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

Contents of /MITgcm/pkg/sphere/sphere.F

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


Revision 1.3 - (show annotations) (download)
Tue Oct 9 00:11:39 2007 UTC (16 years, 6 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint64y, checkpoint64x, checkpoint64z, checkpoint64q, checkpoint64p, checkpoint64s, checkpoint64r, checkpoint64u, checkpoint64t, checkpoint64w, checkpoint64v, checkpoint64i, checkpoint64h, checkpoint64k, checkpoint64j, checkpoint64m, checkpoint64l, checkpoint64o, checkpoint64n, checkpoint64a, checkpoint64c, checkpoint64b, checkpoint64e, checkpoint64d, checkpoint64g, checkpoint64f, checkpoint63p, checkpoint63q, checkpoint63r, checkpoint63s, checkpoint63l, checkpoint63m, checkpoint63n, checkpoint63o, checkpoint63h, checkpoint63i, checkpoint63j, checkpoint63k, checkpoint63d, checkpoint63e, checkpoint63f, checkpoint63g, checkpoint63a, checkpoint63b, checkpoint63c, checkpoint64, checkpoint65, checkpoint60, checkpoint61, checkpoint62, checkpoint63, checkpoint66g, checkpoint66f, checkpoint66e, checkpoint66d, checkpoint66c, checkpoint66b, checkpoint66a, checkpoint66o, checkpoint66n, checkpoint66m, checkpoint66l, checkpoint66k, checkpoint66j, checkpoint66i, checkpoint66h, checkpoint65z, checkpoint65x, checkpoint65y, checkpoint65r, checkpoint65s, checkpoint65p, checkpoint65q, checkpoint65v, checkpoint65w, checkpoint65t, checkpoint65u, checkpoint65j, checkpoint65k, checkpoint65h, checkpoint65i, checkpoint65n, checkpoint65o, checkpoint65l, checkpoint65m, checkpoint65b, checkpoint65c, checkpoint65a, checkpoint65f, checkpoint65g, checkpoint65d, checkpoint65e, checkpoint59q, checkpoint59p, checkpoint59r, checkpoint59m, checkpoint59l, checkpoint59o, checkpoint59n, checkpoint59i, checkpoint59k, checkpoint59j, checkpoint62c, checkpoint62b, checkpoint62a, checkpoint62g, checkpoint62f, checkpoint62e, checkpoint62d, checkpoint62k, checkpoint62j, checkpoint62i, checkpoint62h, checkpoint62o, checkpoint62n, checkpoint62m, checkpoint62l, checkpoint62s, checkpoint62r, checkpoint62q, checkpoint62p, checkpoint62w, checkpoint62v, checkpoint62u, checkpoint62t, checkpoint62z, checkpoint62y, checkpoint62x, checkpoint61f, checkpoint61g, checkpoint61d, checkpoint61e, checkpoint61b, checkpoint61c, checkpoint61a, checkpoint61n, checkpoint61o, checkpoint61l, checkpoint61m, checkpoint61j, checkpoint61k, checkpoint61h, checkpoint61i, checkpoint61v, checkpoint61w, checkpoint61t, checkpoint61u, checkpoint61r, checkpoint61s, checkpoint61p, checkpoint61q, checkpoint61z, checkpoint61x, checkpoint61y, HEAD
Changes since 1.2: +4 -2 lines
add missing cvs $Header:$ or $Name:$

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

  ViewVC Help
Powered by ViewVC 1.1.22