Warning
Page under construction.
The following subroutines can be called by programs that are linked against ecTrans. Here we describe in detail their signatures.
Note
ecTrans is a legacy code with an accumulated 30 years of history. Over this time certain features enabled through optional arguments will have fallen out of use. We are currently reviewing all options to identify those that can be safely deleted, but this takes time. In the mean time, we have tagged below all options we deem to be "potentially deprecatable".
ecTrans in principle follows the coding standard and conventions outlined in the IFS Documentation - Part VI: Technical and Computational Procedures section 1.5. Following these standards, all variable names must begin with a one- or two-character prefix denoting their scope (module level, dummy argument, local variables, loop index, or parameter) and type. These are outlined in Table 1.2. Dummy variables have the following prefixes:
K
- integerP
- real (single or double precision)LD
- logicalCD
- characterYD
- derived typeKIND
parametersAs with the IFS, integer and real variables in ecTrans always have an explicit KIND
specification.
These are defined in the PARKIND1
module which is part of the Fiat library (a dependency of ecTrans). To understand the
subroutines described here, only two must be considered:
INTEGER, PARAMETER :: JPIM = SELECTED_INT_KIND(9)
(i.e. 4-byte integer)INTEGER, PARAMETER :: JPRD = SELECTED_REAL_KIND(13,300)
(i.e. 8-byte float)SETUP_TRANS0
SUBROUTINE SETUP_TRANS0(KOUT, KERR, KPRINTLEV, KMAX_RESOL, KPROMATR, KPRGPNS, KPRGPEW, &
& KPRTRW, KCOMBFLEN, LDMPOFF, LDSYNC_TRANS, KTRANS_SYNC_LEVEL, &
& LDEQ_REGIONS, K_REGIONS_NS, K_REGIONS_EW, K_REGIONS, PRAD, LDALLOPERM, &
& KOPT_MEMORY_TR)
This subroutine initialises all resolution-agnostic aspects of the ecTrans library. It must be called before any other ecTrans subroutines.
OPTIONAL, INTENT(IN)
argumentsINTEGER(KIND=JPIM), OPTIONAL, INTENT(IN) :: KOUT
6
(i.e. STDOUT)INTEGER(KIND=JPIM), OPTIONAL, INTENT(IN) :: KERR
0
(i.e. STDERR)INTEGER(KIND=JPIM), OPTIONAL, INTENT(IN) :: KPRINTLEV
0
INTEGER(KIND=JPIM), OPTIONAL, INTENT(IN) :: KMAX_RESOL
1
INTEGER(KIND=JPIM), OPTIONAL, INTENT(IN) :: KPROMATR
0
INTEGER(KIND=JPIM), OPTIONAL, INTENT(IN) :: KPRGPNS
1
INTEGER(KIND=JPIM), OPTIONAL, INTENT(IN) :: KPRGPEW
1
INTEGER(KIND=JPIM), OPTIONAL, INTENT(IN) :: KPRTRW
1
INTEGER(KIND=JPIM), OPTIONAL, INTENT(IN) :: KCOMBFLEN
LOGICAL, OPTIONAL, INTENT(IN) :: LDMPOFF
.FALSE.
LOGICAL, OPTIONAL, INTENT(IN) :: LDSYNC_TRANS
.FALSE.
INTEGER(KIND=JPIM), OPTIONAL, INTENT(IN) :: KTRANS_SYNC_LEVEL
0
LOGICAL, OPTIONAL, INTENT(IN) :: LDEQ_REGIONS
.FALSE.
LOGICAL, OPTIONAL, INTENT(IN) :: LDALLOPERM
.FALSE.
REAL(KIND=JPRD), OPTIONAL, INTENT(IN) :: PRAD
6371229.0
INTEGER(KIND=JPIM), OPTIONAL, INTENT(IN) :: KOPT_MEMORY_TR
0
OPTIONAL, INTENT(OUT)
argumentsINTEGER(KIND=JPIM), OPTIONAL, INTENT(OUT) :: K_REGIONS(:)
INTEGER(KIND=JPIM), OPTIONAL, INTENT(OUT) :: K_REGIONS_NS
INTEGER(KIND=JPIM), OPTIONAL, INTENT(OUT) :: K_REGIONS_EW
SETUP_TRANS
SUBROUTINE SETUP_TRANS(KSMAX, KDGL, KDLON, KLOEN, LDSPLIT, PSTRET, KTMAX, KRESOL, &
& PWEIGHT, LDGRIDONLY, LDUSERPNM, LDKEEPRPNM, LDUSEFLT, &
& LDSPSETUPONLY, LDPNMONLY, LDUSEFFTW, LD_ALL_FFTW, LDLL, LDSHIFTLL, &
& CDIO_LEGPOL, CDLEGPOLFNAME, KLEGPOLPTR, KLEGPOLPTR_LEN)
This subroutine establishes the environment for performing a spectral transform. Each call to this
subroutine creates work arrays for a certain resolution. One can do this for NMAX_RESOL
different
resolutions (internal parameter corresponding to KMAX_RESOL
argument of SETUP_TRANS0
). One
must call SETUP_TRANS0
before calling this subroutine for the first time.
INTENT(IN)
argumentsINTEGER(KIND=JPIM), INTENT(IN) :: KSMAX
INTEGER(KIND=JPIM), INTENT(IN) :: KDGL
OPTIONAL, INTENT(IN)
argumentsINTEGER(KIND=JPIM), OPTIONAL, INTENT(IN) :: KDLON
2 * KDGL
, unless LDLL
is .TRUE.
in which case2 * KDGL + 2
, unless KLOEN
is provided in which case it will takeMAXVAL(KLOEN)
.INTEGER(KIND=JPIM), OPTIONAL, INTENT(IN) :: KLOEN(:)
KDLON
points.LOGICAL, OPTIONAL, INTENT(IN) :: LDSPLIT
.FALSE.
REAL(KIND=JPRD), OPTIONAL, INTENT(IN) :: PSTRET
1.0
INTEGER(KIND=JPIM), OPTIONAL, INTENT(IN) :: KTMAX
KSMAX
REAL(KIND=JPRD), OPTIONAL, INTENT(IN) :: PWEIGHT(:)
LOGICAL, OPTIONAL, INTENT(IN) :: LDGRIDONLY
.FALSE.
LOGICAL, OPTIONAL, INTENT(IN) :: LDUSEFLT
.FALSE.
LOGICAL, OPTIONAL, INTENT(IN) :: LD_ALL_FFTW
.FALSE.
LOGICAL, OPTIONAL, INTENT(IN) :: LDUSERPNM
.TRUE.
LOGICAL, OPTIONAL, INTENT(IN) :: LDKEEPRPNM
.FALSE.
LOGICAL, OPTIONAL, INTENT(IN) :: LDSPSETUPONLY
.FALSE.
LOGICAL, OPTIONAL, INTENT(IN) :: LDPNMONLY
.FALSE.
LOGICAL, OPTIONAL, INTENT(IN) :: LDLL
.FALSE.
LOGICAL, OPTIONAL, INTENT(IN) :: LDSHIFTLL
.FALSE.
CHARACTER(LEN=*), OPTIONAL, INTENT(IN) :: CDIO_LEGPOL
'readf'
is passed, Legendre polynomials will be read from the file
given to CDLEGPOLFNAME
. If 'writef'
is passed, Legendre polynomials will be written to the
file given to CDLEGPOLFNAME
. If 'membuf'
is passed, Legendre polynomials will be read from
shared memory.CHARACTER(LEN=*), OPTIONAL, INTENT(IN) :: CDLEGPOLFNAME
CDIO_LEGPOL
.TYPE(C_PTR), OPTIONAL, INTENT(IN) :: KLEGPOLPTR
CDIO_LEGPOL
.INTEGER(C_SIZE_T), OPTIONAL, INTENT(IN) :: KLEGPOLPTR_LEN
CDIO_LEGPOL
.OPTIONAL, INTENT(OUT)
argumentsINTEGER(KIND=JPIM), OPTIONAL, INTENT(OUT) :: KRESOL
TRANS_INQ
SUBROUTINE TRANS_INQ(KRESOL, KSPEC, KSPEC2, KSPEC2G, KSPEC2MX, KNUMP, KGPTOT, KGPTOTG, KGPTOTMX, &
& KGPTOTL, KMYMS, KASM0, KUMPP, KPOSSP, KPTRMS, KALLMS, KDIM0G, KFRSTLAT, &
& KLSTLAT, KFRSTLOFF, KPTRLAT, KPTRFRSTLAT, KPTRLSTLAT, KPTRFLOFF, KSTA, &
& KONL, KULTPP, KPTRLS, KNMENG, KPRTRW, KMYSETW, KMYSETV, KMY_REGION_NS, &
& KMY_REGION_EW, LDSPLITLAT, KSMAX, PLAPIN, KNVALUE, KDEF_RESOL, LDLAM, PMU, &
& PGW, PRPNM, KLEI3, KSPOLEGL, KPMS, KDGLU)
This subroutine allows one to obtain information on the transform configuration for a given
KRESOL
handle (or the default handle if this argument is omitted).
OPTIONAL, INTENT(IN)
argumentsINTEGER(KIND=JPIM), OPTIONAL, INTENT(IN) :: KRESOL
OPTIONAL, INTENT(OUT)
argumentsINTEGER(KIND=JPIM), OPTIONAL, INTENT(OUT) :: KSPEC
INTEGER(KIND=JPIM), OPTIONAL, INTENT(OUT) :: KSPEC2
2 * KSPEC
.INTEGER(KIND=JPIM), OPTIONAL, INTENT(OUT) :: KSPEC2G
INTEGER(KIND=JPIM), OPTIONAL, INTENT(OUT) :: KSPEC2MX
INTEGER(KIND=JPIM), OPTIONAL, INTENT(OUT) :: KNUMP
INTEGER(KIND=JPIM), OPTIONAL, INTENT(OUT) :: KGPTOT
INTEGER(KIND=JPIM), OPTIONAL, INTENT(OUT) :: KGPTOTG
INTEGER(KIND=JPIM), OPTIONAL, INTENT(OUT) :: KGPTOTMX
INTEGER(KIND=JPIM), OPTIONAL, INTENT(OUT) :: KGPTOTL(:,:)
INTEGER(KIND=JPIM), OPTIONAL, INTENT(OUT) :: KMYMS(:)
INTEGER(KIND=JPIM), OPTIONAL, INTENT(OUT) :: KASM0(0:)
INTEGER(KIND=JPIM), OPTIONAL, INTENT(OUT) :: KUMPP(:)
INTEGER(KIND=JPIM), OPTIONAL, INTENT(OUT) :: KPOSSP(:)
INTEGER(KIND=JPIM), OPTIONAL, INTENT(OUT) :: KPTRMS(:)
INTEGER(KIND=JPIM), OPTIONAL, INTENT(OUT) :: KALLMS(:)
INTEGER(KIND=JPIM), OPTIONAL, INTENT(OUT) :: KDIM0G(0:)
INTEGER(KIND=JPIM), OPTIONAL, INTENT(OUT) :: KFRSTLAT(:)
INTEGER(KIND=JPIM), OPTIONAL, INTENT(OUT) :: KLSTLAT(:)
INTEGER(KIND=JPIM), OPTIONAL, INTENT(OUT) :: KFRSTLOFF
INTEGER(KIND=JPIM), OPTIONAL, INTENT(OUT) :: KPTRLAT(:)
INTEGER(KIND=JPIM), OPTIONAL, INTENT(OUT) :: KPTRFRSTLAT(:)
INTEGER(KIND=JPIM), OPTIONAL, INTENT(OUT) :: KPTRLSTLAT(:)
INTEGER(KIND=JPIM), OPTIONAL, INTENT(OUT) :: KPTRFLOFF
INTEGER(KIND=JPIM), OPTIONAL, INTENT(OUT) :: KSTA(:,:)
INTEGER(KIND=JPIM), OPTIONAL, INTENT(OUT) :: KONL(:,:)
INTEGER(KIND=JPIM), OPTIONAL, INTENT(OUT) :: KULTPP(:)
INTEGER(KIND=JPIM), OPTIONAL, INTENT(OUT) :: KPTRLS(:)
INTEGER(KIND=JPIM), OPTIONAL, INTENT(OUT) :: KNMENG(:)
INTEGER(KIND=JPIM), OPTIONAL, INTENT(OUT) :: KPRTRW
INTEGER(KIND=JPIM), OPTIONAL, INTENT(OUT) :: KMYSETW
INTEGER(KIND=JPIM), OPTIONAL, INTENT(OUT) :: KMYSETV
INTEGER(KIND=JPIM), OPTIONAL, INTENT(OUT) :: KMY_REGION_NS
INTEGER(KIND=JPIM), OPTIONAL, INTENT(OUT) :: KMY_REGION_EW
LOGICAL, OPTIONAL, INTENT(OUT) :: LDSPLITLAT(:)
INTEGER(KIND=JPIM), OPTIONAL, INTENT(OUT) :: KSMAX
REAL(KIND=JPRB), OPTIONAL, INTENT(OUT) :: PLAPIN(-1:)
INTEGER(KIND=JPIM), OPTIONAL, INTENT(OUT) :: KNVALUE(:)
INTEGER(KIND=JPIM), OPTIONAL, INTENT(OUT) :: KDEF_RESOL
LOGICAL, OPTIONAL, INTENT(OUT) :: LDLAM
REAL(KIND=JPRD), OPTIONAL, INTENT(OUT) :: PMU(:)
REAL(KIND=JPRB), OPTIONAL, INTENT(OUT) :: PGW(:)
REAL(KIND=JPRB), OPTIONAL, INTENT(OUT) :: PRPNM(:,:)
INTEGER(KIND=JPIM), OPTIONAL, INTENT(OUT) :: KLEI3
INTEGER(KIND=JPIM), OPTIONAL, INTENT(OUT) :: KSPOLEGL
INTEGER(KIND=JPIM), OPTIONAL, INTENT(OUT) :: KPMS(0:)
INTEGER(KIND=JPIM), OPTIONAL, INTENT(OUT) :: KDGLU(0:)
DIR_TRANS
SUBROUTINE DIR_TRANS(PSPVOR, PSPDIV, PSPSCALAR, PSPSC3A, PSPSC3B, PSPSC2, LDLATLON, KPROMA, &
& KVSETUV, KVSETSC, KRESOL, KVSETSC3A, KVSETSC3B, KVSETSC2, PGP, PGPUV, &
& PGP3A, PGP3B, PGP2)
This subroutine performs a direct spectral transform.
DIR_TRANS
supports two modes of passing arrays (which we call "call modes" below):
PGP
, receive PSPVOR
, PSPDIV
, PSPSCALAR
.PGPUV
, PGP3A
, PGP3B
, PGP2
, receive PSPVOR
, PSPDIV
, PSPSC3A
, PSPSC3B
, PSPSC2
.Although call mode 1 is simpler, it can often entail having to create array copies in order to
prepare the different elements of PGP
. Hence, with call mode 2, the input and output arrays are
categorized according to whether they are wind related, 3D scalar, or 2D scalar.
Note that PSPSC3A
/PGP3A
and PSPSC3B
/PGP3B
are essentially identical in type and size.
There are cases where it is useful to have the ability to pass two sets of scalar fields to be
transformed independently.
OPTIONAL, INTENT(IN)
argumentsINTEGER(KIND=JPIM), OPTIONAL, INTENT(IN) :: KPROMA
D%NGPTOT
INTEGER(KIND=JPIM), OPTIONAL, INTENT(IN) :: KVSETUV(:)
INTEGER(KIND=JPIM), OPTIONAL, INTENT(IN) :: KVSETSC(:)
INTEGER(KIND=JPIM), OPTIONAL, INTENT(IN) :: KVSETSC3A(:)
INTEGER(KIND=JPIM), OPTIONAL, INTENT(IN) :: KVSETSC3B(:)
INTEGER(KIND=JPIM), OPTIONAL, INTENT(IN) :: KVSETSC2(:)
INTEGER(KIND=JPIM), OPTIONAL, INTENT(IN) :: KRESOL
SETUP_TRANS
.LOGICAL, OPTIONAL, INTENT(IN) :: LDLATLON
REAL(KIND=JPRB), OPTIONAL, INTENT(IN) :: PGP(:,:,:)
REAL(KIND=JPRB), OPTIONAL, INTENT(IN) :: PGPUV(:,:,:,:)
DIR_TRANS
only operates on U/V winds, not vorticity/divergence. Hence, the third dimensionPGPUV
must be equal to 2
(one for U and one for V).REAL(KIND=JPRB), OPTIONAL, INTENT(IN) :: PGP3A(:,:,:,:)
PSPSC3A
.REAL(KIND=JPRB), OPTIONAL, INTENT(IN) :: PGP3B(:,:,:,:)
PSPSC3B
.REAL(KIND=JPRB), OPTIONAL, INTENT(IN) :: PGP2(:,:,:)
PSPSC2
.OPTIONAL, INTENT(OUT)
argumentsREAL(KIND=JPRB), OPTIONAL, INTENT(OUT) :: PSPVOR(:,:)
REAL(KIND=JPRB), OPTIONAL, INTENT(OUT) :: PSPDIV(:,:)
REAL(KIND=JPRB), OPTIONAL, INTENT(OUT) :: PSPSCALAR(:,:)
REAL(KIND=JPRB), OPTIONAL, INTENT(OUT) :: PSPSC3A(:,:,:)
PGP3A
.REAL(KIND=JPRB), OPTIONAL, INTENT(OUT) :: PSPSC3B(:,:,:)
PGP3B
.REAL(KIND=JPRB), OPTIONAL, INTENT(OUT) :: PSPSC2(:,:)
PGP2
.INV_TRANS
SUBROUTINE INV_TRANS(PSPVOR, PSPDIV, PSPSCALAR, PSPSC3A, PSPSC3B, PSPSC2, &
& FSPGL_PROC, LDSCDERS, LDVORGP, LDDIVGP, LDUVDER, LDLATLON, KPROMA, KVSETUV, &
& KVSETSC, KRESOL, KVSETSC3A, KVSETSC3B, KVSETSC2, PGP, PGPUV, PGP3A, PGP3B, &
& PGP2)
This subroutine performs an inverse spectral transform.
INV_TRANS
supports two modes of passing arrays (which we call "call modes" below):
PSPVOR
, PSPDIV
, PSPSCALAR
, receive PGP
.PSPVOR
, PSPDIV
, PSPSC3A
, PSPSC3B
, PSPSC2
, receive PGPUV
, PGP3A
, PGP3B
, PGP2
Although call mode 1 is simpler, it can often entail having to create array copies in order to
subsequently operate on the different elements of PGP
. Hence, with call mode 2, the output
arrays are already categorized according to whether they are wind related, 3D scalar, or 2D
scalar.
Note that PSPSC3A
/PGP3A
and PSPSC3B
/PGP3B
are essentially identical in type and size.
There are cases where it is useful to have the ability to pass two sets of scalar fields to be
transformed independently.
OPTIONAL, INTENT(IN)
argumentsREAL(KIND=JPRB), OPTIONAL, INTENT(IN) :: PSPVOR(:,:)
REAL(KIND=JPRB), OPTIONAL, INTENT(IN) :: PSPDIV(:,:)
REAL(KIND=JPRB), OPTIONAL, INTENT(IN) :: PSPSCALAR(:,:)
REAL(KIND=JPRB), OPTIONAL, INTENT(IN) :: PSPSC3A(:,:,:)
PGP3A
.REAL(KIND=JPRB), OPTIONAL, INTENT(IN) :: PSPSC3B(:,:,:)
PGP3B
.REAL(KIND=JPRB), OPTIONAL, INTENT(IN) :: PSPSC2(:,:)
PGP2
.LOGICAL, OPTIONAL, INTENT(IN) :: LDSCDERS
.FALSE.
LOGICAL, OPTIONAL, INTENT(IN) :: LDVORGP
.FALSE.
LOGICAL, OPTIONAL, INTENT(IN) :: LDDIVGP
.FALSE.
LOGICAL, OPTIONAL, INTENT(IN) :: LDUVDER
.FALSE.
LOGICAL, OPTIONAL, INTENT(IN) :: LDLATLON
.FALSE.
INTEGER(KIND=JPIM), OPTIONAL, INTENT(IN) :: KPROMA
D%NGPTOT
INTEGER(KIND=JPIM), OPTIONAL, INTENT(IN) :: KVSETUV(:)
INTEGER(KIND=JPIM), OPTIONAL, INTENT(IN) :: KVSETSC(:)
INTEGER(KIND=JPIM), OPTIONAL, INTENT(IN) :: KVSETSC3A(:)
INTEGER(KIND=JPIM), OPTIONAL, INTENT(IN) :: KVSETSC3B(:)
INTEGER(KIND=JPIM), OPTIONAL, INTENT(IN) :: KVSETSC2(:)
INTEGER(KIND=JPIM), OPTIONAL, INTENT(IN) :: KRESOL
SETUP_TRANS
.OPTIONAL, EXTERNAL
argumentsOPTIONAL, EXTERNAL :: FSPGL_PROC
OPTIONAL, INTENT(OUT)
argumentsREAL(KIND=JPRB), OPTIONAL, INTENT(OUT) :: PGP(:,:,:)
PSPVOR
is present and LDVORGP
)PSPDIV
is present and LDDIVGP
)PSPVOR
and PSPDIV
are present)PSPVOR
and PSPDIV
are present)PSPSCALAR
is present)PSPSCALAR
is present and LDSCDERS
)PSPVOR
and PSPDIV
are present and LDUVDER
)PSPVOR
and PSPDIV
are present and LDUVDER
)PSPSCALAR
is present and LDSCDERS
)REAL(KIND=JPRB), OPTIONAL, INTENT(OUT) :: PGPUV(:,:,:,:)
PSPVOR
is present and LDVORGP
)PSPDIV
is present and LDDIVGP
)PSPVOR
and PSPDIV
are present)PSPVOR
and PSPDIV
are present)PSPVOR
and PSPDIV
are present and LDUVDER
)PSPVOR
and PSPDIV
are present and LDUVDER
)REAL(KIND=JPRB), OPTIONAL, INTENT(OUT) :: PGP3A(:,:,:,:)
PSPSC3A
.LDSCDERS
)LDSCDERS
)REAL(KIND=JPRB), OPTIONAL, INTENT(OUT) :: PGP3B(:,:,:,:)
PSPSC3B
.LDSCDERS
)LDSCDERS
)REAL(KIND=JPRB), OPTIONAL, INTENT(OUT) :: PGP2(:,:,:)
PSPSC2
.LDSCDERS
)LDSCDERS
)SPECNORM
SUBROUTINE SPECNORM(PNORM, PSPEC, KVSET, KMASTER, KRESOL, PMET)
This subroutine computes the spectral norms of the fields in the given input array and returns the values to a specific MPI task.
OPTIONAL, INTENT(IN)
argumentsREAL(KIND=JPRB), OPTIONAL, INTENT(IN) :: PSPEC(:,:)
KVSET
.INTEGER(KIND=JPIM), OPTIONAL, INTENT(IN) :: KVSET(:)
INTEGER(KIND=JPIM), OPTIONAL, INTENT(IN) :: KMASTER
PNORM
will contain values.1
REAL(KIND=JPRB), OPTIONAL, INTENT(IN) :: PMET(:)
0:R%NSMAX
) (R%NSMAX
is the same as KSMAX
passed to SETUP_TRANS
)1.0
INTEGER(KIND=JPIM), OPTIONAL, INTENT(IN) :: KRESOL
SETUP_TRANS
.1
(i.e. first resolution handle)INTENT(OUT)
argumentsREAL(KIND=JPRB), INTENT(OUT) :: PNORM(:)
KMASTER
needs to
actually allocate its array though.GPNORM_TRANS
SUBROUTINE GPNORM_TRANS(PGP, KFIELDS, KPROMA, PAVE, PMIN, PMAX, LDAVE_ONLY, KRESOL)
This subroutine calculates the global average (and optionally the global minimum and maximum) of the input fields in grid point space. The results are communicated to task 1 only.
INTENT(IN)
argumentsREAL(KIND=JPRB), INTENT(IN) :: PGP(:,:,:)
INTEGER(KIND=JPIM), INTENT(IN) :: KFIELDS
INTEGER(KIND=JPIM), INTENT(IN) :: KPROMA
LOGICAL, INTENT(IN) :: LDAVE_ONLY
PMIN
and PMAX
should already contain the minimum and maximum values of each
field for each MPI task's local region, and afterwards they will contain the global maximum and
minimum (again, only on task 1).OPTIONAL, INTENT(IN)
argumentsINTEGER(KIND=JPIM) ,OPTIONAL, INTENT(IN) :: KRESOL
SETUP_TRANS
.1
(i.e. first resolution handle)INTENT(OUT)
argumentsREAL(KIND=JPRB), INTENT(OUT) :: PAVE(:)
KFIELDS
)INTENT(INOUT)
argumentsREAL(KIND=JPRB), INTENT(INOUT) :: PMIN(:)
KFIELDS
)REAL(KIND=JPRB), INTENT(INOUT) :: PMAX(:)
KFIELDS
)DIST_GRID
SUBROUTINE DIST_GRID(PGPG, KPROMA, KFDISTG, KFROM, KRESOL, PGP, KSORT)
The subroutine distributes one or more global grid point fields resident across one or more MPI
tasks (each field is entirely resident on a single MPI task without decomposition) among all other
tasks, decomposing them along the way. It is the opposite of GATH_GRID
.
The figure below illustrates an example in which four fields are distributed equally across four MPI
tasks. Note that at the end, each task receives part of all four input fields. Note also that the
NPROMA-based blocking is introduced at the end of DIST_GRID
.
INTENT(IN)
argumentsINTEGER(KIND=JPIM), INTENT(IN) :: KFDISTG
INTEGER(KIND=JPIM), INTENT(IN) :: KFROM(:)
OPTIONAL, INTENT(IN)
argumentsREAL(KIND=JPRB), OPTIONAL, INTENT(IN) :: PGPG(:,:)
PGPG
).INTEGER(KIND=JPIM), OPTIONAL, INTENT(IN) :: KPROMA
D%NGPTOT
INTEGER(KIND=JPIM), OPTIONAL, INTENT(IN) :: KRESOL
SETUP_TRANS
.1
(i.e. first resolution handle)INTEGER(KIND=JPIM), OPTIONAL, INTENT(IN) :: KSORT(:)
PGP
)
Default: no sortingINTENT(OUT)
argumentsREAL(KIND=JPRB), OPTIONAL, INTENT(OUT) :: PGP(:,:,:)
GATH_GRID
SUBROUTINE GATH_GRID(PGPG, KPROMA, KFGATHG, KTO, KRESOL, PGP)
The subroutine gathers one or more fields distributed across one or more MPI tasks and sends them
to a particular MPI task. It is the opposite of DIST_GRID
.
The figure below illustrates an example in which four fields which are distributed equally across four MPI tasks are gathered so that each task possesses one field in its entirety. The NPROMA blocking is removed along the way so that at the end, each field only has a single horizontal dimension
INTENT(IN)
argumentsINTEGER(KIND=JPIM), INTENT(IN) :: KFGATHG
INTEGER(KIND=JPIM), INTENT(IN) :: KTO(:)
REAL(KIND=JPRB), OPTIONAL, INTENT(IN) :: PGP(:,:,:)
OPTIONAL, INTENT(IN)
argumentsINTEGER(KIND=JPIM), OPTIONAL, INTENT(IN) :: KPROMA
D%NGPTOT
INTEGER(KIND=JPIM), OPTIONAL, INTENT(IN) :: KRESOL
SETUP_TRANS
.1
(i.e. first resolution handle)OPTIONAL, INTENT(OUT)
argumentsREAL(KIND=JPRB), OPTIONAL, INTENT(OUT) :: PGPG(:,:)
GET_CURRENT
SUBROUTINE GET_CURRENT(KRESOL,LDLAM)
This subroutine allows one to retrieve the current resolution handle, NCUR_RESOL
, and the LAM
parameter, LDLAM
, which indicates if the local area version of ecTrans is being used or not.
Note
LDLAM
will always be .FALSE.
for the moment. ecTrans currently only supports global spectral
transforms, though local (bifourier) transforms may be supported in future.
OPTIONAL, INTENT(OUT)
argumentsINTEGER(KIND=JPIM), OPTIONAL, INTENT(OUT) :: KRESOL
LOGICAL, OPTIONAL, INTENT(OUT) :: LDLAM
TRANS_PNM
SUBROUTINE TRANS_PNM(KRESOL, KM, PRPNM, LDTRANSPOSE, LDCHEAP)
This subroutine computes the Legendre polynomials for a given wavenumber.
INTENT(IN)
argumentsINTEGER(KIND=JPIM), OPTIONAL, INTENT(IN) :: KM
OPTIONAL, INTENT(IN)
argumentsINTEGER(KIND=JPIM), OPTIONAL, INTENT(IN) :: KRESOL
SETUP_TRANS
.1
(i.e. first resolution handle)LOGICAL, OPTIONAL, INTENT(IN) :: LDTRANSPOSE
PRPNM
below)..FALSE.
LOGICAL, OPTIONAL, INTENT(IN) :: LDCHEAP
.FALSE.
INTENT(OUT)
argumentsREAL(KIND=JPRB), INTENT(OUT) :: PRPNM
LDTRANSPOSE
.SETUP_DIMS_MOD
for the defininition of R%NLEI3
.R%NTMAX-KM+3
, R%NLEI3
) if LDTRANSPOSE == .TRUE.
, else (R%NLEI3
, R%NTMAX-KM+3
)DIR_TRANSAD
SUBROUTINE DIR_TRANSAD(PSPVOR, PSPDIV, PSPSCALAR, PSPSC3A, PSPSC3B, PSPSC2, KPROMA, KVSETUV, &
& KVSETSC, KRESOL, KVSETSC3A, KVSETSC3B, KVSETSC2, PGP, PGPUV, PGP3A, &
& PGP3B, PGP2)
This is the adjoint version of the subroutine for performing a direct spectral transform.
Like DIR_TRANS
, DIR_TRANSAD
supports two modes of passing arrays (which we call "call modes"
below):
PGP
, receive PSPVOR
, PSPDIV
, PSPSCALAR
.PGPUV
, PGP3A
, PGP3B
, PGP2
, receive PSPVOR
, PSPDIV
, PSPSC3A
, PSPSC3B
, PSPSC2
.Although call mode 1 is simpler, it can often entail having to create array copies in order to
prepare the different elements of PGP
. Hence, with call mode 2, the input and output arrays are
categorized according to whether they are wind related, 3D scalar, or 2D scalar.
Note that PSPSC3A
/PGP3A
and PSPSC3B
/PGP3B
are essentially identical in type and size.
There are cases where it is useful to have the ability to pass two sets of scalar fields to be
transformed independently.
OPTIONAL, INTENT(IN)
argumentsINTEGER(KIND=JPIM), OPTIONAL, INTENT(IN) :: KPROMA
D%NGPTOT
INTEGER(KIND=JPIM), OPTIONAL, INTENT(IN) :: KVSETUV(:)
INTEGER(KIND=JPIM), OPTIONAL, INTENT(IN) :: KVSETSC(:)
INTEGER(KIND=JPIM), OPTIONAL, INTENT(IN) :: KVSETSC3A(:)
INTEGER(KIND=JPIM), OPTIONAL, INTENT(IN) :: KVSETSC3B(:)
INTEGER(KIND=JPIM), OPTIONAL, INTENT(IN) :: KVSETSC2(:)
INTEGER(KIND=JPIM), OPTIONAL, INTENT(IN) :: KRESOL
SETUP_TRANS
.OPTIONAL, INTENT(OUT)
argumentsREAL(KIND=JPRB), OPTIONAL, INTENT(OUT) :: PGP(:,:,:)
REAL(KIND=JPRB), OPTIONAL, INTENT(OUT) :: PGPUV(:,:,:,:)
DIR_TRANS
only operates on U/V winds, not vorticity/divergence. Hence, the third dimensionPGPUV
must be equal to 2
(one for U and one for V).REAL(KIND=JPRB), OPTIONAL, INTENT(OUT) :: PGP3A(:,:,:,:)
PSPSC3A
.REAL(KIND=JPRB), OPTIONAL, INTENT(OUT) :: PGP3B(:,:,:,:)
PSPSC3B
.REAL(KIND=JPRB), OPTIONAL, INTENT(OUT) :: PGP2(:,:,:)
PSPSC2
.OPTIONAL, INTENT(INOUT)
argumentsREAL(KIND=JPRB), OPTIONAL, INTENT(INOUT) :: PSPVOR(:,:)
REAL(KIND=JPRB), OPTIONAL, INTENT(INOUT) :: PSPDIV(:,:)
REAL(KIND=JPRB), OPTIONAL, INTENT(INOUT) :: PSPSCALAR(:,:)
REAL(KIND=JPRB), OPTIONAL, INTENT(INOUT) :: PSPSC3A(:,:,:)
PGP3A
.REAL(KIND=JPRB), OPTIONAL, INTENT(INOUT) :: PSPSC3B(:,:,:)
PGP3B
.REAL(KIND=JPRB), OPTIONAL, INTENT(INOUT) :: PSPSC2(:,:)
PGP2
.INV_TRANSAD
SUBROUTINE INV_TRANSAD(PSPVOR, PSPDIV, PSPSCALAR, PSPSC3A, PSPSC3B, PSPSC2, FSPGL_PROC, &
& LDSCDERS, LDVORGP, LDDIVGP, LDUVDER, KPROMA, KVSETUV, KVSETSC, KRESOL, &
& KVSETSC3A, KVSETSC3B, KVSETSC2, PGP, PGPUV, PGP3A, PGP3B, &
& PGP2)
This is the adjoint version of the subroutine for performing an inverse spectral transform.
Like INV_TRANS
, INV_TRANSAD
supports two modes of passing arrays (which we call "call modes"
below):
PSPVOR
, PSPDIV
, PSPSCALAR
, receive PGP
.PSPVOR
, PSPDIV
, PSPSC3A
, PSPSC3B
, PSPSC2
, receive PGPUV
, PGP3A
, PGP3B
, PGP2
Although call mode 1 is simpler, it can often entail having to create array copies in order to
subsequently operate on the different elements of PGP
. Hence, with call mode 2, the output
arrays are already categorized according to whether they are wind related, 3D scalar, or 2D
scalar.
Note that PSPSC3A
/PGP3A
and PSPSC3B
/PGP3B
are essentially identical in type and size.
There are cases where it is useful to have the ability to pass two sets of scalar fields to be
transformed independently.
OPTIONAL, INTENT(IN)
argumentsLOGICAL, OPTIONAL, INTENT(IN) :: LDSCDERS
.FALSE.
LOGICAL, OPTIONAL, INTENT(IN) :: LDVORGP
.FALSE.
LOGICAL, OPTIONAL, INTENT(IN) :: LDDIVGP
.FALSE.
LOGICAL, OPTIONAL, INTENT(IN) :: LDUVDER
.FALSE.
INTEGER(KIND=JPIM), OPTIONAL, INTENT(IN) :: KPROMA
D%NGPTOT
INTEGER(KIND=JPIM), OPTIONAL, INTENT(IN) :: KVSETUV(:)
INTEGER(KIND=JPIM), OPTIONAL, INTENT(IN) :: KVSETSC(:)
INTEGER(KIND=JPIM), OPTIONAL, INTENT(IN) :: KVSETSC3A(:)
INTEGER(KIND=JPIM), OPTIONAL, INTENT(IN) :: KVSETSC3B(:)
INTEGER(KIND=JPIM), OPTIONAL, INTENT(IN) :: KVSETSC2(:)
INTEGER(KIND=JPIM), OPTIONAL, INTENT(IN) :: KRESOL
SETUP_TRANS
.REAL(KIND=JPRB), OPTIONAL, INTENT(IN) :: PGP(:,:,:)
PSPVOR
is present and LDVORGP
)PSPDIV
is present and LDDIVGP
)PSPVOR
and PSPDIV
are present)PSPVOR
and PSPDIV
are present)PSPSCALAR
is present)PSPSCALAR
is present and LDSCDERS
)PSPVOR
and PSPDIV
are present and LDUVDER
)PSPVOR
and PSPDIV
are present and LDUVDER
)PSPSCALAR
is present and LDSCDERS
)REAL(KIND=JPRB), OPTIONAL, INTENT(IN) :: PGPUV(:,:,:,:)
PSPVOR
is present and LDVORGP
)PSPDIV
is present and LDDIVGP
)PSPVOR
and PSPDIV
are present)PSPVOR
and PSPDIV
are present)PSPVOR
and PSPDIV
are present and LDUVDER
)PSPVOR
and PSPDIV
are present and LDUVDER
)REAL(KIND=JPRB), OPTIONAL, INTENT(IN) :: PGP3A(:,:,:,:)
PSPSC3A
.LDSCDERS
)LDSCDERS
)REAL(KIND=JPRB), OPTIONAL, INTENT(IN) :: PGP3B(:,:,:,:)
PSPSC3B
.LDSCDERS
)LDSCDERS
)REAL(KIND=JPRB), OPTIONAL, INTENT(IN) :: PGP2(:,:,:)
PSPSC2
.LDSCDERS
)LDSCDERS
)OPTIONAL, EXTERNAL
argumentsOPTIONAL, EXTERNAL :: FSPGL_PROC
OPTIONAL, INTENT(OUT)
argumentsREAL(KIND=JPRB), OPTIONAL, INTENT(OUT) :: PSPVOR(:,:)
REAL(KIND=JPRB), OPTIONAL, INTENT(OUT) :: PSPDIV(:,:)
REAL(KIND=JPRB), OPTIONAL, INTENT(OUT) :: PSPSCALAR(:,:)
REAL(KIND=JPRB), OPTIONAL, INTENT(OUT) :: PSPSC3A(:,:,:)
PGP3A
.REAL(KIND=JPRB), OPTIONAL, INTENT(OUT) :: PSPSC3B(:,:,:)
PGP3B
.REAL(KIND=JPRB), OPTIONAL, INTENT(OUT) :: PSPSC2(:,:)
PGP2
.TRANS_RELEASE
SUBROUTINE TRANS_RELEASE(KRESOL)
This subroutine releases (i.e. deallocates and nullifies) all arrays related to the given
resolution tag KRESOL
.
INTENT(IN)
argumentsINTEGER(KIND=JPIM), INTENT(IN) :: KRESOL
TRANS_END
SUBROUTINE TRANS_END(CDMODE)
This subroutine terminates the transform package by releasing (i.e. deallocating and nullifying) all allocated arrays.
OPTIONAL, INTENT(IN)
argumentsCHARACTER*5, OPTIONAL, INTENT(IN) :: CDMODE
'FINAL'
to finalise'INTER'
just to deallocate N_REGIONS
and NPRCIDS
.'FINAL'
DIST_GRID_32
Note
This subroutine is deprecated and will be removed in a future release.
GATH_GRID_32
Note
This subroutine is deprecated and will be removed in a future release.