NCEPLIBS-sp 2.4.0
sptrunmv.f
Go to the documentation of this file.
1C> @file
2C>
3C> Spectrally interpolate vectors to Mercator
4C> @author IREDELL @date 96-02-29
5
6C> THIS SUBPROGRAM SPECTRALLY TRUNCATES VECTOR FIELDS
7C> ON A GLOBAL CYLINDRICAL GRID, RETURNING THE FIELDS
8C> TO A MERCATOR GRID.
9C> THE WAVE-SPACE CAN BE EITHER TRIANGULAR OR RHOMBOIDAL.
10C> THE GRID-SPACE CAN BE EITHER AN EQUALLY-SPACED GRID
11C> (WITH OR WITHOUT POLE POINTS) OR A GAUSSIAN GRID.
12C> THE GRID FIELDS MAY HAVE GENERAL INDEXING.
13C> THE TRANSFORMS ARE ALL MULTIPROCESSED.
14C> TRANSFORM SEVERAL FIELDS AT A TIME TO IMPROVE VECTORIZATION.
15C> SUBPROGRAM CAN BE CALLED FROM A MULTIPROCESSING ENVIRONMENT.
16C>
17C> PROGRAM HISTORY LOG:
18C> 96-02-29 IREDELL
19C> 1998-12-15 IREDELL OPENMP DIRECTIVES INSERTED
20C>
21C> @param IROMB - INTEGER SPECTRAL DOMAIN SHAPE
22C> (0 FOR TRIANGULAR, 1 FOR RHOMBOIDAL)
23C> @param MAXWV - INTEGER SPECTRAL TRUNCATION
24C> @param IDRTI - INTEGER INPUT GRID IDENTIFIER
25C> (IDRTI=4 FOR GAUSSIAN GRID,
26C> IDRTI=0 FOR EQUALLY-SPACED GRID INCLUDING POLES,
27C> IDRTI=256 FOR EQUALLY-SPACED GRID EXCLUDING POLES)
28C> @param IMAXI - INTEGER EVEN NUMBER OF INPUT LONGITUDES.
29C> @param JMAXI - INTEGER NUMBER OF INPUT LATITUDES.
30C> @param KMAX - INTEGER NUMBER OF FIELDS TO TRANSFORM.
31C> @param MI - INTEGER NUMBER OF POINTS IN THE FASTER ZONAL DIRECTION
32C> @param MJ - INTEGER NUMBER OF POINTS IN THE SLOWER MERID DIRECTION
33C> @param IPRIME - INTEGER INPUT LONGITUDE INDEX FOR THE PRIME MERIDIAN.
34C> (DEFAULTS TO 1 IF IPRIME=0)
35C> (OUTPUT LONGITUDE INDEX FOR PRIME MERIDIAN ASSUMED 1.)
36C> @param ISKIPI - INTEGER SKIP NUMBER BETWEEN INPUT LONGITUDES
37C> (DEFAULTS TO 1 IF ISKIPI=0)
38C> @param JSKIPI - INTEGER SKIP NUMBER BETWEEN INPUT LATITUDES FROM SOUTH
39C> (DEFAULTS TO -IMAXI IF JSKIPI=0)
40C> @param KSKIPI - INTEGER SKIP NUMBER BETWEEN INPUT GRID FIELDS
41C> (DEFAULTS TO IMAXI*JMAXI IF KSKIPI=0)
42C> @param KGSKIP - INTEGER SKIP NUMBER BETWEEN GRID FIELDS
43C> (DEFAULTS TO MI*MJ IF KGSKIP=0)
44C> @param NISKIP - INTEGER SKIP NUMBER BETWEEN GRID I-POINTS
45C> (DEFAULTS TO 1 IF NISKIP=0)
46C> @param NJSKIP - INTEGER SKIP NUMBER BETWEEN GRID J-POINTS
47C> (DEFAULTS TO MI IF NJSKIP=0)
48C> @param JCPU - INTEGER NUMBER OF CPUS OVER WHICH TO MULTIPROCESS
49C> (DEFAULTS TO ENVIRONMENT NCPUS IF JCPU=0)
50C> @param RLAT1 - REAL LATITUDE OF THE FIRST GRID POINT IN DEGREES
51C> @param RLON1 - REAL LONGITUDE OF THE FIRST GRID POINT IN DEGREES
52C> @param DLAT - REAL LATITUDE INCREMENT IN DEGREES SUCH THAT
53C> D(PHI)/D(J)=DLAT*COS(PHI) WHERE J IS MERIDIONAL INDEX.
54C> DLAT IS NEGATIVE FOR GRIDS INDEXED SOUTHWARD.
55C> (IN TERMS OF GRID INCREMENT DY VALID AT LATITUDE RLATI,
56C> THE LATITUDE INCREMENT DLAT IS DETERMINED AS
57C> DLAT=DPR*DY/(RERTH*COS(RLATI/DPR))
58C> WHERE DPR=180/PI AND RERTH IS EARTH'S RADIUS)
59C> @param DLON - REAL LONGITUDE INCREMENT IN DEGREES SUCH THAT
60C> D(LAMBDA)/D(I)=DLON WHERE I IS ZONAL INDEX.
61C> DLON IS NEGATIVE FOR GRIDS INDEXED WESTWARD.
62C> @param GRIDUI - REAL (*) INPUT GRID U-WINDS
63C> @param GRIDVI - REAL (*) INPUT GRID V-WINDS
64C> @param LUV - LOGICAL FLAG WHETHER TO RETURN WINDS
65C> @param LDZ - LOGICAL FLAG WHETHER TO RETURN DIVERGENCE AND VORTICITY
66C> @param LPS - LOGICAL FLAG WHETHER TO RETURN POTENTIAL AND STREAMFCN
67C> @param UM - REAL (*) MERCATOR U-WINDS IF LUV
68C> @param VM - REAL (*) MERCATOR V-WINDS IF LUV
69C> @param DM - REAL (*) MERCATOR DIVERGENCES IF LDZ
70C> @param ZM - REAL (*) MERCATOR VORTICITIES IF LDZ
71C> @param PM - REAL (*) MERCATOR POTENTIALS IF LPS
72C> @param SM - REAL (*) MERCATOR STREAMFCNS IF LPS
73C>
74C> SUBPROGRAMS CALLED:
75C> - SPWGET GET WAVE-SPACE CONSTANTS
76C> - SPLAPLAC COMPUTE LAPLACIAN IN SPECTRAL SPACE
77C> - SPTRANV PERFORM A VECTOR SPHERICAL TRANSFORM
78C> - SPTGPM TRANSFORM SPECTRAL SCALAR TO MERCATOR
79C> - SPTGPMV TRANSFORM SPECTRAL VECTOR TO MERCATOR
80C> - NCPUS GETS ENVIRONMENT NUMBER OF CPUS
81C>
82C> REMARKS: MINIMUM GRID DIMENSIONS FOR UNALIASED TRANSFORMS TO SPECTRAL:
83C> DIMENSION |LINEAR |QUADRATIC
84C> ----------------------- |--------- |-------------
85C> IMAX |2*MAXWV+2 |3*MAXWV/2*2+2
86C> JMAX (IDRT=4,IROMB=0) |1*MAXWV+1 |3*MAXWV/2+1
87C> JMAX (IDRT=4,IROMB=1) |2*MAXWV+1 |5*MAXWV/2+1
88C> JMAX (IDRT=0,IROMB=0) |2*MAXWV+3 |3*MAXWV/2*2+3
89C> JMAX (IDRT=0,IROMB=1) |4*MAXWV+3 |5*MAXWV/2*2+3
90C> JMAX (IDRT=256,IROMB=0) |2*MAXWV+1 |3*MAXWV/2*2+1
91C> JMAX (IDRT=256,IROMB=1) |4*MAXWV+1 |5*MAXWV/2*2+1
92 SUBROUTINE sptrunmv(IROMB,MAXWV,IDRTI,IMAXI,JMAXI,KMAX,MI,MJ,
93 & IPRIME,ISKIPI,JSKIPI,KSKIPI,KGSKIP,
94 & NISKIP,NJSKIP,JCPU,RLAT1,RLON1,DLAT,DLON,
95 & GRIDUI,GRIDVI,LUV,UM,VM,LDZ,DM,ZM,LPS,PM,SM)
96
97 LOGICAL LUV,LDZ,LPS
98 REAL GRIDUI(*),GRIDVI(*)
99 REAL UM(*),VM(*),DM(*),ZM(*),PM(*),SM(*)
100 REAL W((MAXWV+1)*((IROMB+1)*MAXWV+2)/2*2+1,KMAX)
101 REAL EPS((MAXWV+1)*((IROMB+1)*MAXWV+2)/2),EPSTOP(MAXWV+1)
102 REAL ENN1((MAXWV+1)*((IROMB+1)*MAXWV+2)/2)
103 REAL ELONN1((MAXWV+1)*((IROMB+1)*MAXWV+2)/2)
104 REAL EON((MAXWV+1)*((IROMB+1)*MAXWV+2)/2),EONTOP(MAXWV+1)
105 REAL WD((MAXWV+1)*((IROMB+1)*MAXWV+2)/2*2+1,KMAX)
106 REAL WZ((MAXWV+1)*((IROMB+1)*MAXWV+2)/2*2+1,KMAX)
107C - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
108C TRANSFORM INPUT GRID TO WAVE
109 jc=jcpu
110 IF(jc.EQ.0) jc=ncpus()
111 mx=(maxwv+1)*((iromb+1)*maxwv+2)/2
112 mdim=2*mx+1
113 jn=-jskipi
114 IF(jn.EQ.0) jn=imaxi
115 js=-jn
116 inp=(jmaxi-1)*max(0,-jn)+1
117 isp=(jmaxi-1)*max(0,-js)+1
118 CALL sptranv(iromb,maxwv,idrti,imaxi,jmaxi,kmax,
119 & iprime,iskipi,jn,js,mdim,kskipi,0,0,jc,
120 & wd,wz,
121 & gridui(inp),gridui(isp),gridvi(inp),gridvi(isp),-1)
122C - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
123C TRANSFORM WAVE TO OUTPUT WINDS
124 IF(luv) THEN
125 CALL sptgpmv(iromb,maxwv,kmax,mi,mj,mdim,kgskip,niskip,njskip,
126 & rlat1,rlon1,dlat,dlon,wd,wz,um,vm)
127 ENDIF
128C - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
129C TRANSFORM WAVE TO OUTPUT DIVERGENCE AND VORTICITY
130 IF(ldz) THEN
131 CALL sptgpm(iromb,maxwv,kmax,mi,mj,mdim,kgskip,niskip,njskip,
132 & rlat1,rlon1,dlat,dlon,wd,dm)
133 CALL sptgpm(iromb,maxwv,kmax,mi,mj,mdim,kgskip,niskip,njskip,
134 & rlat1,rlon1,dlat,dlon,wz,zm)
135 ENDIF
136C - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
137C TRANSFORM WAVE TO OUTPUT POTENTIAL AND STREAMFUNCTION
138 IF(lps) THEN
139 CALL spwget(iromb,maxwv,eps,epstop,enn1,elonn1,eon,eontop)
140C$OMP PARALLEL DO
141 DO k=1,kmax
142 CALL splaplac(iromb,maxwv,enn1,wd(1,k),wd(1,k),-1)
143 CALL splaplac(iromb,maxwv,enn1,wz(1,k),wz(1,k),-1)
144 wd(1:2,k)=0.
145 wz(1:2,k)=0.
146 ENDDO
147 CALL sptgpm(iromb,maxwv,kmax,mi,mj,mdim,kgskip,niskip,njskip,
148 & rlat1,rlon1,dlat,dlon,wd,pm)
149 CALL sptgpm(iromb,maxwv,kmax,mi,mj,mdim,kgskip,niskip,njskip,
150 & rlat1,rlon1,dlat,dlon,wz,sm)
151 ENDIF
152C - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
153 END
function ncpus()
Set number of CPUs - the number of processors over which to parallelize.
Definition: ncpus.F:24
subroutine splaplac(I, M, ENN1, Q, QD2, IDIR)
Computes the laplacian or the inverse laplacian of a scalar field in spectral space.
Definition: splaplac.f:25
subroutine sptgpm(IROMB, MAXWV, KMAX, MI, MJ, KWSKIP, KGSKIP, NISKIP, NJSKIP, RLAT1, RLON1, DLAT, DLON, WAVE, GM)
This subprogram performs a spherical transform from spectral coefficients of scalar quantities to sca...
Definition: sptgpm.f:56
subroutine sptgpmv(IROMB, MAXWV, KMAX, MI, MJ, KWSKIP, KGSKIP, NISKIP, NJSKIP, RLAT1, RLON1, DLAT, DLON, WAVED, WAVEZ, UM, VM)
This subprogram performs a spherical transform from spectral coefficients of divergences and curls to...
Definition: sptgpmv.f:63
subroutine sptranv(IROMB, MAXWV, IDRT, IMAX, JMAX, KMAX, IPRIME, ISKIP, JNSKIP, JSSKIP, KWSKIP, KGSKIP, JBEG, JEND, JCPU, WAVED, WAVEZ, GRIDUN, GRIDUS, GRIDVN, GRIDVS, IDIR)
This subprogram performs a spherical transform between spectral coefficients of divergences and curls...
Definition: sptranv.f:91
subroutine sptrunmv(IROMB, MAXWV, IDRTI, IMAXI, JMAXI, KMAX, MI, MJ, IPRIME, ISKIPI, JSKIPI, KSKIPI, KGSKIP, NISKIP, NJSKIP, JCPU, RLAT1, RLON1, DLAT, DLON, GRIDUI, GRIDVI, LUV, UM, VM, LDZ, DM, ZM, LPS, PM, SM)
THIS SUBPROGRAM SPECTRALLY TRUNCATES VECTOR FIELDS ON A GLOBAL CYLINDRICAL GRID, RETURNING THE FIELDS...
Definition: sptrunmv.f:96
subroutine spwget(IROMB, MAXWV, EPS, EPSTOP, ENN1, ELONN1, EON, EONTOP)
This subprogram gets wave-space constants.
Definition: spwget.f:18