SCALE-RM
|
module atmosphere / adiabat More...
Functions/Subroutines | |
subroutine, public | atmos_adiabat_setup |
Setup. More... | |
subroutine | atmos_adiabat_cape_1d (KA, KS, KE, Kstr, DENS, TEMP, PRES, QV, QC, Qdry, Rtot, CPtot, CZ, FZ, CAPE, CIN, LCL, LFC, LNB, DENS_p, TEMP_p, BUOY_p, QV_p, converged) |
Calc CAPE and CIN Type of parcel method: Pseudo-adiabatic ascend from lowermost layer of the model Reference: Emanuel(1994) More... | |
subroutine | atmos_adiabat_cape_3d (KA, KS, KE, IA, IS, IE, JA, JS, JE, Kstr, DENS, TEMP, PRES, QV, QC, Qdry, Rtot, CPtot, CZ, FZ, CAPE, CIN, LCL, LFC, LNB, DENS_p, TEMP_p, BUOY_p, QV_p) |
subroutine | atmos_adiabat_liftparcel_1d (KA, KS, KE, Kstr, DENS, TEMP, PRES, QV, QC, QDRY, Rtot, CPtot, DENS_p1D, TEMP_p1D, QV_p1D, kLCL, converged) |
Calc temperature profile with lifting parcel Method: Pseudo-adiabatic ascend from lowermost layer of the model Reference: Emanuel(1994) More... | |
subroutine | atmos_adiabat_liftparcel_3d (KA, KS, KE, IA, IS, IE, JA, JS, JE, Kstr, DENS, TEMP, PRES, QV, QC, QDRY, Rtot, CPtot, DENS_p3D, TEMP_p3D, QV_p3D, kLCL, converged) |
module atmosphere / adiabat
name | description | unit | variable |
---|---|---|---|
BUOY_parcel | buoyancy profile in lifting parcel | m/s2 | P_BUOY |
DENS_parcel | density profile in lifting parcel | kg/m3 | P_DENS |
QV_parcel | humidity profile in lifting parcel | kg/kg | P_QV |
TEMP_parcel | temperature profile in lifting parcel | K | P_TEMP |
subroutine, public scale_atmos_adiabat::atmos_adiabat_setup |
Setup.
Definition at line 58 of file scale_atmos_adiabat.F90.
References scale_atmos_hydrometeor::atmos_hydrometeor_setup(), and scale_atmos_saturation::atmos_saturation_setup().
subroutine scale_atmos_adiabat::atmos_adiabat_cape_1d | ( | integer, intent(in) | KA, |
integer, intent(in) | KS, | ||
integer, intent(in) | KE, | ||
integer, intent(in) | Kstr, | ||
real(rp), dimension (ka), intent(in) | DENS, | ||
real(rp), dimension (ka), intent(in) | TEMP, | ||
real(rp), dimension (ka), intent(in) | PRES, | ||
real(rp), dimension (ka), intent(in) | QV, | ||
real(rp), dimension (ka), intent(in) | QC, | ||
real(rp), dimension (ka), intent(in) | Qdry, | ||
real(rp), dimension (ka), intent(in) | Rtot, | ||
real(rp), dimension(ka), intent(in) | CPtot, | ||
real(rp), dimension ( ka), intent(in) | CZ, | ||
real(rp), dimension (0:ka), intent(in) | FZ, | ||
real(rp), intent(out) | CAPE, | ||
real(rp), intent(out) | CIN, | ||
real(rp), intent(out) | LCL, | ||
real(rp), intent(out) | LFC, | ||
real(rp), intent(out) | LNB, | ||
real(rp), dimension(ka), intent(out) | DENS_p, | ||
real(rp), dimension(ka), intent(out) | TEMP_p, | ||
real(rp), dimension(ka), intent(out) | BUOY_p, | ||
real(rp), dimension (ka), intent(out) | QV_p, | ||
logical, intent(out) | converged | ||
) |
Calc CAPE and CIN Type of parcel method: Pseudo-adiabatic ascend from lowermost layer of the model Reference: Emanuel(1994)
Definition at line 83 of file scale_atmos_adiabat.F90.
References atmos_adiabat_liftparcel_1d(), and scale_const::const_grav.
Referenced by atmos_adiabat_cape_3d().
subroutine scale_atmos_adiabat::atmos_adiabat_cape_3d | ( | integer, intent(in) | KA, |
integer, intent(in) | KS, | ||
integer, intent(in) | KE, | ||
integer, intent(in) | IA, | ||
integer, intent(in) | IS, | ||
integer, intent(in) | IE, | ||
integer, intent(in) | JA, | ||
integer, intent(in) | JS, | ||
integer, intent(in) | JE, | ||
integer, intent(in) | Kstr, | ||
real(rp), dimension (ka,ia,ja), intent(in) | DENS, | ||
real(rp), dimension (ka,ia,ja), intent(in) | TEMP, | ||
real(rp), dimension (ka,ia,ja), intent(in) | PRES, | ||
real(rp), dimension (ka,ia,ja), intent(in) | QV, | ||
real(rp), dimension (ka,ia,ja), intent(in) | QC, | ||
real(rp), dimension (ka,ia,ja), intent(in) | Qdry, | ||
real(rp), dimension (ka,ia,ja), intent(in) | Rtot, | ||
real(rp), dimension(ka,ia,ja), intent(in) | CPtot, | ||
real(rp), dimension ( ka,ia,ja), intent(in) | CZ, | ||
real(rp), dimension (0:ka,ia,ja), intent(in) | FZ, | ||
real(rp), dimension(ia,ja), intent(out) | CAPE, | ||
real(rp), dimension (ia,ja), intent(out) | CIN, | ||
real(rp), dimension (ia,ja), intent(out) | LCL, | ||
real(rp), dimension (ia,ja), intent(out) | LFC, | ||
real(rp), dimension (ia,ja), intent(out) | LNB, | ||
real(rp), dimension(ka,ia,ja), intent(out), optional, target | DENS_p, | ||
real(rp), dimension(ka,ia,ja), intent(out), optional, target | TEMP_p, | ||
real(rp), dimension(ka,ia,ja), intent(out), optional, target | BUOY_p, | ||
real(rp), dimension (ka,ia,ja), intent(out), optional, target | QV_p | ||
) |
Definition at line 196 of file scale_atmos_adiabat.F90.
References atmos_adiabat_cape_1d(), and scale_prc::prc_abort().
subroutine scale_atmos_adiabat::atmos_adiabat_liftparcel_1d | ( | integer, intent(in) | KA, |
integer, intent(in) | KS, | ||
integer, intent(in) | KE, | ||
integer, intent(in) | Kstr, | ||
real(rp), dimension (ka), intent(in) | DENS, | ||
real(rp), dimension (ka), intent(in) | TEMP, | ||
real(rp), dimension (ka), intent(in) | PRES, | ||
real(rp), dimension (ka), intent(in) | QV, | ||
real(rp), dimension (ka), intent(in) | QC, | ||
real(rp), dimension (ka), intent(in) | QDRY, | ||
real(rp), dimension (ka), intent(in) | Rtot, | ||
real(rp), dimension(ka), intent(in) | CPtot, | ||
real(rp), dimension(ka), intent(out) | DENS_p1D, | ||
real(rp), dimension(ka), intent(out) | TEMP_p1D, | ||
real(rp), dimension (ka), intent(out) | QV_p1D, | ||
integer, intent(out) | kLCL, | ||
logical, intent(out) | converged | ||
) |
Calc temperature profile with lifting parcel Method: Pseudo-adiabatic ascend from lowermost layer of the model Reference: Emanuel(1994)
Definition at line 311 of file scale_atmos_adiabat.F90.
References scale_atmos_hydrometeor::cp_water.
Referenced by atmos_adiabat_cape_1d(), and atmos_adiabat_liftparcel_3d().
subroutine scale_atmos_adiabat::atmos_adiabat_liftparcel_3d | ( | integer, intent(in) | KA, |
integer, intent(in) | KS, | ||
integer, intent(in) | KE, | ||
integer, intent(in) | IA, | ||
integer, intent(in) | IS, | ||
integer, intent(in) | IE, | ||
integer, intent(in) | JA, | ||
integer, intent(in) | JS, | ||
integer, intent(in) | JE, | ||
integer, intent(in) | Kstr, | ||
real(rp), dimension (ka,ia,ja), intent(in) | DENS, | ||
real(rp), dimension (ka,ia,ja), intent(in) | TEMP, | ||
real(rp), dimension (ka,ia,ja), intent(in) | PRES, | ||
real(rp), dimension (ka,ia,ja), intent(in) | QV, | ||
real(rp), dimension (ka,ia,ja), intent(in) | QC, | ||
real(rp), dimension (ka,ia,ja), intent(in) | QDRY, | ||
real(rp), dimension (ka,ia,ja), intent(in) | Rtot, | ||
real(rp), dimension(ka,ia,ja), intent(in) | CPtot, | ||
real(rp), dimension(ka,ia,ja), intent(out) | DENS_p3D, | ||
real(rp), dimension(ka,ia,ja), intent(out) | TEMP_p3D, | ||
real(rp), dimension (ka,ia,ja), intent(out) | QV_p3D, | ||
integer, dimension(ia,ja), intent(out) | kLCL, | ||
logical, intent(out) | converged | ||
) |
Definition at line 403 of file scale_atmos_adiabat.F90.
References atmos_adiabat_liftparcel_1d(), and scale_prc::prc_abort().