SCALE-RM
Data Types | Functions/Subroutines
scale_atmos_saturation Module Reference

module atmosphere / saturation More...

Functions/Subroutines

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) 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...
 

Detailed Description

module atmosphere / saturation

Description
Saturation module
Author
Team SCALE
NAMELIST
  • PARAM_ATMOS_SATURATION
    nametypedefault valuecomment
    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

Function/Subroutine Documentation

◆ atmos_saturation_setup()

subroutine, public scale_atmos_saturation::atmos_saturation_setup

Setup.

Definition at line 215 of file scale_atmos_saturation.F90.

215  use scale_prc, only: &
216  prc_abort
217  use scale_const, only: &
218  cpvap => const_cpvap, &
219  cvvap => const_cvvap, &
220  cl => const_cl, &
221  ci => const_ci, &
222  lhv00 => const_lhv00, &
223  lhs00 => const_lhs00, &
224  lhv0 => const_lhv0, &
225  lhs0 => const_lhs0, &
227  use scale_atmos_hydrometeor, only: &
229  implicit none
230 
231  namelist / param_atmos_saturation / &
232  atmos_saturation_ulimit_temp, &
233  atmos_saturation_llimit_temp
234 
235  integer :: ierr
236  !---------------------------------------------------------------------------
237 
238  log_newline
239  log_info("ATMOS_SATURATION_setup",*) 'Setup'
240 
241  !--- read namelist
242  rewind(io_fid_conf)
243  read(io_fid_conf,nml=param_atmos_saturation,iostat=ierr)
244  if( ierr < 0 ) then !--- missing
245  log_info("ATMOS_SATURATION_setup",*) 'Not found namelist. Default used.'
246  elseif( ierr > 0 ) then !--- fatal error
247  log_error("ATMOS_SATURATION_setup",*) 'Not appropriate names in namelist PARAM_ATMOS_SATURATION. Check!'
248  call prc_abort
249  endif
250  log_nml(param_atmos_saturation)
251 
252  rtem00 = 1.0_rp / tem00
253 
254  if ( const_thermodyn_type == 'EXACT' ) then
255 
256  cpovr_liq = ( cpvap - cl ) / rvap
257  cpovr_ice = ( cpvap - ci ) / rvap
258  cvovr_liq = ( cvvap - cl ) / rvap
259  cvovr_ice = ( cvvap - ci ) / rvap
260 
261  lovr_liq = lhv00 / rvap
262  lovr_ice = lhs00 / rvap
263 
264  elseif( const_thermodyn_type == 'SIMPLE' ) then
265 
266  cpovr_liq = 0.0_rp
267  cpovr_ice = 0.0_rp
268  cvovr_liq = 0.0_rp
269  cvovr_ice = 0.0_rp
270 
271  lovr_liq = lhv0 / rvap
272  lovr_ice = lhs0 / rvap
273 
274  endif
275 
276  dalphadt_const = 1.0_rp / ( atmos_saturation_ulimit_temp - atmos_saturation_llimit_temp )
277 
278  log_newline
279  log_info("ATMOS_SATURATION_setup",'(1x,A,F7.2,A,F7.2)') 'Temperature range for liquid/ice mixture : ', &
280  atmos_saturation_llimit_temp, ' - ', &
281  atmos_saturation_ulimit_temp
282 
284 
285  call atmos_saturation_psat_liq( tem_min, psat_min_liq ) ! [IN], [OUT]
286  call atmos_saturation_psat_ice( tem_min, psat_min_ice ) ! [IN], [OUT]
287 
288  return

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().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ 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]temptemperature [K]
[out]alphaliquid/ice separation factor (0-1)

Definition at line 296 of file scale_atmos_saturation.F90.

296  implicit none
297 
298  real(RP), intent(in) :: temp
299  real(RP), intent(out) :: alpha
300  !---------------------------------------------------------------------------
301 
302  alpha = ( temp - atmos_saturation_llimit_temp ) &
303  / ( atmos_saturation_ulimit_temp - atmos_saturation_llimit_temp )
304 
305  alpha = min( max( alpha, 0.0_rp ), 1.0_rp )
306 
307  return

◆ 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]temptemperature [K]
[in]denstemperature [K]
[out]dqsdtem(d qsw/d T)_{rho}

Definition at line 1356 of file scale_atmos_saturation.F90.

1356  implicit none
1357  integer, intent(in) :: KA, KS, KE
1358 
1359  real(RP), intent(in) :: temp (KA)
1360  real(RP), intent(in) :: dens (KA)
1361 
1362  real(RP), intent(out) :: dqsdtem(KA)
1363 
1364  integer :: k
1365  !---------------------------------------------------------------------------
1366 
1367  do k = ks, ke
1368  call atmos_saturation_dqs_dtem_dens_liq_0d( temp(k), dens(k), & ! [IN]
1369  dqsdtem(k) ) ! [OUT]
1370  enddo
1371 
1372  return

◆ 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 1446 of file scale_atmos_saturation.F90.

1446  implicit none
1447  integer, intent(in) :: KA, KS, KE
1448 
1449  real(RP), intent(in) :: temp (KA)
1450  real(RP), intent(in) :: dens (KA)
1451 
1452  real(RP), intent(out) :: dqsdtem(KA)
1453 
1454  integer :: k
1455  !---------------------------------------------------------------------------
1456 
1457  do k = ks, ke
1458  call atmos_saturation_dqs_dtem_dens_ice_0d( temp(k), dens(k), & ! [IN]
1459  dqsdtem(k) ) ! [OUT]
1460  enddo
1461 
1462  return

◆ 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 1588 of file scale_atmos_saturation.F90.

1588  implicit none
1589  integer, intent(in) :: KA, KS, KE
1590 
1591  real(RP), intent(in) :: temp (KA)
1592  real(RP), intent(in) :: pres (KA)
1593  real(RP), intent(in) :: qdry (KA)
1594 
1595  real(RP), intent(out) :: dqsat_dT(KA)
1596  real(RP), intent(out) :: dqsat_dP(KA)
1597 
1598  integer :: k
1599  !---------------------------------------------------------------------------
1600 
1601  do k = ks, ke
1602  call atmos_saturation_dqs_dtem_dpre_liq_0d( temp(k), pres(k), qdry(k), & ! [IN]
1603  dqsat_dt(k), dqsat_dp(k) ) ! [OUT]
1604  enddo
1605 
1606  return

◆ 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 1690 of file scale_atmos_saturation.F90.

1690  implicit none
1691  integer, intent(in) :: KA, KS, KE
1692 
1693  real(RP), intent(in) :: temp (KA)
1694  real(RP), intent(in) :: pres (KA)
1695  real(RP), intent(in) :: qdry (KA)
1696 
1697  real(RP), intent(out) :: dqsat_dT(KA)
1698  real(RP), intent(out) :: dqsat_dP(KA)
1699 
1700  integer :: k
1701  !---------------------------------------------------------------------------
1702 
1703  do k = ks, ke
1704  call atmos_saturation_dqs_dtem_dpre_ice( temp(k), pres(k), qdry(k), & ! [IN]
1705  dqsat_dt(k), dqsat_dp(k) ) ! [OUT]
1706  enddo
1707 
1708  return

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 1868 of file scale_atmos_saturation.F90.

1868  integer, intent(in) :: KA, KS, KE
1869 
1870  real(RP), intent(in) :: DENS(KA)
1871  real(RP), intent(in) :: TEMP(KA)
1872  real(RP), intent(in) :: QV (KA)
1873 
1874  real(RP), intent(out) :: Tdew(KA)
1875  logical, intent(out) :: converged
1876 
1877  integer :: k
1878  !---------------------------------------------------------------------------
1879 
1880  do k = ks, ke
1881  call atmos_saturation_tdew_liq_0d( dens(k), temp(k), qv(k), & ! [IN]
1882  tdew(k), converged ) ! [OUT]
1883  if ( .not. converged ) exit
1884  end do
1885 

References scale_prc::prc_abort().

Here is the call graph for this function:

◆ 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) 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 1940 of file scale_atmos_saturation.F90.

1940  use scale_const, only: &
1941  rdry => const_rdry, &
1942  cpdry => const_cpdry
1943  use scale_atmos_hydrometeor, only: &
1944  atmos_hydrometeor_lhv
1945  real(RP), intent(in) :: DENS
1946  real(RP), intent(in) :: POTT
1947  real(RP), intent(in) :: TEMP
1948  real(RP), intent(in) :: QV
1949 
1950  real(RP), intent(out) :: POTE
1951 
1952  real(RP) :: TL
1953  real(RP) :: Pv
1954  real(RP) :: LHV
1955 
1956  pv = dens * qv * rvap * temp
1957  tl = 55.0_rp + 2840.0_rp / ( cpdry / rdry * log(temp) - log(pv) - 4.805_rp )
1958  call atmos_hydrometeor_lhv( temp, lhv ) ! [IN], [OUT]
1959 
1960  pote = pott * exp( lhv * qv / ( cpdry * temp ) &
1961  * 1.0784_rp * ( 1.0_rp + 0.810_rp * qv ) )
1962 
1963  return

References scale_const::const_cpdry, and scale_const::const_rdry.

Referenced by atmos_saturation_pote_1d().

Here is the caller graph for this function:

◆ 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 1971 of file scale_atmos_saturation.F90.

1971  integer, intent(in) :: KA, KS, KE
1972 
1973  real(RP), intent(in) :: DENS(KA)
1974  real(RP), intent(in) :: POTT(KA)
1975  real(RP), intent(in) :: TEMP(KA)
1976  real(RP), intent(in) :: QV (KA)
1977 
1978  real(RP), intent(out) :: POTE(KA)
1979 
1980  integer :: k
1981 
1982  do k = ks, ke
1983  call atmos_saturation_pote_0d( dens(k), pott(k), temp(k), qv(k), & ! [IN]
1984  pote(k) ) ! [OUT]
1985  end do
1986 
1987  return

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, scale_atmos_hydrometeor::lhv, and scale_prc::prc_abort().

Here is the call graph for this function:

◆ 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 2132 of file scale_atmos_saturation.F90.

2132  use scale_prc, only: &
2133  prc_abort
2134  use scale_atmos_hydrometeor, only: &
2135  cp_vapor, &
2136  cp_water, &
2137  cp_ice, &
2138  cv_vapor, &
2139  cv_water, &
2140  cv_ice, &
2141  lhv, &
2142  lhf
2143  implicit none
2144 
2145  real(RP), intent(in) :: DENS
2146  real(RP), intent(in) :: Emoist0
2147 
2148  real(RP), intent(inout) :: TEMP
2149  real(RP), intent(inout) :: QV
2150  real(RP), intent(inout) :: QC
2151  real(RP), intent(inout) :: QI
2152  real(RP), intent(inout) :: CPtot
2153  real(RP), intent(inout) :: CVtot
2154 
2155  logical, intent(out) :: converged
2156 
2157  ! working
2158  real(RP) :: QSUM
2159  real(RP) :: TEMP0
2160  real(RP) :: QV0, QC0, QI0
2161  real(RP) :: CVtot0
2162  real(RP) :: alpha
2163  real(RP) :: qsat
2164  real(RP) :: Emoist ! moist internal energy
2165 
2166  ! d(X)/dT
2167  real(RP) :: dalpha_dT
2168  real(RP) :: dqsat_dT
2169  real(RP) :: dqc_dT, dqi_dT
2170  real(RP) :: dCVtot_dT
2171  real(RP) :: dEmoist_dT
2172  real(RP) :: dtemp
2173 
2174  integer, parameter :: itelim = 100
2175  real(RP), parameter :: dtemp_criteria = 0.1_rp**(2+rp/2)
2176  integer :: ite
2177  !---------------------------------------------------------------------------
2178 
2179  temp0 = temp
2180  cvtot0 = cvtot
2181  qv0 = qv
2182  qc0 = qc
2183  qi0 = qi
2184  qsum = qv + qc + qi
2185 
2186  call atmos_saturation_dens2qsat_all( temp, dens, & ! [IN]
2187  qsat ) ! [OUT]
2188 
2189  if ( qsum <= qsat ) then
2190  qv = qsum
2191  cptot = cptot + qc * ( cp_vapor - cp_water ) + qi * ( cp_vapor - cp_ice )
2192  cvtot = cvtot + qc * ( cv_vapor - cv_water ) + qi * ( cv_vapor - cv_ice )
2193  qc = 0.0_rp
2194  qi = 0.0_rp
2195  temp = ( emoist0 - lhv * qv ) / cvtot
2196  converged = .true.
2197  return
2198  end if
2199 
2200  converged = .false.
2201  do ite = 1, itelim
2202 
2203  ! dX/dT
2204  call atmos_saturation_dqs_dtem_dens_all_0d( temp, dens, & ! [IN]
2205  dqsat_dt, & ! [OUT]
2206  qsat=qsat, alpha=alpha ) ! [OUT]
2207  call atmos_saturation_dalphadt( temp, dalpha_dt ) ! [IN], [OUT]
2208 
2209  dqc_dt = ( qsum - qv ) * dalpha_dt - dqsat_dt * ( alpha )
2210  dqi_dt = -( qsum - qv ) * dalpha_dt - dqsat_dt * ( 1.0_rp-alpha )
2211 
2212  dcvtot_dt = dqsat_dt * cv_vapor &
2213  + dqc_dt * cv_water &
2214  + dqi_dt * cv_ice
2215 
2216  demoist_dt = temp * dcvtot_dt + cvtot + dqsat_dt * lhv - dqi_dt * lhf
2217 
2218  ! Saturation
2219  qv = qsat
2220  qc = ( qsum - qsat ) * ( alpha )
2221  qi = ( qsum - qsat ) * ( 1.0_rp - alpha )
2222 
2223  cvtot = cvtot0 &
2224  + cv_vapor * ( qv - qv0 ) &
2225  + cv_water * ( qc - qc0 ) &
2226  + cv_ice * ( qi - qi0 )
2227 
2228  emoist = temp * cvtot + qv * lhv - qi * lhf
2229 
2230  ! update temp by the newtonian method
2231  dtemp = ( emoist - emoist0 ) / demoist_dt
2232  temp = temp - dtemp
2233 
2234  if ( abs(dtemp) < dtemp_criteria ) then
2235  converged = .true.
2236  exit
2237  endif
2238 
2239  if( temp*0.0_rp /= 0.0_rp ) exit
2240  enddo
2241 
2242  cptot = cptot &
2243  + cp_vapor * ( qv - qv0 ) &
2244  + cp_water * ( qc - qc0 ) &
2245  + cp_ice * ( qi - qi0 )
2246 
2247  return

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, scale_atmos_hydrometeor::lhv, and scale_prc::prc_abort().

Here is the call graph for this function:

◆ 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 2259 of file scale_atmos_saturation.F90.

2259  use scale_const, only: &
2260  eps => const_eps, &
2261  lhv0 => const_lhv0
2262  use scale_atmos_hydrometeor, only: &
2263  atmos_hydrometeor_entr, &
2264  atmos_hydrometeor_entr2temp, &
2265  cp_vapor, &
2266  cp_water
2267 
2268  real(RP), intent(in) :: PRES
2269  real(RP), intent(in) :: ENTR
2270  real(RP), intent(in) :: Qdry
2271 
2272  real(RP), intent(inout) :: QV
2273  real(RP), intent(inout) :: QC
2274  real(RP), intent(inout) :: CPtot
2275  real(RP), intent(inout) :: Rtot
2276 
2277  real(RP), intent(out) :: TEMP
2278  logical, intent(out) :: converged
2279 
2280  real(RP), parameter :: TEMMIN = 0.1_rp
2281  real(RP), parameter :: criteria = 1.e-8_rp
2282  integer, parameter :: itelim = 100
2283  integer :: ite
2284 
2285  real(RP) :: Qsum
2286  real(RP) :: qsat, psat
2287 
2288  real(RP) :: TEMP_ite
2289  real(RP) :: QV_ite
2290  real(RP) :: ENTR_ite
2291  real(RP) :: Rtot_ite
2292  real(RP) :: CPtot_ite
2293 
2294  real(RP) :: dqsat_dT, dqsat_dP
2295 
2296  real(RP) :: TEMP_prev
2297  real(RP) :: dENTR_dT
2298  !---------------------------------------------------------------------------
2299 
2300  qsum = qv + qc
2301 
2302  ! all liquid evaporates
2303  rtot = rtot + qc * rvap
2304  cptot = cptot + qc * ( cp_vapor - cp_water )
2305  qv = qsum
2306  qc = 0.0_rp
2307 
2308  call atmos_hydrometeor_entr2temp( entr, pres, qsum, 0.0_rp, qdry, & ! [IN]
2309  rtot, cptot, & ! [IN]
2310  temp ) ! [OUT]
2311 
2312  call atmos_saturation_pres2qsat_liq( temp, pres, qdry, & ! [IN]
2313  qsat ) ! [OUT]
2314  if ( qsum <= qsat ) then
2315  ! unsaturated
2316 
2317  converged = .true.
2318 
2319  return
2320 
2321  else
2322  ! saturated
2323 
2324  temp_ite = temp
2325 
2326  do ite = 1, itelim
2327 
2328  call atmos_saturation_dqs_dtem_dpre_liq_0d( temp_ite, pres, qdry, & ! [IN]
2329  dqsat_dt, dqsat_dp, & ! [OUT]
2330  qsat=qsat, psat=psat ) ! [OUT]
2331 
2332  qv_ite = min( qsum, qsat )
2333 
2334  rtot_ite = rtot - ( qsum - qv_ite ) * rvap
2335  cptot_ite = cptot - ( qsum - qv_ite ) * ( cp_vapor - cp_water )
2336 
2337  dentr_dt = cptot_ite / temp_ite &
2338  + ( ( cp_vapor - cp_water ) * log( temp_ite / tem00 ) &
2339  - rvap * log( psat/ psat0 ) &
2340  + lhv0 / tem00 &
2341  ) * dqsat_dt
2342 
2343  call atmos_hydrometeor_entr( temp_ite, pres, & ! [IN]
2344  qv_ite, 0.0_rp, qdry, & ! [IN] ! QI = 0
2345  rtot_ite, cptot_ite, & ! [IN]
2346  entr_ite ) ! [OUT]
2347 
2348  temp_prev = temp_ite
2349  temp_ite = temp_ite - ( entr_ite - entr ) / max( dentr_dt, eps )
2350  temp_ite = max( temp_ite, temmin )
2351 
2352  if( abs(temp_ite-temp_prev) < criteria ) then
2353  converged = .true.
2354  exit
2355  end if
2356 
2357  enddo
2358 
2359  end if
2360 
2361  qv = qv_ite
2362  qc = qsum - qv
2363  cptot = cptot_ite
2364  rtot = rtot_ite
2365  temp = temp_ite
2366 
2367  return

References scale_const::const_eps, scale_const::const_lhv0, scale_atmos_hydrometeor::cp_vapor, and scale_atmos_hydrometeor::cp_water.

scale_const::const_lhv0
real(rp), parameter, public const_lhv0
latent heat of vaporizaion at 0C [J/kg]
Definition: scale_const.F90:75
scale_prc::prc_abort
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:342
scale_atmos_hydrometeor::cp_water
real(rp), public cp_water
CP for water [J/kg/K].
Definition: scale_atmos_hydrometeor.F90:133
scale_const::const_lhs00
real(rp), public const_lhs00
latent heat of sublimation at 0K [J/kg]
Definition: scale_const.F90:78
scale_const::const_lhv00
real(rp), public const_lhv00
latent heat of vaporizaion at 0K [J/kg]
Definition: scale_const.F90:76
scale_const::const_eps
real(rp), public const_eps
small number
Definition: scale_const.F90:33
scale_atmos_hydrometeor
module atmosphere / hydrometeor
Definition: scale_atmos_hydrometeor.F90:12
scale_const::const_thermodyn_type
character(len=h_short), public const_thermodyn_type
internal energy type
Definition: scale_const.F90:99
scale_const::const_cpvap
real(rp), parameter, public const_cpvap
specific heat (water vapor, constant pressure) [J/kg/K]
Definition: scale_const.F90:64
scale_prc
module PROCESS
Definition: scale_prc.F90:11
scale_atmos_hydrometeor::cv_vapor
real(rp), public cv_vapor
CV for vapor [J/kg/K].
Definition: scale_atmos_hydrometeor.F90:130
scale_const
module CONSTANT
Definition: scale_const.F90:11
scale_const::const_cpdry
real(rp), public const_cpdry
specific heat (dry air,constant pressure) [J/kg/K]
Definition: scale_const.F90:56
scale_const::const_ci
real(rp), parameter, public const_ci
specific heat (ice) [J/kg/K]
Definition: scale_const.F90:67
scale_const::const_cl
real(rp), parameter, public const_cl
specific heat (liquid water) [J/kg/K]
Definition: scale_const.F90:66
scale_const::const_cvvap
real(rp), public const_cvvap
specific heat (water vapor, constant volume) [J/kg/K]
Definition: scale_const.F90:65
scale_atmos_hydrometeor::lhf
real(rp), public lhf
latent heat of fusion for use [J/kg]
Definition: scale_atmos_hydrometeor.F90:128
scale_const::const_rdry
real(rp), public const_rdry
specific gas constant (dry air) [J/kg/K]
Definition: scale_const.F90:55
scale_const::const_lhs0
real(rp), parameter, public const_lhs0
latent heat of sublimation at 0C [J/kg]
Definition: scale_const.F90:77
scale_atmos_hydrometeor::atmos_hydrometeor_setup
subroutine, public atmos_hydrometeor_setup
Setup.
Definition: scale_atmos_hydrometeor.F90:154
scale_atmos_hydrometeor::cp_vapor
real(rp), public cp_vapor
CP for vapor [J/kg/K].
Definition: scale_atmos_hydrometeor.F90:131
scale_atmos_hydrometeor::cv_water
real(rp), public cv_water
CV for water [J/kg/K].
Definition: scale_atmos_hydrometeor.F90:132
scale_atmos_hydrometeor::cv_ice
real(rp), public cv_ice
CV for ice [J/kg/K].
Definition: scale_atmos_hydrometeor.F90:134
scale_atmos_hydrometeor::cp_ice
real(rp), public cp_ice
CP for ice [J/kg/K].
Definition: scale_atmos_hydrometeor.F90:135