forked from GomezGroup/1DSeaLevelModel_FWTW
-
Notifications
You must be signed in to change notification settings - Fork 3
/
Copy pathspharmt.f90
2050 lines (1771 loc) · 61.8 KB
/
spharmt.f90
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
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
module spharmt
!
! fortran95 spherical harmonic module.
! For gaussian grids and triangular truncation only.
!
! version 1.1 9/30/2003 (second version - now all fortran using fft99)
! Jeff Whitaker <[email protected]>
! Modified by Sam Goldberg to add subroutines that perform SH transforms using
! 2-dimensional matrices of coefficients (U(l,m)) instead of 1-dimensional vectors
! - spat2spec, spec2spat
!-----------------------------------------------------------------------------
! to use this module, put "use spharmt" at top of main program
! (just after program statement).
! to initialize arrays needed to compute spherical harmonic transforms
! at T(ntrunc) resolution on a nlon x nlat gaussian grid,
! do something like this:
!
! type (sphere) :: sphere_dat
! parameter (nlon=192,nlat=nlon/2)
! parameter (ntrunc=62)
! rsphere = 6.3712e6
! call spharmt_init(sphere_dat,nlon,nlat,ntrunc,re)
! derived data type "sphere" contains stuff needed for legendre
! transforms and fft. By creating multiple instances of sphere
! data types, spherical harmonic transforms on multiple grids can
! be done easily in the same code.
!
! Components of sphere derived data type are:
!
! NLON (integer) - number of longitudes
! NLAT (integer) - number of latitudes
! NTRUNC (integer) - triangular truncation limit
! RE (real) - radius of sphere (m).
! PNM (real pointer array dimension ((ntrunc+1)*(ntrunc+2)/2, nlat) -
! Associated legendre polynomials.
! HNM (real pointer array, same size as pnm) = -(1 - x*x) * d(pnm)/dx
! at x = sin(latitude).
! GAULATS (real pointer array dimension nlat) - sin(gaussian latitudes).
! WEIGHTS (real pointer array dimension nlat) - gaussian weights.
! GWRC (real pointer array dimension nlat) - gaussian weights divided by
! rsphere*(1-x**2).
! INDXM (integer pointer array dimension (ntrunc+1)*(ntrunc+2)/2) - value of
! zonal wavenumber m corresponding to spectral array index
! nm=1,(ntrunc+1)*(ntrunc+2)/2)
! INDXN (integer pointer array same size as indxm) - value of
! spherical harmonic degree n corresponding to spectral array index
! nm=1,(ntrunc+1)*(ntrunc+2)/2)
! LAP (real pointer array dimension (ntrunc+1)*(ntrunc+2)/2) -
! lapacian operator in spectral space =
! -n*(n+1)/rsphere**2, where n is degree of spherical harmonic.
! ILAP (real pointer array same size as lap) - inverse laplacian operator in
! spectral space.
! TRIGS, IFAX: arrays needed by Temperton FFT.
! ISINITIALIZED (logical) - true if derived data type has been
! initialized by call to spharm_init, false if not initialized.
!
! public routines:
!
! SPHARMT_INIT(sphere_dat,nlon,nlat,ntrunc,re): initialize a sphere object
! (sphere_dat - derived data type "sphere"). inputs are nlon (number of unique
! longitudes), nlat (number of gaussian latitudes), and re
! (radius of sphere in m). Must be called before anything else.
!
! SPHARMT_DESTROY(sphere_dat): cleans up pointer arrays allocated by
! spharmt_init.
!
! SPHARM(sphere_dat, ugrid, anm, idir): spherical harmonic transform
! (forward, i.e. grid to spectral, for idir=1 and backward for idir=-1).
! Arguments are sphere derived data type (sphere_dat), gridded data (ugrid),
! complex spectral coefficients (anm), and flag specifying direction of
! transform (idir). See "Import Details" below for information
! about indexing of grid and spectral arrays.
! GAULW(gaulats,weights): compute sin(gaussian latitudes) and gaussian
! weights. Number of latitudes determined by size of input arrays.
! LEGEND(x,pnm,hnm): compute associated legendre functions (pnm) and
! their derivates hnm = (X**2-1)*D(PNM)/DX at x = sin(latitude). The
! input arrays pnm and hnm should have dimension (ntrunc+1)*(ntrunc+2)/2,
! where ntrunc is triangular truncation limit (see Import Details below
! for a description of the spectral indexing).
! GETUV(sphere_dat,vrtnm,divnm,ug,vg): get U,V (u*cos(lat),v*cos(lat) on grid)
! from spectral coefficients of vorticity and divergence.
! Input: sphere_dat, spectral coeffs. of vort and div (vrtnm,divnm).
! Output: gridded U,V (ug,vg).
!
! GETVRTDIV(sphere_dat,vrtnm,divnm,ug,vg): get spectral coefficients of
! vorticity and divergence (vrtnm,divnm) from gridded U,V (ug,vg).
!
! COSGRAD(sphere_dat,chinm,uchig,vchig): compute cos(lat)*vector gradient.
! Inputs: sphere_dat, spectral coefficient array (chinm).
! Outputs: longitudinal and latitudinal components of gradient on the gaussian
! grid (uchig,vchig).
! SPECSMOOTH(sphere_dat,datagrid,smooth): isotropic spectral smoothing.
! Inputs: sphere_dat, smooth(ntrunc+1) (smoothing factor as a function
! of spherical harmonic degree), datagrid (gridded data to be smoothed).
! Outputs: datagrid (smoothed gridded data).
!
! SUMNM(sphere_dat,am,bm,anm,isign1,isign2):
! given the arrays of fourier coeffs, am and bm,
! compute the complex spherical harmonic coeffs (anm) of:
! isign1*( (1./rsphere*(1.-x**2))*d(ag)/d(lon) + (isign2/rsphere)*d(bg)/dx )
! where ag and bg are the grid pt counterparts of am,bm,
! isign1,isign2 are +1 or -1, rsphere is radius of sphere, x=sin(lat))
! am,bm can be computed from gridded data (ag,bg) using RFFT.
!
! RFFT(sphere_dat, data, coeff, idir): computes fourier harmonics in zonal
! direction of a gridded array. idir=+1 for forward (grid
! to fourier) and -1 for backward (fourier to grid) transform.
! data (nlon,nlat) contains gridded data, coeff (ntrunc+1,nlat) contains
! complex zonal fourier harmonics.
!
! Important Details:
!
! The gridded data is assumed to be oriented such that i=1 is the Greenwich
! meridian and j=1 is the northernmost point. Grid indices increase eastward
! and southward. The grid increment in longitude is 2*pi/nlon radians.
! For example nlon = 72 for a five degree grid. nlon must be greater than or
! equal to 4. The efficiency of the computation is improved when nlon is a
! product of small prime numbers.
! The spectral data is assumed to be in a complex array of dimension
! (NTRUNC+1)*(NTRUNC+2)/2. NTRUNC is the triangular truncation limit (NTRUNC=42
! for T42). NTRUNC must be <= nlon/2. Coefficients are ordered so that first
! (nm=1) is m=0,n=0, second is m=0,n=1, nm=mtrunc is m=0,n=mtrunc, nm=mtrunc+1
! is m=1,n=1, etc. In Fortran90 syntax, values of m (degree) and n (order)
! corresponding as a function of the index nm are:
! indxm = (/((m,n=m,mtrunc),m=0,mtrunc)/)
! indxn = (/((n,n=m,mtrunc),m=0,mtrunc)/)
! Conversely, the index nm as a function of m and n is:
! nm = sum((/(i,i=mtrunc+1,mtrunc-m+2,-1)/))+n-m+1
! The associated legendre polynomials are normalized so that the integral
! (pbar(n,m,theta)**2)*sin(theta) on the interval theta=0 to theta=pi is 1,
! where: pbar(m,n,theta) = sqrt((2*n+1)*factorial(n-m)/(2*factorial(n+m)))
! *sin(theta)**m/(2**n*factorial(n)) times the (n+m)th derivative of
! (x**2-1)**n with respect to x=cos(theta).
! note: theta = 0.5*pi - phi, where phi is latitude and theta is colatitude,
! Therefore, cos(theta) = sin(phi) and sin(theta) = cos(phi).
! Note that pbar(0,0,theta)=sqrt(2)/2,
! and pbar(1,0,theta)=0.5*sqrt(6.)*sin(lat).
!----------------------------------------------------------------------------
! explicit typing
implicit none
! everything private to module, unless otherwise specified.
private
public :: sphere,spharmt_init,spharmt_destroy,spharm,rfft,cosgrad,&
specsmooth,getuv,getvrtdiv,sumnm,gaulw,legend,spat2spec,spec2spat
type sphere
! rsphere is radius of sphere in m.
real :: rsphere = 0.
! nlons is number of longitudes (not including cyclic point).
integer :: nlons = 0
! nlats is number of gaussian latitudes.
integer :: nlats = 0
! for fft.
real, dimension(:), pointer :: trigs
integer, dimension(:), pointer :: ifax
! ntrunc is triangular truncation limit (e.g. 42 for T42 truncation)
integer :: ntrunc = 0
! gaulats is sin(gaussian lats), weights are gaussian weights,
! gwrc is gaussian weights divided by re*cos(lat)**2, lap is
! the laplacian operatore -n*(n+1)/re**2 (where n is the degree
! of the spherical harmonic),
! ilap is the inverse laplacian (1/lap, set to zero for n=0).
real, dimension(:), pointer :: gaulats, weights, gwrc, lap, ilap
! pnm is associated legendre polynomials, hnm is -(1-x**2)d(pnm)/dx,
! where x = gaulats.
real, dimension(:,:), pointer :: pnm, hnm
! indxm is zonal wavenumber index, indxn is spherical harmonic degree index.
integer, dimension(:), pointer :: indxm, indxn
logical :: isinitialized=.FALSE.
end type sphere
contains
subroutine spharmt_init(sphere_dat,nlon,nlat,ntrunc,re)
! initialize a sphere object.
integer, intent(in) :: nlon,nlat,ntrunc
real, intent(in) :: re
type (sphere), intent(inout) :: sphere_dat
real, dimension(:), allocatable :: pnm_tmp,hnm_tmp
double precision, dimension(:), allocatable :: gaulats_tmp,weights_tmp
integer :: nmdim,m,n,j
sphere_dat%nlons = nlon
sphere_dat%nlats = nlat
sphere_dat%ntrunc = ntrunc
sphere_dat%rsphere = re
nmdim = (ntrunc+1)*(ntrunc+2)/2
allocate(gaulats_tmp(nlat))
allocate(weights_tmp(nlat))
call gaulw(gaulats_tmp,weights_tmp)
allocate(sphere_dat%gwrc(nlat))
sphere_dat%gwrc = weights_tmp/(dble(re)*(1.d0-gaulats_tmp**2))
allocate(sphere_dat%gaulats(nlat))
allocate(sphere_dat%weights(nlat))
sphere_dat%gaulats = gaulats_tmp
sphere_dat%weights = weights_tmp
deallocate(weights_tmp)
allocate(sphere_dat%indxm(nmdim))
allocate(sphere_dat%indxn(nmdim))
allocate(sphere_dat%lap(nmdim))
allocate(sphere_dat%ilap(nmdim))
sphere_dat%indxm = (/((m,n=m,ntrunc),m=0,ntrunc)/)
sphere_dat%indxn = (/((n,n=m,ntrunc),m=0,ntrunc)/)
sphere_dat%lap(:)=-real(sphere_dat%indxn(:))*real(sphere_dat%indxn(:)+1)/re**2
sphere_dat%ilap(1) = 0.
sphere_dat%ilap(2:nmdim) = 1./sphere_dat%lap(2:nmdim)
allocate(sphere_dat%pnm(nmdim,nlat))
allocate(sphere_dat%hnm(nmdim,nlat))
allocate(pnm_tmp(nmdim))
allocate(hnm_tmp(nmdim))
do j=1,nlat
call legend(gaulats_tmp(j),pnm_tmp,hnm_tmp)
sphere_dat%pnm(:,j) = pnm_tmp(:)
sphere_dat%hnm(:,j) = hnm_tmp(:)
enddo
deallocate(gaulats_tmp)
deallocate(pnm_tmp)
deallocate(hnm_tmp)
allocate(sphere_dat%trigs((3*nlon/2)+1))
allocate(sphere_dat%ifax(13))
call set99(sphere_dat%trigs,sphere_dat%ifax,nlon)
sphere_dat%isinitialized = .TRUE.
end subroutine spharmt_init
subroutine spharmt_destroy(sphere_dat)
! deallocate pointers in sphere object.
type (sphere), intent(inout) :: sphere_dat
if (.not. sphere_dat%isinitialized) return
deallocate(sphere_dat%gaulats)
deallocate(sphere_dat%weights)
deallocate(sphere_dat%gwrc)
deallocate(sphere_dat%pnm)
deallocate(sphere_dat%hnm)
deallocate(sphere_dat%indxm)
deallocate(sphere_dat%indxn)
deallocate(sphere_dat%lap)
deallocate(sphere_dat%ilap)
deallocate(sphere_dat%trigs)
deallocate(sphere_dat%ifax)
end subroutine spharmt_destroy
subroutine gaulw(sinlats, wts)
! compute sin of gaussian latitudes and gaussian weights.
! uses the iterative method presented in Numerical Recipes.
double precision, intent (inout), dimension(:) :: sinlats, wts
integer :: itermax
integer :: i, iter, j, nlat, nprec
double precision :: pi, pp, p1, p2, p3, z, z1, converg, ten
ten = 10.d0
converg = ten * EPSILON(ten)
nprec = precision(converg)
converg = .1**nprec
itermax = 10
pi = 4.0*ATAN(1.0)
nlat = size(sinlats)
if (size(sinlats) .ne. size(wts)) then
print *, 'sinlats and wts must be same size in gaulw!'
stop
end if
do i=1,nlat
z = cos(pi*(i - 0.25)/(nlat + 0.5))
do iter=1,itermax
p1 = 1.0
p2 = 0.0
do j=1,nlat
p3 = p2
p2 = p1
p1 = ((2.0*j - 1.0)*z*p2 - (j - 1.0)*p3)/j
end do
pp = nlat*(z*p1 - p2)/(z*z - 1.0E+00)
z1 = z
z = z1 - p1/pp
if(ABS(z - z1) .LT. converg) go to 10
end do
print *, 'abscissas failed to converge in itermax iterations'
print *, 'stopping in gaulw!'
stop
10 continue
sinlats (i) = z
wts (i) = 2.0/((1.0 - z*z)*pp*pp)
end do
return
end subroutine gaulw
SUBROUTINE LEGEND(X,PMN,HMN)
!
! THIS SUBROUTINE COMPUTES ASSOCIATED LEGENDRE
! FUNCTIONS, PMN, AND THE DERIVATIVE QUANTITY
! HMN = -(1 - X*X) * D(PMN)/DX AT X = COS( COLATITUDE )
! GIVEN SOME COLATITUDE IN THE DOMAIN.
!
! ACCURACY:
!
! THE RECURSION RELATION EMPLOYED IS LESS ACCURATE
! NEAR THE POLES FOR M + N .GE. ROUGHLY 75 BECAUSE THE RECURSION
! INVOLVES THE DIFFERENCE OF NEARLY EQUAL NUMBERS, SO
! THE ASSOCIATED LEGENDRE FUNCTIONS ARE COMPUTED USING DOUBLE
! PRECISION ACCUMULATORS.
! SEE BELOUSOV, TABLES OF NORMALIZED ASSOCIATED LEGENDRE
! POLYNOMIALS, C. 1962, FOR INFORMATION ON THE ACCURATE
! COMPUTATION OF ASSOCIATED LEGENDRE FUNCTIONS.
real, dimension(:), intent(inout) :: PMN, HMN
double precision, intent(in) :: x
integer :: m,n,nm,i,nmax,np1,nmstrt,j,nmdim,ntrunc
double precision :: A, B, PROD, SINSQ, &
eps,epsnmp1,EPSPMN, PMNJ, PMNJM1, PMNJM2
!**** SET PARAMETERS FOR ENTRY INTO THE RECURSIVE FORMULAE.
SINSQ = 1.D0 - X * X
A = 1.D0
B = 0.D0
PROD = 1.D0
nmdim = size(pmn)
ntrunc = nint((-1.+sqrt(1+8*float(nmdim)))/2.)-1
if (size(pmn) .ne. size(hmn)) then
print *, 'pnm and hnm must be same size in subroutine legend!'
stop
end if
!**** LOOP FOR THE 'M' INDEX.
NMSTRT = 0
DO I = 1, NTRUNC+1
M = (I - 1)
NMAX = NTRUNC+1-M
! WHEN M=0 (I=1), STANDARD LEGENDRE POLYNOMIALS ARE
! GENERATED.
IF (M .NE. 0) THEN
A = A + 2.D0
B = B + 2.D0
PROD = PROD * SINSQ * A / B
END IF
!**** GENERATE PMN AND HMN FOR J = 1 AND 2.
PMNJM2 = SQRT(0.5D0 * PROD)
NM = NMSTRT + 1
PMN(NM) = PMNJM2
PMNJM1 = SQRT( DBLE(2 * M + 3) ) * X * PMNJM2
IF (NM .NE. NMDIM) PMN(NM+1) = PMNJM1
NP1 = M + 1
EPSNMP1 = SQRT( DBLE(NP1*NP1 - M*M) / DBLE(4*NP1*NP1 - 1) )
EPSPMN = X * PMNJM1 - EPSNMP1 * PMNJM2
HMN(NM) = DBLE(M) * EPSNMP1 * PMNJM1
IF (NM .NE. NMDIM) &
HMN(NM+1) = DBLE(M+1) * EPSPMN - &
DBLE(M+2) * EPSNMP1 * PMNJM2
!**** LOOP FOR THE 'N' INDEX.
! NOW APPLY THE RECURSION FORMULAE FOR J .GE. 3.
DO J = 3, NMAX
N = M + J - 1
NM = NMSTRT + J
EPS = SQRT( DBLE(N*N - M*M) / DBLE(4*N*N - 1) )
PMNJ = EPSPMN / EPS
PMN(NM) = PMNJ
! COMPUTE EPS * PMN FOR J+1.
EPSPMN = X * PMNJ - EPS * PMNJM1
! COMPUTE THE DERIVATIVE.
HMN(NM) = DBLE(N) * EPSPMN - &
DBLE(N+1) * EPS * PMNJM1
PMNJM2 = PMNJM1
PMNJM1 = PMNJ
ENDDO
NMSTRT = NMSTRT + NMAX
ENDDO
END SUBROUTINE LEGEND
subroutine rfft(sphere_dat, data, coeff, idir)
! real multiple fft (uses temperton fft991)
type (sphere), intent(in) :: sphere_dat
real, dimension(sphere_dat%nlons,sphere_dat%nlats), intent(inout) :: data
real, dimension(sphere_dat%nlats*(sphere_dat%nlons+2)) :: wrk1
real, dimension(sphere_dat%nlats*(sphere_dat%nlons+1)) :: wrk2
complex, dimension(sphere_dat%ntrunc+1,sphere_dat%nlats), intent(inout) :: coeff
integer :: nlons,nlats,ntrunc,mwaves,i,j,m,n
integer, intent(in) :: idir
if (.not. sphere_dat%isinitialized) then
print *, 'uninitialized sphere object in rfft!'
stop
end if
nlons = sphere_dat%nlons
nlats = sphere_dat%nlats
ntrunc = sphere_dat%ntrunc
if (ntrunc .gt. nlons/2) then
print *, 'ntrunc must be less than or equal to nlons in rfft'
stop
end if
mwaves = ntrunc+1
!==> forward transform.
if (idir .eq. +1) then
!==> copy the data into the work array.
! transforms are computed pairwise using a complex fft.
n = 0
wrk1 = 0.
do j=1,nlats
do i=1,nlons+2
n = n + 1
wrk1(n) = 0.0
if (i .le. nlons) then
wrk1(n) = data(i,j)
end if
enddo
enddo
call fft991(wrk1,wrk2,sphere_dat%trigs,sphere_dat%ifax,1,nlons+2,nlons,nlats,-1)
n = -1
do j=1,nlats
do m=1,(nlons/2)+1
n = n + 2
if (m .le. mwaves) then
coeff(m,j) = cmplx(wrk1(n),wrk1(n+1))
end if
enddo
enddo
!==> inverse transform.
else if (idir .eq. -1) then
wrk1 = 0.
n = -1
do j=1,nlats
do m=1,(nlons/2)+1
n = n + 2
if (m .le. mwaves) then
wrk1(n) = real(coeff(m,j))
wrk1(n+1) = aimag(coeff(m,j))
end if
enddo
enddo
call fft991(wrk1,wrk2,sphere_dat%trigs,sphere_dat%ifax,1,nlons+2,nlons,nlats,1)
n = 0
do j=1,nlats
do i=1,nlons+2
n = n + 1
if (i .le. nlons) then
data(i,j) = wrk1(n)
end if
enddo
enddo
else
write(6,*) ' idir must be +1 or -1 in RFFT!'
write(6,*) ' execution terminated.'
stop
end if
end subroutine rfft
subroutine spharm(sphere_dat, ugrid, anm, idir)
! spherical harmonic transform
type (sphere), intent(in) :: sphere_dat
real, dimension(sphere_dat%nlons,sphere_dat%nlats), intent(inout) :: ugrid
complex, dimension((sphere_dat%ntrunc+1)*(sphere_dat%ntrunc+2)/2), intent(inout) :: anm
complex, dimension(sphere_dat%ntrunc+1,sphere_dat%nlats) :: am
integer :: nlats,ntrunc,mwaves,nmstrt,nm,m,n,j
integer, intent(in) :: idir
if (.not. sphere_dat%isinitialized) then
print *, 'uninitialized sphere object in spharm!'
stop
end if
nlats = sphere_dat%nlats
ntrunc = sphere_dat%ntrunc
mwaves = ntrunc+1
IF (IDIR .eq. +1) then
!==> GRID SPACE TO SPECTRAL SPACE TRANSFORMATION
! FIRST, INITIALIZE ARRAY.
anm = 0.
!==> perform ffts on each latitude.
call rfft(sphere_dat, ugrid, am, 1)
!==> SUM OVER ALL GAUSSIAN LATITUDES FOR EACH MODE AND EACH WAVE TO
! OBTAIN THE TRANSFORMED VARIABLE IN SPECTRAL SPACE.
do j=1,nlats
NMSTRT = 0
DO m = 1, mwaves
DO n = 1, mwaves-m+1
NM = NMSTRT + n
anm(NM)=anm(NM)+sphere_dat%pnm(nm,j)*sphere_dat%weights(j)*am(m,j)
enddo
nmstrt = nmstrt + mwaves-m+1
enddo
enddo
!==> SPECTRAL SPACE TO GRID SPACE TRANSFORMATION.
else if (idir .eq. -1) then
DO J = 1, NLATS
!==> INVERSE LEGENDRE TRANSFORM TO GET VALUES OF THE ZONAL FOURIER
! TRANSFORM AT LATITUDE j.
!==> SUM THE VARIOUS MERIDIONAL MODES TO BUILD THE FOURIER SERIES
! COEFFICIENT FOR ZONAL WAVENUMBER M=I-1 AT the GIVEN LATITUDE.
NMSTRT = 0
do m = 1, MWAVES
am(m,j) = cmplx(0., 0.)
DO n = 1, mwaves-m+1
NM = NMSTRT + n
am(m,j) = am(m,j) + anm(NM) * sphere_dat%pnm(NM,j)
enddo
NMSTRT = NMSTRT + mwaves-m+1
enddo
ENDDO
!==> FOURIER TRANSFORM TO COMPUTE THE VALUES OF THE VARIABLE IN GRID
! SPACE at THE J-TH LATITUDE.
call rfft(sphere_dat, ugrid, am, -1)
else
print *, 'error in spharm: idir must be -1 or +1!'
print *, 'execution terminated in subroutine spharm'
end if
end subroutine spharm
subroutine getuv(sphere_dat,vrtnm,divnm,ug,vg)
! compute U,V (winds times cos(lat)) from vrtnm,divnm
! (spectral coeffs of vorticity and divergence).
type (sphere), intent(in) :: sphere_dat
real, dimension(sphere_dat%nlons,sphere_dat%nlats), intent(out) :: ug,vg
complex, dimension((sphere_dat%ntrunc+1)*(sphere_dat%ntrunc+2)/2), intent(in) :: vrtnm,divnm
complex, dimension(sphere_dat%ntrunc+1,sphere_dat%nlats) :: um,vm
integer :: nlats,ntrunc,mwaves,m,j,n,nm,nmstrt
real :: rm
if (.not. sphere_dat%isinitialized) then
print *, 'uninitialized sphere object in getuv!'
stop
end if
nlats = sphere_dat%nlats
ntrunc = sphere_dat%ntrunc
mwaves = ntrunc+1
do j=1,nlats
nmstrt = 0
do m=1,mwaves
rm = m-1
um(m,j) = cmplx(0.,0.)
vm(m,j) = cmplx(0.,0.)
do n=1,mwaves-m+1
nm = nmstrt + n
um(m,j) = um(m,j) + (sphere_dat%ilap(nm)/sphere_dat%rsphere)*( &
cmplx(0.,rm)*divnm(nm)*sphere_dat%pnm(nm,j) + &
vrtnm(nm)*sphere_dat%hnm(nm,j) )
vm(m,j) = vm(m,j) + (sphere_dat%ilap(nm)/sphere_dat%rsphere)*( &
cmplx(0.,rm)*vrtnm(nm)*sphere_dat%pnm(nm,j) - &
divnm(nm)*sphere_dat%hnm(nm,j) )
enddo
nmstrt = nmstrt + mwaves-m+1
enddo
enddo
call rfft(sphere_dat, ug, um, -1)
call rfft(sphere_dat, vg, vm, -1)
end subroutine getuv
subroutine getvrtdiv(sphere_dat,vrtnm,divnm,ug,vg)
! compute vrtnm,divnm (spectral coeffs of vorticity and
! divergence) from U,V (winds time cos(lat)).
type (sphere), intent(in) :: sphere_dat
real, dimension(sphere_dat%nlons,sphere_dat%nlats), intent(inout) :: ug,vg
complex, dimension((sphere_dat%ntrunc+1)*(sphere_dat%ntrunc+2)/2), intent(in) :: vrtnm,divnm
complex, dimension(sphere_dat%ntrunc+1,sphere_dat%nlats) :: um,vm
if (.not. sphere_dat%isinitialized) then
print *, 'uninitialized sphere object in getvrtdiv!'
stop
end if
call rfft(sphere_dat, ug, um, 1)
call rfft(sphere_dat, vg, vm, 1)
call sumnm(sphere_dat,um,vm,divnm,1,1)
call sumnm(sphere_dat,vm,um,vrtnm,1,-1)
end subroutine getvrtdiv
subroutine sumnm(sphere_dat,am,bm,anm,isign1,isign2)
!
! given the arrays of fourier coeffs, am and bm,
! compute the complex spherical harmonic coeffs of:
!
! isign1*( (1./re*(1.-x**2))*d(ag)/d(lon) + (isign2/re)*d(bg)/dx )
!
! where x = sin(lat), isign1 and isign2 are either +1 or -1,
! ag and bg are the physical space values of am,bm and re
! is the radius of the sphere.
!
! the result is returned in anm.
!
! for example on how to use this routine, see subroutine getvrtdiv.
!
!
type (sphere), intent(in) :: sphere_dat
integer, intent(in) :: isign1,isign2
complex, dimension((sphere_dat%ntrunc+1)*(sphere_dat%ntrunc+2)/2) :: anm
complex, dimension(sphere_dat%ntrunc+1,sphere_dat%nlats), intent(in) :: am,bm
integer :: nlats,ntrunc,mwaves,j,m,n,nm,nmstrt
real :: sign1,sign2,rm
if (.not. sphere_dat%isinitialized) then
print *, 'uninitialized sphere object in sumnm!'
stop
end if
sign1 = float(isign1)
sign2 = float(isign2)
if (isign2 .ne. 1 .and. isign2 .ne. -1) then
print *, ' isign2 must be +1 or -1 in sumnm!'
print *, ' execution terminated in sumnm'
end if
if (isign1 .ne. 1 .and. isign1 .ne. -1) then
print *, ' isign1 must be +1 or -1 in sumnm!'
print *, ' execution terminated in sumnm'
stop
end if
nlats = sphere_dat%nlats
ntrunc = sphere_dat%ntrunc
mwaves = ntrunc+1
anm = 0.
DO J=1,NLATS
NMSTRT = 0
DO m = 1, MWAVES
RM = M-1
DO n = 1, mwaves-m+1
NM = NMSTRT + n
aNM(NM) = aNM(NM) + sign1*sphere_dat%GWrC(j)*(CMPLX(0.,RM) &
* sphere_dat%PNM(NM,J) * am(m,j) &
+ sign2 * sphere_dat%HNM(NM,J) * bm(m,j))
enddo
nmstrt = nmstrt + mwaves - m +1
enddo
enddo
end subroutine sumnm
subroutine cosgrad(sphere_dat,divnm,ug,vg)
! compute coslat * gradient of spectral coefficients (divnm)
! vector gradient returned on grid as (ug,vg)
type (sphere), intent(in) :: sphere_dat
real, dimension(sphere_dat%nlons,sphere_dat%nlats), intent(out) :: ug,vg
complex, dimension((sphere_dat%ntrunc+1)*(sphere_dat%ntrunc+2)/2), intent(in) :: divnm
complex, dimension(sphere_dat%ntrunc+1,sphere_dat%nlats) :: um,vm
integer :: nlats,ntrunc,mwaves,j,m,n,nm,nmstrt
real :: rm
if (.not. sphere_dat%isinitialized) then
print *, 'uninitialized sphere object in cosgrad!'
stop
end if
nlats = sphere_dat%nlats
ntrunc = sphere_dat%ntrunc
mwaves = ntrunc+1
do j=1,nlats
nmstrt = 0
do m=1,mwaves
rm = (m-1)
um(m,j) = cmplx(0.,0.)
vm(m,j) = cmplx(0.,0.)
do n=1,mwaves-m+1
nm = nmstrt + n
um(m,j) = um(m,j) + (1./sphere_dat%rsphere)* &
cmplx(0.,rm)*divnm(nm)*sphere_dat%pnm(nm,j)
vm(m,j) = vm(m,j) - (1./sphere_dat%rsphere)* &
divnm(nm)*sphere_dat%hnm(nm,j)
enddo
nmstrt = nmstrt + mwaves - m +1
enddo
enddo
call rfft(sphere_dat, ug, um, -1)
call rfft(sphere_dat, vg, vm, -1)
end subroutine cosgrad
subroutine specsmooth(sphere_dat,datagrid,smooth)
! isotropic spectral smoothing of datagrid.
! input: smooth(sphere_dat%ntrunc+1) - smoothing factor as a
! function of degree (sphere_dat%indxn).
type (sphere), intent(in) :: sphere_dat
real, dimension(sphere_dat%nlons,sphere_dat%nlats), intent(inout) :: datagrid
real, dimension(sphere_dat%ntrunc+1), intent(in) :: smooth
complex, dimension((sphere_dat%ntrunc+1)*(sphere_dat%ntrunc+2)/2) :: dataspec
integer :: n,nm,nmdim
if (.not. sphere_dat%isinitialized) then
print *, 'uninitialized sphere object in specsmooth!'
stop
end if
nmdim = (sphere_dat%ntrunc+1)*(sphere_dat%ntrunc+2)/2
call spharm(sphere_dat, datagrid, dataspec, 1)
do nm=1,nmdim
n = sphere_dat%indxn(nm)
dataspec(nm) = dataspec(nm)*smooth(n+1)
enddo
call spharm(sphere_dat, datagrid, dataspec, -1)
end subroutine specsmooth
subroutine fft99 (a,work,trigs,ifax,inc,jump,n,lot,isign)
! purpose performs multiple fast fourier transforms. this package
! will perform a number of simultaneous real/half-complex
! periodic fourier transforms or corresponding inverse
! transforms, i.e. given a set of real data vectors, the
! package returns a set of 'half-complex' fourier
! coefficient vectors, or vice versa. the length of the
! transforms must be an even number greater than 4 that has
! no other factors except possibly powers of 2, 3, and 5.
! this is an all-fortran version of a optimized routine
! fft99 written for xmp/ymps by dr. clive temperton of
! ecmwf.
!
! the package fft99f contains several user-level routines:
!
! subroutine set99
! an initialization routine that must be called once
! before a sequence of calls to the fft routines
! (provided that n is not changed).
!
! subroutines fft99 and fft991
! two fft routines that return slightly different
! arrangements of the data in gridpoint space.
!
! usage let n be of the form 2**p * 3**q * 5**r, where p .ge. 1,
! q .ge. 0, and r .ge. 0. then a typical sequence of
! calls to transform a given set of real vectors of length
! n to a set of 'half-complex' fourier coefficient vectors
! of length n is
!
! dimension ifax(13),trigs(3*n/2+1),a(m*(n+2)),
! + work(m*(n+1))
!
! call set99 (trigs, ifax, n)
! call fft99 (a,work,trigs,ifax,inc,jump,n,m,isign)
!
! see the individual write-ups for set99, fft99, and
! fft991 below, for a detailed description of the
! arguments.
!
! history the package was written by clive temperton at ecmwf in
! november, 1978. it was modified, documented, and tested
! for ncar by russ rew in september, 1980.
!
!-----------------------------------------------------------------------
!
! subroutine set99 (trigs, ifax, n)
!
! purpose a set-up routine for fft99 and fft991. it need only be
! called once before a sequence of calls to the fft
! routines (provided that n is not changed).
!
! argument ifax(13),trigs(3*n/2+1)
! dimensions
!
! arguments
!
! on input trigs
! a floating point array of dimension 3*n/2 if n/2 is
! even, or 3*n/2+1 if n/2 is odd.
!
! ifax
! an integer array. the number of elements actually used
! will depend on the factorization of n. dimensioning
! ifax for 13 suffices for all n less than a million.
!
! n
! an even number greater than 4 that has no prime factor
! greater than 5. n is the length of the transforms (see
! the documentation for fft99 and fft991 for the
! definitions of the transforms).
!
! on output ifax
! contains the factorization of n/2. ifax(1) is the
! number of factors, and the factors themselves are stored
! in ifax(2),ifax(3),... if set99 is called with n odd,
! or if n has any prime factors greater than 5, ifax(1)
! is set to -99.
!
! trigs
! an array of trigonometric function values subsequently
! used by the fft routines.
!
!-----------------------------------------------------------------------
!
! subroutine fft991 (a,work,trigs,ifax,inc,jump,n,m,isign)
! and
! subroutine fft99 (a,work,trigs,ifax,inc,jump,n,m,isign)
!
! purpose perform a number of simultaneous real/half-complex
! periodic fourier transforms or corresponding inverse
! transforms, using ordinary spatial order of gridpoint
! values (fft991) or explicit cyclic continuity in the
! gridpoint values (fft99). given a set
! of real data vectors, the package returns a set of
! 'half-complex' fourier coefficient vectors, or vice
! versa. the length of the transforms must be an even
! number that has no other factors except possibly powers
! of 2, 3, and 5. this is an all-fortran version of
! optimized routine fft991 written for xmp/ymps by
! dr. clive temperton of ecmwf.
!
! argument a(m*(n+2)), work(m*(n+1)), trigs(3*n/2+1), ifax(13)
! dimensions
!
! arguments
!
! on input a
! an array of length m*(n+2) containing the input data
! or coefficient vectors. this array is overwritten by
! the results.
!
! work
! a work array of dimension m*(n+1)
!
! trigs
! an array set up by set99, which must be called first.
!
! ifax
! an array set up by set99, which must be called first.
!
! inc
! the increment (in words) between successive elements of
! each data or coefficient vector (e.g. inc=1 for
! consecutively stored data).
!
! jump
! the increment (in words) between the first elements of
! successive data or coefficient vectors. on crays,
! try to arrange data so that jump is not a multiple of 8
! (to avoid memory bank conflicts). for clarification of
! inc and jump, see the examples below.
!
! n
! the length of each transform (see definition of
! transforms, below).
!
! m
! the number of transforms to be done simultaneously.
!
! isign
! = +1 for a transform from fourier coefficients to
! gridpoint values.
! = -1 for a transform from gridpoint values to fourier
! coefficients.
!
! on output a
! if isign = +1, and m coefficient vectors are supplied
! each containing the sequence:
!
! a(0),b(0),a(1),b(1),...,a(n/2),b(n/2) (n+2 values)
!
! then the result consists of m data vectors each
! containing the corresponding n+2 gridpoint values:
!
! for fft991, x(0), x(1), x(2),...,x(n-1),0,0.
! for fft99, x(n-1),x(0),x(1),x(2),...,x(n-1),x(0).
! (explicit cyclic continuity)
!
! when isign = +1, the transform is defined by:
! x(j)=sum(k=0,...,n-1)(c(k)*exp(2*i*j*k*pi/n))
! where c(k)=a(k)+i*b(k) and c(n-k)=a(k)-i*b(k)
! and i=sqrt (-1)
!
! if isign = -1, and m data vectors are supplied each
! containing a sequence of gridpoint values x(j) as