UPP (develop)
|
pvetc() computes potential vorticity, etc. More...
Go to the source code of this file.
Functions/Subroutines | |
subroutine | mptgen (mpirank, mpisize, nd, jt1, jt2, j1, j2, jx, jm, jn) |
mptgen() generates grid decomposition dimensions. | |
subroutine | mptranr4 (mpicomm, mpisize, im, ida, idb, jm, jma, jmb, jda, km, kma, kmb, kdb, a, b, ta, tb) |
mptranr4() transposes grid decompositions. | |
subroutine | mxwind (km, p, u, v, t, h, pmw, umw, vmw, tmw, hmw) |
mxwind() computes maximum wind level fields. | |
subroutine | p2pv (km, pvu, h, t, p, u, v, kpv, pv, pvpt, pvpb, lpv, upv, vpv, hpv, tpv, ppv, spv) |
p2pv() interpolates to potential vorticity level. | |
subroutine | p2th (km, theta, u, v, h, t, pvu, sigma, rh, omga, kth, th, lth, uth, vth, hth, tth, zth, sigmath, rhth, oth) |
p2th() interpolates from pressure level to isentropic (theta) level. | |
subroutine | pvetc (km, p, px, py, t, tx, ty, h, u, v, av, hm, s, bvf2, pvn, theta, sigma, pvu) |
pvetc() computes potential vorticity, etc. | |
subroutine | rsearch1 (km1, z1, km2, z2, l2) |
rsearch1() searches for a surrounding real interval. | |
subroutine | tpause (km, p, u, v, t, h, ptp, utp, vtp, ttp, htp, shrtp) |
tpause() computes tropopause level fields. | |
pvetc() computes potential vorticity, etc.
This subprogram computes
computation | equation |
---|---|
Montgomery streamfunction | hm=cp*t+g*z |
Specific entropy | s=cp*log(t/t0)-r*log(p/p0) |
Brunt-Vaisala frequency squared | bvf2=g/cp*ds/dz |
Potential vorticity | pvn=(av*ds/dz-dv/dz*ds/dx+du/dz*ds/dy)/rho/cp |
Potential temperature | theta=t0*exp(s/cp) |
Static stability | sigma=t/g*bvf2 |
Potential vorticity in PV units | pvu=10**-6*theta*pvn |
[in] | km | integer number of levels. |
[in] | p | real (km) pressure (Pa). |
[in] | px | real (km) pressure x-gradient (Pa/m). |
[in] | py | real (km) pressure y-gradient (Pa/m). |
[in] | t | real (km) (virtual) temperature (K). |
[in] | tx | real (km) (virtual) temperature x-gradient (K/m). |
[in] | ty | real (km) (virtual) temperature y-gradient (K/m). |
[in] | h | real (km) height (m). |
[in] | u | real (km) x-component wind (m/s). |
[in] | v | real (km) y-component wind (m/s). |
[in] | av | real (km) absolute vorticity (1/s). |
[out] | hm | real (km) Montgomery streamfunction (m**2/s**2). |
[out] | s | real (km) specific entropy (J/K/kg). |
[out] | bvf2 | real (km) Brunt-Vaisala frequency squared (1/s**2). |
[out] | pvn | real (km) potential vorticity (m**2/kg/s). |
[out] | theta | real (km) (virtual) potential temperature (K). |
[out] | sigma | real (km) static stability (K/m). |
[out] | pvu | real (km) potential vorticity (10**-6*K*m**2/kg/s). |
Date | Programmer | Comments |
---|---|---|
1999-10-18 | Mark Iredell | Initial |
Definition in file GFSPOST.F.
subroutine mptgen | ( | integer(kint_mpi), intent(in) | mpirank, |
integer(kint_mpi), intent(in) | mpisize, | ||
integer(kint_mpi), intent(in) | nd, | ||
integer(kint_mpi), dimension(nd), intent(in) | jt1, | ||
integer(kint_mpi), dimension(nd), intent(in) | jt2, | ||
integer(kint_mpi), dimension(nd), intent(out) | j1, | ||
integer(kint_mpi), dimension(nd), intent(out) | j2, | ||
integer(kint_mpi), dimension(nd), intent(out) | jx, | ||
integer(kint_mpi), dimension(nd), intent(out) | jm, | ||
integer(kint_mpi), dimension(nd), intent(out) | jn | ||
) |
mptgen() generates grid decomposition dimensions.
This subprogram decomposes total dimensions of a problem into smaller domains to be managed on a distributed memory system. The last dimension given is decomposed first. If more decompositions are possible, the next to last dimension is decomposed next, and so on. The transpositions between decompositions should be done by mptran*.
[in] | mpirank | integer(kint_mpi) rank of the process (from mpi_comm_rank). |
[in] | mpisize | integer(kint_mpi) size of the process (from mpi_comm_size). |
[in] | nd | integer(kint_mpi) number of dimensions to decompose. |
[in] | jt1 | integer(kint_mpi) (nd) lower bounds of total dimensions. |
[in] | jt2 | integer(kint_mpi) (nd) upper bounds of total dimensions. |
[out] | j1 | integer(kint_mpi) (nd) lower bounds of local decompositions. |
[out] | j2 | integer(kint_mpi) (nd) upper bounds of local decompositions. |
[out] | jx | integer(kint_mpi) (nd) local size of decompositions. |
[out] | jm | integer(kint_mpi) (nd) maximum size of decompositions. |
[out] | jn | integer(kint_mpi) (nd) number of decompositions. |
Date | Programmer | Comments |
---|---|---|
1999-02-12 | Mark Iredell | Initial |
subroutine mptranr4 | ( | integer(kint_mpi), intent(in) | mpicomm, |
integer(kint_mpi), intent(in) | mpisize, | ||
integer(kint_mpi), intent(in) | im, | ||
integer(kint_mpi), intent(in) | ida, | ||
integer(kint_mpi), intent(in) | idb, | ||
integer(kint_mpi), intent(in) | jm, | ||
integer(kint_mpi), intent(in) | jma, | ||
integer(kint_mpi), intent(in) | jmb, | ||
integer(kint_mpi), intent(in) | jda, | ||
integer(kint_mpi), intent(in) | km, | ||
integer(kint_mpi), intent(in) | kma, | ||
integer(kint_mpi), intent(in) | kmb, | ||
integer(kint_mpi), intent(in) | kdb, | ||
real(4), dimension(ida,jda,kma), intent(in) | a, | ||
real(4), dimension(idb,kdb,jmb), intent(out) | b, | ||
real(4), dimension(im,jm,km,mpisize), intent(inout) | ta, | ||
real(4), dimension(im,jm,km,mpisize), intent(inout) | tb | ||
) |
mptranr4() transposes grid decompositions.
This subprogram transposes an array of data from one grid decomposition to another by using message passing. The grid decompositions should be generated by mptgen.
[in] | mpicomm | integer(kint_mpi) mpi communicator. |
[in] | mpisize | integer(kint_mpi) size of the process (from mpi_comm_size). |
[in] | im | integer(kint_mpi) undecomposed range. |
[in] | ida | integer(kint_mpi) undecomposed input dimension. |
[in] | idb | integer(kint_mpi) undecomposed output dimension. |
[in] | jm | integer(kint_mpi) output grid decomposition size. |
[in] | jma | integer(kint_mpi) input grid undecomposed range. |
[in] | jmb | integer(kint_mpi) output grid decomposed range. |
[in] | jda | integer(kint_mpi) input grid undecomposed dimension. |
[in] | km | integer(kint_mpi) input grid decomposition size. |
[in] | kma | integer(kint_mpi) input grid decomposed range. |
[in] | kmb | integer(kint_mpi) output grid undecomposed range. |
[in] | kdb | integer(kint_mpi) output grid undecomposed dimension. |
[in] | a | real(4) (ida,jda,kma) input array. |
[out] | b | real(4) (idb,kdb,jmb) output array. |
[out] | ta,tb | real(4) (im,jm,km,mpisize) work arrays. |
It does not do any communication hiding.
This subprogram must be used rather than mpi_alltoall in any of the following cases:
The output grid range is not the full extent (either kmb<mpisize or kmb<kda or jma<mpisize or jma<jda)
If none of these conditions apply, mpi_alltoall could be used directly rather than this subprogram and would be more efficient.
include 'mpif.h' ! use mpi parameter(jt=1000,kt=10000) ! set problem size real,allocatable:: a(:,:),b(:,:) ! declare arrays call mpi_init(ierr) ! initialize mpi call mpi_comm_rank(MPI_COMM_WORLD,mpirank,ierr) ! get mpi rank call mpi_comm_size(MPI_COMM_WORLD,mpisize,ierr) ! get mpi size call mptgen(mpirank,mpisize,1,1,jt,j1,j2,jx,jm,jn) ! decompose output call mptgen(mpirank,mpisize,1,1,kt,k1,k2,kx,km,kn) ! decompose input allocate(a(jt,k1:k2),b(kt,j1:j2)) ! allocate arrays a=reshape((/((j+k,j=1,jt),k=k1,k2)/),(/jt,k2-k1+1/)) ! initialize input call mptranr4(MPI_COMM_WORLD,mpisize,1,1,1, ! transpose arrays & jm,jt,j2-j1+1,jt,km,k2-k1+1,kt,kt,a,b) print '(2i8,f16.1)',((k,j,b(k,j),k=2000,kt,2000), ! print some values & j=((j1-1)/200+1)*200,j2,200) call mpi_finalize(ierr) ! finalize mpi endThis transpose took 0.6 seconds on 4 2-way winterhawk nodes.
Date | Programmer | Comments |
---|---|---|
1999-02-12 | Mark Iredell | Initial |
subroutine mxwind | ( | integer, intent(in) | km, |
real, dimension(km), intent(in) | p, | ||
real, dimension(km), intent(in) | u, | ||
real, dimension(km), intent(in) | v, | ||
real, dimension(km), intent(in) | t, | ||
real, dimension(km), intent(in) | h, | ||
real, intent(out) | pmw, | ||
real, intent(out) | umw, | ||
real, intent(out) | vmw, | ||
real, intent(out) | tmw, | ||
real, intent(out) | hmw | ||
) |
mxwind() computes maximum wind level fields.
This subprogram finds the maximum wind level and computes fields at the maximum wind level. The maximum wind level is searched for between 500 mb and 100 mb. The height and wind speed at the maximum wind speed level is calculated by assuming the wind speed varies quadratically in height in the neighborhood of the maximum wind level. The other fields are interpolated linearly in height to the maximum wind level. The maximum wind level pressure is found hydrostatically.
[in] | km | integer number of levels. |
[in] | p | real (km) pressure (Pa). |
[in] | u | real (km) x-component wind (m/s). |
[in] | v | real (km) y-component wind (m/s). |
[in] | t | real (km) temperature (K). |
[in] | h | real (km) height (m). |
[out] | pmw | real maximum wind level pressure (Pa). |
[out] | umw | real maximum wind level x-component wind (m/s). |
[out] | vmw | real maximum wind level y-component wind (m/s). |
[out] | tmw | maximum wind level temperature (K). |
[out] | hmw | real maximum wind level height (m). |
Date | Programmer | Comments |
---|---|---|
1999-10-18 | Mark Iredell | Initial |
2005-02-02 | Mark Iredell | Changed upper limit to 100 mb |
Definition at line 459 of file GFSPOST.F.
References rsearch1().
Referenced by miscln().
subroutine p2pv | ( | integer, intent(in) | km, |
real, dimension(km), intent(in) | pvu, | ||
real, dimension(km), intent(in) | h, | ||
real, dimension(km), intent(in) | t, | ||
real, dimension(km), intent(in) | p, | ||
real, dimension(km), intent(in) | u, | ||
real, dimension(km), intent(in) | v, | ||
integer, intent(in) | kpv, | ||
real, dimension(kpv), intent(in) | pv, | ||
real, dimension(kpv), intent(in) | pvpt, | ||
real, dimension(kpv), intent(in) | pvpb, | ||
logical*1, dimension(kpv), intent(out) | lpv, | ||
real, dimension(kpv), intent(out) | upv, | ||
real, dimension(kpv), intent(out) | vpv, | ||
real, dimension(kpv), intent(out) | hpv, | ||
real, dimension(kpv), intent(out) | tpv, | ||
real, dimension(kpv), intent(out) | ppv, | ||
real, dimension(kpv), intent(out) | spv | ||
) |
p2pv() interpolates to potential vorticity level.
This subprogram interpolates fields to given potential vorticity levels within given pressure limits. The output level is the first encountered from the top pressure limit. If the given potential vorticity level is not found, the outputs are zero and the bitmap is false. The interpolation is linear in potential vorticity.
[in] | km | integer number of levels. |
[in] | pvu | real (km) potential vorticity in PV units (10**-6*K*m**2/kg/s). |
[in] | h | real (km) height (m). |
[in] | t | real (km) temperature (K). |
[in] | p | real (km) pressure (Pa). |
[in] | u | real (km) x-component wind (m/s). |
[in] | v | real (km) y-component wind (m/s). |
[in] | kpv | integer number of potential vorticity levels. |
[in] | pv | real (kpv) potential vorticity levels (10**-6*K*m**2/kg/s). |
[in] | pvpt | real (kpv) top pressures for PV search (Pa). |
[in] | pvpb | real (kpv) bottom pressures for PV search (Pa). |
[out] | lpv | logical*1 (kpv) bitmap. |
[out] | upv | real (kpv) x-component wind (m/s). |
[out] | vpv | real (kpv) y-component wind (m/s). |
[out] | hpv | real (kpv) temperature (K). |
[out] | tpv | real (kpv) temperature (K). |
[out] | ppv | real (kpv) pressure (Pa). |
[out] | spv | real (kpv) wind speed shear (1/s). |
Date | Programmer | Comments |
---|---|---|
1999-10-18 | Mark Iredell | Initial |
2021-08-31 | Hui-ya Chuang | Increase depth criteria for identifying PV layer from 25 to 50 to avoid finding shallow high level PV layer in high latitudes |
Definition at line 203 of file GFSPOST.F.
References rsearch1().
Referenced by mdl2thandpv().
subroutine p2th | ( | integer, intent(in) | km, |
real, dimension(km), intent(in) | theta, | ||
real, dimension(km), intent(in) | u, | ||
real, dimension(km), intent(in) | v, | ||
real, dimension(km), intent(in) | h, | ||
real, dimension(km), intent(in) | t, | ||
real, dimension(km), intent(in) | pvu, | ||
real, dimension(km), intent(in) | sigma, | ||
real, dimension(km), intent(in) | rh, | ||
real, dimension(km), intent(in) | omga, | ||
integer, intent(in) | kth, | ||
real, dimension(kth), intent(in) | th, | ||
logical*1, dimension(kth), intent(out) | lth, | ||
real, dimension(kth), intent(out) | uth, | ||
real, dimension(kth), intent(out) | vth, | ||
real, dimension(kth), intent(out) | hth, | ||
real, dimension(kth), intent(out) | tth, | ||
real, dimension(kth), intent(out) | zth, | ||
real, dimension(kth), intent(out) | sigmath, | ||
real, dimension(kth), intent(out) | rhth, | ||
real, dimension(kth), intent(out) | oth | ||
) |
p2th() interpolates from pressure level to isentropic (theta) level.
This subprogram interpolates fields to given isentropic levels. The interpolation is linear in entropy. Outside the domain the bitmap is set to false.
[in] | km | integer number of levels. |
[in] | theta | real (km) potential temperature (K). |
[in] | u | real (km) x-component wind (m/s). |
[in] | v | real (km) y-component wind (m/s). |
[in] | h | real (km) height (m). |
[in] | t | real (km) temperature (K). |
[in] | pvu | real (km) potential vorticity in PV units (10**-6*K*m**2/kg/s). |
[in] | sigma | real (km) static stability (K/m). |
[in] | rh | real (km) relative humidity. |
[in] | omga | real (km) vertical velocity in pressure coordinates (Pa/s). |
[in] | kth | integer number of isentropic levels. |
[in] | th | real (kth) isentropic levels (K). |
[in] | lth | logical*1 (kth) bitmap. |
[out] | uth | real (kth) x-component wind (m/s). |
[out] | vth | real (kth) y-component wind (m/s). |
[out] | hth | real (kth) height (m). |
[out] | tth | real (kth) temperature (K). |
[out] | zth | real (kth) potential vorticity in PV units (10**-6*K*m**2/kg/s). |
[out] | sigmath | (kth) real static stability (K/m). |
[out] | rhth | (kth) real relative humidity. |
[out] | oth | (kth) real vertical velocity in pressure coordinates (Pa/s). |
Date | Programmer | Comments |
---|---|---|
1999-10-18 | Mark Iredell | Initial |
Definition at line 137 of file GFSPOST.F.
References rsearch1().
Referenced by mdl2thandpv().
subroutine pvetc | ( | integer, intent(in) | km, |
real, dimension(km), intent(in) | p, | ||
real, dimension(km), intent(in) | px, | ||
real, dimension(km), intent(in) | py, | ||
real, dimension(km), intent(in) | t, | ||
real, dimension(km), intent(in) | tx, | ||
real, dimension(km), intent(in) | ty, | ||
real, dimension(km), intent(in) | h, | ||
real, dimension(km), intent(in) | u, | ||
real, dimension(km), intent(in) | v, | ||
real, dimension(km), intent(in) | av, | ||
real, dimension(km), intent(out) | hm, | ||
real, dimension(km), intent(out) | s, | ||
real, dimension(km), intent(out) | bvf2, | ||
real, dimension(km), intent(out) | pvn, | ||
real, dimension(km), intent(out) | theta, | ||
real, dimension(km), intent(out) | sigma, | ||
real, dimension(km), intent(out) | pvu | ||
) |
pvetc() computes potential vorticity, etc.
[in] | km | integer number of levels. |
[in] | p | real (km) pressure (Pa). |
[in] | px | real (km) pressure x-gradient (Pa/m). |
[in] | py | real (km) pressure y-gradient (Pa/m). |
[in] | t | real (km) (virtual) temperature (K). |
[in] | tx | real (km) (virtual) temperature x-gradient (K/m). |
[in] | ty | real (km) (virtual) temperature y-gradient (K/m). |
[in] | h | real (km) height (m). |
[in] | u | real (km) x-component wind (m/s). |
[in] | v | real (km) y-component wind (m/s). |
[in] | av | real (km) absolute vorticity (1/s). |
[out] | hm | real (km) Montgomery streamfunction (m**2/s**2). |
[out] | s | real (km) specific entropy (J/K/kg). |
[out] | bvf2 | real (km) Brunt-Vaisala frequency squared (1/s**2). |
[out] | pvn | real (km) potential vorticity (m**2/kg/s). |
[out] | theta | real (km) (virtual) potential temperature (K). |
[out] | sigma | real (km) static stability (K/m). |
[out] | pvu | real (km) potential vorticity (10**-6*K*m**2/kg/s). |
Definition at line 62 of file GFSPOST.F.
References rsearch1().
Referenced by mdl2thandpv().
subroutine rsearch1 | ( | integer, intent(in) | km1, |
real, dimension(km1), intent(in) | z1, | ||
integer, intent(in) | km2, | ||
real, dimension(km2), intent(in) | z2, | ||
integer, dimension(km2), intent(out) | l2 | ||
) |
rsearch1() searches for a surrounding real interval.
This subprogram searches a monotonic sequences of real numbers for intervals that surround a given search set of real numbers. the sequences may be monotonic in either direction; the real numbers may be single or double precision.
[in] | km1 | integer number of points in the sequence. |
[in] | z1 | real (km1) sequence values to search. (z1 must be monotonic in either direction) |
[in] | km2 | integer number of points to search for. |
[in] | z2 | real (km2) set of values to search for. (z2 need not be monotonic) |
[out] | l2 | integer (km2) interval locations from 0 to km1. (z2 will be between z1(l2) and z1(l2+1)) |
Date | Programmer | Comments |
---|---|---|
1998-05-01 | Mark Iredell | Initial |
Definition at line 297 of file GFSPOST.F.
Referenced by mxwind(), p2pv(), p2th(), pvetc(), and tpause().
subroutine tpause | ( | integer, intent(in) | km, |
real, dimension(km), intent(in) | p, | ||
real, dimension(km), intent(in) | u, | ||
real, dimension(km), intent(in) | v, | ||
real, dimension(km), intent(in) | t, | ||
real, dimension(km), intent(in) | h, | ||
real, intent(out) | ptp, | ||
real, intent(out) | utp, | ||
real, intent(out) | vtp, | ||
real, intent(out) | ttp, | ||
real, intent(out) | htp, | ||
real, intent(out) | shrtp | ||
) |
tpause() computes tropopause level fields.
This subprogram finds the tropopause level and computes fields at the tropopause level. The tropopause is defined as the lowest level above 500 mb which has a temperature lapse rate of less than 2 K/km. The lapse rate must average less than 2 K/km over a 2 km depth. If no such level is found below 50 mb, the tropopause is set to 50 mb. The tropopause fields are interpolated linearly in lapse rate. The tropopause pressure is found hydrostatically. The tropopause wind shear is computed as the partial derivative of wind speed with respect to height at the tropopause level.
[in] | km | integer number of levels. |
[in] | p | real (km) pressure (Pa). |
[in] | u | real (km) x-component wind (m/s). |
[in] | v | real (km) y-component wind (m/s). |
[in] | t | real (km) temperature (K). |
[in] | h | real (km) height (m). |
[out] | ptp | real tropopause pressure (Pa). |
[out] | utp | real tropopause x-component wind (m/s). |
[out] | vtp | real tropopause y-component wind (m/s). |
[out] | ttp | real tropopause temperature (K). |
[out] | htp | real tropopause height (m). |
[out] | shrtp | real tropopause wind shear (1/s). |
Date | Programmer | Comments |
---|---|---|
1999-10-18 | Mark Iredell | Initial |
Definition at line 372 of file GFSPOST.F.
References rsearch1().
Referenced by miscln().