Processing math: 100%

derivylm Subroutine

public subroutine derivylm(v1, v2, v3, lmax, rabs, ylm, dydth, dydfi, d2ydth2, d2ydfi2, d2ydthdfi)

Uses

Calculate the 1st and 2nd derivatives of real spherical harmonics with respect to θ, ϕ. Use recursion relations for the assoc. Legendre functions P[l,m] to generate the derivatives. These are (taken from Abramowitz and Stegun, Handbook of Mathematical Functions, chapt. 8.): P[l,m+1]=(x21)12((lm)xP[l,m](l+m)P[l1,m]) (x21)ddxP[l,m]=(l+m)(lm+1)(x21)12P[l,m1]mxP[l,m] (x21)ddxP[l,m]=lxP[l,m](l+m)P[l1,m] where x=cosθ, (x21)12=sinθ, ddθ=sinθddx Adding these equations: ddθP[l,m](cosθ)=12((l+m)(lm+1)P[l,m1]+P[l,m+1]) It is implied that P[l,m]=0 if m>l or m<l. Also, the term (x21)12 is ambiguous for real x, 0<x<1; here it is interpreted as (x21)12=sinθ but (x21)12=1sinθ otherwise the result from Eq.4 (which is cross-checked and correct) does not follow. For the 2nd derivative apply Eq.4 twice. Result: d2dθ2P[l,m](cosθ)=14((l+m)(lm+1)(l+m1)(lm+2)P[l,m2]((lm)(l+m+1)+(l+m)(lm+1))P[l,m]+P[l,m+2]) The ϕ-derivatives act on \cos{\phi},\sin{\phi} and are trivial. For the associated Legendre functions use the recursion formulas: (lm+1)P[l+1,m]=(2l+1)cosθP[l,m](l+m)P[l1,m] P[l+1,m]=P[l1,m](2l+1)sinθP[l,m1] ( with x=cosθ. Recursion algorithm for the calculation of P[l,m] and calculation of Yml taken over from subr. ymy of KKR program (implemented there by M. Weinert, B. Drittler). For m<0, use P[l,m]=P[l,m](lm)!(l+m)! Taking into account the lm-prefactors of the spherical harmonics, we construct and use the functions Q[l,m]=2l+14π(lm)!(l+m)!P[l,m] whence Eq.4 and Eq.7 become ddθQ[l,m]=12((l+m)(lm+1)Q[l,m1]+(l+m+1)(lm)Q[l,m+1]) d2dθ2Q[l,m]=14((l+m)(l+m1)(lm+1)(lm+2)Q[l,m2]+((lm)(l+m+1)+(l+m)(lm+1))Q[l,m]+(lm)(lm1)(l+m+1)(l+m+2)Q[l,m+2]) Note on sign convension: For the needs of GGA PW91 as implemented here, ylm and derivatives come with a different sign convention compared to the usual in the program: sinϕm(1)msinϕm. Thus some signs change.

Arguments

Type IntentOptional Attributes Name
real(kind=dp) :: v1
real(kind=dp) :: v2
real(kind=dp) :: v3
integer :: lmax
real(kind=dp) :: rabs
real(kind=dp) :: ylm(*)
real(kind=dp) :: dydth(*)
real(kind=dp) :: dydfi(*)
real(kind=dp) :: d2ydth2(*)
real(kind=dp) :: d2ydfi2(*)
real(kind=dp) :: d2ydthdfi(*)