-
Notifications
You must be signed in to change notification settings - Fork 106
/
MakeGrid2D.f95
416 lines (356 loc) · 14.5 KB
/
MakeGrid2D.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
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
subroutine MakeGrid2D(grid, cilm, lmax, interval, nlat, nlong, norm, csphase, &
f, a, north, south, east, west, dealloc, exitstatus)
!------------------------------------------------------------------------------
!
! Given the Spherical Harmonic coefficients cilm, this subroutine
! will compute a 2D grid with equal latitude and longitude spacings.
! Note that since this is NOT done using FFTs, this routine is therefore
! relatively SLOW! You might want to instead use MakeGridDH.
!
! The value at grid(1,1) correspons to 90 degrees latitude
! and 0 degrees longitude, and the longitude spacing goes from 0 to 360,
! with both points being calculated. If the optional
! parameters NORTH, SOUTH, EAST and WEST are specified, the upper-left and
! lower right coordinates of the output grid are (NORTH, WEST) and
! (SOUTH, EAST), respectively.
!
! Calling Parameters
!
! IN
! cilm Spherical harmonic coefficients.
! lmax Maximum degree of expansions to be performed.
! interval Spacing of output grid in DEGREES.
!
! OUT
! grid Gridded expansion of spherical harmonic coefficients.
! nlat Number of latitude points for the grid.
! nlong Number of longitude points for the grid.
!
! OPTIONAL (IN)
! norm Spherical harmonic normalization:
! (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.
! f Flattening of the function. If included, an ellipsoid
! with these parameters will be subtracted from the data.
! a Semimajor axis of the function. If included,
! an ellipsoid with these parameters will be subtracted
! from the data.
! north Maximum latitude to compute, in degrees.
! south Minimum latitude to compute, in degrees.
! east Maximum longitude to compute, in degrees.
! west Minimum latitude to compute, in degrees.
! dealloc 0 (default): Do not deallocate saved variables in the
! Legendre function routines. (1) Dellocate this memory
! at the end of the routine.
!
! 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.
!
! Notes:
! 1. If lmax is greater than the the maximum spherical harmonic
! degree of the input file, then this file will be ZERO PADDED!
! (i.e., those degrees after lmax are assumed to be zero).
!
! Dependencies: PlmBar, PlBar, PlmSchmidt, PlSchmidt, PLegendreA, &
! PLegendre, PlmON, CSPHASE_DEFAULT
!
! Copyright (c) 2016, SHTOOLS
! All rights reserved.
!
!------------------------------------------------------------------------------
use SHTOOLS, only: PlmBar, PlBar, PlmSchmidt, PlSchmidt, PLegendreA, &
PLegendre, PlmON, PlON, CSPHASE_DEFAULT
implicit none
real*8, intent(in) :: cilm(:,:,:), interval
real*8, intent(out) :: grid(:,:)
integer, intent(in) :: lmax
integer, intent(out) :: nlat, nlong
integer, intent(in), optional :: norm, csphase, dealloc
real*8, intent(in), optional :: f, a, north, south, east, west
integer, intent(out), optional :: exitstatus
integer :: l, m, j, k, index, l1, m1, lmax_comp, phase, lnorm, temp, &
astat(4)
real*8 :: pi, latmax, latmin, longmin, longmax, lat, longitude, &
x, intervalrad, r_ref
real*8, allocatable :: pl(:), cosm(:, :), sinm(:, :), cilm2(:,:, :)
if (present(exitstatus)) exitstatus = 0
temp = 0
if (present(north)) temp = temp + 1
if (present(south)) temp = temp + 1
if (present(east)) temp = temp + 1
if (present(west)) temp = temp + 1
if (temp /=0 .and. temp /=4) then
print*, "Error --- MakeGrid2d"
print*, "The optional parameters NORTH, SOUTH, EAST, and WEST " // &
"must all be specified", present(north), present(south), &
present(east), present(west)
if (present(exitstatus)) then
exitstatus = 5
return
else
stop
endif
end if
if (temp == 4) then
latmax = north
latmin = south
longmin = west
longmax = east
if (latmax < latmin) then
print*, "Error --- MakeGrid2d"
print*, "NORTH must be larger than SOUTH."
print*, "NORTH = ", latmax
print*, "SOUTH = ", latmin
if (present(exitstatus)) then
exitstatus = 2
return
else
stop
endif
end if
if (latmax > 90.0d0 .or. latmin < -90.0d0) then
print*, "Error --- MakeGrid2d"
print*, "NORTH and SOUTH must lie between 90 and -90."
print*, "NORTH = ", latmax
print*, "SOUTH = ", latmin
if (present(exitstatus)) then
exitstatus = 2
return
else
stop
endif
end if
if (longmin > longmax) longmax = longmax + 360.0d0
else
latmax = 90.0d0
latmin = -90.0d0
longmin = 0.0d0
longmax = 360.d0
end if
nlat = int((latmax - latmin) / interval + 1)
nlong = int((longmax - longmin) / interval + 1)
if (present(norm)) then
if (norm > 4 .or. norm < 1) then
print*, "Error --- MakeGrid2d"
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
endif
end if
lnorm = norm
else
lnorm = 1
end if
if (present(csphase)) then
if (csphase /= -1 .and. csphase /= 1) then
print*, "Error --- MakeGrid2D"
print*, "CSPHASE must be 1 (exclude) or -1 (include)"
print*, "Input valuse is ", csphase
if (present(exitstatus)) then
exitstatus = 2
return
else
stop
endif
else
phase = csphase
end if
else
phase = CSPHASE_DEFAULT
end if
if (size(grid(:,1)) < nlat .or. size(grid(1,:)) < nlong ) then
print*, "Error --- MakeGrid2D"
print*, "GRID must be dimensioned ( (LATMAX-LATMIN)/INTERVAL+1, " // &
"(LONGMAX-LONGMIN)/INTERVAL+1 ) where"
print*, "INTERVAL = ", interval
print*, "LATMAX = ", latmax
print*, "LATMIN = ", latmin
print*, "LONGMIN = ", longmin
print*, "LONGMAX = ", longmax
print*, "Input array is dimensioned ", size(grid(:,1)), size(grid(1,:))
if (present(exitstatus)) then
exitstatus = 1
return
else
stop
endif
elseif (size(cilm(:,1,1)) < 2 .or. size(cilm(1,:,1)) < lmax+1 .or. &
size(cilm(1,1,:)) < lmax+1) then
print*, "Error --- MakeGrid2D"
print*, "CILM must be dimensioned as (2, LMAX+1, LMAX+1) " // &
"where LMAX is ", lmax
print*, "Input dimension is ", size(cilm(:,1,1)), size(cilm(1,:,1)), &
size(cilm(1,1,:))
if (present(exitstatus)) then
exitstatus = 1
return
else
stop
endif
end if
if ((present(f) .and. .not. present(a)) .or. (present(a) .and. &
.not. present(f)) ) then
print*, "Error --- MakeGrid2D"
print*, "Both F and A must be present"
print*, "F ", present(f)
print*, "A ", present(a)
if (present(exitstatus)) then
exitstatus = 5
return
else
stop
endif
end if
allocate (pl(((lmax+1) * (lmax+2)) / 2), stat = astat(1))
allocate (cosm(lmax+1, int(360./interval + 1)), stat = astat(2))
allocate (sinm(lmax+1, int(360./interval + 1)), stat = astat(3))
allocate (cilm2(2,lmax+1,lmax+1), stat = astat(4))
if (astat(1) /= 0 .or. astat(2) /= 0 .or. astat(3) /=0 &
.or. astat(4) /= 0) then
print*, "Error --- MakeGrid2D"
print*, "Problem allocating arrays PL, COSM, SINM, and CILM2", &
astat(1), astat(2), astat(3), astat(4)
if (present(exitstatus)) then
exitstatus = 3
return
else
stop
endif
end if
pi = acos(-1.0d0)
grid = 0.0d0
lmax_comp = min(lmax, size(cilm(1,1,:))-1)
intervalrad = interval*pi/180.0d0
cilm2(1:2,1:lmax+1,1:lmax+1) = cilm(1:2,1:lmax+1,1:lmax+1)
!--------------------------------------------------------------------------
!
! Precomputing sines and cosines leads to an increase in speed by a
! factor of almost 4 with no optimization, and by a factor of about 15
! with normal optimizations.
!
!--------------------------------------------------------------------------
do k = 1, nlong
longitude = longmin * pi / 180.0d0 + dble(k-1) * intervalrad
sinm(1,k) = 0.0d0
cosm(1,k) = 1.0d0
if (lmax > 0) then
sinm(2,k) = sin(longitude)
cosm(2,k) = cos(longitude)
end if
do m = 2, lmax, 1
sinm(m+1,k) = 2 * sinm(m,k)*cosm(2,k) - sinm(m-1,k)
cosm(m+1,k) = 2 * cosm(m,k)*cosm(2,k) - cosm(m-1,k)
end do
end do
do j = 1, nlat
lat = latmax - (j-1)*interval
x = sin(lat*pi/180.0d0)
if (lat == 90.0d0 .or. lat == -90.0d0) then
if (present(exitstatus)) then
select case (lnorm)
case (1); call PlBar(pl, lmax_comp, x, &
exitstatus = exitstatus)
case (2); call PlSchmidt(pl, lmax_comp, x, &
exitstatus = exitstatus)
case (3); call PLegendre(pl, lmax_comp, x, &
exitstatus = exitstatus)
case (4); call PlON(pl, lmax_comp, x, &
exitstatus = exitstatus)
end select
if (exitstatus /= 0) return
else
select case (lnorm)
case (1); call PlBar(pl, lmax_comp, x)
case (2); call PlSchmidt(pl, lmax_comp, x)
case (3); call PLegendre(pl, lmax_comp, x)
case (4); call PlON(pl, lmax_comp, x)
end select
endif
do l = lmax_comp, 0, -1
l1 = l + 1
grid(j,1) = grid(j,1) + cilm2(1,l1,1) * pl(l1)
end do
grid(j, 2:nlong) = grid(j,1)
if (present(f)) grid(j,1:nlong) = grid(j,1:nlong) - a * (1.0d0 - f)
else
if (present(exitstatus)) then
select case (lnorm)
case (1); call PlmBar(pl, lmax_comp, x, csphase = phase, &
exitstatus = exitstatus)
case (2); call PlmSchmidt(pl, lmax_comp, x, &
csphase = phase, &
exitstatus = exitstatus)
case (3); call PLegendreA(pl, lmax_comp, x, &
csphase = phase, &
exitstatus = exitstatus)
case (4); call PlmON(pl, lmax_comp, x, csphase = phase, &
exitstatus = exitstatus)
end select
if (exitstatus /= 0) return
else
select case (lnorm)
case (1); call PlmBar(pl, lmax_comp, x, csphase = phase)
case (2); call PlmSchmidt(pl, lmax_comp, x, csphase = phase)
case (3); call PLegendreA(pl, lmax_comp, x, csphase = phase)
case (4); call PlmON(pl, lmax_comp, x, csphase = phase)
end select
endif
do k = 1, nlong
! do m = 0 term first
m1 = 1
do l = 0, lmax_comp, 1
l1 = l + 1
index = ((l+1) * l) / 2 + 1
grid(j,k) = grid(j,k) + cilm2(1,l1,1) * cosm(1,k) * &
pl(index)
end do
do m = 1, lmax_comp, 1
m1 = m + 1
do l = m, lmax_comp, 1
l1 = l + 1
index = ((l+1) * l) / 2 + m + 1
grid(j,k) = grid(j,k) + (cilm2(1,l1,m1) * cosm(m1,k) +&
cilm2(2,l1,m1) * sinm(m1,k)) &
* pl(index)
end do
end do
end do
if (present(f)) then
r_ref = a**2 * (1.0d0 + tan(lat * pi / 180.0d0)**2) / &
(1.0d0 + tan(lat * pi / 180.0d0)**2 / (1.0d0 - f)**2)
r_ref = sqrt(r_ref)
grid(j,1:nlong) = grid(j,1:nlong) - r_ref
end if
end if
end do
! deallocate memory
if (present(dealloc)) then
if (dealloc == 1) then
select case (lnorm)
case (1); call PlmBar(pl, -1, x, csphase = phase)
case (2); call PlmSchmidt(pl, -1, x, csphase = phase)
case (4); call PlmON(pl, -1, x, csphase = phase)
end select
end if
end if
deallocate (pl)
deallocate (cosm)
deallocate (sinm)
deallocate (cilm2)
end subroutine MakeGrid2D