grid_tools 1.14.0
Loading...
Searching...
No Matches
pfun.f90
Go to the documentation of this file.
1
4
11module pfun
12!=============================================================================
13use pkind, only: sp,dp
14implicit none
15private
18
19interface gd; module procedure gd_s, gd_d; end interface
20interface gdi; module procedure gdi_s, gdi_d; end interface
21interface hav; module procedure hav_s, hav_d; end interface
22interface havh; module procedure havh_s, havh_d; end interface
23interface ahav; module procedure ahav_s, ahav_d; end interface
24interface ahavh; module procedure ahavh_s, ahavh_d; end interface
25interface atanh; module procedure atanh_s, atanh_d; end interface
26interface sech; module procedure sech_s, sech_d; end interface
27interface sechs; module procedure sechs_s, sechs_d; end interface
28interface sinoxm; module procedure sinoxm_d; end interface
29interface sinox; module procedure sinox_d; end interface
30interface sinhoxm; module procedure sinhoxm_d;end interface
31interface sinhox; module procedure sinhox_d; end interface
32
33contains
34
40function gd_s(x) result(y)! [gd]
41! Gudermannian function
42implicit none
43real(sp),intent(in ):: x
44real(sp) :: y
45y=atan(sinh(x))
46end function gd_s
47
53function gd_d(x) result(y)! [gd]
54implicit none
55real(dp),intent(in ):: x
56real(dp) :: y
57y=atan(sinh(x))
58end function gd_d
59
65function gdi_s(y) result(x)! [gdi]
66implicit none
67real(sp),intent(in ):: y
68real(sp) :: x
69x=atanh(sin(y))
70end function gdi_s
71
77function gdi_d(y) result(x)! [gdi]
78implicit none
79real(dp),intent(in ):: y
80real(dp) :: x
81x=atanh(sin(y))
82end function gdi_d
83
89function hav_s(t) result(a)! [hav]
90use pietc_s, only: o2
91implicit none
92real(sp),intent(in ):: t
93real(sp) :: a
94a=(sin(t*o2))**2
95end function hav_s
96
102function hav_d(t) result(a)! [hav]
103use pietc, only: o2
104implicit none
105real(dp),intent(in ):: t
106real(dp) :: a
107a=(sin(t*o2))**2
108end function hav_d
109
117function havh_s(t) result(a)! [havh]
118use pietc_s, only: o2
119implicit none
120real(sp),intent(in ):: t
121real(sp) :: a
122a=-(sinh(t*o2))**2
123end function havh_s
124
130function havh_d(t) result(a)! [havh]
131use pietc, only: o2
132implicit none
133real(dp),intent(in ):: t
134real(dp) :: a
135a=-(sinh(t*o2))**2
136end function havh_d
137
143function ahav_s(a) result(t)! [ahav]
144use pietc_s, only: u2
145implicit none
146real(sp),intent(in ):: a
147real(sp) :: t
148t=u2*asin(sqrt(a))
149end function ahav_s
150
156function ahav_d(a) result(t)! [ahav]
157use pietc, only: u2
158implicit none
159real(dp),intent(in ):: a
160real(dp) :: t
161t=u2*asin(sqrt(a))
162end function ahav_d
163
171function ahavh_s(a) result(t)! [ahavh]
172use pietc_s, only: u2
173implicit none
174real(sp),intent(in ):: a
175real(sp) :: t
176t=u2*asinh(sqrt(-a))
177end function ahavh_s
178
184function ahavh_d(a) result(t)! [ahavh]
185use pietc, only: u2
186implicit none
187real(dp),intent(in ):: a
188real(dp) :: t
189t=u2*asinh(sqrt(-a))
190end function ahavh_d
191
197function atanh_s(t) result(a)! [atanh]
198use pietc_s, only: u1,o2,o3,o5
199implicit none
200real(sp),intent(IN ):: t
201real(sp) :: a,tt
202real(sp),parameter :: o7=u1/7_sp,o9=u1/9_sp
203!=============================================================================
204if(abs(t)>=u1)stop 'In atanh; no solution'
205if(abs(t)>1.e-3_sp)then; a=log((u1+t)/(u1-t))*o2
206else; tt=t*t; a=t*(u1+tt*(o3+tt*(o5+tt*(o7+tt*o9))))
207endif
208end function atanh_s
209
215function atanh_d(t) result(a)! [atanh]
216use pietc, only: u1,o2,o3,o5
217implicit none
218real(dp),intent(IN ):: t
219real(dp) :: a,tt
220real(dp),parameter :: o7=u1/7_dp,o9=u1/9_dp
221!=============================================================================
222if(abs(t)>=u1)stop 'In atanh; no solution'
223if(abs(t)>1.e-3_dp)then; a=log((u1+t)/(u1-t))*o2
224else; tt=t*t; a=t*(u1+tt*(o3+tt*(o5+tt*(o7+tt*o9))))
225endif
226end function atanh_d
227
233function sech_s(x)result(r)! [sech]
234! This indirect way of computing 1/cosh(x) avoids overflows at large x
235use pietc_s, only: u1,u2
236implicit none
237real(sp),intent(in ):: x
238real(sp) :: r
239real(sp) :: e,ax
240ax=abs(x)
241e=exp(-ax)
242r=e*u2/(u1+e*e)
243end function sech_s
244
250function sech_d(x)result(r)! [sech]
251use pietc, only: u1,u2
252implicit none
253real(dp),intent(in ):: x
254real(dp) :: r
255real(dp) :: e,ax
256ax=abs(x)
257e=exp(-ax)
258r=e*u2/(u1+e*e)
259end function sech_d
260
266function sechs_s(x)result(r)! [sechs]
267implicit none
268real(sp),intent(in ):: x
269real(sp) :: r
270r=sech(x)**2
271end function sechs_s
272
278function sechs_d(x)result(r)! [sechs]
279implicit none
280real(dp),intent(in ):: x
281real(dp) :: r
282r=sech(x)**2
283end function sechs_d
284
290function sinoxm_d(x) result(r)! [sinoxm]
291use pietc, only: u1
292implicit none
293real(dp),intent(in ):: x
294real(dp) :: r
295!-----------------------------------------------------------------------------
296real(dp):: xx
297!=============================================================================
298xx=x*x
299if(xx > .05_dp)then; r=sin(x)/x-u1
300else ; r=-xx*(u1-xx*(u1-xx*(u1-xx*(u1-xx*(u1-xx/&
301 156._dp)/110._dp)/72._dp)/42._dp)/20._dp)/6._dp
302endif
303end function sinoxm_d
304
310function sinox_d(x) result(r)! [sinox]
311use pietc, only: u1
312implicit none
313real(dp),intent(in ):: x
314real(dp) :: r
315!=============================================================================
316r=sinoxm(x)+u1
317end function sinox_d
318
324function sinhoxm_d(x) result(r)! [sinhoxm]
325use pietc, only: u1
326implicit none
327real(dp),intent(in ):: x
328real(dp) :: r
329!-----------------------------------------------------------------------------
330real(dp):: xx
331!=============================================================================
332xx=x*x
333if(xx > .05_dp)then; r=sinh(x)/x-u1
334else; r=xx*(u1+xx*(u1+xx*(u1+xx*(u1+xx*(u1+xx/&
335 156._dp)/110._dp)/72._dp)/42._dp)/20._dp)/6._dp
336endif
337end function sinhoxm_d
338
344function sinhox_d(x) result(r)! [sinhox]
345use pietc, only: u1
346implicit none
347real(dp),intent(in ):: x
348real(dp) :: r
349!=============================================================================
350r=sinhoxm(x)+u1
351end function sinhox_d
352
353end module pfun
This module is for evaluating several useful real-valued functions that are not always available in F...
Definition pfun.f90:11
real(dp) function sinhox_d(x)
Evaluate the symmetric real function sinh(x)/x.
Definition pfun.f90:345
real(dp) function atanh_d(t)
Hyperbolic arc-tangent for double precision real.
Definition pfun.f90:216
real(dp) function sech_d(x)
Hyperbolic secant for double precision real.
Definition pfun.f90:251
real(dp) function sinhoxm_d(x)
Evaluate the symmetric real function sinh(x)/x-1.
Definition pfun.f90:325
real(sp) function hav_s(t)
Haversine function for single precision real (for geometry on the sphere).
Definition pfun.f90:90
real(dp) function gd_d(x)
Gudermannian function (related to Mercator map projections)
Definition pfun.f90:54
real(sp) function gd_s(x)
Gudermannian function (related to Mercator map projections)
Definition pfun.f90:41
real(sp) function ahav_s(a)
Arc-haversine function for single precision real.
Definition pfun.f90:144
real(dp) function gdi_d(y)
Inverse Gudermannian function for double precision real.
Definition pfun.f90:78
real(dp) function havh_d(t)
Hyperbolic-haversine for double precision real (for pseudosphere geometry).
Definition pfun.f90:131
real(sp) function gdi_s(y)
Inverse Gudermannian function for single precision real.
Definition pfun.f90:66
real(dp) function ahav_d(a)
Arc-haversine function for double precision real.
Definition pfun.f90:157
real(dp) function sinox_d(x)
Evaluate the symmetric real function sin(x)/x.
Definition pfun.f90:311
real(sp) function atanh_s(t)
Hyperbolic arc-tangent for single precision real.
Definition pfun.f90:198
real(dp) function sechs_d(x)
Hyperbolic secant-squared function (logistic distribution).
Definition pfun.f90:279
real(dp) function hav_d(t)
Haversine function for double precision real (for geometry on the sphere).
Definition pfun.f90:103
real(sp) function havh_s(t)
Hyperbolic-haversine for single precision real (for pseudosphere geometry).
Definition pfun.f90:118
real(dp) function ahavh_d(a)
Hyperbolic arc-haversine for double precision real.
Definition pfun.f90:185
real(dp) function sinoxm_d(x)
Evaluate the symmetric real function sin(x)/x-1, still accurate near x=0.
Definition pfun.f90:291
real(sp) function ahavh_s(a)
Hyperbolic arc-haversine for single precision real.
Definition pfun.f90:172
real(sp) function sechs_s(x)
Hyperbolic secant-squared function (logistic distribution).
Definition pfun.f90:267
real(sp) function sech_s(x)
Hyperbolic secant for single precision real.
Definition pfun.f90:234
Some of the commonly used constants (pi etc) mainly for double-precision subroutines.
Definition pietc.f90:14
real(dp), parameter o5
fifth
Definition pietc.f90:35
real(dp), parameter u1
one
Definition pietc.f90:20
real(dp), parameter o2
half
Definition pietc.f90:32
real(dp), parameter o3
third
Definition pietc.f90:33
logical, parameter t
for pain-relief in logical ops
Definition pietc.f90:17
real(dp), parameter u2
two
Definition pietc.f90:22
Standard integer, real, and complex single and double precision kinds.
Definition pkind.f90:7
integer, parameter dp
Double precision real kind.
Definition pkind.f90:11
integer, parameter sp
Single precision real kind.
Definition pkind.f90:10