SCALE-RM
Data Types | Functions/Subroutines | Variables
scale_atmos_hydrometeor Module Reference

module atmosphere / hydrometeor More...

Functions/Subroutines

subroutine, public atmos_hydrometeor_setup
 Setup. More...
 
subroutine, public atmos_hydrometeor_regist (NL, NI, NAME, DESC, UNIT, Q0, ADVC)
 ATMOS_HYDROMETEOR_regist Regist tracer. More...
 
subroutine atmos_hydrometeor_lhv_0d (temp, lhv)
 
subroutine atmos_hydrometeor_lhv_1d (KA, KS, KE, temp, lhv)
 
subroutine atmos_hydrometeor_lhs_1d (KA, KS, KE, temp, lhs)
 
subroutine atmos_hydrometeor_lhf_1d (KA, KS, KE, temp, lhf)
 
subroutine atmos_hydrometeor_entr_2d (IA, IS, IE, JA, JS, JE, TEMP, PRES, QV, QI, Qdry, Rtot, CPtot, entr)
 calc temp, pres, q -> entropy (2D) More...
 

Variables

logical, public atmos_hydrometeor_ice_phase
 
integer, public i_qv = -1
 
integer, parameter, public n_hyd = 6
 
integer, parameter, public i_hc = 1
 liquid water cloud More...
 
integer, parameter, public i_hr = 2
 liquid water rain More...
 
integer, parameter, public i_hi = 3
 ice water cloud More...
 
integer, parameter, public i_hs = 4
 snow More...
 
integer, parameter, public i_hg = 5
 graupel More...
 
integer, parameter, public i_hh = 6
 hail More...
 
character(len=h_short), dimension(n_hyd), parameter, public hyd_name = (/ "QC", "QR", "QI", "QS", "QG", "QH" /)
 
character(len=h_mid), dimension(n_hyd), parameter, public hyd_desc = (/ "cloud ", "rain ", "ice water", "snow ", "graupel ", "hail " /)
 
character(len=h_short), dimension(n_hyd), parameter, public num_name = (/ "NC", "NR", "NI", "NS", "NG", "NH" /)
 
real(rp), dimension(n_hyd), public hyd_dens
 
logical, public atmos_hydrometeor_dry = .true.
 
integer, public qha = 0
 
integer, public qhs = -1
 
integer, public qhe = -2
 
integer, public qla = 0
 
integer, public qls = -1
 
integer, public qle = -2
 
integer, public qia = 0
 
integer, public qis = -1
 
integer, public qie = -2
 
real(rp), public lhv
 latent heat of vaporization for use [J/kg] More...
 
real(rp), public lhs
 latent heat of sublimation for use [J/kg] More...
 
real(rp), public lhf
 latent heat of fusion for use [J/kg] More...
 
real(rp), public cv_vapor
 CV for vapor [J/kg/K]. More...
 
real(rp), public cp_vapor
 CP for vapor [J/kg/K]. More...
 
real(rp), public cv_water
 CV for water [J/kg/K]. More...
 
real(rp), public cp_water
 CP for water [J/kg/K]. More...
 
real(rp), public cv_ice
 CV for ice [J/kg/K]. More...
 
real(rp), public cp_ice
 CP for ice [J/kg/K]. More...
 

Detailed Description

module atmosphere / hydrometeor

Description
Hydrometeor module
Author
Team SCALE

Function/Subroutine Documentation

◆ atmos_hydrometeor_setup()

subroutine, public scale_atmos_hydrometeor::atmos_hydrometeor_setup

Setup.

Definition at line 154 of file scale_atmos_hydrometeor.F90.

154  use scale_const, only: &
155  const_setup, &
156  cpvap => const_cpvap, &
157  cvvap => const_cvvap, &
158  cl => const_cl, &
159  ci => const_ci, &
160  lhv00 => const_lhv00, &
161  lhs00 => const_lhs00, &
162  lhf00 => const_lhf00, &
163  lhv0 => const_lhv0, &
164  lhs0 => const_lhs0, &
165  lhf0 => const_lhf0, &
166  dwatr => const_dwatr, &
167  dice => const_dice, &
168  thermodyn_type => const_thermodyn_type
169  use scale_prc, only: &
170  prc_abort
171  implicit none
172  !---------------------------------------------------------------------------
173 
174  if ( initialized ) return
175  initialized = .true.
176 
177  call const_setup
178 
179  log_newline
180  log_info("ATMOS_HYDROMETEOR_setup",*) 'Setup'
181 
182  if ( thermodyn_type == 'EXACT' ) then
183 
184  cv_vapor = cvvap
185  cp_vapor = cpvap
186  cv_water = cl
187  cp_water = cv_water
188  cv_ice = ci
189  cp_ice = cv_ice
190 
191  lhv = lhv00
192  lhs = lhs00
193  lhf = lhf00
194  thermodyn_emask = 1.0_rp
195 
196  elseif( thermodyn_type == 'SIMPLE' ) then
197 
198  cv_vapor = cvvap
199  cp_vapor = cpvap
200  cv_water = cvvap
201  cp_water = cv_water
202  cv_ice = cvvap
203  cp_ice = cv_ice
204 
205  lhv = lhv0
206  lhs = lhs0
207  lhf = lhf0
208  thermodyn_emask = 0.0_rp
209 
210  else
211  log_error("ATMOS_HYDROMETEOR_setup",*) 'Not appropriate ATMOS_THERMODYN_ENERGY_TYPE. Check!', trim(thermodyn_type)
212  call prc_abort
213  endif
214 
215  hyd_dens(:) = (/ dwatr, & ! HC
216  dwatr, & ! HR
217  dice, & ! HI
218  dice, & ! HS
219  dice, & ! HG
220  dice /) ! HH
221 
222  return

References scale_const::const_ci, scale_const::const_cl, scale_const::const_cpvap, scale_const::const_cvvap, scale_const::const_dice, scale_const::const_dwatr, scale_const::const_lhf0, scale_const::const_lhf00, scale_const::const_lhs0, scale_const::const_lhs00, scale_const::const_lhv0, scale_const::const_lhv00, scale_const::const_setup(), scale_const::const_thermodyn_type, cp_ice, cp_vapor, cp_water, cv_ice, cv_vapor, cv_water, hyd_dens, lhf, lhs, lhv, and scale_prc::prc_abort().

Referenced by scale_atmos_adiabat::atmos_adiabat_setup(), scale_atmos_saturation::atmos_saturation_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_hydrometeor_regist()

subroutine, public scale_atmos_hydrometeor::atmos_hydrometeor_regist ( integer, intent(in)  NL,
integer, intent(in)  NI,
character(len=*), dimension(1+nl+ni), intent(in)  NAME,
character(len=*), dimension(1+nl+ni), intent(in)  DESC,
character(len=*), dimension(1+nl+ni), intent(in)  UNIT,
integer, intent(out)  Q0,
logical, dimension(1+nl+ni), intent(in), optional  ADVC 
)

ATMOS_HYDROMETEOR_regist Regist tracer.

Parameters
[in]nlnumber of liquid water tracers
[in]ninumber of ice water tracers

Definition at line 237 of file scale_atmos_hydrometeor.F90.

237  use scale_prc, only: &
238  prc_abort
239  use scale_tracer, only: &
241  use scale_const, only: &
242  rvap => const_rvap
243  implicit none
244  integer, intent(in) :: NL
245  integer, intent(in) :: NI
246  character(len=*), intent(in) :: NAME(1+NL+NI)
247  character(len=*), intent(in) :: DESC(1+NL+NI)
248  character(len=*), intent(in) :: UNIT(1+NL+NI)
249 
250  integer, intent(out) :: Q0
251 
252  logical, intent(in), optional :: ADVC(1+NL+NI)
253 
254  real(RP) :: CV (1+NL+NI)
255  real(RP) :: CP (1+NL+NI)
256  real(RP) :: R (1+NL+NI)
257  real(RP) :: EI0 (1+NL+NI)
258  logical :: MASS (1+NL+NI)
259  logical :: ADVC_(1+NL+NI)
260 
261  integer :: NQ
262  integer :: n
263  !---------------------------------------------------------------------------
264 
265  if ( .not. atmos_hydrometeor_dry ) then
266  log_error("ATMOS_HYDROMETEOR_regist",*) 'tracer for hydrometeor is already registerd'
267  call prc_abort
268  endif
269 
270  atmos_hydrometeor_dry = .false.
271 
272  nq = 0
273 
274  ! vapor
275  nq = nq + 1
276  cv(nq) = cv_vapor
277  cp(nq) = cp_vapor
278  r(nq) = rvap
279  ei0(nq) = lhv
280 
281  ! liquid
282  do n = 1, nl
283  nq = nq + 1
284  cv(nq) = cv_water
285  cp(nq) = cp_water
286  r(nq) = 0.0_rp
287  ei0(nq) = 0.0_rp
288  end do
289 
290  ! ice
291  do n = 1, ni
292  nq = nq + 1
293  cv(nq) = cv_ice
294  cp(nq) = cp_ice
295  r(nq) = 0.0_rp
296  ei0(nq) = - lhf
297  end do
298 
299  ! NQ = 1 + NL + NI, vapor + liqid + ice
300 
301  if ( present(advc) ) then
302  advc_(:) = advc(:)
303  else
304  advc_(:) = .true.
305  endif
306 
307  do n = 1, nq
308  mass(n) = .true.
309  end do
310 
311  call tracer_regist( q0, & ! [OUT]
312  nq, & ! [IN]
313  name, & ! [IN]
314  desc, & ! [IN]
315  unit, & ! [IN]
316  cv, cp, r, ei0, & ! [IN], optional
317  advc_, mass ) ! [IN], optional
318 
319  i_qv = q0
320 
321  if ( nq > 1 ) then
322  qhs = q0 + 1
323  qha = nl + ni
324  qhe = qhs + qha - 1
325  endif
326 
327  if ( nl > 0 ) then
328  qls = q0 + 1
329  qla = nl
330  qle = qls + nl - 1
331  endif
332 
333  if ( ni > 0 ) then
334  qis = qle + 1
335  qia = ni
336  qie = qis + ni - 1
337  endif
338 
339  atmos_hydrometeor_ice_phase = qia > 0
340 
341  return

References atmos_hydrometeor_dry, atmos_hydrometeor_ice_phase, scale_const::const_rvap, cp_ice, cp_vapor, cp_water, cv_ice, cv_vapor, cv_water, i_qv, lhf, lhv, scale_prc::prc_abort(), qha, qhe, qhs, qia, qie, qis, qla, qle, qls, and scale_tracer::tracer_regist().

Referenced by mod_atmos_driver::atmos_driver_tracer_setup(), and mod_atmos_phy_mp_driver::atmos_phy_mp_driver_tracer_setup().

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

◆ atmos_hydrometeor_lhv_0d()

subroutine scale_atmos_hydrometeor::atmos_hydrometeor_lhv_0d ( real(rp), intent(in)  temp,
real(rp), intent(out)  lhv 
)
Parameters
[in]temptemperature [K]
[out]lhvlatent heat of vaporization [J/kg]

Definition at line 348 of file scale_atmos_hydrometeor.F90.

348  use scale_const, only: &
349  tem00 => const_tem00, &
350  lhv0 => const_lhv0
351  implicit none
352  real(RP), intent(in) :: temp
353  real(RP), intent(out) :: lhv
354  !---------------------------------------------------------------------------
355 
356  lhv = lhv0 + ( cp_vapor - cp_water ) * ( temp - tem00 ) * thermodyn_emask
357 
358  return

References scale_const::const_lhv0, scale_const::const_tem00, cp_vapor, and cp_water.

Referenced by atmos_hydrometeor_lhv_1d().

Here is the caller graph for this function:

◆ atmos_hydrometeor_lhv_1d()

subroutine scale_atmos_hydrometeor::atmos_hydrometeor_lhv_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(out)  lhv 
)
Parameters
[in]temptemperature [K]
[out]lhvlatent heat of vaporization [J/kg]

Definition at line 367 of file scale_atmos_hydrometeor.F90.

367  implicit none
368  integer, intent(in) :: KA, KS, KE
369 
370  real(RP), intent(in) :: temp(KA)
371 
372  real(RP), intent(out) :: lhv (KA)
373 
374  integer :: k
375  !---------------------------------------------------------------------------
376 
377  do k = ks, ke
378  call atmos_hydrometeor_lhv_0d( temp(k), lhv(k) )
379  enddo
380 
381  return

References atmos_hydrometeor_lhv_0d(), scale_const::const_lhs0, scale_const::const_tem00, cp_ice, and cp_vapor.

Here is the call graph for this function:

◆ atmos_hydrometeor_lhs_1d()

subroutine scale_atmos_hydrometeor::atmos_hydrometeor_lhs_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(out)  lhs 
)
Parameters
[in]temptemperature [K]
[out]lhslatent heat of sublimation [J/kg]

Definition at line 461 of file scale_atmos_hydrometeor.F90.

461  implicit none
462  integer, intent(in) :: KA, KS, KE
463 
464  real(RP), intent(in) :: temp(KA)
465 
466  real(RP), intent(out) :: lhs (KA)
467 
468  integer :: k
469  !---------------------------------------------------------------------------
470 
471  do k = ks, ke
472  call atmos_hydrometeor_lhs( temp(k), lhs(k) )
473  enddo
474 
475  return

References scale_const::const_lhf0, scale_const::const_tem00, cp_ice, and cp_water.

◆ atmos_hydrometeor_lhf_1d()

subroutine scale_atmos_hydrometeor::atmos_hydrometeor_lhf_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(out)  lhf 
)
Parameters
[in]temptemperature [K]
[out]lhflatent heat of fusion [J/kg]

Definition at line 555 of file scale_atmos_hydrometeor.F90.

555  implicit none
556  integer, intent(in) :: KA, KS, KE
557 
558  real(RP), intent(in) :: temp(KA)
559  real(RP), intent(out) :: lhf (KA)
560 
561  integer :: k
562  !---------------------------------------------------------------------------
563 
564  do k = ks, ke
565  call atmos_hydrometeor_lhf( temp(k), lhf(k) )
566  enddo
567 
568  return

References scale_const::const_lhf0, scale_const::const_lhv0, scale_const::const_pre00, scale_const::const_psat0, scale_const::const_rdry, scale_const::const_rvap, and scale_const::const_tem00.

◆ atmos_hydrometeor_entr_2d()

subroutine scale_atmos_hydrometeor::atmos_hydrometeor_entr_2d ( integer, intent(in)  IA,
integer, intent(in)  IS,
integer, intent(in)  IE,
integer, intent(in)  JA,
integer, intent(in)  JS,
integer, intent(in)  JE,
real(rp), dimension (ia,ja), intent(in)  TEMP,
real(rp), dimension (ia,ja), intent(in)  PRES,
real(rp), dimension (ia,ja), intent(in)  QV,
real(rp), dimension (ia,ja), intent(in)  QI,
real(rp), dimension (ia,ja), intent(in)  Qdry,
real(rp), dimension (ia,ja), intent(in)  Rtot,
real(rp), dimension(ia,ja), intent(in)  CPtot,
real(rp), dimension(ia,ja), intent(out)  entr 
)

calc temp, pres, q -> entropy (2D)

Parameters
[in]temptemperature [K]
[in]prespressure [Pa]
[in]qvwater vapor mass concentration [kg/kg]
[in]qiice water mass concentration [kg/kg]
[in]qdrydry air mass concentration [kg/kg]
[in]rtotgas constant
[in]cptotspecific heat
[out]entrentropy [J/K]

Definition at line 675 of file scale_atmos_hydrometeor.F90.

675  implicit none
676  integer, intent(in) :: IA, IS, IE
677  integer, intent(in) :: JA, JS, JE
678 
679  real(RP), intent(in) :: TEMP (IA,JA)
680  real(RP), intent(in) :: PRES (IA,JA)
681  real(RP), intent(in) :: QV (IA,JA)
682  real(RP), intent(in) :: QI (IA,JA)
683  real(RP), intent(in) :: Qdry (IA,JA)
684  real(RP), intent(in) :: Rtot (IA,JA)
685  real(RP), intent(in) :: CPtot(IA,JA)
686 
687  real(RP), intent(out) :: entr(IA,JA)
688 
689  integer :: i, j
690  !---------------------------------------------------------------------------
691 
692  !$omp parallel do OMP_SCHEDULE_ collapse(2)
693  do j = js, je
694  do i = is, ie
695  call atmos_hydrometeor_entr_0d( &
696  temp(i,j), pres(i,j), & ! [IN]
697  qv(i,j), qi(i,j), qdry(i,j), & ! [IN]
698  rtot(i,j), cptot(i,j), & ! [IN]
699  entr(i,j) ) ! [OUT]
700  enddo
701  enddo
702 
703  return

References scale_const::const_lhf0, scale_const::const_lhv0, scale_const::const_pre00, scale_const::const_psat0, scale_const::const_rdry, scale_const::const_rvap, and scale_const::const_tem00.

Variable Documentation

◆ atmos_hydrometeor_ice_phase

logical, public scale_atmos_hydrometeor::atmos_hydrometeor_ice_phase

Definition at line 74 of file scale_atmos_hydrometeor.F90.

74  logical, public :: ATMOS_HYDROMETEOR_ice_phase

Referenced by atmos_hydrometeor_regist(), and mod_atmos_phy_cp_driver::atmos_phy_cp_driver_setup().

◆ i_qv

integer, public scale_atmos_hydrometeor::i_qv = -1

◆ n_hyd

integer, parameter, public scale_atmos_hydrometeor::n_hyd = 6

Definition at line 79 of file scale_atmos_hydrometeor.F90.

79  integer, public, parameter :: N_HYD = 6

Referenced by mod_atmos_bnd_driver::atmos_boundary_driver_send(), mod_atmos_bnd_driver::atmos_boundary_set_online(), mod_atmos_phy_cp_driver::atmos_phy_cp_driver_calc_tendency(), scale_atmos_phy_cp_kf::atmos_phy_cp_kf_tendency(), mod_atmos_phy_cp_vars::atmos_phy_cp_vars_check(), mod_atmos_phy_cp_vars::atmos_phy_cp_vars_fillhalo(), mod_atmos_phy_cp_vars::atmos_phy_cp_vars_restart_def_var(), mod_atmos_phy_cp_vars::atmos_phy_cp_vars_restart_read(), mod_atmos_phy_cp_vars::atmos_phy_cp_vars_restart_write(), mod_atmos_phy_cp_vars::atmos_phy_cp_vars_setup(), mod_atmos_phy_mp_driver::atmos_phy_mp_driver_calc_tendency(), mod_atmos_phy_mp_driver::atmos_phy_mp_driver_qhyd2qtrc(), mod_atmos_phy_mp_driver::atmos_phy_mp_driver_qhyd2qtrc_onlyqv(), mod_atmos_phy_mp_driver::atmos_phy_mp_driver_setup(), scale_atmos_phy_mp_kessler::atmos_phy_mp_kessler_effective_radius(), scale_atmos_phy_mp_kessler::atmos_phy_mp_kessler_qhyd2qtrc(), scale_atmos_phy_mp_kessler::atmos_phy_mp_kessler_qtrc2qhyd(), scale_atmos_phy_mp_sn14::atmos_phy_mp_sn14_effective_radius(), scale_atmos_phy_mp_sn14::atmos_phy_mp_sn14_qhyd2qtrc(), scale_atmos_phy_mp_sn14::atmos_phy_mp_sn14_qtrc2nhyd(), scale_atmos_phy_mp_sn14::atmos_phy_mp_sn14_qtrc2qhyd(), scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_crg_qtrc2qhyd(), scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_effective_radius(), scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_qhyd2qtrc(), scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_qtrc2nhyd(), scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_qtrc2qhyd(), scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08_adjustment(), scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08_effective_radius(), scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08_qhyd2qtrc(), scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08_qtrc2qhyd(), mod_atmos_phy_mp_vars::atmos_phy_mp_vars_get_diagnostic(), mod_atmos_phy_mp_vars::atmos_phy_mp_vars_history(), mod_atmos_phy_mp_vars::atmos_phy_mp_vars_setup(), mod_atmos_phy_rd_driver::atmos_phy_rd_driver_calc_tendency(), scale_atmos_phy_rd_mstrnx::atmos_phy_rd_mstrnx_flux(), scale_atmos_phy_rd_mstrnx::atmos_phy_rd_mstrnx_setup(), mod_atmos_vars::atmos_vars_calc_diagnostics(), mod_atmos_vars::atmos_vars_setup(), scale_comm_cartesc_nest::comm_cartesc_nest_domain_shape(), mod_mkinit::mkinit(), mod_realinput_grads::parentatmosinputgrads(), mod_realinput_scale::parentatmosinputscale(), mod_realinput_wrfarw::parentatmosinputwrfarw(), and mod_realinput::realinput_surface().

◆ i_hc

integer, parameter, public scale_atmos_hydrometeor::i_hc = 1

liquid water cloud

Definition at line 81 of file scale_atmos_hydrometeor.F90.

81  integer, public, parameter :: I_HC = 1

Referenced by scale_atmos_phy_cp_kf::atmos_phy_cp_kf_tendency(), scale_atmos_phy_mp_kessler::atmos_phy_mp_kessler_effective_radius(), scale_atmos_phy_mp_kessler::atmos_phy_mp_kessler_qhyd2qtrc(), scale_atmos_phy_mp_kessler::atmos_phy_mp_kessler_qtrc2qhyd(), scale_atmos_phy_mp_sn14::atmos_phy_mp_sn14_effective_radius(), scale_atmos_phy_mp_sn14::atmos_phy_mp_sn14_qhyd2qtrc(), scale_atmos_phy_mp_sn14::atmos_phy_mp_sn14_qtrc2nhyd(), scale_atmos_phy_mp_sn14::atmos_phy_mp_sn14_qtrc2qhyd(), scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_crg_qtrc2qhyd(), scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_effective_radius(), scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_qhyd2qtrc(), scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_qtrc2nhyd(), scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_qtrc2qhyd(), scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_setup(), scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08_effective_radius(), scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08_qhyd2qtrc(), scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08_qtrc2qhyd(), mod_atmos_phy_mp_vars::atmos_phy_mp_vars_get_diagnostic(), mod_atmos_phy_rd_driver::atmos_phy_rd_driver_calc_tendency(), scale_atmos_phy_rd_mstrnx::atmos_phy_rd_mstrnx_flux(), mod_atmos_vars::atmos_vars_setup(), mod_mkinit::mkinit(), mod_realinput_grads::parentatmosinputgrads(), and mod_realinput_wrfarw::parentatmosinputwrfarw().

◆ i_hr

integer, parameter, public scale_atmos_hydrometeor::i_hr = 2

◆ i_hi

integer, parameter, public scale_atmos_hydrometeor::i_hi = 3

◆ i_hs

integer, parameter, public scale_atmos_hydrometeor::i_hs = 4

◆ i_hg

integer, parameter, public scale_atmos_hydrometeor::i_hg = 5

◆ i_hh

integer, parameter, public scale_atmos_hydrometeor::i_hh = 6

◆ hyd_name

character(len=h_short), dimension(n_hyd), parameter, public scale_atmos_hydrometeor::hyd_name = (/ "QC", "QR", "QI", "QS", "QG", "QH" /)

Definition at line 88 of file scale_atmos_hydrometeor.F90.

88  character(len=H_SHORT), public, parameter :: HYD_NAME(N_HYD) = &
89  (/ "QC", "QR", "QI", "QS", "QG", "QH" /)

Referenced by mod_atmos_phy_cp_driver::atmos_phy_cp_driver_calc_tendency(), mod_atmos_phy_cp_vars::atmos_phy_cp_vars_setup(), mod_atmos_phy_mp_vars::atmos_phy_mp_vars_setup(), and mod_realinput_scale::parentatmosinputscale().

◆ hyd_desc

character(len=h_mid), dimension(n_hyd), parameter, public scale_atmos_hydrometeor::hyd_desc = (/ "cloud ", "rain ", "ice water", "snow ", "graupel ", "hail " /)

Definition at line 90 of file scale_atmos_hydrometeor.F90.

90  character(len=H_MID), public, parameter :: HYD_DESC(N_HYD) = &
91  (/ "cloud ", "rain ", "ice water", "snow ", "graupel ", "hail " /)

Referenced by mod_atmos_phy_mp_vars::atmos_phy_mp_vars_setup().

◆ num_name

character(len=h_short), dimension(n_hyd), parameter, public scale_atmos_hydrometeor::num_name = (/ "NC", "NR", "NI", "NS", "NG", "NH" /)

Definition at line 92 of file scale_atmos_hydrometeor.F90.

92  character(len=H_SHORT), public, parameter :: NUM_NAME(N_HYD) = &
93  (/ "NC", "NR", "NI", "NS", "NG", "NH" /)

Referenced by mod_atmos_phy_mp_vars::atmos_phy_mp_vars_setup(), and mod_realinput_scale::parentatmosinputscale().

◆ hyd_dens

real(rp), dimension(n_hyd), public scale_atmos_hydrometeor::hyd_dens

Definition at line 95 of file scale_atmos_hydrometeor.F90.

95  real(RP), public :: HYD_DENS(N_HYD)

Referenced by atmos_hydrometeor_setup(), and scale_atmos_phy_rd_mstrnx::atmos_phy_rd_mstrnx_flux().

◆ atmos_hydrometeor_dry

logical, public scale_atmos_hydrometeor::atmos_hydrometeor_dry = .true.

◆ qha

integer, public scale_atmos_hydrometeor::qha = 0

◆ qhs

integer, public scale_atmos_hydrometeor::qhs = -1

◆ qhe

integer, public scale_atmos_hydrometeor::qhe = -2

◆ qla

integer, public scale_atmos_hydrometeor::qla = 0

◆ qls

integer, public scale_atmos_hydrometeor::qls = -1

◆ qle

integer, public scale_atmos_hydrometeor::qle = -2

◆ qia

integer, public scale_atmos_hydrometeor::qia = 0

◆ qis

integer, public scale_atmos_hydrometeor::qis = -1

Definition at line 123 of file scale_atmos_hydrometeor.F90.

123  integer, public :: QIS = -1

Referenced by atmos_hydrometeor_regist().

◆ qie

integer, public scale_atmos_hydrometeor::qie = -2

Definition at line 124 of file scale_atmos_hydrometeor.F90.

124  integer, public :: QIE = -2

Referenced by atmos_hydrometeor_regist().

◆ lhv

real(rp), public scale_atmos_hydrometeor::lhv

◆ lhs

real(rp), public scale_atmos_hydrometeor::lhs

latent heat of sublimation for use [J/kg]

Definition at line 127 of file scale_atmos_hydrometeor.F90.

127  real(RP), public :: LHS

Referenced by atmos_hydrometeor_setup(), and scale_atmos_phy_cp_kf::atmos_phy_cp_kf_tendency().

◆ lhf

real(rp), public scale_atmos_hydrometeor::lhf

◆ cv_vapor

real(rp), public scale_atmos_hydrometeor::cv_vapor

◆ cp_vapor

real(rp), public scale_atmos_hydrometeor::cp_vapor

◆ cv_water

real(rp), public scale_atmos_hydrometeor::cv_water

CV for water [J/kg/K].

Definition at line 132 of file scale_atmos_hydrometeor.F90.

132  real(RP), public :: CV_WATER

Referenced by atmos_hydrometeor_regist(), atmos_hydrometeor_setup(), scale_atmos_hydrostatic::atmos_hydrostatic_buildrho_1d(), scale_atmos_hydrostatic::atmos_hydrostatic_buildrho_atmos_0d(), scale_atmos_hydrostatic::atmos_hydrostatic_buildrho_atmos_1d(), scale_atmos_hydrostatic::atmos_hydrostatic_buildrho_atmos_rev_1d(), scale_atmos_hydrostatic::atmos_hydrostatic_buildrho_bytemp_3d(), scale_atmos_hydrostatic::atmos_hydrostatic_buildrho_bytemp_atmos_rev_1d(), scale_atmos_hydrostatic::atmos_hydrostatic_buildrho_real_3d(), mod_atmos_phy_cp_driver::atmos_phy_cp_driver_calc_tendency(), scale_atmos_phy_cp_kf::atmos_phy_cp_kf_tendency(), scale_atmos_phy_mp_kessler::atmos_phy_mp_kessler_adjustment(), scale_atmos_phy_mp_common::atmos_phy_mp_negative_fixer(), scale_atmos_phy_mp_common::atmos_phy_mp_precipitation_semilag(), scale_atmos_phy_mp_common::atmos_phy_mp_precipitation_upwind(), scale_atmos_phy_mp_common::atmos_phy_mp_saturation_adjustment_3d(), scale_atmos_phy_mp_sn14::atmos_phy_mp_sn14_terminal_velocity(), scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_crg_qtrc2qhyd(), scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08_adjustment(), scale_atmos_saturation::atmos_saturation_moist_conversion_dens_all_0d(), scale_atmos_saturation::atmos_saturation_pote_1d(), scale_cpl_phy_sfc_skin::cpl_phy_sfc_skin(), mod_cpl_vars::cpl_putatm(), scale_atmos_phy_mp_sn14::freezing_water(), mod_land_driver::land_driver_calc_tendency(), scale_land_dyn_bucket::land_dyn_bucket(), scale_land_dyn_bucket::land_dyn_bucket_setup(), mod_land_vars::land_vars_monitor(), mod_ocean_driver::ocean_driver_calc_tendency(), scale_ocean_dyn_slab::ocean_dyn_slab(), scale_ocean_dyn_slab::ocean_dyn_slab_setup(), scale_ocean_phy_ice_simple::ocean_phy_ice_adjustment(), scale_ocean_phy_ice_simple::ocean_phy_ice_simple(), mod_ocean_vars::ocean_vars_monitor(), and mod_urban_vars::urban_vars_monitor().

◆ cp_water

real(rp), public scale_atmos_hydrometeor::cp_water

CP for water [J/kg/K].

Definition at line 133 of file scale_atmos_hydrometeor.F90.

133  real(RP), public :: CP_WATER

Referenced by scale_atmos_adiabat::atmos_adiabat_liftparcel_1d(), atmos_hydrometeor_lhs_1d(), atmos_hydrometeor_lhv_0d(), atmos_hydrometeor_regist(), atmos_hydrometeor_setup(), scale_atmos_hydrostatic::atmos_hydrostatic_buildrho_1d(), scale_atmos_hydrostatic::atmos_hydrostatic_buildrho_atmos_0d(), scale_atmos_hydrostatic::atmos_hydrostatic_buildrho_atmos_1d(), scale_atmos_hydrostatic::atmos_hydrostatic_buildrho_atmos_rev_1d(), scale_atmos_hydrostatic::atmos_hydrostatic_buildrho_bytemp_3d(), scale_atmos_hydrostatic::atmos_hydrostatic_buildrho_bytemp_atmos_rev_1d(), scale_atmos_hydrostatic::atmos_hydrostatic_buildrho_real_3d(), scale_atmos_phy_cp_kf::atmos_phy_cp_kf_tendency(), scale_atmos_phy_mp_kessler::atmos_phy_mp_kessler_adjustment(), scale_atmos_phy_mp_common::atmos_phy_mp_negative_fixer(), scale_atmos_phy_mp_common::atmos_phy_mp_precipitation_semilag(), scale_atmos_phy_mp_common::atmos_phy_mp_precipitation_upwind(), scale_atmos_phy_mp_common::atmos_phy_mp_saturation_adjustment_3d(), scale_atmos_phy_mp_sn14::atmos_phy_mp_sn14_terminal_velocity(), scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_crg_qtrc2qhyd(), scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08_adjustment(), scale_atmos_saturation::atmos_saturation_moist_conversion_dens_all_0d(), scale_atmos_saturation::atmos_saturation_moist_conversion_pres_liq_0d(), scale_atmos_saturation::atmos_saturation_pote_1d(), and scale_atmos_phy_mp_sn14::freezing_water().

◆ cv_ice

real(rp), public scale_atmos_hydrometeor::cv_ice

◆ cp_ice

real(rp), public scale_atmos_hydrometeor::cp_ice
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_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_rvap
real(rp), parameter, public const_rvap
specific gas constant (water vapor) [J/kg/K]
Definition: scale_const.F90:63
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_const
module CONSTANT
Definition: scale_const.F90:11
scale_tracer::tracer_regist
subroutine, public tracer_regist(QS, NQ, NAME, DESC, UNIT, CV, CP, R, ENGI0, ADVC, MASS)
Regist tracer.
Definition: scale_tracer.F90:65
scale_tracer
module TRACER
Definition: scale_tracer.F90:12
scale_const::const_ci
real(rp), parameter, public const_ci
specific heat (ice) [J/kg/K]
Definition: scale_const.F90:67
scale_const::const_dwatr
real(rp), parameter, public const_dwatr
density of water [kg/m3]
Definition: scale_const.F90:82
scale_const::const_tem00
real(rp), parameter, public const_tem00
temperature reference (0C) [K]
Definition: scale_const.F90:90
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_const::const_lhf0
real(rp), public const_lhf0
latent heat of fusion at 0C [J/kg]
Definition: scale_const.F90:79
scale_const::const_lhs0
real(rp), parameter, public const_lhs0
latent heat of sublimation at 0C [J/kg]
Definition: scale_const.F90:77
scale_const::const_dice
real(rp), parameter, public const_dice
density of ice [kg/m3]
Definition: scale_const.F90:83
scale_const::const_setup
subroutine, public const_setup
Setup.
Definition: scale_const.F90:115
scale_const::const_lhf00
real(rp), public const_lhf00
latent heat of fusion at 0K [J/kg]
Definition: scale_const.F90:80