module atmosphere / saturation
More...
|
subroutine, public | atmos_saturation_setup |
| Setup. More...
|
|
subroutine | atmos_saturation_alpha_0d (temp, alpha) |
| calc liquid/ice separation factor (0D) More...
|
|
subroutine | atmos_saturation_dqs_dtem_dens_liq_1d (KA, KS, KE, temp, dens, dqsdtem) |
|
subroutine | atmos_saturation_dqs_dtem_dens_ice_1d (KA, KS, KE, temp, dens, dqsdtem) |
|
subroutine | atmos_saturation_dqs_dtem_dpre_liq_1d (KA, KS, KE, temp, pres, qdry, dqsat_dT, dqsat_dP) |
|
subroutine | atmos_saturation_dqs_dtem_dpre_ice_1d (KA, KS, KE, temp, pres, qdry, dqsat_dT, dqsat_dP) |
|
subroutine | atmos_saturation_tdew_liq_1d (KA, KS, KE, DENS, TEMP, QV, Tdew, converged) |
|
subroutine | atmos_saturation_pote_0d (DENS, POTT, TEMP, QV, POTE) |
| calculate equivalent potential temperature Bolton, D., 1980: The computation of equivalent potential temperature. Monthly Weather Rev., 108, 1046-1053. PT_E = PT exp( L QV / (CPdry T_L) f ) f ~ 1.0784 ( 1 + 0.810 QV ) Here T_L is temperature at the lifting condensation level and T_L ~ 55 + 2840 / ( CPdry/Rdry log(T) - log(P_v) - 4.805 ) More...
|
|
subroutine | atmos_saturation_pote_1d (KA, KS, KE, DENS, POTT, TEMP, QV, POTE) |
|
subroutine | atmos_saturation_moist_conversion_dens_all_0d (DENS, Emoist0, TEMP, QV, QC, QI, CPtot, CVtot, converged) |
| Iterative moist conversion (liquid/ice mixture) at constant density (volume) More...
|
|
subroutine | atmos_saturation_moist_conversion_pres_liq_0d (PRES, Entr, Qdry, QV, QC, Rtot, CPtot, TEMP, converged) |
| Iterative moist conversion for liquid water at constant pressure. More...
|
|
module atmosphere / saturation
- Description
- Saturation module
- Author
- Team SCALE
- NAMELIST
-
PARAM_ATMOS_SATURATION
name | type | default value | comment |
ATMOS_SATURATION_ULIMIT_TEMP | real(RP) | 273.15_RP | upper limit temperature |
ATMOS_SATURATION_LLIMIT_TEMP | real(RP) | 233.15_RP | lower limit temperature |
- History Output
- No history output
◆ atmos_saturation_setup()
subroutine, public scale_atmos_saturation::atmos_saturation_setup |
Setup.
Definition at line 244 of file scale_atmos_saturation.F90.
260 namelist / param_atmos_saturation / &
261 atmos_saturation_ulimit_temp, &
262 atmos_saturation_llimit_temp
268 log_info(
"ATMOS_SATURATION_setup",*)
'Setup'
272 read(io_fid_conf,nml=param_atmos_saturation,iostat=ierr)
274 log_info(
"ATMOS_SATURATION_setup",*)
'Not found namelist. Default used.'
275 elseif( ierr > 0 )
then
276 log_error(
"ATMOS_SATURATION_setup",*)
'Not appropriate names in namelist PARAM_ATMOS_SATURATION. Check!'
279 log_nml(param_atmos_saturation)
281 rtem00 = 1.0_rp / tem00
285 cpovr_liq = ( cpvap - cl ) / rvap
286 cpovr_ice = ( cpvap - ci ) / rvap
287 cvovr_liq = ( cvvap - cl ) / rvap
288 cvovr_ice = ( cvvap - ci ) / rvap
290 lovr_liq = lhv00 / rvap
291 lovr_ice = lhs00 / rvap
300 lovr_liq = lhv0 / rvap
301 lovr_ice = lhs0 / rvap
305 dalphadt_const = 1.0_rp / ( atmos_saturation_ulimit_temp - atmos_saturation_llimit_temp )
308 log_info(
"ATMOS_SATURATION_setup",
'(1x,A,F7.2,A,F7.2)')
'Temperature range for liquid/ice mixture : ', &
309 atmos_saturation_llimit_temp,
' - ', &
310 atmos_saturation_ulimit_temp
318 call atmos_saturation_psat_liq( tem_min, psat_min_liq )
319 call atmos_saturation_psat_ice( tem_min, psat_min_ice )
References scale_atmos_hydrometeor::atmos_hydrometeor_setup(), scale_const::const_ci, scale_const::const_cl, scale_const::const_cpvap, scale_const::const_cvvap, scale_const::const_lhs0, scale_const::const_lhs00, scale_const::const_lhv0, scale_const::const_lhv00, scale_const::const_thermodyn_type, scale_io::io_fid_conf, and scale_prc::prc_abort().
Referenced by scale_atmos_adiabat::atmos_adiabat_setup(), mod_rm_driver::rm_driver(), and mod_rm_prep::rm_prep().
◆ atmos_saturation_alpha_0d()
subroutine scale_atmos_saturation::atmos_saturation_alpha_0d |
( |
real(rp), intent(in) |
temp, |
|
|
real(rp), intent(out) |
alpha |
|
) |
| |
calc liquid/ice separation factor (0D)
- Parameters
-
[in] | temp | temperature [K] |
[out] | alpha | liquid/ice separation factor (0-1) |
Definition at line 331 of file scale_atmos_saturation.F90.
334 real(RP),
intent(in) :: temp
335 real(RP),
intent(out) :: alpha
338 alpha = ( temp - atmos_saturation_llimit_temp ) &
339 / ( atmos_saturation_ulimit_temp - atmos_saturation_llimit_temp )
341 alpha = min( max( alpha, 0.0_rp ), 1.0_rp )
◆ atmos_saturation_dqs_dtem_dens_liq_1d()
subroutine scale_atmos_saturation::atmos_saturation_dqs_dtem_dens_liq_1d |
( |
integer, intent(in) |
KA, |
|
|
integer, intent(in) |
KS, |
|
|
integer, intent(in) |
KE, |
|
|
real(rp), dimension (ka), intent(in) |
temp, |
|
|
real(rp), dimension (ka), intent(in) |
dens, |
|
|
real(rp), dimension(ka), intent(out) |
dqsdtem |
|
) |
| |
- Parameters
-
[in] | temp | temperature [K] |
[in] | dens | temperature [K] |
[out] | dqsdtem | (d qsw/d T)_{rho} |
Definition at line 1949 of file scale_atmos_saturation.F90.
1951 integer,
intent(in) :: KA, KS, KE
1953 real(RP),
intent(in) :: temp (KA)
1954 real(RP),
intent(in) :: dens (KA)
1956 real(RP),
intent(out) :: dqsdtem(KA)
1962 call atmos_saturation_dqs_dtem_dens_liq_0d( temp(k), dens(k), &
◆ atmos_saturation_dqs_dtem_dens_ice_1d()
subroutine scale_atmos_saturation::atmos_saturation_dqs_dtem_dens_ice_1d |
( |
integer, intent(in) |
KA, |
|
|
integer, intent(in) |
KS, |
|
|
integer, intent(in) |
KE, |
|
|
real(rp), dimension (ka), intent(in) |
temp, |
|
|
real(rp), dimension (ka), intent(in) |
dens, |
|
|
real(rp), dimension(ka), intent(out) |
dqsdtem |
|
) |
| |
Definition at line 2042 of file scale_atmos_saturation.F90.
2044 integer,
intent(in) :: KA, KS, KE
2046 real(RP),
intent(in) :: temp (KA)
2047 real(RP),
intent(in) :: dens (KA)
2049 real(RP),
intent(out) :: dqsdtem(KA)
2055 call atmos_saturation_dqs_dtem_dens_ice_0d( temp(k), dens(k), &
◆ atmos_saturation_dqs_dtem_dpre_liq_1d()
subroutine scale_atmos_saturation::atmos_saturation_dqs_dtem_dpre_liq_1d |
( |
integer, intent(in) |
KA, |
|
|
integer, intent(in) |
KS, |
|
|
integer, intent(in) |
KE, |
|
|
real(rp), dimension (ka), intent(in) |
temp, |
|
|
real(rp), dimension (ka), intent(in) |
pres, |
|
|
real(rp), dimension (ka), intent(in) |
qdry, |
|
|
real(rp), dimension(ka), intent(out) |
dqsat_dT, |
|
|
real(rp), dimension(ka), intent(out) |
dqsat_dP |
|
) |
| |
Definition at line 2189 of file scale_atmos_saturation.F90.
2191 integer,
intent(in) :: KA, KS, KE
2193 real(RP),
intent(in) :: temp (KA)
2194 real(RP),
intent(in) :: pres (KA)
2195 real(RP),
intent(in) :: qdry (KA)
2197 real(RP),
intent(out) :: dqsat_dT(KA)
2198 real(RP),
intent(out) :: dqsat_dP(KA)
2204 call atmos_saturation_dqs_dtem_dpre_liq_0d( temp(k), pres(k), qdry(k), &
2205 dqsat_dt(k), dqsat_dp(k) )
◆ atmos_saturation_dqs_dtem_dpre_ice_1d()
subroutine scale_atmos_saturation::atmos_saturation_dqs_dtem_dpre_ice_1d |
( |
integer, intent(in) |
KA, |
|
|
integer, intent(in) |
KS, |
|
|
integer, intent(in) |
KE, |
|
|
real(rp), dimension (ka), intent(in) |
temp, |
|
|
real(rp), dimension (ka), intent(in) |
pres, |
|
|
real(rp), dimension (ka), intent(in) |
qdry, |
|
|
real(rp), dimension(ka), intent(out) |
dqsat_dT, |
|
|
real(rp), dimension(ka), intent(out) |
dqsat_dP |
|
) |
| |
Definition at line 2295 of file scale_atmos_saturation.F90.
2297 integer,
intent(in) :: KA, KS, KE
2299 real(RP),
intent(in) :: temp (KA)
2300 real(RP),
intent(in) :: pres (KA)
2301 real(RP),
intent(in) :: qdry (KA)
2303 real(RP),
intent(out) :: dqsat_dT(KA)
2304 real(RP),
intent(out) :: dqsat_dP(KA)
2310 call atmos_saturation_dqs_dtem_dpre_ice( temp(k), pres(k), qdry(k), &
2311 dqsat_dt(k), dqsat_dp(k) )
References scale_const::const_undef.
◆ atmos_saturation_tdew_liq_1d()
subroutine scale_atmos_saturation::atmos_saturation_tdew_liq_1d |
( |
integer, intent(in) |
KA, |
|
|
integer, intent(in) |
KS, |
|
|
integer, intent(in) |
KE, |
|
|
real(rp), dimension(ka), intent(in) |
DENS, |
|
|
real(rp), dimension(ka), intent(in) |
TEMP, |
|
|
real(rp), dimension (ka), intent(in) |
QV, |
|
|
real(rp), dimension(ka), intent(out) |
Tdew, |
|
|
logical, intent(out) |
converged |
|
) |
| |
Definition at line 2481 of file scale_atmos_saturation.F90.
2482 integer,
intent(in) :: KA, KS, KE
2484 real(RP),
intent(in) :: DENS(KA)
2485 real(RP),
intent(in) :: TEMP(KA)
2486 real(RP),
intent(in) :: QV (KA)
2488 real(RP),
intent(out) :: Tdew(KA)
2489 logical,
intent(out) :: converged
2495 call atmos_saturation_tdew_liq_0d( dens(k), temp(k), qv(k), &
2496 tdew(k), converged )
2497 if ( .not. converged )
exit
References scale_prc::prc_abort().
◆ atmos_saturation_pote_0d()
subroutine scale_atmos_saturation::atmos_saturation_pote_0d |
( |
real(rp), intent(in) |
DENS, |
|
|
real(rp), intent(in) |
POTT, |
|
|
real(rp), intent(in) |
TEMP, |
|
|
real(rp), intent(in) |
QV, |
|
|
real(rp), intent(out) |
POTE |
|
) |
| |
calculate equivalent potential temperature Bolton, D., 1980: The computation of equivalent potential temperature. Monthly Weather Rev., 108, 1046-1053. PT_E = PT exp( L QV / (CPdry T_L) f ) f ~ 1.0784 ( 1 + 0.810 QV ) Here T_L is temperature at the lifting condensation level and T_L ~ 55 + 2840 / ( CPdry/Rdry log(T) - log(P_v) - 4.805 )
Definition at line 2561 of file scale_atmos_saturation.F90.
2566 atmos_hydrometeor_lhv
2567 real(RP),
intent(in) :: DENS
2568 real(RP),
intent(in) :: POTT
2569 real(RP),
intent(in) :: TEMP
2570 real(RP),
intent(in) :: QV
2572 real(RP),
intent(out) :: POTE
2578 pv = dens * qv * rvap * temp
2579 tl = 55.0_rp + 2840.0_rp / ( cpdry / rdry * log(temp) - log(pv) - 4.805_rp )
2580 call atmos_hydrometeor_lhv( temp, lhv )
2582 pote = pott * exp( lhv * qv / ( cpdry * tl ) &
2583 * 1.0784_rp * ( 1.0_rp + 0.810_rp * qv ) )
References scale_const::const_cpdry, and scale_const::const_rdry.
Referenced by atmos_saturation_pote_1d().
◆ atmos_saturation_pote_1d()
subroutine scale_atmos_saturation::atmos_saturation_pote_1d |
( |
integer, intent(in) |
KA, |
|
|
integer, intent(in) |
KS, |
|
|
integer, intent(in) |
KE, |
|
|
real(rp), dimension(ka), intent(in) |
DENS, |
|
|
real(rp), dimension(ka), intent(in) |
POTT, |
|
|
real(rp), dimension(ka), intent(in) |
TEMP, |
|
|
real(rp), dimension (ka), intent(in) |
QV, |
|
|
real(rp), dimension(ka), intent(out) |
POTE |
|
) |
| |
Definition at line 2593 of file scale_atmos_saturation.F90.
2594 integer,
intent(in) :: KA, KS, KE
2596 real(RP),
intent(in) :: DENS(KA)
2597 real(RP),
intent(in) :: POTT(KA)
2598 real(RP),
intent(in) :: TEMP(KA)
2599 real(RP),
intent(in) :: QV (KA)
2601 real(RP),
intent(out) :: POTE(KA)
2606 call atmos_saturation_pote_0d( dens(k), pott(k), temp(k), qv(k), &
References atmos_saturation_pote_0d(), scale_atmos_hydrometeor::cp_vapor, scale_atmos_hydrometeor::cp_water, scale_atmos_hydrometeor::cv_vapor, scale_atmos_hydrometeor::cv_water, and scale_atmos_hydrometeor::lhv.
◆ atmos_saturation_moist_conversion_dens_all_0d()
subroutine scale_atmos_saturation::atmos_saturation_moist_conversion_dens_all_0d |
( |
real(rp), intent(in) |
DENS, |
|
|
real(rp), intent(in) |
Emoist0, |
|
|
real(rp), intent(inout) |
TEMP, |
|
|
real(rp), intent(inout) |
QV, |
|
|
real(rp), intent(inout) |
QC, |
|
|
real(rp), intent(inout) |
QI, |
|
|
real(rp), intent(inout) |
CPtot, |
|
|
real(rp), intent(inout) |
CVtot, |
|
|
logical, intent(out) |
converged |
|
) |
| |
Iterative moist conversion (liquid/ice mixture) at constant density (volume)
Definition at line 2770 of file scale_atmos_saturation.F90.
2781 real(RP),
intent(in) :: DENS
2782 real(RP),
intent(in) :: Emoist0
2784 real(RP),
intent(inout) :: TEMP
2785 real(RP),
intent(inout) :: QV
2786 real(RP),
intent(inout) :: QC
2787 real(RP),
intent(inout) :: QI
2788 real(RP),
intent(inout) :: CPtot
2789 real(RP),
intent(inout) :: CVtot
2791 logical,
intent(out) :: converged
2796 real(RP) :: QV0, QC0, QI0
2797 real(RP) :: CPtot0, CVtot0
2803 real(RP) :: dalpha_dT
2804 real(RP) :: dqsat_dT
2805 real(RP) :: dqc_dT, dqi_dT
2806 real(RP) :: dCVtot_dT
2807 real(RP) :: dEmoist_dT
2810 integer,
parameter :: itelim = 100
2811 real(RP),
parameter :: dtemp_criteria = 0.1_rp**(2+rp/2)
2829 temp = ( emoist0 - lhv * qv ) / cvtot
2831 call atmos_saturation_dens2qsat_all( temp, dens, &
2834 if ( qsum <= qsat )
then
2852 call atmos_saturation_dqs_dtem_dens_all_0d( temp, dens, &
2854 qsat=qsat, alpha=alpha )
2855 call atmos_saturation_dalphadt( temp, dalpha_dt )
2857 dqc_dt = ( qsum - qv ) * dalpha_dt - dqsat_dt * ( alpha )
2858 dqi_dt = -( qsum - qv ) * dalpha_dt - dqsat_dt * ( 1.0_rp-alpha )
2864 demoist_dt = temp * dcvtot_dt + cvtot + dqsat_dt * lhv - dqi_dt *
lhf
2868 qc = ( qsum - qsat ) * ( alpha )
2869 qi = ( qsum - qsat ) * ( 1.0_rp - alpha )
2878 emoist = temp * cvtot + qsat * lhv - qi *
lhf
2881 dtemp = ( emoist - emoist0 ) / demoist_dt
2884 if ( abs(dtemp) < dtemp_criteria )
then
2889 if( temp*0.0_rp /= 0.0_rp )
exit
References scale_atmos_hydrometeor::cp_ice, scale_atmos_hydrometeor::cp_vapor, scale_atmos_hydrometeor::cp_water, scale_atmos_hydrometeor::cv_ice, scale_atmos_hydrometeor::cv_vapor, scale_atmos_hydrometeor::cv_water, scale_atmos_hydrometeor::lhf, and scale_atmos_hydrometeor::lhv.
◆ atmos_saturation_moist_conversion_pres_liq_0d()
subroutine scale_atmos_saturation::atmos_saturation_moist_conversion_pres_liq_0d |
( |
real(rp), intent(in) |
PRES, |
|
|
real(rp), intent(in) |
Entr, |
|
|
real(rp), intent(in) |
Qdry, |
|
|
real(rp), intent(inout) |
QV, |
|
|
real(rp), intent(inout) |
QC, |
|
|
real(rp), intent(inout) |
Rtot, |
|
|
real(rp), intent(inout) |
CPtot, |
|
|
real(rp), intent(out) |
TEMP, |
|
|
logical, intent(out) |
converged |
|
) |
| |
Iterative moist conversion for liquid water at constant pressure.
Definition at line 2909 of file scale_atmos_saturation.F90.
2914 atmos_hydrometeor_entr, &
2915 atmos_hydrometeor_entr2temp, &
2919 real(RP),
intent(in) :: PRES
2920 real(RP),
intent(in) :: ENTR
2921 real(RP),
intent(in) :: Qdry
2923 real(RP),
intent(inout) :: QV
2924 real(RP),
intent(inout) :: QC
2925 real(RP),
intent(inout) :: CPtot
2926 real(RP),
intent(inout) :: Rtot
2928 real(RP),
intent(out) :: TEMP
2929 logical,
intent(out) :: converged
2931 real(RP),
parameter :: TEMMIN = 0.1_rp
2932 real(RP),
parameter :: criteria = 1.e-8_rp
2933 integer,
parameter :: itelim = 100
2937 real(RP) :: qsat, psat
2939 real(RP) :: TEMP_ite
2941 real(RP) :: ENTR_ite
2942 real(RP) :: Rtot_ite
2943 real(RP) :: CPtot_ite
2945 real(RP) :: dqsat_dT, dqsat_dP
2947 real(RP) :: TEMP_prev
2948 real(RP) :: dENTR_dT
2954 rtot = rtot + qc * rvap
2959 call atmos_hydrometeor_entr2temp( entr, pres, qsum, 0.0_rp, qdry, &
2963 call atmos_saturation_pres2qsat_liq( temp, pres, qdry, &
2965 if ( qsum <= qsat )
then
2980 call atmos_saturation_dqs_dtem_dpre_liq_0d( temp_ite, pres, qdry, &
2981 dqsat_dt, dqsat_dp, &
2982 qsat=qsat, psat=psat )
2984 qv_ite = min( qsum, qsat )
2986 rtot_ite = rtot - ( qsum - qv_ite ) * rvap
2989 dentr_dt = cptot_ite / temp_ite &
2991 - rvap * log( psat/ psat0 ) &
2995 call atmos_hydrometeor_entr( temp_ite, pres, &
2996 qv_ite, 0.0_rp, qdry, &
2997 rtot_ite, cptot_ite, &
3000 temp_prev = temp_ite
3001 temp_ite = temp_ite - ( entr_ite - entr ) / max( dentr_dt, eps )
3002 temp_ite = max( temp_ite, temmin )
3004 if( abs(temp_ite-temp_prev) < criteria )
then
References scale_const::const_eps, scale_const::const_lhv0, scale_atmos_hydrometeor::cp_vapor, and scale_atmos_hydrometeor::cp_water.
real(rp), parameter, public const_lhv0
latent heat of vaporizaion at 0C [J/kg]
subroutine, public prc_abort
Abort Process.
real(rp), public cp_water
CP for water [J/kg/K].
real(rp), public const_lhs00
latent heat of sublimation at 0K [J/kg]
real(rp), public const_lhv00
latent heat of vaporizaion at 0K [J/kg]
real(rp), public const_eps
small number
module atmosphere / hydrometeor
character(len=h_short), public const_thermodyn_type
internal energy type
real(rp), parameter, public const_cpvap
specific heat (water vapor, constant pressure) [J/kg/K]
real(rp), public cv_vapor
CV for vapor [J/kg/K].
real(rp), public const_cpdry
specific heat (dry air,constant pressure) [J/kg/K]
real(rp), parameter, public const_ci
specific heat (ice) [J/kg/K]
real(rp), parameter, public const_cl
specific heat (liquid water) [J/kg/K]
real(rp), public const_cvvap
specific heat (water vapor, constant volume) [J/kg/K]
real(rp), public lhf
latent heat of fusion for use [J/kg]
real(rp), public const_rdry
specific gas constant (dry air) [J/kg/K]
real(rp), parameter, public const_lhs0
latent heat of sublimation at 0C [J/kg]
subroutine, public atmos_hydrometeor_setup
Setup.
real(rp), public cp_vapor
CP for vapor [J/kg/K].
real(rp), public cv_water
CV for water [J/kg/K].
real(rp), public cv_ice
CV for ice [J/kg/K].
real(rp), public cp_ice
CP for ice [J/kg/K].