use the logger everywhere

This converts the last remaining warnings to use the logging system.
Also drops `catand` and replace it with the intrinsic `atan`, which
supports complex as well as real numbers.

Note: before 3eab989d the `catand` function was actually incorrent!
The definition of arctan(z) can be obtained starting from the identity

  d/dz arctan(z) = 1/(1 + z²) = ½ [1/(1 + iz) + 1/(1 - iz)],

integrating and using the definition log(z) = ∫₁^z dz/z,

  arctan(z) = -i/2 [log(1 + iz) - log(1 - iz)].

If log is the principal branch, log(z) = log|z| + i arg(z), then

  arctan(z) = -i/2 log(w) = 1/2 arg(w) -i/2 log|w|

where w = (1 + iz)/(1 - iz). Finally, the real part is

  Re arctan(z) = 1/2 atan2(2Re(z), 1 - |z|²).

The term -|z|² is missing from the `catand` definition of GRAY,
but is present in the original Fortran 77 code from [SLATEC]:
it has probably been lost in the translation.

[SLATEC]: https://people.math.sc.edu/Burkardt/f_src/slatec/slatec.f90
This commit is contained in:
Michele Guerini Rocco 2022-07-15 00:09:36 +02:00
parent 0a63a20e73
commit cc889bb5a0
Signed by: rnhmjoj
GPG Key ID: BFBAF4C975F76450
4 changed files with 25 additions and 90 deletions

View File

@ -140,7 +140,8 @@ contains
subroutine setcdcoeff_ncl(zeff,rbx,fc,amu,rhop,cst2,eccdpar)
use magsurf_data, only : ch,tjp,tlm,njpt,nlmt
use dierckx, only : profil
use dierckx, only : profil
use logger, only : log_warning
implicit none
integer, parameter :: ksp=3
real(wp_), intent(in) :: zeff,rbx,fc,amu,rhop
@ -154,7 +155,14 @@ contains
call SpitzFuncCoeff(amu,Zeff,fc)
nlm=nlmt
call profil(0,tjp,njpt,tlm,nlmt,ch,ksp,ksp,rhop,nlm,chlm,ierr)
if(ierr>0) print*,' Hlambda profil =',ierr
if (ierr > 0) then
block
character(256) :: msg
write(msg, '(a, g0)') &
'when computing H(λ) `profil` returned ier=', ierr
call log_warning(msg, mod='eccd', proc='setcdcoeff_ncl')
end block
end if
npar=3+2*nlm
allocate(eccdpar(npar))
eccdpar(1)=zeff

View File

@ -652,7 +652,6 @@ contains
! beam tracing initial conditions igrad=1
! !!!!!! check ray tracing initial conditions igrad=0 !!!!!!
use const_and_precisions, only : zero,one,pi,half,two,degree,ui=>im
use math, only : catand
use gray_params, only : gray_parameters
use beamdata, only : nray,nrayr,nrayth,rwmax
@ -730,7 +729,7 @@ contains
dw = wwcsi - wweta
ts = -(dk*sin(2*phirrad) - ui*dw*sin(2*phiwrad))
tc = (dk*cos(2*phirrad) - ui*dw*cos(2*phiwrad))
phic = half*catand(ts/tc)
phic = half*atan(ts/tc)
ddd = dk*cos(2*(phirrad+phic)) - ui*dw*cos(2*(phiwrad+phic))
sss = sk - ui*sw
qi1 = half*(sss + ddd)

View File

@ -469,6 +469,7 @@ contains
subroutine contours_psi(h,rcn,zcn,rup,zup,rlw,zlw)
use const_and_precisions, only : wp_,pi
use gray_params, only : iequil
use logger, only : log_warning
use dierckx, only : profil,sproota
use equilibrium, only : rmaxis,zmaxis,aminor,frhotor,tr,nsr,tz,nsz,cceq, &
kspl,psiant,psinop,points_tgo
@ -517,7 +518,14 @@ contains
zc=zlw+(zup-zlw)*(1.0_wp_-cos(th*(ic-1)))/2.0_wp_
iopt=1
call profil(iopt,tr,nsr,tz,nsz,cceq,kspl,kspl,zc,nsr,czc,ier)
if(ier.gt.0) print*,' profil =',ier
if (ier > 0) then
block
character(256) :: msg
write(msg, '(a, a, g0)') &
'when computing ψ(R,z) contour `profil` returned ier=', ier
call log_warning(msg, mod='magsurf_data', proc='contours_psi')
end block
end if
val=h*psiant+psinop
call sproota(val,tr,nsr,czc,zeroc,mest,m,ier)
if (zeroc(1).gt.rwallm) then

View File

@ -1,100 +1,20 @@
! Miscellaneous mathematical functions
module math
use const_and_precisions, only : wp_, zero, one
implicit none
implicit none
contains
function catand(z)
!***begin prologue catan
!***purpose compute the complex arc tangent.
!***library slatec (fnlib)
!***category c4a
!***type complex (catan-c)
!***keywords arc tangent, elementary functions, fnlib, trigonometric
!***author fullerton, w., (lanl)
!***description
!
! catan(z) calculates the complex trigonometric arc tangent of z.
! the result is in units of radians, and the real part is in the first
! or fourth quadrant.
!
!***references (none)
!***routines called (none)
!***revision history (yymmdd)
! 770801 date written
! 890531 changed all specific intrinsics to generic. (wrb)
! 890531 revision date from version 3.2
! 891214 prologue converted to version 4.0 format. (bab)
! 900315 calls to xerror changed to calls to xermsg. (thj)
! 900326 removed duplicate information from description section.
! (wrb)
!***end prologue catan
use const_and_precisions, only : comp_eps, pi2=>pihalf, czero, cunit
implicit none
complex(wp_) :: catand
complex(wp_), intent(in) :: z
complex(wp_) :: z2
real(wp_) :: r,x,y,r2,xans,yans,twoi
integer :: i
logical, save :: first=.true.
integer, save :: nterms
real(wp_), save :: rmin, rmax, sqeps
!***first executable statement catan
if (first) then
! nterms = log(eps)/log(rbnd) where rbnd = 0.1
nterms = int(-0.4343_wp_*log(0.5_wp_*comp_eps) + 1.0_wp_)
sqeps = sqrt(comp_eps)
rmin = sqrt (1.5_wp_*comp_eps)
rmax = 2.0_wp_/comp_eps
endif
first = .false.
!
r = abs(z)
if (r<=0.1_wp_) then
!
catand = z
if (r<rmin) return
!
catand = czero
z2 = z*z
do i=1,nterms
twoi = 2*(nterms-i) + 1
catand = 1.0_wp_/twoi - z2*catand
end do
catand = z*catand
!
else if (r<=rmax) then
x = real(z)
y = aimag(z)
r2 = r*r
if (r2==one.and.x==zero) print*,'catand, z is +i or -i'
if (abs(r2-one)<=sqeps) then
if (abs(cunit+z*z) < sqeps) &
print*,'catand, answer lt half precision, z**2 close to -1'
!
end if
xans = 0.5_wp_*atan2(2.0_wp_*x, one-r2)
yans = 0.25_wp_*log((r2+2.0_wp_*y+one)/(r2-2.0_wp_*y+one))
catand = cmplx(xans, yans, wp_)
!
else
catand = cmplx(pi2, zero, wp_)
if (real(z)<zero) catand = cmplx(-pi2, zero, wp_)
end if
end function catand
function fact(k)
pure function fact(k)
! Factorial function k!
! Note: the result is real
! Note: the result is a real number
use const_and_precisions, only : wp_
implicit none
integer, intent(in) :: k
real(wp_) :: fact
fact = gamma(real(1 + k, kind=wp_))
end function fact
end module math