-
Notifications
You must be signed in to change notification settings - Fork 106
/
PLegendreA_d1.f95
185 lines (154 loc) · 5.36 KB
/
PLegendreA_d1.f95
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
subroutine PLegendreA_d1(p, dp, lmax, z, csphase)
!-------------------------------------------------------------------------------
!
! This function evalutates all of the unnormalized associated legendre
! polynomials and their first derivatives up to degree lmax.
!
! Calling Parameters
!
! IN
! lmax Maximum spherical harmonic degree to compute.
! z [-1, 1], cos(colatitude) or sin(latitude).
!
! OPTIONAL (IN)
! csphase 1: Do not include the phase factor of (-1)^m
! -1: Apply the phase factor of (-1)^m.
!
! OUT
! p A vector of all associated Legendgre polynomials
! evaluated at z up to lmax. The length must by greater
! or equal to (lmax+1)*(lmax+2)/2.
!
! Notes:
!
! 1. The integral of plm**2 over (-1,1) is 2 * (l+m)! / (2l+1) / (l-m)!.
! 2. The index of the array p corresponds to l*(l+1)/2 + m + 1.
! 3. The index of the array p corresponds to l*(l+1)/2 + m + 1.
! 4. The derivative is evaluated with respecte to z, and NOT cos(colatitude)
! or sin(latitude).
! 5. Derivatives are calculated using the unnormalized identities
! P'l,m = ( (l+m) Pl-1,m - l z Plm ) / (1-z**2) (for l>m), and
! P'll = - l z Pll / (1-z**2) (for l=m).
! 6. The derivative is not defined at z=+-1 for all m>0, and is therefore not
! calculated here.
! 7. The default is to exlude the Condon-Shortley phase of (-1)^m.
!
! Dependencies: CSPHASE_DEFAULT
!
! Copyright (c) 2005-2012, Mark A. Wieczorek
! All rights reserved.
!
!-------------------------------------------------------------------------------
use SHTOOLS, only: CSPHASE_DEFAULT
implicit none
integer, intent(in) :: lmax
real*8, intent(out) :: p(:), dp(:)
real*8, intent(in) :: z
integer, intent(in), optional :: csphase
real*8 :: pm2, pm1, pmm, sinsq, sinsqr, fact, plm
integer :: k, kstart, m, l, sdim
integer*1 :: phase
sdim = (lmax+1)*(lmax+2)/2
if (size(p) < sdim) then
print*, "Error --- PLegendreA_d1"
print*, "P must be dimensioned as (LMAX+1)*(LMAX+2)/2 " // &
"where LMAX is ", lmax
print*, "Input array is dimensioned ", size(p)
stop
else if (size(dp) < sdim) then
print*, "Error --- PLegendreA_d1"
print*, "DP must be dimensioned as (LMAX+1)*(LMAX+2)/2 " // &
"where LMAX is ", lmax
print*, "Input array is dimensioned ", size(dp)
stop
else if (lmax < 0) then
print*, "Error --- PLegendreA_d1"
print*, "LMAX must be greater than or equal to 0."
print*, "Input value is ", lmax
stop
elseif(abs(z) > 1.0d0) then
print*, "Error --- PLegendreA_d1"
print*, "ABS(Z) must be less than or equal to 1."
print*, "Input value is ", z
stop
else if (abs(z) == 1.0d0) then
print*, "Error --- PLegendreA_d1"
print*, "Derivative can not be calculated at Z = 1 or -1."
print*, "Input value is ", z
stop
end if
if (present(csphase)) then
if (csphase == -1) then
phase = -1
else if (csphase == 1) then
phase = 1
else
print*, "Error --- PLegendreA_d1"
print*, "CSPHASE must be 1 (exclude) or -1 (include)."
print*, "Input value is ", csphase
stop
end if
else
phase = CSPHASE_DEFAULT
end if
!---------------------------------------------------------------------------
!
! Calculate P(l,0)
!
!---------------------------------------------------------------------------
sinsq = (1.0d0-z) * (1.0d0+z)
sinsqr = sqrt(sinsq)
pm2 = 1.d0
p(1) = 1.d0
dp(1) = 0.0d0
if (lmax == 0) return
pm1 = z
p(2) = pm1
dp(2) = 1.0d0
k = 2
do l = 2, lmax, 1
k = k+l
plm = ( z*(2*l-1)*pm1-(l-1)*pm2 ) / dble(l)
p(k) = plm
dp(k) = l * (pm1 -z*plm) / sinsq
pm2 = pm1
pm1 = plm
end do
!---------------------------------------------------------------------------
!
! Calculate P(m,m), P(m+1,m), and P(l,m)
!
!---------------------------------------------------------------------------
pmm = 1.0d0
fact = -1.0d0
kstart = 1
do m = 1, lmax - 1, 1
! Calculate P(m,m)
kstart = kstart+m+1
fact = fact+2.0d0
pmm = phase * pmm * sinsqr * fact
p(kstart) = pmm
dp(kstart) = -m*z*pmm / sinsq
pm2 = pmm
! Calculate P(m+1,m)
k = kstart+m+1
pm1 = z*pmm*(2*m+1)
p(k) = pm1
dp(k) = ( (2*m+1)*pmm - (m+1)*z*pm1) / sinsq
! Calculate P(l,m)
do l = m+2, lmax, 1
k = k+l
plm = ( z*(2*l-1)*pm1-(l+m-1)*pm2 ) / dble(l-m)
p(k) = plm
dp(k) = ( (l+m) * pm1 - l*z*plm) / sinsq
pm2 = pm1
pm1 = plm
end do
end do
! P(lmax, lmax)
kstart = kstart+m+1
fact = fact+2.0d0
pmm = phase * pmm * sinsqr * fact
p(kstart) = pmm
dp(kstart) = -lmax*z*pmm / sinsq
end subroutine PLegendreA_d1