Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

doxygen updates in sorc/grid_tools.fd - applied comments from R. Purser #409

Merged
merged 3 commits into from
Mar 11, 2021
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
123 changes: 63 additions & 60 deletions sorc/grid_tools.fd/regional_esg_grid.fd/pfun.f90
Original file line number Diff line number Diff line change
@@ -1,8 +1,11 @@
!> @file
!! ???
!! @brief Routine evaluating useful functions not always available in Fortran.
!! @author R. J. Purser

!> This module is for ???
!> This module is for evaluating several useful real-valued functions
!! that are not always available in Fortran compilers. These have applications
!! in map transformations, spherical and pseudo-spherical surface geometry,
!! probability distributions, and splines.
!!
!! @author R. J. Purser
module pfun
Expand All @@ -29,10 +32,10 @@ module pfun

contains

!> ???
!> Gudermannian function (related to Mercator map projections)
!!
!! @param[in] x ???
!! @return y ???
!! @param[in] x single precision real argument of function
!! @return y Gudermannian function of x
!! @author R. J. Purser
function gd_s(x) result(y)! [gd]
! Gudermannian function
Expand All @@ -42,10 +45,10 @@ function gd_s(x) result(y)! [gd]
y=atan(sinh(x))
end function gd_s

!> ???
!> Gudermannian function (related to Mercator map projections)
!!
!! @param[in] x ???
!! @return y ???
!! @param[in] x double precision real argument of function
!! @return y Gudermannian function of x
!! @author R. J. Purser
function gd_d(x) result(y)! [gd]
implicit none
Expand All @@ -56,8 +59,8 @@ end function gd_d

!> Inverse Gudermannian function for single precision real.
!!
!! @param[in] y ???
!! @return x ???
!! @param[in] y single precision real argument
!! @return x inverse Gudermannian function of y
!! @author R. J. Purser
function gdi_s(y) result(x)! [gdi]
implicit none
Expand All @@ -68,8 +71,8 @@ end function gdi_s

!> Inverse Gudermannian function for double precision real.
!!
!! @param[in] y ???
!! @return x ???
!! @param[in] y double precision real argument
!! @return x inverse Gudermannian function of y
!! @author R. J. Purser
function gdi_d(y) result(x)! [gdi]
implicit none
Expand All @@ -78,10 +81,10 @@ function gdi_d(y) result(x)! [gdi]
x=atanh(sin(y))
end function gdi_d

!> Haversine function for single precision real.
!> Haversine function for single precision real (for geometry on the sphere).
!!
!! @param[in] t ???
!! @return a ???
!! @param[in] t single precision real argument
!! @return a haversine function of t
!! @author R. J. Purser
function hav_s(t) result(a)! [hav]
use pietc_s, only: o2
Expand All @@ -91,10 +94,10 @@ function hav_s(t) result(a)! [hav]
a=(sin(t*o2))**2
end function hav_s

!> Haversine function for double precision real.
!> Haversine function for double precision real (for geometry on the sphere).
!!
!! @param[in] t ???
!! @return a ???
!! @param[in] t double precision real argument
!! @return a haversine function of t
!! @author R. J. Purser
function hav_d(t) result(a)! [hav]
use pietc, only: o2
Expand All @@ -104,12 +107,12 @@ function hav_d(t) result(a)! [hav]
a=(sin(t*o2))**2
end function hav_d

!> Hyperbolic-haversine for single precision real.
!> Hyperbolic-haversine for single precision real (for pseudosphere geometry).
!!
!! @note The minus sign in the hyperbolic-haversine definition.
!!
!! @param[in] t ???
!! @return a ???
!! @param[in] t single precision real argument
!! @return a hyperbolic-haversine function of t
!! @author R. J. Purser
function havh_s(t) result(a)! [havh]
use pietc_s, only: o2
Expand All @@ -119,10 +122,10 @@ function havh_s(t) result(a)! [havh]
a=-(sinh(t*o2))**2
end function havh_s

!> Hyperbolic-haversine for double precision real.
!> Hyperbolic-haversine for double precision real (for pseudosphere geometry).
!!
!! @param[in] t ???
!! @return a ???
!! @param[in] t double precision real argument
!! @return a hyperbolic-haversine function of t
!! @author R. J. Purser
function havh_d(t) result(a)! [havh]
use pietc, only: o2
Expand All @@ -134,8 +137,8 @@ end function havh_d

!> Arc-haversine function for single precision real.
!!
!! @param[in] a ???
!! @return t ???
!! @param[in] a single precision real argument
!! @return t arc-haversine function of a
!! @author R. J. Purser
function ahav_s(a) result(t)! [ahav]
use pietc_s, only: u2
Expand All @@ -147,8 +150,8 @@ end function ahav_s

!> Arc-haversine function for double precision real.
!!
!! @param[in] a ???
!! @return t ???
!! @param[in] a double precision real argument
!! @return t arc-haversine function of a
!! @author R. J. Purser
function ahav_d(a) result(t)! [ahav]
use pietc, only: u2
Expand All @@ -162,8 +165,8 @@ end function ahav_d
!!
!! @note The minus sign in the hyperbolic arc-haversine definition.
!!
!! @param[in] a ???
!! @return t ???
!! @param[in] a single precision real argument
!! @return t hyperbolic arc-haversine of a
!! @author R. J. Purser
function ahavh_s(a) result(t)! [ahavh]
use pietc_s, only: u2
Expand All @@ -175,8 +178,8 @@ end function ahavh_s

!> Hyperbolic arc-haversine for double precision real.
!!
!! @param[in] a ???
!! @return t ???
!! @param[in] a double precision real argument
!! @return t hyperbolic arc-haversine of a
!! @author R. J. Purser
function ahavh_d(a) result(t)! [ahavh]
use pietc, only: u2
Expand All @@ -186,10 +189,10 @@ function ahavh_d(a) result(t)! [ahavh]
t=u2*asinh(sqrt(-a))
end function ahavh_d

!> ???
!> Hyperbolic arc-tangent for single precision real. (compilers now have this)
!!
!! @param[in] t ???
!! @return a ???
!! @param[in] t single precision real argument
!! @return a hyperbolic arc-tangent of t
!! @author R. J. Purser
function atanh_s(t) result(a)! [atanh]
use pietc_s, only: u1,o2,o3,o5
Expand All @@ -204,10 +207,10 @@ function atanh_s(t) result(a)! [atanh]
endif
end function atanh_s

!> ???
!> Hyperbolic arc-tangent for double precision real.
!!
!! @param[in] t ???
!! @return a ???
!! @param[in] t double precision real argument
!! @return a hyperbolic arc-tangent of t
!! @author R. J. Purser
function atanh_d(t) result(a)! [atanh]
use pietc, only: u1,o2,o3,o5
Expand All @@ -222,10 +225,10 @@ function atanh_d(t) result(a)! [atanh]
endif
end function atanh_d

!> ???
!> Hyperbolic secant for single precision real.
!!
!! @param[in] x ???
!! @return t ???
!! @param[in] x single precision real argument
!! @return r hyperbolic secant of x
!! @author R. J. Purser
function sech_s(x)result(r)! [sech]
! This indirect way of computing 1/cosh(x) avoids overflows at large x
Expand All @@ -239,10 +242,10 @@ function sech_s(x)result(r)! [sech]
r=e*u2/(u1+e*e)
end function sech_s

!> ???
!> Hyperbolic secant for double precision real.
!!
!! @param[in] x ???
!! @return r ???
!! @param[in] x double precision real argument
!! @return r hyperbolic secant of x
!! @author R. J. Purser
function sech_d(x)result(r)! [sech]
use pietc, only: u1,u2
Expand All @@ -255,10 +258,10 @@ function sech_d(x)result(r)! [sech]
r=e*u2/(u1+e*e)
end function sech_d

!> ???
!> Hyperbolic secant-squared function (logistic distribution).
!!
!! @param[in] x ???
!! @return r ???
!! @param[in] x single precision real argument
!! @return r sech-squared of x
!! @author R. J. Purser
function sechs_s(x)result(r)! [sechs]
implicit none
Expand All @@ -267,10 +270,10 @@ function sechs_s(x)result(r)! [sechs]
r=sech(x)**2
end function sechs_s

!> ???
!> Hyperbolic secant-squared function (logistic distribution).
!!
!! @param[in] x ???
!! @return r ???
!! @param[in] x double precision real argument
!! @return r sech-squared of x
!! @author R. J. Purser
function sechs_d(x)result(r)! [sechs]
implicit none
Expand All @@ -279,10 +282,10 @@ function sechs_d(x)result(r)! [sechs]
r=sech(x)**2
end function sechs_d

!> Evaluate the symmetric real function sin(x)/x-1.
!> Evaluate the symmetric real function sin(x)/x-1, still accurate near x=0.
!!
!! @param[in] x ???
!! @return r ???
!! @param[in] x double precision real argument
!! @return r sin(x)/x-1
!! @author R. J. Purser
function sinoxm_d(x) result(r)! [sinoxm]
use pietc, only: u1
Expand All @@ -301,8 +304,8 @@ end function sinoxm_d

!> Evaluate the symmetric real function sin(x)/x.
!!
!! @param[in] x ???
!! @return r ???
!! @param[in] x double precision real argument
!! @return r sin(x)/x
!! @author R. J. Purser
function sinox_d(x) result(r)! [sinox]
use pietc, only: u1
Expand All @@ -313,10 +316,10 @@ function sinox_d(x) result(r)! [sinox]
r=sinoxm(x)+u1
end function sinox_d

!> Evaluate the symmetric real function sinh(x)/x-1.
!> Evaluate the symmetric real function sinh(x)/x-1. still accurate near x=0.
!!
!! @param[in] x ???
!! @return r ???
!! @param[in] x double precision real argument
!! @return r sinh(x)-1
!! @author R. J. Purser
function sinhoxm_d(x) result(r)! [sinhoxm]
use pietc, only: u1
Expand All @@ -335,8 +338,8 @@ end function sinhoxm_d

!> Evaluate the symmetric real function sinh(x)/x.
!!
!! @param[in] x ???
!! @return r ???
!! @param[in] x double precision real argument
!! @return r sinh(x)/x
!! @author R. J. Purser
function sinhox_d(x) result(r)! [sinhox]
use pietc, only: u1
Expand Down
5 changes: 5 additions & 0 deletions sorc/grid_tools.fd/regional_esg_grid.fd/pkind.f90
Original file line number Diff line number Diff line change
@@ -1,4 +1,9 @@
!> @file
!! @brief Standard integer, real, and complex single and double precision kinds.
!! @author R. J. Purser

!> Standard integer, real, and complex single and double precision kinds.
!! @author R. J. Purser
module pkind
integer,parameter:: spi=selected_int_kind(6) !< Single precision integer kind.
integer,parameter:: dpi=selected_int_kind(12) !< Double precision integer kind.
Expand Down