/[MITgcm]/MITgcm/pkg/sphere/sphere.F
ViewVC logotype

Annotation of /MITgcm/pkg/sphere/sphere.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.1 - (hide annotations) (download)
Thu Feb 17 18:55:10 2005 UTC (19 years, 3 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint57g_post, checkpoint57g_pre, checkpoint57f_pre, checkpoint57f_post
Adding sphere package to do spherical harmonic decomposition

1 heimbach 1.1
2     #include "CPP_OPTIONS.h"
3    
4     c ==================================================================
5     c
6     c shpere.F: Routines that handle the projection onto sherical
7     c harmonics.
8     c
9     c Routines:
10     c
11     c o adfsc4dat - Adjoint routine of fsc4dat.
12     c o shc2grid - Evaluate a spherical harmonics model on a
13     c regular grid.
14     c o shc4grid - Evaluate a spherical harmonics model on a
15     c regular grid.
16     c
17     c o shcrotate - s/r used for the sph analysis ...
18     c o shc2zone - s/r used for the sph analysis ...
19     c o shc4zone - s/r used for the sph analysis ...
20     c o fsc2dat - s/r used for the sph analysis ...
21     c o frsbase - s/r used for the sph analysis ...
22     c o helmholtz - s/r used for the sph analysis ...
23     c o recur_z - s/r used for the sph analysis ...
24     c
25     c o shcError - Print error message and stop program (see below).
26     c
27     c IMSL Routines:
28     c
29     c o fftrb - Compute the real periodic sequence from its Fourier
30     c coefficients.
31     c --> replace by NAGLIB routine. C06...
32     c
33     c Platform-specific:
34     c
35     c M abort - FORTRAN intrinsic on SUN and, presumably, on CRAY to
36     c exit the program if an error condition is met.
37     c
38     c --> Replaced ABORT by shcError ( Christian Eckert MIT 22-Mar-2000 )
39     c
40     c Note:
41     c =====
42     c
43     c Where code is written in lower case letters, changes were
44     c introduced. The changes are mainly related to F90 syntax.
45     c Additionally, I replaced the call to the intrinsic *AMOD*
46     c by its generic name *MOD*. ( Ch.E. 25-May-1999 )
47     c
48     c
49     c Documentation:
50     c
51     c
52     c ==================================================================
53    
54     SUBROUTINE FSC4DAT(N,FSC)
55    
56     IMPLICIT NONE
57    
58     INTEGER N
59     REAL FSC(N)
60     REAL SCALE
61     integer i
62    
63     #ifdef USE_SPH_PROJECTION
64     CALL FFTRF(N,FSC,FSC)
65     #else
66     do i=1,n
67     fsc(i) = 0.0
68     enddo
69     #endif
70    
71     SCALE = 2.0/FLOAT(N)
72     FSC(1) = FSC(1)/FLOAT(N)
73    
74     CE FSC(2:N-1:2) = FSC(2:N-1:2)*SCALE
75     CE FSC(3:N :2) =-FSC(3:N: 2)*SCALE ! change sign of SINE coeffs.
76    
77     do i = 2,n-1,2
78     fsc( i ) = fsc( i )*scale
79     enddo
80     do i = 3,n,2
81     fsc( i ) = -fsc( i )*scale ! change sign of sine coeffs.
82     enddo
83    
84     IF (MOD(N,2).EQ.0) THEN
85     FSC(N) = FSC(N)/FLOAT(N)
86     ENDIF
87    
88     RETURN
89     END
90    
91     c ==================================================================
92    
93     subroutine adfsc4dat(n,adfsc)
94    
95     implicit none
96    
97     integer n
98     real adfsc(n)
99    
100     integer i
101    
102     #ifdef USE_SPH_PROJECTION
103     call fftrb(n,adfsc,adfsc)
104     #else
105     do i=1,n
106     adfsc(i) = 0.0
107     enddo
108     #endif
109    
110     do i=1,n
111     adfsc(i) = adfsc(i)/float(n)
112     enddo
113    
114     end
115    
116     c ==================================================================
117    
118     SUBROUTINE SHC2GRID( LMAX, SHC, NLON, NLAT, GRID, P )
119     C
120     C --- Evaluate a spherical harmonics model on a regular grid.
121     C
122     C SHC ! spherical harmonic coefficients in the order of
123     C !
124     C ! C00 : SHC(1)
125     C ! S11 C10 C11 : SHC(2) SHC(3) SHC(4)
126     C ! S22 S21 C20 C21 C22 : SHC(5) SHC(6) SHC(7) SHC(8) SHC(9)
127     C !
128     C ! i.e., C(L,M) = SHC(L*L+L+1+M)
129     C ! S(L,M) = SHC(L*L+L+1-M)
130     C NLAT ! number of latitude (zonal) lines.
131     C NLON ! number of longitude (meridional) lines.
132     C !
133     C GRID ! grid values in the order of
134     C ! ((west to east),south to north)
135     C ! ((LON=1,NLON),LAT=1,NLAT)
136     C ! grid values are defined on the "centers"
137     C ! of geographically regular equi-angular blocks.
138     C
139     C --- Coded by Dr. Myung-Chan Kim, Center for Space Research, 1997
140     C University of Texas at Austin, Austin, Texas, 78712
141     C
142     IMPLICIT NONE
143    
144     INTEGER LMAX ! INPUT max. degree
145     REAL SHC((1+LMAX)*(1+LMAX)) ! INPUT spherical harmonics
146     INTEGER NLAT ! INPUT number of latitude lines.
147     INTEGER NLON ! INPUT number of longitude lines.
148     REAL GRID(NLON,NLAT) ! OUTPUT grid values
149     REAL P((LMAX+1)*(LMAX+2)/2) ! work space
150    
151     INTEGER NLONLMT
152     PARAMETER (NLONLMT=10000)
153     REAL HS (NLONLMT)
154     REAL HN (NLONLMT)
155    
156     REAL DLAT, ANGLE, XLAT1, XLAT2
157     INTEGER LATS, LATN
158    
159     ce integer js, jn
160    
161     integer i
162    
163     ce IF(NLON.GT.NLONLMT) CALL ABORT('NLON.GT.NLONLMT')
164     ce IF(NLON.LT.LMAX*2 ) CALL ABORT('NLON.LT.LMAX*2 ')
165     IF(NLON.GT.NLONLMT) CALL shcError('NLON.GT.NLONLMT',1)
166     IF(NLON.LT.LMAX*2 ) CALL shcError('NLON.LT.LMAX*2 ',1)
167    
168     DLAT = 180.0/FLOAT(NLAT)
169     ANGLE = 180.0/FLOAT(NLON)
170    
171     CALL SHCROTATE(LMAX,SHC,ANGLE)
172     DO LATS = 1,(NLAT+1)/2
173     LATN = NLAT-(LATS-1)
174     XLAT2 =-90.0 + FLOAT(LATS)*DLAT
175     XLAT1 = XLAT2-DLAT
176     do i=1,nlon
177     hs(i) = 0.0
178     hn(i) = 0.0
179     enddo
180     CALL SHC2ZONE( LMAX, SHC, XLAT1, XLAT2, HS, HN, P )
181     CALL FSC2DAT( NLON, HS )
182     CALL FSC2DAT( NLON, HN )
183     do i=1,nlon
184     grid(i,lats) = hs(i)
185     grid(i,latn) = hn(i)
186     enddo
187     ENDDO
188     CALL SHCROTATE(LMAX,SHC,-ANGLE)
189     RETURN
190     END
191    
192     c ==================================================================
193    
194     SUBROUTINE SHC4GRID( LMAX, SHC, NLON, NLAT, GRID, P )
195    
196     IMPLICIT NONE
197    
198     INTEGER LMAX ! INPUT max. degree
199     REAL SHC((1+LMAX)*(1+LMAX)) ! OUTPUT spherical harmonics
200     INTEGER NLAT ! INPUT number of latitude lines.
201     INTEGER NLON ! INPUT number of longitude lines.
202     REAL GRID(NLON,NLAT) ! INPUT grid values
203     REAL P((LMAX+1)*(LMAX+2)/2) ! work space
204    
205     INTEGER NLONLMT
206     PARAMETER (NLONLMT=10000)
207     REAL HS (NLONLMT)
208     REAL HN (NLONLMT)
209    
210     REAL DLAT, ANGLE, XLAT1, XLAT2
211     INTEGER LATS, LATN
212    
213     ce integer js, jn
214    
215     integer i
216    
217     IF(NLON.GT.NLONLMT) THEN
218     PRINT *, 'NLON = ', NLON
219     PRINT *, 'NLONLMT = ', NLONLMT
220     ce CALL ABORT('NLON.GT.NLONLMT')
221     CALL shcError('NLON.GT.NLONLMT',1)
222     END IF
223    
224     IF(NLON.LT.LMAX*2 ) THEN
225     PRINT *, 'NLON = ', NLON
226     PRINT *, 'LMAX = ', LMAX
227     ce CALL ABORT('NLON.LT.LMAX*2')
228     CALL shcError('NLON.LT.LMAX*2',1)
229     END IF
230    
231     DLAT = 180.0/FLOAT(NLAT)
232     ANGLE = 180.0/FLOAT(NLON)
233     DO LATS = 1,(NLAT+1)/2
234     LATN = NLAT-(LATS-1)
235     do i = 1,nlon
236     hs(i) = grid(i,lats)
237     hn(i) = grid(i,latn)
238     enddo
239     CALL FSC4DAT(NLON,HS)
240     CALL FSC4DAT(NLON,HN)
241     XLAT2 =-90.0 + FLOAT(LATS)*DLAT
242     XLAT1 = XLAT2-DLAT
243     CALL SHC4ZONE(LMAX,SHC,XLAT1,XLAT2,HS,HN,P)
244     ENDDO
245     CALL SHCROTATE(LMAX,SHC,-ANGLE)
246    
247     RETURN
248     END
249    
250     c ==================================================================
251    
252     SUBROUTINE SHCROTATE(LMAX,SHC,ANGLE)
253    
254     IMPLICIT NONE
255    
256     INTEGER LMAX ! max. degree of spherical harmonics.
257     REAL SHC((1+LMAX)*(1+LMAX)) ! spherical harmonic coeffs.
258     REAL ANGLE ! in degree.
259    
260     INTEGER LMAXLMT
261     PARAMETER (LMAXLMT=10000)
262    
263     REAL H(LMAXLMT*2+1)
264     INTEGER K, L, M
265     REAL SINA, COSA, C, S
266    
267     if (mod(angle,360.0) .ne. 0.0) then
268    
269     CALL FRSBASE(ANGLE,H,1,1+LMAX+LMAX)
270    
271     DO M = 0,LMAX
272     IF(M.EQ.0) THEN
273     SINA = 0.0
274     COSA = 1.0
275     ELSE
276     COSA = H(M+M)
277     SINA = H(M+M+1)
278     ENDIF
279     DO L = M,LMAX
280     K = L*L+L+1
281     C = SHC(K+M)
282     S = SHC(K-M)
283     SHC(K+M) = COSA*C + SINA*S
284     SHC(K-M) =-SINA*C + COSA*S
285     ENDDO
286     ENDDO
287    
288     ENDIF
289    
290     RETURN
291     END
292    
293     c ==================================================================
294    
295     SUBROUTINE SHC2ZONE(LMAX,SHC,XLAT1,XLAT2,HS,HN,P)
296    
297     IMPLICIT NONE
298    
299     INTEGER LMAX ! INPUT max. degree of spher. harmonics.
300     REAL SHC((1+LMAX)*(1+LMAX)) ! INPUT spherical harmonic coeffs.
301     REAL XLAT1 ! INPUT latitude.
302     REAL XLAT2 ! INPUT latitude.
303     REAL HS(1+LMAX+LMAX) ! OUTPUT
304     REAL HN(1+LMAX+LMAX) ! OUTPUT
305     REAL P((LMAX+1)*(LMAX+2)/2) ! INPUT
306    
307     INTEGER LMAXLMT
308     PARAMETER (LMAXLMT=5000)
309    
310     REAL FACT(0:LMAXLMT) !
311     INTEGER LMAX1, J, K, L, M
312     REAL DEG2RAD, SLAT
313     REAL A, B, E, F, Q, DA, DB, XLAT
314    
315     ce IF (LMAX.GT.LMAXLMT) CALL ABORT('LMAX.GT.LMAXLMT')
316     IF (LMAX.GT.LMAXLMT) CALL shcError('LMAX.GT.LMAXLMT',1)
317    
318     DEG2RAD = ACOS(-1.0)/180.0
319     XLAT = 0.5*(XLAT1+XLAT2)
320     do k = 1,lmax
321     fact(k) = 1.0
322     enddo
323     SLAT = SIN(XLAT*DEG2RAD)
324    
325     CALL HELMHOLTZ(LMAX,SLAT,P)
326     LMAX1 = LMAX+1
327     K = 0
328     DO M = 0,LMAX
329     A = 0.0
330     B = 0.0
331     E = 0.0
332     F = 0.0
333     DO L = M,LMAX-1,2
334     K = K+1
335     J = L*L+L+1
336     Q = FACT(L)*P(K)
337     DA = Q*SHC(J+M)
338     DB = Q*SHC(J-M)
339     A = A+DA
340     B = B+DB
341     E = E+DA
342     F = F+DB
343     K = K+1
344     J = J+L+L+2
345     Q = FACT(L+1)*P(K)
346     DA = Q*SHC(J+M)
347     DB = Q*SHC(J-M)
348     A = A+DA
349     B = B+DB
350     E = E-DA
351     F = F-DB
352     ENDDO
353     IF(MOD(LMAX-L,2).EQ.0) THEN
354     K = K+1
355     J = LMAX*LMAX+LMAX+1
356     Q = FACT(LMAX)*P(K)
357     DA = Q*SHC(J+M)
358     DB = Q*SHC(J-M)
359     A = A+DA
360     B = B+DB
361     E = E+DA
362     F = F+DB
363     ENDIF
364     IF(M.EQ.0) THEN
365     HS(1) = A
366     HN(1) = E
367     ELSE
368     HS(M+M ) = A
369     HS(M+M+1) = B
370     HN(M+M ) = E
371     HN(M+M+1) = F
372     ENDIF
373     ENDDO
374    
375     RETURN
376     END
377    
378     c ==================================================================
379    
380     SUBROUTINE SHC4ZONE(LMAX,SHC,XLAT1,XLAT2,HS,HN,P)
381    
382     IMPLICIT NONE
383    
384     INTEGER LMAX
385     REAL SHC((1+LMAX)*(1+LMAX)) ! spherical harmonic coeffs.
386     REAL XLAT1 ! latitude.
387     REAL XLAT2 ! latitude.
388     REAL P((LMAX+1)*(LMAX+2)/2)
389     C
390     C - output -
391     C
392     REAL HS(1+LMAX+LMAX)
393     REAL HN(1+LMAX+LMAX)
394     C
395     C - work space -
396     C
397     INTEGER LMAXLMT
398     PARAMETER (LMAXLMT=5000)
399    
400     REAL XLAT, SIN0, SIN1, SIN2, SCALE
401     REAL DEG2RAD, CE, CO, SE, SO
402     INTEGER I, J, K, L, M
403     INTEGER JP, I1, I2
404    
405     ce IF(LMAX.GT.LMAXLMT) CALL ABORT('LMAX.GT.LMAXLMT')
406     IF(LMAX.GT.LMAXLMT) CALL shcError('LMAX.GT.LMAXLMT',1)
407    
408     IF(XLAT1 .LT. -90.0+1.E-10) THEN
409     do i = 1,(1+lmax)*(1+lmax)
410     shc(i) = 0.0
411     enddo
412     ENDIF
413    
414     XLAT = 0.5*(XLAT1+XLAT2)
415     DEG2RAD = ACOS(-1.0)/180.0
416     SIN0 = SIN(XLAT *DEG2RAD)
417     SIN1 = SIN(AMAX1(-90.0,XLAT1)*DEG2RAD)
418     SIN2 = SIN(AMIN1( 0.0,XLAT2)*DEG2RAD)
419     SCALE = (SIN2 - SIN1)*0.25
420     do i = 1,lmax+lmax+1
421     hs(i) = hs(i)*scale
422     hn(i) = hn(i)*scale
423     enddo
424    
425     CALL HELMHOLTZ(LMAX,SIN0,P)
426     I = 0
427     cadj loop = parallel
428     DO M = 0,LMAX
429     IF (M .EQ. 0) THEN
430     J = 1
431     JP = 1
432     ELSE
433     J = 2*M
434     JP = J+1
435     ENDIF
436     CE = HS(J)+HN(J)
437     CO = HS(J)-HN(J)
438     SE = HS(J)+HN(JP)
439     SO = HS(J)-HN(JP)
440     I1 = I+1
441     I2 = I+2
442     cadj loop = parallel
443     DO L = M,LMAX-1,2
444     K = L*L+L+1
445     SHC(K+M) = SHC(K+M) + P(I1) * CE
446     SHC(K-M) = SHC(K-M) + P(I1) * SE
447     K = K+L+L+2
448     SHC(K+M) = SHC(K+M) + P(I2) * CO
449     SHC(K-M) = SHC(K-M) + P(I2) * SO
450     ENDDO
451     I = I + 2
452     IF(MOD(LMAX-M,2).NE.1) THEN
453     I = I+1
454     K = LMAX*LMAX+LMAX+1
455     SHC(K+M) = SHC(K+M) + P(I) * CE
456     SHC(K-M) = SHC(K-M) + P(I) * SE
457     ENDIF
458     ENDDO
459     RETURN
460     END
461    
462     c ==================================================================
463    
464     subroutine fsc2dat(n,fsc)
465     c
466     c --- this routine are coded to clarify the imsl's fftrf/fftrb.
467     c
468     implicit none
469    
470     integer n
471     real fsc(n)
472    
473     integer i
474    
475     do i = 2,n-1,2
476     fsc(i ) = fsc(i )*0.5
477     fsc(i+1) = -fsc(i+1)*0.5 ! change sign of sine coeffs.
478     enddo
479    
480     #ifdef USE_SPH_PROJECTION
481     call fftrb(n,fsc,fsc)
482     #else
483     do i = 1,n
484     fsc( i ) = 0.0
485     enddo
486     #endif
487    
488     return
489     end
490    
491     c ==================================================================
492    
493     SUBROUTINE FRSBASE(A,H,I,J)
494     C
495     C --- Given A (in degree), return H(I:J) = 1,COS(A),SIN(A).....
496     C
497     IMPLICIT NONE
498    
499     REAL A
500     REAL H(1)
501     INTEGER I, J
502    
503     REAL DEG2RAD, T, ARG
504     INTEGER N, L
505    
506     DEG2RAD = ACOS(-1.0)/180.0
507    
508     N = J-I+1
509     IF(N.LE.0) RETURN
510     H(I) = 1.0
511     IF(N.EQ.1) RETURN
512     ARG = A * DEG2RAD
513     H(I+1) = COS(ARG)
514     IF(N.EQ.2) RETURN
515     H(I+2) = SIN(ARG)
516     IF(N.EQ.3) RETURN
517     T = H(I+1)+H(I+1)
518     H(I+3) = T*H(I+1) - 1.0
519     IF(N.EQ.4) RETURN
520     H(I+4) = T*H(I+2)
521     IF(N.EQ.5) RETURN
522     DO L = I+5,J
523     H(L) = T*H(L-2)-H(L-4)
524     END DO
525    
526     RETURN
527     END
528    
529     c ==================================================================
530    
531     SUBROUTINE HELMHOLTZ(LMAX,S,P)
532     C
533     C --- compute the fully normalized associated legendre polynomials.
534     C
535     IMPLICIT NONE
536    
537     INTEGER LMAX ! INPUT max. degree.
538     REAL S ! INPUT sin(latitude).
539     REAL P(1) ! OUTPUT assoc. legendre polynomials.
540    
541     INTEGER LMAXLMT
542     PARAMETER(LMAXLMT=10800)
543    
544     REAL A(0:LMAXLMT)
545     REAL ISECT(0:LMAXLMT)
546     REAL X(LMAXLMT)
547     REAL Y(LMAXLMT+LMAXLMT)
548     REAL Z(LMAXLMT)
549    
550     INTEGER LMAXOLD
551     DATA LMAXOLD/-1/
552     SAVE LMAXOLD, ISECT, X, Y, Z
553    
554     REAL C, CM, AK
555     INTEGER K, L, N, M
556    
557     IF(LMAX.NE.LMAXOLD) THEN
558     ISECT(0) = 1
559     DO L = 1,LMAX
560     X(L) = SNGL(1.0D0/DSQRT(DBLE(FLOAT(4*L*L-1))))
561     Y(L) = SNGL(DSQRT(DBLE(FLOAT(L))))
562     Y(L+LMAX) = SNGL(DSQRT(DBLE(FLOAT(L+LMAX))))
563     ISECT(L) = L*LMAX-L*(L-3)/2+1
564     ENDDO
565     CALL RECUR_Z(Z,LMAX)
566     LMAXOLD = LMAX
567     ENDIF
568    
569     C = SNGL(DSQRT(1.0D0-DBLE(S)*DBLE(S)))
570     CM = 1.0
571     P(1) = 1.0
572     DO M = 1,LMAX
573     K = ISECT(M)
574     CM = C * CM
575     P(K) = Z(M) * CM
576     ENDDO
577     N = 1
578     DO M = 0,LMAX-N
579     K = ISECT(M)+N
580     L = N+M
581     AK = X(L)*Y(L+M)*Y(N)
582     P(K) = S*P(K-1)/AK
583     A(M) = AK
584     ENDDO
585     DO N = 2,LMAX
586     DO M = 0,LMAX-N
587     K = ISECT(M)+N
588     L = N+M
589     AK = X(L)*Y(L+M)*Y(N)
590     P(K) =(S*P(K-1)-P(K-2)*A(M))/AK
591     A(M) = AK
592     ENDDO
593     ENDDO
594     RETURN
595     END
596    
597     c ==================================================================
598    
599     SUBROUTINE RECUR_Z(Z,LMAX)
600     C
601     C --- Coefficients for fully normalized sectorial Legendre polynomials.
602     C
603     IMPLICIT NONE
604    
605     INTEGER LMAX
606     REAL Z(LMAX)
607    
608     DOUBLE PRECISION ZZ
609     INTEGER L
610    
611     ZZ = 2.0D0
612     DO L = 1,LMAX
613     ZZ = ZZ * DBLE(FLOAT(L+L+1))/DBLE(FLOAT(L+L))
614     Z(L) = SNGL(DSQRT(ZZ))
615     ENDDO
616    
617     RETURN
618     END
619    
620     c ==================================================================
621    
622     subroutine shcError( errstring, ierr )
623     c
624     c --- Print an error message and exit. Written in order to replace
625     c the machine specific routine ABORT.
626     c ( Christian Eckert MIT 22-Mar-2000 )
627     c
628     implicit none
629    
630     integer ierr
631     character*(*) errstring
632    
633     if ( ierr .ne. 0 ) then
634     print*
635     print*,' sphere: ',errstring
636     print*
637     stop ' ... program stopped.'
638     endif
639    
640     return
641     end
642    
643     c ==================================================================

  ViewVC Help
Powered by ViewVC 1.1.22