-
Notifications
You must be signed in to change notification settings - Fork 106
/
SHGLQ.f95
278 lines (233 loc) · 9.16 KB
/
SHGLQ.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
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
subroutine SHGLQ(lmax, zero, w, plx, norm, csphase, cnorm, exitstatus)
!------------------------------------------------------------------------------
!
! This routine will precompute the zeros and weights that are
! necessary for the Gauss-Legendre quadrature sherical harmonic
! expansion routines. This routine can also calculate a matrix of
! legendre functions on the grid nodes to speed up computations
! in the case where more than one expansion is performed.
! Note, if memory is an issue because of large spherical harmonic
! degrees, DO NOT PRECOMPUTE THE ARRAY PLX!
!
! Calling Parameters
!
! IN
! lmax Maximum spherical harmonic degree to
! be used in all following expansions.
!
! OUT
! zero An array of length (lmax + 1) which contains the Gauss
! points (i.e., zeros of P(l+1,m=0)).
! w An array of length lmax + 1 which contains the weights
! used in the Gauss-Legendre quadrature routines.
!
! OPTIONAL
! plx Array of normalized associated Legendre Polnomials
! computed at the Gauss points. The size of this array
! must be (1:lmax+1, 1:(lmax+1)*(lmax+2)/2).
! norm Normalization to be used when calculating Legendre
! functions
! (1) "geodesy" (default)
! (2) Schmidt
! (3) unnormalized
! (4) orthonormalized
! csphase 1: Do not include the phase factor of (-1)^m
! -1: Apply the phase factor of (-1)^m.
! cnorm: 1: compute complex normalized Legendre functions.
! 0 (default): compute real normalized Legendre functions.
!
! OPTIONAL (OUT)
! exitstatus If present, instead of executing a STOP when an error
! is encountered, the variable exitstatus will be
! returned describing the error.
! 0 = No errors;
! 1 = Improper dimensions of input array;
! 2 = Improper bounds for input variable;
! 3 = Error allocating memory;
! 4 = File IO error.
!
! Dependencies: PreGLQ, PlmBar, PlmSchmidt, PLegendreA, PlmON,
! CSPHASE_DEFAULT, PlmIndex
!
! Copyright (c) 2016, SHTOOLS
! All rights reserved.
!
!------------------------------------------------------------------------------
use SHTOOLS, only: PreGLQ, PlmBar, PlmSchmidt, PLegendreA, PlmON, &
CSPHASE_DEFAULT, PlmIndex
implicit none
integer, intent(in) :: lmax
real*8, intent(out) :: zero(:), w(:)
real*8, intent(out), optional :: plx(:,:)
integer, intent(in), optional :: norm, csphase, cnorm
integer, intent(out), optional :: exitstatus
integer :: n, i, astat, phase, l, m, i_s, cnormin, lnorm
real*8 :: upper, lower, pi
real*8, allocatable :: pl(:)
if (present(exitstatus)) exitstatus = 0
if (size(zero) < lmax+1) then
print*, "Error --- SHGLQ"
print*, "ZERO must be dimensioned as (LMAX+1) where LMAX is ", lmax
print*, "Input array is dimensioned ", size(zero)
if (present(exitstatus)) then
exitstatus = 1
return
else
stop
end if
else if (size(w) < lmax+1) then
print*, "Error --- SHGLQ"
print*, "W must be dimensioned as (LMAX+1) where LMAX is ", lmax
print*, "Input array is dimensioned ", size(w)
if (present(exitstatus)) then
exitstatus = 1
return
else
stop
end if
end if
if (present(plx)) then
if (size(plx(:,1)) < lmax +1 .or. size(plx(1,:)) &
< (lmax+1)*(lmax+2)/2) then
print*, "Error --- SHGLQ"
print*, "PLX must be dimensioned as (LMAX+1, " // &
"(LMAX+1)*(LMAX+2)/2) where LMAX is ", lmax
print*, "Input array is dimensioned as ", size(plx(:,1)), &
size(plx(1,:))
if (present(exitstatus)) then
exitstatus = 1
return
else
stop
end if
end if
end if
if (present(norm)) then
if (norm > 4 .or. norm < 1) then
print*, "Error --- SHGLQ"
print*, "Parameter NORM must be 1 (geodesy), 2 (Schmidt), " // &
"3 (unnormalized), or 4 (orthonormalized)."
print*, "Input value is ", norm
if (present(exitstatus)) then
exitstatus = 2
return
else
stop
end if
end if
lnorm = norm
else
lnorm = 1
end if
if (present(csphase)) then
if (csphase /= -1 .and. csphase /= 1) then
print*, "Error --- SHGLQ"
print*, "CSPHASE must be 1 (exclude) or -1 (include)."
print*, "Input value is ", csphase
if (present(exitstatus)) then
exitstatus = 2
return
else
stop
end if
else
phase = csphase
end if
else
phase = CSPHASE_DEFAULT
end if
allocate (pl(((lmax+2)*(lmax+1))/2), stat = astat)
if (astat /= 0) then
print*, "Error --- SHGLQ"
print*, "Problem allocating array PL", astat
if (present(exitstatus)) then
exitstatus = 3
return
else
stop
end if
end if
if (present(cnorm)) then
cnormin = cnorm
else
cnormin = 0
end if
pi = acos(-1.0d0)
upper = 1.0d0
lower = -1.0d0
n = lmax + 1
! lmax is the maxium degree of the equation we
! are going to integrate.
! Determine Gauss Points and Weights.
if (present(exitstatus)) then
call PreGLQ(lower, upper, n, zero, w, exitstatus = exitstatus)
if (exitstatus /= 0) return
else
call PreGLQ(lower, upper, n, zero, w)
endif
!--------------------------------------------------------------------------
!
! Next, fill the array plx with all of the associated legendre polynomials
! evaluated at the Gauss points. The first index of this array corresponds
! to the Gauss point, whereas the (lmax+2)*(lmax+1)/2 elements in the
! second column correspond to all of the assiated legendre polynomials,
! ordered according to index = (l+1)*l/2 + m + 1.
!
!--------------------------------------------------------------------------
if (present(plx)) then
do i = 1, (n + 1) / 2
if (present(exitstatus)) then
select case(lnorm)
case(1)
call PlmBar(pl, lmax, zero(i), csphase = phase, &
cnorm=cnormin, exitstatus = exitstatus)
case(2)
call PlmSchmidt(pl, lmax, zero(i), csphase = phase, &
cnorm=cnormin, exitstatus = exitstatus)
case(3)
call PLegendreA(pl, lmax, zero(i), csphase = phase, &
exitstatus = exitstatus)
case(4)
call PlmON(pl, lmax, zero(i), csphase = phase, &
cnorm=cnormin, exitstatus = exitstatus)
if (exitstatus /= 0) return
end select
else
select case(lnorm)
case(1)
call PlmBar(pl, lmax, zero(i), csphase = phase, &
cnorm=cnormin)
case(2)
call PlmSchmidt(pl, lmax, zero(i), csphase = phase, &
cnorm=cnormin)
case(3)
call PLegendreA(pl, lmax, zero(i), csphase = phase)
case(4)
call PlmON(pl, lmax, zero(i), csphase = phase, &
cnorm=cnormin)
end select
end if
if (i == (n + 1) / 2 .and. mod(n,2) /= 0) then
plx(i, 1:((lmax+1)*(lmax+2))/2) = pl(1:((lmax+1)*(lmax+2))/2)
else
i_s = n + 1 - i
plx(i, 1:((lmax+1)*(lmax+2))/2) = pl(1:((lmax+1)*(lmax+2))/2)
do l = 0, lmax
do m = 0, l
plx(i_s,PlmIndex(l,m)) = plx(i,PlmIndex(l,m)) &
* (-1.0d0)**(l-m)
end do
end do
end if
end do
select case(lnorm)
case(1)
call PlmBar(pl, -1, zero(1), csphase = phase, cnorm=cnormin)
case(2)
call PlmSchmidt(pl, -1, zero(1), csphase = phase, cnorm=cnormin)
case(4)
call PlmON(pl, -1, zero(1), csphase = phase, cnorm=cnormin)
end select
end if
deallocate (pl)
end subroutine SHGLQ