-
Notifications
You must be signed in to change notification settings - Fork 106
/
CilmMinusRhoH.f95
628 lines (556 loc) · 23.7 KB
/
CilmMinusRhoH.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
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
subroutine CilmMinusRhoH(cilm, gridin, lmax, nmax, mass, d, rho, gridtype, w, &
zero, plx, n, dref, exitstatus)
!------------------------------------------------------------------------------
!
! This routine will compute the potential coefficients associated
! with the input gridded relief using the method of Wieczorek and
! Phillips (1998). The input grid must contain the degree-0 term
! and the computed coefficients will be referenced to the corresonding
! degree-0 radius. Note that the array plx is optional, and should not
! be precomputed when memory is an issue (i.e., lmax>360). In contrast
! to CilmPlus, takes into account lateral variations in density by
! inputing both a relief grid (GRIDIN) and density grid (RHO).
!
! Calling Parameters
!
! IN
! gridin Input grid to be transformed to spherical harmonics.
! lmax Maximum spherical harmonic degree to compute. For
! gridtype 3 and 4, this must be less than or equal
! to N/2 - 1.
! nmax Order of expansion.
! mass Mass of planet.
! rho Input grid of the density with same dimensions as
! GRIDIN.
! gridtype 1 = Gauss-Legendre quadrature grid corresponding
! to LMAX.
! 2 = N by N Driscoll and Healy grid corresponding
! to LMAX.
! 3 = N by 2N Driscoll and Healy grid corresponding
! to LMAX.
! (4 = 2D Cartesian using MakeGrid2D is not implemented).
!
! OUT
! cilm Output spherical harmonic coefficients with
! dimensions (2, lmax+1, lmax+1).
! d The radius that the coefficients are referenced
! to. This parameter corresponds to the degree zero term
! of the data.
!
! OPTIONAL
! w Gauss-Legendre points used in the integrations of
! dimension lmax+1.
! zero Array of dimension lmax+1 that contains the latitudinal
! gridpoints used in the Gauss-Legendre quadrature
! integration scheme. Only needed when PLX is not
! included (determined from a call to SHGLQ).
! plx Input array of Associated Legendre Polnomials computed
! at the Gauss points (determined from a call to
! SHGLQ). If this is not included, then the optional
! array zero MUST be inlcuded.
! N Number of latitude points in the Driscoll and Healy
! grid. Required for Gridtype 2 and 3.
! dref The reference radius used to be used when calculating
! the spherical harmonic coefficients. If not specified,
! this will be set equal to the mean radius of GRIDIN.
!
! 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.
!
! All units assumed to be SI.
!
! Copyright (c) 2005-2019, SHTOOLS
! All rights reserved.
!
!------------------------------------------------------------------------------
use SHTOOLS, only: SHExpandGLQ, SHExpandDH
use ftypes
implicit none
real(dp), intent(in) :: gridin(:,:), mass, rho(:,:)
real(dp), intent(in), optional :: w(:), zero(:), plx(:,:), dref
real(dp), intent(out) :: cilm(:,:,:), d
integer(int32), intent(in) :: lmax, nmax, gridtype
integer(int32), intent(in), optional :: n
integer(int32), intent(out), optional :: exitstatus
real(dp) :: prod, pi, scalef
real(dp), allocatable :: cilmn(:, :, :), grid(:, :)
integer(int32) :: j, l, k, nlat, nlong, astat(2), lmax_dh
if (present(exitstatus)) exitstatus = 0
pi = acos(-1.0_dp)
if (size(cilm(:,1,1)) < 2 .or. size(cilm(1,:,1)) < lmax+1 .or. &
size(cilm(1,1,:)) < lmax+1) then
print*, "Error --- CilmMinusRhoH"
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
end if
end if
if (gridtype == 4) then
print*, "Error --- CilmMinusRhoH"
print*, "GRIDTYPE 4 (Cartesian obtained from MakeGrid2D) is " // &
"not allowed."
if (present(exitstatus)) then
exitstatus = 2
return
else
stop
end if
else if (gridtype == 1) then
if (present(n)) then
print*, "Error --- CilmMinusRhoH"
print*, "N can not be present when using GLQ grids."
if (present(exitstatus)) then
exitstatus = 5
return
else
stop
end if
end if
if (present(w) .and. present(zero) .and. present(plx)) then
print*, "Error --- CilmMinusRhoH"
print*, "For GLQ grids, either W and ZERO or W and PLX " // &
"must be present, but not all three."
if (present(exitstatus)) then
exitstatus = 5
return
else
stop
end if
else if (present(w) .and. present(zero)) then
if (size(w) < lmax + 1) then
print*, "Error --- CilmMinusRhoH"
print*, "W must be dimensioned as (LMAX+1) where LMAX is ", lmax
print*, "Input dimension is ", size(w)
if (present(exitstatus)) then
exitstatus = 1
return
else
stop
end if
end if
if (size(zero) < lmax + 1) then
print*, "Error --- CilmMinusRhoH"
print*, "ZERO must be dimensioned as (LMAX+1) " // &
"where LMAX is ", lmax
print*, "Input dimension is ", size(zero)
if (present(exitstatus)) then
exitstatus = 1
return
else
stop
end if
end if
else if (present(plx) .and. present(w)) then
if (size(w) < lmax + 1) then
print*, "Error --- CilmMinusRhoH"
print*, "W must be dimensioned as (LMAX+1) where LMAX is ", lmax
print*, "Input dimension is ", size(w)
if (present(exitstatus)) then
exitstatus = 1
return
else
stop
end if
end if
if (size(plx(:,1)) < lmax+1 .or. &
size(plx(1,:)) < (lmax+1)*(lmax+2)/2) then
print*, "Error --- CilmMinusRhoH"
print*, "PLX must be dimensioned as (LMAX+1, " // &
"(LMAX+1)*(LMAX+2)/2) where LMAX is ", lmax
print*, "Input dimension is ", size(plx(:,1)), size(plx(1,:))
stop
end if
else
print*, "Error --- CilmMinusRhoH"
print*, "For GLQ grids, either W and ZERO or W and " // &
"PLX must be present"
if (present(exitstatus)) then
exitstatus = 5
return
else
stop
end if
end if
else if (gridtype == 2) then
if (present(w) .or. present(zero) .or. present(plx)) then
print*, "Error --- CilmMinusRhoH"
print*, "W, ZERO and PLX can not be present for " // &
"Driscoll-Healy grids."
if (present(exitstatus)) then
exitstatus = 5
return
else
stop
end if
else if (.not. present(N)) then
print*, "Error --- CilmMinusRhoH"
print*, "N must be present when GRIDTYPE is 2 or 3."
if (present(exitstatus)) then
exitstatus = 5
return
else
stop
end if
end if
else if (gridtype == 3) then
if (present(w) .or. present(zero) .or. present(plx)) then
print*, "Error --- CilmMinusRhoH"
print*, "W, ZERO and PLX can not be present for " // &
"Driscoll-Healy grids."
if (present(exitstatus)) then
exitstatus = 5
return
else
stop
end if
else if (.not. present(N)) then
print*, "Error --- CilmMinusRhoH"
print*, "N must be present when GRIDTYPE is 2 or 3."
if (present(exitstatus)) then
exitstatus = 5
return
else
stop
end if
end if
else
print*, "Error --- CilmMinusRhoH"
print*, "GRIDTYPE must be 2 (GLQ), 3 (NxN) or 4 (Nx2N)"
print*, "Input value is ", gridtype
if (present(exitstatus)) then
exitstatus = 2
return
else
stop
end if
end if
if (gridtype == 1) then
nlat = lmax+1
nlong = 2 * lmax + 1
else if (gridtype == 2) then
nlat = N
nlong = N
lmax_dh = N / 2 - 1
if (lmax > lmax_dh) then
print*, "Error --- CilmMinusRhoH"
print*, "For Driscoll-Healy grids, LMAX must be less " // &
"than or equal to N/2 -1, where N is ", N
print*, "Input value of LMAX is ", lmax
if (present(exitstatus)) then
exitstatus = 2
return
else
stop
end if
end if
else if (gridtype == 3) then
nlat = N
nlong = 2 * N
lmax_dh = N / 2 - 1
if (lmax > lmax_dh) then
print*, "Error --- CilmMinusRhoH"
print*, "For Driscoll-Healy grids, LMAX must be less " // &
"than or equal to N/2 -1, where N is ", N
print*, "Input value of LMAX is ", lmax
if (present(exitstatus)) then
exitstatus = 2
return
else
stop
end if
end if
end if
if (size(gridin(1,:)) < nlong .or. size(gridin(:,1)) < nlat .or. &
size(rho(1,:)) < nlong .or. size(rho(:,1)) < nlat) then
print*, "Error --- CilmMinusRhoH"
if (gridtype == 2) then
print*, "GRIDIN and RHO must be dimensioned as " // &
"(LMAX+1, 2*LMAX+1) where LMAX is ", lmax
print*, "Input dimension of GRIDIN is ", size(gridin(1,:)), &
size(gridin(:,1))
print*, "Input dimension of RHO is ", size(rho(1,:)), size(rho(:,1))
if (present(exitstatus)) then
exitstatus = 1
return
else
stop
end if
else if (gridtype == 3) then
print*, "GRIDIN and RHO must be dimensioned as " // &
"(N, N) where N is ", n
print*, "Input dimension of GRIDIN is ", size(gridin(1,:)), &
size(gridin(:,1))
print*, "Input dimension of RHO is ", size(rho(1,:)), size(rho(:,1))
if (present(exitstatus)) then
exitstatus = 1
return
else
stop
end if
else if (gridtype == 4) then
print*, "GRIDIN and RHO must be dimensioned as " // &
"(N, 2N) where N is ", n
print*, "Input dimension of GRIDIN is ", size(gridin(1,:)), &
size(gridin(:,1))
print*, "Input dimension of RHO is ", size(gridin(1,:)), &
size(gridin(:,1))
if (present(exitstatus)) then
exitstatus = 1
return
else
stop
end if
end if
end if
allocate (cilmn(2, lmax+1, lmax+1), stat=astat(1))
allocate (grid(nlat, nlong), stat=astat(2))
if (astat(1) /= 0 .or. astat(2) /= 0) then
print*, "Error --- CilmMinusRhoH"
print*, "Problem allocating arrays CILMN and GRID", astat(1), astat(2)
if (present(exitstatus)) then
exitstatus = 3
return
else
stop
end if
end if
!--------------------------------------------------------------------------
!
! Do the expansion.
!
!--------------------------------------------------------------------------
grid(1:nlat, 1:nlong) = gridin(1:nlat,1:nlong)
cilm = 0.0_dp
cilmn = 0.0_dp
! Determine reference radius, if necessary
if (present(dref)) then
d = dref
else
select case (gridtype)
case (1)
if (present(plx)) then
if (present(exitstatus)) then
call SHExpandGLQ(cilmn(1:2,1:lmax+1,1:lmax+1), lmax, &
grid(1:nlat,1:nlong), w, plx = plx, &
norm = 1, csphase = 1, lmax_calc = 0,&
exitstatus = exitstatus)
if (exitstatus /= 0) return
else
call SHExpandGLQ(cilmn(1:2,1:lmax+1,1:lmax+1), lmax, &
grid(1:nlat,1:nlong), w, plx = plx, &
norm = 1, csphase = 1, lmax_calc = 0)
end if
else
if (present(exitstatus)) then
call SHExpandGLQ(cilmn(1:2,1:lmax+1,1:lmax+1), lmax, &
grid(1:nlat,1:nlong), w, zero=zero, &
norm = 1, csphase = 1, lmax_calc = 0,&
exitstatus = exitstatus)
if (exitstatus /= 0) return
else
call SHExpandGLQ(cilmn(1:2,1:lmax+1,1:lmax+1), lmax, &
grid(1:nlat,1:nlong), w, zero=zero, &
norm = 1, csphase = 1, lmax_calc = 0)
end if
end if
case (2)
if (present(exitstatus)) then
call SHExpandDH(grid(1:nlat,1:nlong), n, &
cilmn(1:2,1:lmax+1,1:lmax+1), lmax_dh, &
norm = 1, sampling = 1, csphase = 1, &
lmax_calc = 0, exitstatus = exitstatus)
if (exitstatus /= 0) return
else
call SHExpandDH(grid(1:nlat,1:nlong), n, &
cilmn(1:2,1:lmax+1,1:lmax+1), lmax_dh, &
norm = 1, sampling = 1, csphase = 1, &
lmax_calc = 0)
end if
case (3)
if (present(exitstatus)) then
call SHExpandDH(grid(1:nlat,1:nlong), n, &
cilmn(1:2,1:lmax+1,1:lmax+1), lmax_dh, &
norm = 1, sampling = 2, csphase = 1, &
lmax_calc = 0, exitstatus = exitstatus)
if (exitstatus /= 0) return
else
call SHExpandDH(grid(1:nlat,1:nlong), n, &
cilmn(1:2,1:lmax+1,1:lmax+1), lmax_dh, &
norm = 1, sampling = 2, csphase = 1, &
lmax_calc = 0)
end if
end select
d = cilmn(1,1,1) ! use mean radius of relief for reference sphere
end if
! Do k = 1 terms first
grid(1:nlat, 1:nlong) = (gridin(1:nlat,1:nlong) - d) &
* rho(1:nlat,1:nlong) / d
select case (gridtype)
case (1)
if (present(plx)) then
if (present(exitstatus)) then
call SHExpandGLQ(cilmn(1:2,1:lmax+1,1:lmax+1), lmax, &
grid(1:nlat,1:nlong), w, plx = plx, &
norm = 1, csphase = 1, &
exitstatus = exitstatus)
if (exitstatus /= 0) return
else
call SHExpandGLQ(cilmn(1:2,1:lmax+1,1:lmax+1), lmax, &
grid(1:nlat,1:nlong), w, plx = plx, &
norm = 1, csphase = 1)
end if
else
if (present(exitstatus)) then
call SHExpandGLQ(cilmn(1:2,1:lmax+1,1:lmax+1), lmax, &
grid(1:nlat,1:nlong), w, zero=zero, &
norm = 1, csphase = 1, &
exitstatus = exitstatus)
if (exitstatus /= 0) return
else
call SHExpandGLQ(cilmn(1:2,1:lmax+1,1:lmax+1), lmax, &
grid(1:nlat,1:nlong), w, zero=zero, &
norm = 1, csphase = 1)
end if
end if
case (2)
if (present(exitstatus)) then
call SHExpandDH(grid(1:nlat,1:nlong), n, &
cilmn(1:2,1:lmax+1,1:lmax+1), lmax_dh, &
norm = 1, sampling = 1, csphase = 1, &
lmax_calc = lmax, exitstatus = exitstatus)
if (exitstatus /= 0) return
else
call SHExpandDH(grid(1:nlat,1:nlong), n, &
cilmn(1:2,1:lmax+1,1:lmax+1), lmax_dh, &
norm = 1, sampling = 1, csphase = 1, &
lmax_calc = lmax)
end if
case (3)
if (present(exitstatus)) then
call SHExpandDH(grid(1:nlat,1:nlong), n, &
cilmn(1:2,1:lmax+1,1:lmax+1), lmax_dh, &
norm = 1, sampling = 2, csphase = 1, &
lmax_calc = lmax, exitstatus = exitstatus)
if (exitstatus /= 0) return
else
call SHExpandDH(grid(1:nlat,1:nlong), n, &
cilmn(1:2,1:lmax+1,1:lmax+1), lmax_dh, &
norm = 1, sampling = 2, csphase = 1, &
lmax_calc = lmax)
end if
end select
do l = 0, lmax
cilm(1:2,l+1,1:l+1) = 4.0_dp * pi * (d**3) * cilmn(1:2,l+1,1:l+1) &
/ mass / dble(2*l+1)
end do
scalef = maxval(abs(gridin(1:nlat,1:nlong) - d))
do k = 2, nmax
grid(1:nlat,1:nlong) = rho(1:nlat,1:nlong) &
* ((gridin(1:nlat,1:nlong) - d) / scalef)**k
select case (gridtype)
case (1)
if (present(plx)) then
if (present(exitstatus)) then
call SHExpandGLQ(cilmn(1:2,1:lmax+1,1:lmax+1), lmax, &
grid(1:nlat,1:nlong), w, plx=plx, &
norm = 1, csphase = 1, &
exitstatus = exitstatus)
if (exitstatus /= 0) return
else
call SHExpandGLQ(cilmn(1:2,1:lmax+1,1:lmax+1), lmax, &
grid(1:nlat,1:nlong), w, plx=plx, &
norm = 1, csphase = 1)
end if
else
if (present(exitstatus)) then
call SHExpandGLQ(cilmn(1:2,1:lmax+1,1:lmax+1), lmax, &
grid(1:nlat,1:nlong), w, zero=zero, &
norm = 1, csphase = 1, &
exitstatus = exitstatus)
if (exitstatus /= 0) return
else
call SHExpandGLQ(cilmn(1:2,1:lmax+1,1:lmax+1), lmax, &
grid(1:nlat,1:nlong), w, zero=zero, &
norm = 1, csphase = 1)
end if
end if
case(2)
if (present(exitstatus)) then
call SHExpandDH(grid(1:nlat,1:nlong), n, &
cilmn(1:2,1:lmax+1,1:lmax+1), lmax_dh, &
norm = 1, sampling = 1, csphase = 1, &
lmax_calc = lmax, exitstatus = exitstatus)
if (exitstatus /= 0) return
else
call SHExpandDH(grid(1:nlat,1:nlong), n, &
cilmn(1:2,1:lmax+1,1:lmax+1), lmax_dh, &
norm = 1, sampling = 1, csphase = 1, &
lmax_calc = lmax)
end if
case(3)
if (present(exitstatus)) then
call SHExpandDH(grid(1:nlat,1:nlong), n, &
cilmn(1:2,1:lmax+1,1:lmax+1), lmax_dh, &
norm = 1, sampling = 2, csphase = 1, &
lmax_calc = lmax, exitstatus = exitstatus)
if (exitstatus /= 0) return
else
call SHExpandDH(grid(1:nlat,1:nlong), n, &
cilmn(1:2,1:lmax+1,1:lmax+1), lmax_dh, &
norm = 1, sampling = 2, csphase = 1, &
lmax_calc = lmax)
end if
end select
do l = 0, lmax
prod = 4.0_dp * pi * (d**3) / mass * (scalef / d)**k
do j = 2, k, 1
prod = prod * dble(l+j-3)
end do
prod = prod / (dble(2*l+1) * dble(fact(k)))
cilm(1:2,l+1,1:l+1) = cilm(1:2,l+1,1:l+1) &
+ cilmn(1:2,l+1,1:l+1) * prod
end do
end do
deallocate(cilmn)
deallocate(grid)
CONTAINS
function fact(i)
!----------------------------------------------------------------------
!
! This function computes the factorial of an integer.
!
!----------------------------------------------------------------------
implicit none
integer(int32) :: i, j
real(dp) :: fact
if (i == 0) then
fact = 1.0_dp
else if (i < 0) then
print*, "Argument to FACT must be positive"
if (present(exitstatus)) then
exitstatus = 2
return
else
stop
end if
else
fact = 1.0_dp
do j = 1, i
fact = fact * j
end do
end if
end function fact
end subroutine CilmMinusRhoH