UPP  001
 All Data Structures Files Functions Pages
upp_physics Module Reference

calcape() computes CAPE/CINS and other storm related variables. More...

Public Member Functions

subroutine, public CALCAPE (ITYPE, DPBND, P1D, T1D, Q1D, L1D, CAPE, CINS, PPARC, ZEQL, THUND)
 calcape() computes CAPE and CINS. More...
 
subroutine, public CALCAPE2 (ITYPE, DPBND, P1D, T1D, Q1D, L1D, CAPE, CINS, LFC, ESRHL, ESRHH, DCAPE, DGLD, ESP)
 calcape2() computes CAPE and CINS. More...
 
subroutine, public CALRH (P1, T1, Q1, RH)
 
subroutine, public CALRH_GFS (P1, T1, Q1, RH)
 calrh_gfs() computes relative humidity. More...
 
subroutine, public CALRH_GSD (P1, T1, Q1, RHB)
 
subroutine, public CALRH_NAM (P1, T1, Q1, RH)
 calrh_nam() computes relative humidity. More...
 
subroutine, public CALRH_PW (RHPW)
 
elemental real function, public fpvsnew (t)
 
elemental real function, public TVIRTUAL (T, Q)
 

Detailed Description

calcape() computes CAPE/CINS and other storm related variables.

calcape2() computes additional storm related variables.

calrh(), calrh_nam(), calrh_gfs(), calrh_gsd() compute RH using various algorithms.

The NAM v4.1.18 algorithm (calrh_nam()) is selected as default for NMMB and FV3GFS, FV3GEFS, and FV3R for the UPP 2020 unification.

calrh_pw() algorithm use at GSD for RUC and Rapid Refresh.

fpvsnew() computes saturation vapor pressure.

tvirtual() computes virtual temperature.

Program history log:

Date Programmer Comments
2020-05-20 Jesse Meng Initial
Author
Jesse Meng
Date
2020-05-20

Definition at line 27 of file UPP_PHYSICS.f.

Member Function/Subroutine Documentation

subroutine, public upp_physics::CALCAPE ( integer, intent(in)  ITYPE,
real, intent(in)  DPBND,
real, dimension(im,jsta:jend), intent(in)  P1D,
real, dimension(im,jsta:jend), intent(in)  T1D,
real, dimension(im,jsta:jend), intent(inout)  Q1D,
integer, dimension(im,jsta:jend), intent(in)  L1D,
real, dimension(im,jsta:jend), intent(inout)  CAPE,
real, dimension(im,jsta:jend), intent(inout)  CINS,
real, dimension(im,jsta:jend), intent(inout)  PPARC,
real, dimension(im,jsta:jend), intent(inout)  ZEQL,
real, dimension(im,jsta:jend)  THUND 
)

calcape() computes CAPE and CINS.

This routine computes CAPE and CINS given temperature, pressure, and specific humidty. In "storm and cloud dynamics" (1989, academic press) cotton and anthes define CAPE (equation 9.16, p501) as

el
cape = sum g * ln(thetap/thetaa) dz
lcl
Where,
el = equilibrium level,
lcl = lifting condenstation level,
g = gravitational acceleration,
thetap = lifted parcel potential temperature,
thetaa = ambient potential temperature.
Note that the integrand ln(THETAP/THETAA) approximately
equals (THETAP-THETAA)/THETAA.  This ratio is often used
in the definition of CAPE/CINS.

Two types of CAPE/CINS can be computed by this routine.  The
summation process is the same For both cases.  What differs
is the definition of the parcel to lift.  FOR ITYPE=1 the
parcel with the warmest THETA-E in A DPBND pascal layer above
the model surface is lifted.  the arrays P1D, T1D, and Q1D
are not used.  For itype=2 the arrays P1D, T1D, and Q1D
define the parcel to lift in each column.  Both types of
CAPE/CINS may be computed in a single execution of the post
processor.

This algorithm proceeds as follows.
For each column, 
   (1)  Initialize running CAPE and CINS SUM TO 0.0
   (2)  Compute temperature and pressure at the LCL using
        look up table (PTBL).  Use either parcel that gives
        max THETAE in lowest DPBND above ground (ITYPE=1)
        or given parcel from t1D,Q1D,...(ITYPE=2).
   (3)  Compute the temp of a parcel lifted from the LCL.
        We know that the parcel's
        equivalent potential temperature (THESP) remains
        constant through this process.  we can
        compute tpar using this knowledge using look
        up table (subroutine TTBLEX).
   (4)  Find the equilibrium level.  This is defined as the
        highest positively buoyant layer.
        (If there is no positively buoyant layer, CAPE/CINS
         will be zero)
   (5)  Compute CAPE/CINS.  
        (A) Compute THETAP.  We know TPAR and P.
        (B) Compute THETAA.  We know T and P.  
   (6)  Add G*(THETAP-THETAA)*DZ to the running CAPE or CINS sum.
        (A) If THETAP > THETAA, add to the CAPE sum.
        (B) If THETAP < THETAA, add to the CINS sum.
   (7)  Are we at equilibrium level? 
        (A) If yes, stop the summation.
        (b) if no, contiunue the summation.
   (8)  Enforce limits on CAPE and CINS (i.e. no negative CAPE)
Parameters
[in]ITYPEINTEGER Flag specifying how parcel to lift is identified. See comments above.
[in]DPBNDDepth over which one searches for most unstable parcel.
[in]P1DArray of pressure of parcels to lift.
[in]T1DArray of temperature of parcels to lift.
[in]Q1DArray of specific humidity of parcels to lift.
[in]L1DArray of model level of parcels to lift.
[out]CAPEConvective available potential energy (J/kg).
[out]CINSConvective inhibition (J/kg).
[out]PPARCPressure level of parcel lifted when one searches over a particular depth to compute CAPE/CIN.

Program history log:

Date Programmer Comments
1993-02-10 Russ Treadon Initial
1993-06-19 Russ Treadon Generalized routine to allow for type 2 CAPE/CINS calculations
1994-09-23 Mike Baldwin Modified to use look up tables instead of complicated equations
1994-10-13 Mike Baldwin Modified to continue CAPE/CINS calc up to at highest buoyant layer
1998-06-12 T Black Conversion from 1-D TO 2-D
1998-08-18 T Black Compute APE internally
2000-01-04 Jim Tuccillo MPI Version
2002-01-15 Mike Baldwin WRF Version
2003-08-24 G Manikin Added level of parcel being lifted as output from the routine and added the depth over which one searches for the most unstable parcel as input
2010-09-09 G Manikin Changed computation to use virtual temp added eq lvl hght and thunder parameter
2015-??-?? S Moorthi Optimization and threading
2021-07-28 W Meng Restrict computation from undefined grids
2021-09-01 E Colon Equivalent level height index for RTMA
Author
Russ Treadon W/NP2
Date
1993-02-10

Definition at line 523 of file UPP_PHYSICS.f.

References fpvsnew().

subroutine, public upp_physics::CALCAPE2 ( integer, intent(in)  ITYPE,
real, intent(in)  DPBND,
real, dimension(im,jsta:jend), intent(in)  P1D,
real, dimension(im,jsta:jend), intent(in)  T1D,
real, dimension(im,jsta:jend), intent(inout)  Q1D,
integer, dimension(im,jsta:jend), intent(in)  L1D,
real, dimension(im,jsta:jend), intent(inout)  CAPE,
real, dimension(im,jsta:jend), intent(inout)  CINS,
real, dimension(im,jsta:jend), intent(inout)  LFC,
real, dimension(im,jsta:jend), intent(inout)  ESRHL,
real, dimension(im,jsta:jend), intent(inout)  ESRHH,
real, dimension(im,jsta:jend), intent(inout)  DCAPE,
real, dimension(im,jsta:jend), intent(inout)  DGLD,
real, dimension(im,jsta:jend), intent(inout)  ESP 
)

calcape2() computes CAPE and CINS.

This routine computes CAPE and CINS given temperature, pressure, and specific humidty. In "storm and cloud dynamics" (1989, academic press) cotton and anthes define CAPE (equation 9.16, p501) as

el
cape = sum g * ln(thetap/thetaa) dz
lcl
Where,
el = equilibrium level,
lcl = lifting condenstation level,
g = gravitational acceleration,
thetap = lifted parcel potential temperature,
thetaa = ambient potential temperature.
Note that the integrand ln(THETAP/THETAA) approximately
equals (THETAP-THETAA)/THETAA.  This ratio is often used
in the definition of CAPE/CINS.

Two types of CAPE/CINS can be computed by this routine.  The
summation process is the same For both cases.  What differs
is the definition of the parcel to lift.  FOR ITYPE=1 the
parcel with the warmest THETA-E in A DPBND pascal layer above
the model surface is lifted.  the arrays P1D, T1D, and Q1D
are not used.  For itype=2 the arrays P1D, T1D, and Q1D
define the parcel to lift in each column.  Both types of
CAPE/CINS may be computed in a single execution of the post
processor.

This algorithm proceeds as follows.
For each column, 
   (1)  Initialize running CAPE and CINS SUM TO 0.0
   (2)  Compute temperature and pressure at the LCL using
        look up table (PTBL).  Use either parcel that gives
        max THETAE in lowest DPBND above ground (ITYPE=1)
        or given parcel from t1D,Q1D,...(ITYPE=2).
   (3)  Compute the temp of a parcel lifted from the LCL.
        We know that the parcel's
        equivalent potential temperature (THESP) remains
        constant through this process.  we can
        compute tpar using this knowledge using look
        up table (subroutine TTBLEX).
   (4)  Find the equilibrium level.  This is defined as the
        highest positively buoyant layer.
        (If there is no positively buoyant layer, CAPE/CINS
         will be zero)
   (5)  Compute CAPE/CINS.  
        (A) Compute THETAP.  We know TPAR and P.
        (B) Compute THETAA.  We know T and P.  
   (6)  Add G*(THETAP-THETAA)*DZ to the running CAPE or CINS sum.
        (A) If THETAP > THETAA, add to the CAPE sum.
        (B) If THETAP < THETAA, add to the CINS sum.
   (7)  Are we at equilibrium level? 
        (A) If yes, stop the summation.
        (b) if no, contiunue the summation.
   (8)  Enforce limits on CAPE and CINS (i.e. no negative CAPE)
Parameters
[in]ITYPEINTEGER Flag specifying how parcel to lift is identified. See comments above.
[in]DPBNDDepth over which one searches for most unstable parcel.
[in]P1DArray of pressure of parcels to lift.
[in]T1DArray of temperature of parcels to lift.
[in]Q1DArray of specific humidity of parcels to lift.
[in]L1DArray of model level of parcels to lift.
[out]CAPEConvective available potential energy (J/kg).
[out]CINSConvective inhibition (J/kg).
[out]LFClevel of free convection (m).
[out]ESRHLLower bound to account for effective helicity calculation.
[out]ESRHHUpper bound to account for effective helicity calculation.
[out]DCAPEdowndraft CAPE (J/KG).
[out]DGLDDendritic growth layer depth (m).
[out]ESPEnhanced stretching potential.

Program history log:

Date Programmer Comments
1993-02-10 Russ Treadon Initial
1993-06-19 Russ Treadon Generalized routine to allow for type 2 CAPE/CINS calculations
1994-09-23 Mike Baldwin Modified to use look up tables instead of complicated equations
1994-10-13 Mike Baldwin Modified to continue CAPE/CINS calc up to at highest buoyant layer
1998-06-12 T Black Conversion from 1-D TO 2-D
1998-08-18 T Black Compute APE internally
2000-01-04 Jim Tuccillo MPI Version
2002-01-15 Mike Baldwin WRF Version
2003-08-24 G Manikin Added level of parcel being lifted as output from the routine and added the depth over which one searches for the most unstable parcel as input
2010-09-09 G Manikin Changed computation to use virtual temp added eq lvl hght and thunder parameter
2015-??-?? S Moorthi Optimization and threading
2021-09-03 J Meng Modified to add 0-3km CAPE/CINS, LFC, effective helicity, downdraft CAPE, dendritic growth layer depth, ESP
2021-09-01 E Colon Equivalent level height index for RTMA
Author
Russ Treadon W/NP2
Date
1993-02-10

Definition at line 998 of file UPP_PHYSICS.f.

References fpvsnew().

subroutine, public upp_physics::CALRH_GFS ( real, dimension(im,jsta:jend), intent(in)  P1,
real, dimension(im,jsta:jend), intent(in)  T1,
real, dimension(im,jsta:jend), intent(inout)  Q1,
real, dimension(im,jsta:jend), intent(inout)  RH 
)

calrh_gfs() computes relative humidity.

This routine computes relative humidity given pressure, temperature, specific humidity. an upper and lower bound of 100 and 1 percent relative humidity is enforced. When these bounds are applied the passed specific humidity array is adjusted as necessary to produce the set relative humidity.

Parameters
[in]P1Pressure (pa)
[in]T1Temperature (K)
[in]Q1Specific humidity (kg/kg)
[out]RHRelative humidity (decimal form)
[out]Q1Specific humidity (kg/kg)

Program History Log

Date Programmer Comments
????-??-?? DENNIS DEAVEN Initial
1992-12-22 Russ Treadon Modified as described above
1998-06-08 T Black Conversion from 1-D to 2-D
1998-08-18 Mike Baldwin Modify to compute RH over ice as in model
1998-12-16 Geoff Manikin undo RH computation over ice
2000-01-04 Jim Tuccillo MPI Version
2002-06-11 Mike Baldwin WRF Version
2013-08-13 S. Moorthi Threading
2006-03-19 Wen Meng Modify top pressure to 1 pa
Author
Russ Treadon W/NP2
Date
1992-12-22

Definition at line 171 of file UPP_PHYSICS.f.

References fpvsnew().

subroutine, public upp_physics::CALRH_NAM ( real, dimension(im,jsta:jend), intent(in)  P1,
real, dimension(im,jsta:jend), intent(in)  T1,
real, dimension(im,jsta:jend), intent(inout)  Q1,
real, dimension(im,jsta:jend), intent(out)  RH 
)

calrh_nam() computes relative humidity.

This routine computes relative humidity given pressure, temperature, specific humidity. an upper and lower bound of 100 and 1 percent relative humidity is enforced. When these bounds are applied the passed specific humidity array is adjusted as necessary to produce the set relative humidity.

Parameters
[in]P1Pressure (pa)
[in]T1Temperature (K)
[in]Q1Specific humidity (kg/kg)
[out]RHRelative humidity (decimal form)
[out]Q1Specific humidity (kg/kg)

Program History Log

Date Programmer Comments
????-??-?? DENNIS DEAVEN Initial
1992-12-22 Russ Treadon Modified as described above
1998-06-08 T Black Conversion from 1-D to 2-D
1998-08-18 Mike Baldwin Modify to compute RH over ice as in model
1998-12-16 Geoff Manikin undo RH computation over ice
2000-01-04 Jim Tuccillo MPI Version
2002-06-11 Mike Baldwin WRF Version
2006-03-19 Wen Meng Modify top pressure to 1 pa
Author
Russ Treadon W/NP2
Date
1992-12-22

Definition at line 91 of file UPP_PHYSICS.f.

elemental real function, public upp_physics::fpvsnew ( real, intent(in)  t)

fpvsnew() computes saturation vapor pressure.

Compute saturation vapor pressure from the temperature. A linear interpolation is done between values in a lookup table computed in gpvs. See documentation for fpvsx for details. Input values outside table range are reset to table extrema. The interpolation accuracy is almost 6 decimal places. On the Cray, fpvs is about 4 times faster than exact calculation. This function should be expanded inline in the calling routine.

Parameters
[in]tReal(krealfp) Temperature in Kelvin.
[out]fpvsnewReal(krealfp) Saturation vapor pressure in Pascals.

Program history log:

Date Programmer Comments
1991-05-07 Iredell Initial. Made into inlinable function
1994-12-30 Iredell Expand table
1999-03-01 Iredell F90 module
2001-02-26 Iredell Ice phase
Author
N Phillips w/NMC2X2
Date
1982-12-30

Definition at line 341 of file UPP_PHYSICS.f.

Referenced by BNDLYR(), CALCAPE(), CALCAPE2(), CALRH_GFS(), FDLVL_MASS(), INITPOST_GFS_NEMS_MPIIO(), INITPOST_NETCDF(), MISCLN(), and SURFCE().


The documentation for this module was generated from the following file: