Class sncndnK_kennedy_mod
In: ZolotarevModule/sncndnK_kennedy_mod.F90
sncndnK_kennedy_mod dot/f_93.png

Compute Jacobi elliptic function sn,cn,dn and complete elliptic integral K

Refs.

1
A.D.Kennedy, help-lat/0402038.
2
B.Laszkiewicza and K.Zietak, Appl.Math.Comput. 206, 298 (2008)
3
E.L.Wachspress, Comp.Math.Appl. 39 (2000) 131.
4
NEC, ASL/SX User Guide.

Methods

Public Instance methods

PI
Constant :
PI = 4*ATAN(1.0_QP) :real(QP), parameter
QP
Constant :
QP =SELECTED_REAL_KIND(33) :integer, parameter
Function :recursive
sn :real(QP)
z :real(QP), intent(in)
a :real(QP), intent(in)
b :real(QP), intent(in)
agm :real(QP), intent(inout)

[Source]

recursive function arithgeom(z,a,b,agm) result(sn)
  implicit none
  real(QP), intent(in)    :: z,a,b
  real(QP), intent(inout) :: agm
  real(QP), save :: pb = -1.0_QP
  real(QP) :: sn
  real(QP) :: xi
  if (b <= pb) then
    pb = -1.0_QP
    agm = a
    sn = sin(z*a)
    return
  endif
  pb = b
  xi = arithgeom(z,(a+b)/2.0_QP,sqrt(a*b),agm)
  sn = 2*a*xi/((a+b)+(a-b)*xi**2)
  return
end function
Function :
ceK :real(QP)
kpm :real(QP), intent(in)

Return complete elliptic integral 1st kind K with modulus k

K(k)

 kpm : complementary modulus kappa'  = sqrt(1-kappa^2)
 ceK : K(kappa) = K(sqrt(1-kappa'^2))

[Source]

function comp_elliptic_integ_K(kpm) result(ceK)
!
! Return complete elliptic integral 1st kind K with modulus k
!
! K(k)
!
!  kpm : complementary modulus kappa'  = sqrt(1-kappa^2)
!  ceK : K(kappa) = K(sqrt(1-kappa'^2))
!
  implicit none
  real(QP), intent(in) :: kpm
  real(QP) :: ceK
  real(QP) :: a,b,a0,b0
  a = 1.0_QP
  b = kpm
  do
    a0 = a
    b0 = b
    a = (a0+b0)/2.0_QP
    b = sqrt(a0*b0)
    if (b <= b0) exit
  enddo
  ceK = 2*atan(1.0_QP)/a
  return
end function
Function :
et :real(QP)
v :real(QP), intent(in)
q :real(QP), intent(in)

Elliptic Theta Function theta(4,v,q) (Double prec.)

theta(4,v,q) = 1 + 2 sum_{n=0,inf}(-1)^n q^(n^2) cos(2 n v)

Converted form Algorithm used section on Ref.[4].

 (Note. typo on the formula exsits in Ref.[4])

[Source]

function detf4(v,q) result(et)
!
! Elliptic Theta Function theta(4,v,q) (Double prec.)
!
! theta(4,v,q) = 1 + 2 sum_{n=0,inf}(-1)^n q^(n^2) cos(2 n v)
!
! Converted form Algorithm used section on Ref.[4].
! 
!  (Note. typo on the formula exsits in Ref.[4])
!  
!
  implicit none
  real(QP), intent(in) :: v,q
  real(QP) :: et
  real(QP), parameter :: th=0.125_QP
  real(QP) :: u,cw,c2w,c3w,c4w
  real(QP) :: PLQ,wq,w

  et = abs(v/PI)

  u = et-INT(et)

  if (q == 0.0_QP) then
    et = 1.0_QP
    return
  endif
  if ( q <= th ) then
    cw = COS(2*PI*u)
    c2w =  2*cw**2-1.0_QP
    c3w = (4*cw**2-3.0_QP)*cw
    c4w = (8*cw**2-8.0_QP)*cw**2+1.0_QP
    et = 2*(((c4w*q**7-c3w)*q**5+c2w)*q**3-cw)*q+1.0_QP
  else
    if (u > 0.5_QP) u=1.0_QP-u
    PLQ = PI**2*(1.0_QP/log(q))
    wq = exp(PLQ)
    w = exp(2*u*PLQ)
    et = sqrt(PI*(-1.0_QP/log(q)))*exp((u**2-u+0.25_QP)*PLQ)* ((((wq**12*w+wq**6)*w+wq**2)*w+1.0_QP)*w+1.0_QP+((wq**8/w+wq**2)*(wq**2/w)+1.0_QP)*(wq**2/w))
  endif

  return
end function
Function :
f :real(QP)
v :real(QP), intent(in)
q :real(QP), intent(in)

theta function theta(4,v,q) = alp theta(2, i v’, q’)

via Jacobi identity. (q’ is small for q near 1)

approximation is good for q near 1.

[Source]

function etf2i(v,q) result(f)
!
! theta function  theta(4,v,q) = alp theta(2, i v', q') 
!
! via Jacobi identity.  (q' is small for q near 1)
!
! approximation is good for q near 1.
!
  implicit none
  real(QP), intent(in) :: q,v
  real(QP) :: f,fa,fb,u
  real(QP) :: c,wp,wm,w,a,b,wq,PLQ
  integer :: i,ns

  u = v/PI
  u = u-INT(u)

  if (u > 0.5_QP) u=1.0_QP-u

  PLQ = PI**2*(1.0_QP/log(q))
  wq = exp(PLQ)
  w = exp(2*u*PLQ)

  ns = INT(sqrt( (log(EPSILON(1.0_QP)*2)-PI*u)/log(q)+(PI*u/log(q)+0.5_QP)**2) -(PI*u/log(q)+0.5_QP))+1
  
  i=ns
  fa = wq**(i*(i+1))*w
  do i=ns,1,-1
    fa=(fa + wq**(i*(i-1)))*w
  enddo

  i=ns
  fb = wq**(i*(i-1))*(wq**2/w)
  do i=ns,2,-1
    fb=(fb + wq**((i-1)*(i-2)))*(wq**2/w)
  enddo
  fb=fb + 1.0_QP
     
  f = sqrt(PI*(-1.0_QP/log(q)))*exp((u**2-u+0.25_QP)*PLQ)*(fa+fb)

  return
end function
Function :
f :real(QP)
v :real(QP), intent(in)
q :real(QP), intent(in)

Elliptic Theta Function theta(4,v,q)

theta(4,v,q) = 1 + 2 sum_{n=0,inf}(-1)^n q^(n^2) cos(2 n v)

[Source]

function etf4(v,q) result(f)
!
! Elliptic Theta Function theta(4,v,q)
!
! theta(4,v,q) = 1 + 2 sum_{n=0,inf}(-1)^n q^(n^2) cos(2 n v)
!
  implicit none
  real(QP), intent(in) :: q,v
  real(QP), parameter :: QUAD=0.125_QP
  real(QP) :: f
  real(QP) :: vp,qp,Lq
  if ( q < QUAD ) then
    f = etf4r(v,q)
  else
    f = etf2i(v,q)
  endif
  return
end function
Function :
f :real(QP)
v :real(QP), intent(in)
q :real(QP), intent(in)

theta function theta(4,v,q)

via the definitnion series

approximation is good for small q.

[Source]

function etf4r(v,q) result(f)
!
! theta function theta(4,v,q)
!
! via the definitnion series
!
! approximation is good for small q.
!
  implicit none
  real(QP), intent(in) :: q,v
  real(QP) :: f
  real(QP) :: u,c,wp,wm,w
  integer :: i,ns

  if (q==0.0_QP) then
    f = 1.0_QP
    return
  endif

  u = abs(v/PI)
  u = u - INT(u)

  ns = sqrt(log(EPSILON(1.0_QP))/log(abs(q)))+2
#ifdef _GGG_
  wp = 0.0_QP
  wm = 0.0_QP
  do i=ns,1,-1
    w = (-1)**i*q**(i**2)*cos(2*PI*i*u)
    if (w > 0.0_QP) then
      wp = wp + w
    else
      wm = wm + w
    endif
  enddo
  f=(wp+wm)*2+1.0_QP
#else
  f=1.0_QP
  do i=ns,2,-1
    f=1.0_QP-q**(2*i-1)*cos(2*PI*i*u)/cos(2*PI*(i-1)*u)*f
  enddo
  f=1.0_QP-q*cos(2*PI*u)*2*f
#endif
  return
end function
Function :
e :real(QP)
x :real(QP), intent(in)
kpm :real(QP), intent(in)

Return imcomplete elliptic integral 1st kind with modulus k

  x : argument

kpm : complementary modulus kappa’ = sqrt(1-kappa^2)

  e : F(x;kappa) imcomplete elliptic intagral of 1st kind

[Source]

function imcomp_elliptic_integ1(x,kpm) result(e)
!
! Return imcomplete elliptic integral 1st kind with modulus k
!
!   x : argument
! kpm : complementary modulus kappa' = sqrt(1-kappa^2)
!   e : F(x;kappa) imcomplete elliptic intagral of 1st kind
!
  real(QP), intent(in) :: x,kpm
  real(QP) :: e
  real(QP) :: a,b,d,t,a0,b0,d0,t0
  integer :: s,s0
  integer :: j

  a = 1.0_QP
  b = kpm
  d = b/a
  t = x/sqrt((1.0_QP-x)*(1.0_QP+x))
  s = 0
  j = 0
  do
    a0 = a
    b0 = b
    d0 = d
    t0 = t
    s0 = s
    j = j + 1
    a = (a0 + b0)/2.0_QP
    b = sqrt(a0*b0)
    d = b/a
    t = ((1.0_QP + d0)*t0)/(1.0_QP-d0*t0**2)
    if (t*t0 > 0.0_QP) then
      s = 2*s0
    else if ( t < 0.0_QP .and. t0 > 0.0_QP) then
      s = 2*s0+1
    else if ( t > 0.0_QP .and. t0 < 0.0_QP) then
      s = 2*s0-1
    endif
    if (b <= b0) exit
  enddo
  e = (atan(t)+s*4*atan(1.0_QP))/2**j/a

  return
end function
Subroutine :
z :real(QP), intent(in)
kap :real(QP), intent(in)
sn :real(QP), intent(out)
cn :real(QP), intent(out)
dn :real(QP), intent(out)
K :real(QP), intent(out)

return Jacobi elliptic function

sn(z,kappa) cn(z,kappa) dn(z,kappa)

and complete elliptic integral

   K(kappa)

  z : real argument for sn,cn,dn

kap : modulus

[Source]

subroutine sncndnK(z,kap,sn,cn,dn,K)
!
! return Jacobi elliptic function
!
! sn(z,kappa)
! cn(z,kappa)
! dn(z,kappa)
!
! and complete elliptic integral
!
!    K(kappa)
!
!   z : real argument for sn,cn,dn
! kap : modulus
!
!
  implicit none
  real(QP), intent(in)  :: z,kap
  real(QP), intent(out) :: sn,cn,dn,K
  real(QP) :: agm,kappm
  integer :: sgn
  kappm = sqrt((1.0_QP-kap)*(1.0_QP+kap))
  sn = arithgeom(z,1.0_QP,kappm,agm)
  K = 2*atan(1.0_QP)/agm
  sgn = mod(INT(abs(z)/K),4)
  sgn = ieor(sgn,ishft(sgn,-1))
  sgn = 1 - ishft(iand(sgn,1),1)
  cn = sgn*sqrt((1.0_QP-    sn)*(1.0_QP+    sn))
  dn =     sqrt((1.0_QP-kap*sn)*(1.0_QP+kap*sn))
  return
end subroutine
Subroutine :
z :real(QP), intent(in)
kappm :real(QP), intent(in)
sn :real(QP), intent(out)
cn :real(QP), intent(out)
dn :real(QP), intent(out)
K :real(QP), intent(out)

return Jacobi elliptic function

sn(z,kappa) cn(z,kappa) dn(z,kappa)

and complete elliptic integral

   K(kappa)

    z : real argument for sn,cn,dn

kappm : complementary modulus kappa’ = sqrt(1-kappa**2)

[Source]

subroutine sncndnK2(z,kappm,sn,cn,dn,K)
!
! return Jacobi elliptic function
!
! sn(z,kappa)
! cn(z,kappa)
! dn(z,kappa)
!
! and complete elliptic integral
!
!    K(kappa)
!
!     z : real argument for sn,cn,dn
! kappm : complementary modulus kappa' = sqrt(1-kappa**2)
!
  implicit none
  real(QP), intent(in)  :: z,kappm
  real(QP), intent(out) :: sn,cn,dn,K
  real(QP) :: agm
  integer :: sgn
  sn = arithgeom(z,1.0_QP,kappm,agm)
  K = 2*atan(1.0_QP)/agm
  sgn = mod(INT(abs(z)/K),4)
  sgn = ieor(sgn,ishft(sgn,-1))
  sgn = 1 - ishft(iand(sgn,1),1)
  cn = sgn*sqrt((1.0_QP-sn)*(1.0_QP+sn))
  dn =     sqrt((1.0_QP-sn)*(1.0_QP+sn)+(kappm*sn)**2)
  return
end subroutine