NCEPLIBS-sp 2.4.0
sptrunsv.f
Go to the documentation of this file.
1C> @file
2C> @brief Spectrally interpolate vectors to polar stereo.
3C>
4C> 96-02-29 | Iredell | Initial.
5C> 1998-12-15 | Iredell | Openmp directives inserted.
6C>
7C> @author Iredell @date 96-02-29
8
9C> This subprogram spectrally truncates vector fields
10C> on a global cylindrical grid, returning the fields
11C> to specific pairs of polar stereographic scalar fields.
12C>
13C> The wave-space can be either triangular or rhomboidal.
14C>
15C> The grid-space can be either an equally-spaced grid
16C> (with or without pole points) or a gaussian grid.
17C>
18C> The grid fields may have general indexing.
19C>
20C> The transforms are all multiprocessed.
21C>
22C> Transform several fields at a time to improve vectorization.
23C>
24C> Subprogram can be called from a multiprocessing environment.
25C>
26C> Minimum grid dimensions for unaliased transforms to spectral:
27C> Dimension |Linear |Quadratic
28C> ----------------------- |--------- |-------------
29C> IMAX |2*MAXWV+2 |3*MAXWV/2*2+2
30C> JMAX (IDRT=4,IROMB=0) |1*MAXWV+1 |3*MAXWV/2+1
31C> JMAX (IDRT=4,IROMB=1) |2*MAXWV+1 |5*MAXWV/2+1
32C> JMAX (IDRT=0,IROMB=0) |2*MAXWV+3 |3*MAXWV/2*2+3
33C> JMAX (IDRT=0,IROMB=1) |4*MAXWV+3 |5*MAXWV/2*2+3
34C> JMAX (IDRT=256,IROMB=0) |2*MAXWV+1 |3*MAXWV/2*2+1
35C> JMAX (IDRT=256,IROMB=1) |4*MAXWV+1 |5*MAXWV/2*2+1
36C>
37C> @param IROMB integer spectral domain shape
38C> (0 for triangular, 1 for rhomboidal)
39C> @param MAXWV integer spectral truncation
40C> @param IDRTI integer input grid identifier
41C> - IDRTI=4 for Gaussian grid
42C> - IDRTI=0 for equally-spaced grid including poles
43C> - IDRTI=256 for equally-spaced grid excluding poles
44C> @param IMAXI integer even number of input longitudes.
45C> @param JMAXI integer number of input latitudes.
46C> @param KMAX integer number of fields to transform.
47C> @param NPS integer odd order of the polar stereographic grids
48C> @param IPRIME integer input longitude index for the prime meridian.
49C> (defaults to 1 if IPRIME=0)
50C> (output longitude index for prime meridian assumed 1.)
51C> @param ISKIPI integer skip number between input longitudes
52C> (defaults to 1 if ISKIPI=0)
53C> @param JSKIPI integer skip number between input latitudes from south
54C> (defaults to -IMAXI if JSKIPI=0)
55C> @param KSKIPI integer skip number between input grid fields
56C> (defaults to IMAXI*JMAXI if KSKIPI=0)
57C> @param KGSKIP integer skip number between grid fields
58C> (defaults to NPS*NPS if KGSKIP=0)
59C> @param NISKIP integer skip number between grid i-points
60C> (defaults to 1 if NISKIP=0)
61C> @param NJSKIP integer skip number between grid j-points
62C> (defaults to NPS if NJSKIP=0)
63C> @param JCPU integer number of cpus over which to multiprocess
64C> (defaults to environment NCPUS if JCPU=0)
65C> @param TRUE real latitude at which ps grid is true (usually 60.)
66C> @param XMESH real grid length at true latitude (m)
67C> @param ORIENT real longitude at bottom of Northern PS grid
68C> (Southern PS grid will have opposite orientation.)
69C> @param GRIDUI real input grid u-winds
70C> @param GRIDVI real input grid v-winds
71C> @param LUV logical flag whether to return winds
72C> @param LDZ logical flag whether to return divergence and vorticity
73C> @param LPS logical flag whether to return potential and streamfcn
74C> @param UN real northern ps u-winds if luv
75C> @param VN real northern ps v-winds if luv
76C> @param US real southern ps u-winds if luv
77C> @param VS real southern ps v-winds if luv
78C> @param DN real northern divergences if ldz
79C> @param ZN real northern vorticities if ldz
80C> @param DS real southern divergences if ldz
81C> @param ZS real southern vorticities if ldz
82C> @param PN real northern potentials if lps
83C> @param SN real northern streamfcns if lps
84C> @param PS real southern potentials if lps
85C> @param SS real southern streamfcns if lps
86C>
87C> @author Iredell @date 96-02-29
88 SUBROUTINE sptrunsv(IROMB,MAXWV,IDRTI,IMAXI,JMAXI,KMAX,NPS,
89 & IPRIME,ISKIPI,JSKIPI,KSKIPI,KGSKIP,
90 & NISKIP,NJSKIP,JCPU,TRUE,XMESH,ORIENT,
91 & GRIDUI,GRIDVI,
92 & LUV,UN,VN,US,VS,LDZ,DN,ZN,DS,ZS,
93 & LPS,PN,SN,PS,SS)
94 LOGICAL LUV,LDZ,LPS
95 REAL GRIDUI(*),GRIDVI(*)
96 REAL UN(*),VN(*),US(*),VS(*),DN(*),ZN(*),DS(*),ZS(*)
97 REAL PN(*),SN(*),PS(*),SS(*)
98 REAL EPS((MAXWV+1)*((IROMB+1)*MAXWV+2)/2),EPSTOP(MAXWV+1)
99 REAL ENN1((MAXWV+1)*((IROMB+1)*MAXWV+2)/2)
100 REAL ELONN1((MAXWV+1)*((IROMB+1)*MAXWV+2)/2)
101 REAL EON((MAXWV+1)*((IROMB+1)*MAXWV+2)/2),EONTOP(MAXWV+1)
102 REAL WD((MAXWV+1)*((IROMB+1)*MAXWV+2)/2*2+1,KMAX)
103 REAL WZ((MAXWV+1)*((IROMB+1)*MAXWV+2)/2*2+1,KMAX)
104
105C TRANSFORM INPUT GRID TO WAVE
106 jc=jcpu
107 IF(jc.EQ.0) jc=ncpus()
108 mx=(maxwv+1)*((iromb+1)*maxwv+2)/2
109 mdim=2*mx+1
110 jn=-jskipi
111 IF(jn.EQ.0) jn=imaxi
112 js=-jn
113 inp=(jmaxi-1)*max(0,-jn)+1
114 isp=(jmaxi-1)*max(0,-js)+1
115 CALL sptranv(iromb,maxwv,idrti,imaxi,jmaxi,kmax,
116 & iprime,iskipi,jn,js,mdim,kskipi,0,0,jc,
117 & wd,wz,
118 & gridui(inp),gridui(isp),gridvi(inp),gridvi(isp),-1)
119
120C TRANSFORM WAVE TO OUTPUT WINDS
121 IF(luv) THEN
122 CALL sptgpsv(iromb,maxwv,kmax,nps,mdim,kgskip,niskip,njskip,
123 & true,xmesh,orient,wd,wz,un,vn,us,vs)
124 ENDIF
125
126C TRANSFORM WAVE TO OUTPUT DIVERGENCE AND VORTICITY
127 IF(ldz) THEN
128 CALL sptgps(iromb,maxwv,kmax,nps,mdim,kgskip,niskip,njskip,
129 & true,xmesh,orient,wd,dn,ds)
130 CALL sptgps(iromb,maxwv,kmax,nps,mdim,kgskip,niskip,njskip,
131 & true,xmesh,orient,wz,zn,zs)
132 ENDIF
133
134C TRANSFORM WAVE TO OUTPUT POTENTIAL AND STREAMFUNCTION
135 IF(lps) THEN
136 CALL spwget(iromb,maxwv,eps,epstop,enn1,elonn1,eon,eontop)
137C$OMP PARALLEL DO
138 DO k=1,kmax
139 CALL splaplac(iromb,maxwv,enn1,wd(1,k),wd(1,k),-1)
140 CALL splaplac(iromb,maxwv,enn1,wz(1,k),wz(1,k),-1)
141 wd(1:2,k)=0.
142 wz(1:2,k)=0.
143 ENDDO
144 CALL sptgps(iromb,maxwv,kmax,nps,mdim,kgskip,niskip,njskip,
145 & true,xmesh,orient,wd,pn,ps)
146 CALL sptgps(iromb,maxwv,kmax,nps,mdim,kgskip,niskip,njskip,
147 & true,xmesh,orient,wz,sn,ss)
148 ENDIF
149 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 sptgps(IROMB, MAXWV, KMAX, NPS, KWSKIP, KGSKIP, NISKIP, NJSKIP, TRUE, XMESH, ORIENT, WAVE, GN, GS)
This subprogram performs a spherical transform from spectral coefficients of scalar quantities to sca...
Definition: sptgps.f:81
subroutine sptgpsv(IROMB, MAXWV, KMAX, NPS, KWSKIP, KGSKIP, NISKIP, NJSKIP, TRUE, XMESH, ORIENT, WAVED, WAVEZ, UN, VN, US, VS)
This subprogram performs a spherical transform from spectral coefficients of divergences and curls to...
Definition: sptgpsv.f:83
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 sptrunsv(IROMB, MAXWV, IDRTI, IMAXI, JMAXI, KMAX, NPS, IPRIME, ISKIPI, JSKIPI, KSKIPI, KGSKIP, NISKIP, NJSKIP, JCPU, TRUE, XMESH, ORIENT, GRIDUI, GRIDVI, LUV, UN, VN, US, VS, LDZ, DN, ZN, DS, ZS, LPS, PN, SN, PS, SS)
This subprogram spectrally truncates vector fields on a global cylindrical grid, returning the fields...
Definition: sptrunsv.f:94
subroutine spwget(IROMB, MAXWV, EPS, EPSTOP, ENN1, ELONN1, EON, EONTOP)
This subprogram gets wave-space constants.
Definition: spwget.f:18