/
SHMultiply.f95
374 lines (324 loc) · 13.9 KB
/
SHMultiply.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
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
subroutine SHMultiply(shout, sh1, lmax1, sh2, lmax2, precomp, norm, &
csphase, exitstatus)
!------------------------------------------------------------------------------
!
! This subroutine will multiply two spherical harmonic fields which are
! expressed up to maximum spherical harmonic degrees lmax1 and lmax2. The
! output spherical harmonic coefficients will have a maximum spherical
! harmonic degree equal to lmax1 + lmax2.
!
! Calling Parameters:
!
! IN
! sh1 Spherical harmonic field with maximum spherical harmonic
! degree lmax1.
! sh2 Spherical harmonic field with maximun spherical harmonic
! degree lmax2.
! lmax1 Maximum spherical harmonic degree of sh1.
! lmax2 Maximum spherical harmonic degree of sh2.
!
! OUT
! shout Spherical harmonic expansion of spatial multiplication
! of sh1 and sh2, with a maximum spherical harmonic
! degree of lmax1 + lmax2.
!
! OPTIONAL (IN)
! precomp If 1, the array plx will be precomputed when calling
! the subroutine SHGLQ. If 0 (default), then this array
! will not be precomputed.
! csphase 1: Do not include the phase factor of (-1)^m
! -1: Apply the phase factor of (-1)^m.
! norm: Normalization to be used when calculating Legendre
! functions
! (1) "geodesy" (default)
! (2) Schmidt
! (3) unnormalized
! (4) orthonormalized
!
! 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.
!
! Copyright (c) 2005-2019, SHTOOLS
! All rights reserved.
!
!------------------------------------------------------------------------------
use SHTOOLS, only: SHGLQ, MakeGridGLQ, SHExpandGLQ, CSPHASE_DEFAULT
use ftypes
implicit none
real(dp), intent(out) :: shout(:,:,:)
real(dp), intent(in) :: sh1(:,:,:), sh2(:,:,:)
integer, intent(in) :: lmax1, lmax2
integer, intent(in), optional :: precomp, norm, csphase
integer, intent(out), optional :: exitstatus
integer :: lmaxout, phase, mnorm, astat(2), nlat, nlong
real(dp), allocatable, save :: zero(:), w(:)
integer, save :: first = 1, lmaxout_last = -1
real(dp), allocatable :: grid1glq(:,:), grid2glq(:,:), plx(:,:)
!$OMP threadprivate(zero, w, first, lmaxout_last)
if (present(exitstatus)) exitstatus = 0
if (size(sh1(:,1,1)) < 2 .or. size(sh1(1,:,1)) < lmax1+1 .or. &
size(sh1(1,1,:)) < lmax1+1) then
print*, "Error --- SHMultiply"
print*, "SH1 must be dimensioned as (2, LMAX1+1, LMAX1+1) " // &
"where LMAX1 is", lmax1
print*, "Input array is dimensioned ", size(sh1(:,1,1)), &
size(sh1(1,:,1)), size(sh1(1,1,:))
if (present(exitstatus)) then
exitstatus = 1
return
else
stop
end if
else if (size(sh2(:,1,1)) < 2 .or. size(sh2(1,:,1)) < lmax2+1 &
.or. size(sh2(1,1,:)) < lmax2+1) then
print*, "Error --- SHMultiply"
print*, "SH2 must be dimensioned as (2,LMAX2+1, LMAX2+1) " // &
"where LMAX2 is", lmax2
print*, "Input array is dimensioned ", size(sh2(:,1,1)), &
size(sh2(1,:,1)), size(sh2(1,1,:))
if (present(exitstatus)) then
exitstatus = 1
return
else
stop
end if
else if (size(shout(:,1,1)) < 2 .or. size(shout(1,:,1)) < lmax1+lmax2+1 &
.or. size(shout(1,1,:)) < lmax1+lmax2+1) then
print*, "Error --- SHMultiply"
print*, "SHOUT must be dimensioned as (2, LMAX1+LMAX2+1, " // &
"LMAX1+LMAX2+1) where LMAX1 and LMAX2 are", lmax1, lmax2
print*, "Input array is dimensioned ", size(shout(:,1,1)), &
size(shout(1,:,1)), size(shout(1,1,:))
if (present(exitstatus)) then
exitstatus = 1
return
else
stop
end if
end if
if (present(csphase)) then
if (csphase /= -1 .and. csphase /= 1) then
print*, "Error --- SHMultiply"
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
if (present(precomp)) then
if (precomp /= 1 .and. precomp /= 0) then
print*, "Error --- SHMultiply"
print*, "PRECOMP must be either 0 (do not precompute PLX) " // &
"or 1 (precompute PLX)."
print*, "Input value is ", precomp
if (present(exitstatus)) then
exitstatus = 2
return
else
stop
end if
end if
end if
if (present(norm)) then
if (norm > 4 .or. norm < 1) then
print*, "Error - SHMultiply"
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
mnorm = norm
else
mnorm = 1
end if
lmaxout = lmax1 + lmax2
nlat = lmax1 + lmax2 + 1
nlong = 2 * (lmax1 + lmax2) + 1
if (first == 1) then
first = 0
lmaxout_last = lmaxout
allocate (zero(lmaxout+1), stat = astat(1))
allocate (w(lmaxout+1), stat = astat(2))
if (sum(astat(1:2)) /= 0) then
print*, "Error --- SHMultiply"
print*, "Problem allocating arrays ZERO and W", astat(1), astat(2)
if (present(exitstatus)) then
exitstatus = 3
return
else
stop
end if
end if
if (present(exitstatus)) then
call SHGLQ(lmaxout, zero, w, csphase = phase, norm = mnorm, &
exitstatus = exitstatus)
if (exitstatus /= 0) return
else
call SHGLQ(lmaxout, zero, w, csphase = phase, norm = mnorm)
end if
end if
if (lmaxout /= lmaxout_last) then
lmaxout_last = lmaxout
deallocate (zero)
deallocate (w)
allocate (zero(lmaxout+1), stat = astat(1))
allocate (w(lmaxout+1), stat = astat(2))
if (sum(astat(1:2)) /= 0) then
print*, "Error --- SHMultiply"
print*, "Problem allocating arrays ZERO and W", astat(1), astat(2)
if (present(exitstatus)) then
exitstatus = 3
return
else
stop
end if
end if
if (present(exitstatus)) then
call SHGLQ(lmaxout, zero, w, csphase = phase, norm = mnorm, &
exitstatus = exitstatus)
if (exitstatus /= 0) return
else
call SHGLQ(lmaxout, zero, w, csphase = phase, norm = mnorm)
end if
end if
allocate (grid1glq(nlat, nlong), stat = astat(1))
allocate (grid2glq(nlat, nlong), stat = astat(2))
if (sum(astat(1:2)) /= 0) then
print*, "Error --- SHMultiply"
print*, "Problem allocating arrays GRID1GLQ and GRID2GLQ", &
astat(1), astat(2)
if (present(exitstatus)) then
exitstatus = 3
return
else
stop
end if
end if
if (present(precomp)) then
if (precomp == 0) then
if (present(exitstatus)) then
call MakeGridGLQ(grid1glq, sh1(1:2,1:lmax1+1, 1:lmax1+1), &
lmaxout, zero = zero, csphase = phase, &
norm = mnorm, exitstatus = exitstatus)
if (exitstatus /= 0) return
call MakeGridGLQ(grid2glq, sh2(1:2,1:lmax2+1, 1:lmax2+1), &
lmaxout, zero = zero, csphase = phase, &
norm = mnorm, exitstatus = exitstatus)
if (exitstatus /= 0) return
grid1glq(1:nlat,1:nlong) = grid1glq(1:nlat,1:nlong) &
* grid2glq(1:nlat,1:nlong)
call SHExpandGLQ(shout, lmaxout, grid1glq, w, zero = zero, &
csphase = phase, norm = mnorm, &
exitstatus = exitstatus)
if (exitstatus /= 0) return
else
call MakeGridGLQ(grid1glq, sh1(1:2,1:lmax1+1, 1:lmax1+1), &
lmaxout, zero = zero, csphase = phase, &
norm = mnorm)
call MakeGridGLQ(grid2glq, sh2(1:2,1:lmax2+1, 1:lmax2+1), &
lmaxout, zero = zero, csphase = phase, &
norm = mnorm)
grid1glq(1:nlat,1:nlong) = grid1glq(1:nlat,1:nlong) &
* grid2glq(1:nlat,1:nlong)
call SHExpandGLQ(shout, lmaxout, grid1glq, w, zero = zero, &
csphase = phase, norm = mnorm)
end if
else
allocate (plx(lmax1+lmax2+1, (lmax1+lmax2+1)*(lmax1+lmax2+2)/2 ), &
stat=astat(1))
if (astat(1) /= 0) then
print*, "Error --- SHMultiply"
print*, "Problem allocating array PLX", astat(1)
if (present(exitstatus)) then
exitstatus = 3
return
else
stop
end if
end if
if (present(exitstatus)) then
call SHGLQ(lmaxout, zero, w, plx = plx, csphase = phase, &
norm = mnorm, exitstatus = exitstatus)
if (exitstatus /= 0) return
call MakeGridGLQ(grid1glq, sh1(1:2,1:lmax1+1, 1:lmax1+1), &
lmaxout, plx = plx, csphase = phase, &
norm = mnorm, exitstatus = exitstatus)
if (exitstatus /= 0) return
call MakeGridGLQ(grid2glq, sh2(1:2,1:lmax2+1, 1:lmax2+1), &
lmaxout, plx = plx, csphase = phase, &
norm = mnorm, exitstatus = exitstatus)
if (exitstatus /= 0) return
grid1glq(1:nlat,1:nlong) = grid1glq(1:nlat,1:nlong) &
* grid2glq(1:nlat,1:nlong)
call SHExpandGLQ(shout, lmaxout, grid1glq, w, plx = plx, &
csphase = phase, norm = mnorm, &
exitstatus = exitstatus)
if (exitstatus /= 0) return
else
call SHGLQ(lmaxout, zero, w, plx = plx, csphase = phase, &
norm = mnorm)
call MakeGridGLQ(grid1glq, sh1(1:2,1:lmax1+1, 1:lmax1+1), &
lmaxout, plx = plx, csphase = phase, &
norm = mnorm)
call MakeGridGLQ(grid2glq, sh2(1:2,1:lmax2+1, 1:lmax2+1), &
lmaxout, plx = plx, csphase = phase, &
norm = mnorm)
grid1glq(1:nlat,1:nlong) = grid1glq(1:nlat,1:nlong) &
* grid2glq(1:nlat,1:nlong)
call SHExpandGLQ(shout, lmaxout, grid1glq, w, plx = plx, &
csphase = phase, norm = mnorm)
end if
deallocate (plx)
end if
else
if (present(exitstatus)) then
call MakeGridGLQ(grid1glq, sh1(1:2,1:lmax1+1, 1:lmax1+1), &
lmaxout, zero = zero, csphase = phase, &
norm = mnorm, exitstatus = exitstatus)
if (exitstatus /= 0) return
call MakeGridGLQ(grid2glq, sh2(1:2,1:lmax2+1, 1:lmax2+1), &
lmaxout, zero = zero, csphase = phase, &
norm = mnorm, exitstatus = exitstatus)
if (exitstatus /= 0) return
grid1glq(1:nlat,1:nlong) = grid1glq(1:nlat,1:nlong) &
* grid2glq(1:nlat,1:nlong)
call SHExpandGLQ(shout, lmaxout, grid1glq, w, zero = zero, &
csphase = phase, norm = mnorm, &
exitstatus = exitstatus)
if (exitstatus /= 0) return
else
call MakeGridGLQ(grid1glq, sh1(1:2,1:lmax1+1, 1:lmax1+1), &
lmaxout, zero = zero, csphase = phase, &
norm = mnorm)
call MakeGridGLQ(grid2glq, sh2(1:2,1:lmax2+1, 1:lmax2+1), &
lmaxout, zero = zero, csphase = phase, &
norm = mnorm)
grid1glq(1:nlat,1:nlong) = grid1glq(1:nlat,1:nlong) &
* grid2glq(1:nlat,1:nlong)
call SHExpandGLQ(shout, lmaxout, grid1glq, w, zero = zero, &
csphase = phase, norm = mnorm)
end if
end if
deallocate(grid1glq)
deallocate(grid2glq)
end subroutine SHMultiply