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

module Hydrometeor More...

Functions/Subroutines

subroutine, public atmos_hydrometeor_setup
 Setup. More...
 
subroutine, public atmos_hydrometeor_regist (Q0, NV, NL, NI, NAME, DESC, UNIT, ADVC)
 Regist tracer. More...
 
subroutine atmos_hydrometeor_lhv_0d (lhv, temp)
 
subroutine atmos_hydrometeor_lhv_1d (lhv, temp)
 
subroutine atmos_hydrometeor_lhv_2d (lhv, temp)
 
subroutine atmos_hydrometeor_lhv_3d (lhv, temp)
 
subroutine atmos_hydrometeor_lhs_0d (lhs, temp)
 
subroutine atmos_hydrometeor_lhs_1d (lhs, temp)
 
subroutine atmos_hydrometeor_lhs_2d (lhs, temp)
 
subroutine atmos_hydrometeor_lhs_3d (lhs, temp)
 
subroutine atmos_hydrometeor_lhf_0d (lhf, temp)
 
subroutine atmos_hydrometeor_lhf_1d (lhf, temp)
 
subroutine atmos_hydrometeor_lhf_2d (lhf, temp)
 
subroutine atmos_hydrometeor_lhf_3d (lhf, temp)
 
subroutine atmos_hydrometeor_entr_0d (entr, temp, pres, q, Rq)
 calc temp, pres, q -> entropy (0D) More...
 
subroutine atmos_hydrometeor_entr_2d (entr, temp, pres, q, Rq)
 calc temp, pres, q -> entropy (2D) More...
 
subroutine atmos_hydrometeor_entr_3d (entr, temp, pres, q, Rq)
 calc temp, pres, q -> entropy (3D) More...
 
subroutine, public atmos_hydrometeor_diagnose_number_concentration (QTRC)
 

Variables

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...
 
integer, public i_qc = -1
 
integer, public i_qr = -1
 
integer, public i_qi = -1
 
integer, public i_qs = -1
 
integer, public i_qg = -1
 
integer, public i_qh = -1
 
integer, public i_nc = -1
 
integer, public i_nr = -1
 
integer, public i_ni = -1
 
integer, public i_ns = -1
 
integer, public i_ng = -1
 
integer, public i_nh = -1
 
integer, public qhs = -1
 
integer, public qhe = -2
 
integer, public qls = -1
 
integer, public qle = -2
 
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...
 

Detailed Description

module Hydrometeor

Description
Hydrometeor module
Author
Team SCALE
History
  • 2016-08-06 (S.Nishizawa) [new]

Function/Subroutine Documentation

◆ atmos_hydrometeor_setup()

subroutine, public scale_atmos_hydrometeor::atmos_hydrometeor_setup ( )

Setup.

Definition at line 129 of file scale_atmos_sub_hydrometeor.F90.

References scale_const::const_ci, scale_const::const_cl, scale_const::const_cpvap, scale_const::const_cvvap, 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_thermodyn_type, scale_stdio::io_fid_log, scale_stdio::io_l, lhf, lhs, lhv, and scale_process::prc_mpistop().

Referenced by mod_rm_driver::scalerm(), and mod_rm_prep::scalerm_prep().

129  use scale_const, only: &
130  cpvap => const_cpvap, &
131  cvvap => const_cvvap, &
132  cl => const_cl, &
133  ci => const_ci, &
134  lhv00 => const_lhv00, &
135  lhs00 => const_lhs00, &
136  lhf00 => const_lhf00, &
137  lhv0 => const_lhv0, &
138  lhs0 => const_lhs0, &
139  lhf0 => const_lhf0, &
140  thermodyn_type => const_thermodyn_type
141  use scale_process, only: &
143  implicit none
144  !---------------------------------------------------------------------------
145 
146  if( io_l ) write(io_fid_log,*)
147  if( io_l ) write(io_fid_log,*) '++++++ Module[HYDEROMETER] / Categ[ATMOS SHARE] / Origin[SCALElib]'
148 
149  if ( thermodyn_type == 'EXACT' ) then
150 
151  cv_vapor = cvvap
152  cp_vapor = cpvap
153  cv_water = cl
154  cp_water = cv_water
155  cv_ice = ci
156  cp_ice = cv_ice
157 
158  lhv = lhv00
159  lhs = lhs00
160  lhf = lhf00
161  thermodyn_emask = 1.0_rp
162 
163  elseif( thermodyn_type == 'SIMPLE' ) then
164 
165  cv_vapor = cvvap
166  cp_vapor = cpvap
167  cv_water = cvvap
168  cp_water = cv_water
169  cv_ice = cvvap
170  cp_ice = cv_ice
171 
172  lhv = lhv0
173  lhs = lhs0
174  lhf = lhf0
175  thermodyn_emask = 0.0_rp
176 
177  else
178  write(*,*) 'xxx Not appropriate ATMOS_THERMODYN_ENERGY_TYPE. Check!', trim(thermodyn_type)
179  call prc_mpistop
180  endif
181 
182  return
subroutine, public prc_mpistop
Abort MPI.
real(rp), parameter, public const_ci
specific heat (ice) [J/kg/K]
Definition: scale_const.F90:69
real(rp), parameter, public const_cl
specific heat (liquid water) [J/kg/K]
Definition: scale_const.F90:68
real(rp), public const_cvvap
specific heat (water vapor, constant volume) [J/kg/K]
Definition: scale_const.F90:67
real(rp), public const_lhf0
latent heat of fusion at 0C [J/kg]
Definition: scale_const.F90:81
real(rp), parameter, public const_lhs0
latent heat of sublimation at 0C [J/kg]
Definition: scale_const.F90:79
real(rp), parameter, public const_lhv0
latent heat of vaporizaion at 0C [J/kg]
Definition: scale_const.F90:77
real(rp), public const_lhf00
latent heat of fusion at 0K [J/kg]
Definition: scale_const.F90:82
module PROCESS
real(rp), public const_lhs00
latent heat of sublimation at 0K [J/kg]
Definition: scale_const.F90:80
real(rp), public const_lhv00
latent heat of vaporizaion at 0K [J/kg]
Definition: scale_const.F90:78
module CONSTANT
Definition: scale_const.F90:14
real(rp), parameter, public const_cpvap
specific heat (water vapor, constant pressure) [J/kg/K]
Definition: scale_const.F90:66
character(len=h_short), public const_thermodyn_type
internal energy type
Definition: scale_const.F90:98
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(out)  Q0,
integer, intent(in)  NV,
integer, intent(in)  NL,
integer, intent(in)  NI,
character(len=*), dimension(nv+nl+ni), intent(in)  NAME,
character(len=*), dimension(nv+nl+ni), intent(in)  DESC,
character(len=*), dimension(nv+nl+ni), intent(in)  UNIT,
logical, dimension(nv+nl+ni), intent(in), optional  ADVC 
)

Regist tracer.

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

Definition at line 196 of file scale_atmos_sub_hydrometeor.F90.

References scale_const::const_rvap, i_qv, scale_process::prc_mpistop(), qhe, qhs, qie, qis, qle, qls, and scale_tracer::tracer_regist().

Referenced by scale_atmos_phy_mp_kessler::atmos_phy_mp_kessler_config(), scale_atmos_phy_mp_sn14::atmos_phy_mp_sn14_config(), scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_config(), and scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08_config().

196  use scale_process, only: &
198  use scale_tracer, only: &
200  use scale_const, only: &
201  rvap => const_rvap
202  implicit none
203 
204  integer, intent(out) :: Q0
205  integer, intent(in) :: NV
206  integer, intent(in) :: NL
207  integer, intent(in) :: NI
208  character(len=*), intent(in) :: NAME(NV+NL+NI)
209  character(len=*), intent(in) :: DESC(NV+NL+NI)
210  character(len=*), intent(in) :: UNIT(NV+NL+NI)
211 
212  logical, intent(in), optional :: ADVC(NV+NL+NI)
213 
214  real(RP) :: CV (NV+NL+NI)
215  real(RP) :: CP (NV+NL+NI)
216  real(RP) :: R (NV+NL+NI)
217  logical :: MASS (NV+NL+NI)
218  logical :: ADVC_(NV+NL+NI)
219 
220  integer :: NQ
221  integer :: n
222  !---------------------------------------------------------------------------
223 
224  if ( i_qv > 0 ) then
225  write(*,*) 'xxx tracer for hydrometeor is already registerd'
226  call prc_mpistop
227  endif
228 
229  if ( nv /= 1 ) then
230  write(*,*) 'xxx number of vapor must be 1 at this moment'
231  call prc_mpistop
232  endif
233 
234  nq = 0
235 
236  do n = 1, nv
237  nq = nq + 1
238  cv(nq) = cv_vapor
239  cp(nq) = cp_vapor
240  r(nq) = rvap
241  end do
242 
243  do n = 1, nl
244  nq = nq + 1
245  cv(nq) = cv_water
246  cp(nq) = cp_water
247  r(nq) = 0.0_rp
248  end do
249 
250  do n = 1, ni
251  nq = nq + 1
252  cv(nq) = cv_ice
253  cp(nq) = cp_ice
254  r(nq) = 0.0_rp
255  end do
256 
257  ! NQ = 1 + NL + NI, vapor + liqid + ice
258 
259  if ( present(advc) ) then
260  advc_(:) = advc(:)
261  else
262  advc_(:) = .true.
263  endif
264 
265  do n = 1, nq
266  mass(n) = .true.
267  end do
268 
269  call tracer_regist( q0, & ! [OUT]
270  nq, & ! [IN]
271  name, & ! [IN]
272  desc, & ! [IN]
273  unit, & ! [IN]
274  cv, cp, r, & ! [IN], optional
275  advc_, mass ) ! [IN], optional
276 
277  i_qv = q0
278 
279  if ( nq > 1 ) then
280  qhs = i_qv + 1
281  qhe = qhs + nl + ni - 1
282  endif
283 
284  if ( nl > 0 ) then
285  qls = i_qv + 1
286  qle = qls + nl - 1
287  endif
288 
289  if ( ni > 0 ) then
290  qis = qle + 1
291  qie = qis + ni - 1
292  endif
293 
294  return
subroutine, public prc_mpistop
Abort MPI.
module TRACER
module PROCESS
real(rp), parameter, public const_rvap
specific gas constant (water vapor) [J/kg/K]
Definition: scale_const.F90:65
module CONSTANT
Definition: scale_const.F90:14
subroutine, public tracer_regist(QS, NQ, NAME, DESC, UNIT, CV, CP, R, ADVC, MASS)
Regist tracer.
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(out)  lhv,
real(rp), intent(in)  temp 
)
Parameters
[out]lhvlatent heat of vaporization [J/kg]
[in]temptemperature [K]

Definition at line 301 of file scale_atmos_sub_hydrometeor.F90.

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

301  use scale_const, only: &
302  tem00 => const_tem00, &
303  lhv0 => const_lhv0
304  implicit none
305 
306  real(RP), intent(out) :: lhv
307  real(RP), intent(in) :: temp
308  !---------------------------------------------------------------------------
309 
310  lhv = lhv0 + ( cp_vapor - cp_water ) * ( temp - tem00 ) * thermodyn_emask
311 
312  return
real(rp), parameter, public const_tem00
temperature reference (0C) [K]
Definition: scale_const.F90:92
real(rp), parameter, public const_lhv0
latent heat of vaporizaion at 0C [J/kg]
Definition: scale_const.F90:77
module CONSTANT
Definition: scale_const.F90:14

◆ atmos_hydrometeor_lhv_1d()

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

Definition at line 319 of file scale_atmos_sub_hydrometeor.F90.

References scale_const::const_lhv0, scale_const::const_tem00, scale_grid_index::ke, and scale_grid_index::ks.

319  use scale_const, only: &
320  tem00 => const_tem00, &
321  lhv0 => const_lhv0
322  implicit none
323 
324  real(RP), intent(out) :: lhv (KA)
325  real(RP), intent(in) :: temp(KA)
326 
327  integer :: k
328  !---------------------------------------------------------------------------
329 
330  do k = ks, ke
331  lhv(k) = lhv0 + ( cp_vapor - cp_water ) * ( temp(k) - tem00 ) * thermodyn_emask
332  enddo
333 
334  return
integer, public ke
end point of inner domain: z, local
real(rp), parameter, public const_tem00
temperature reference (0C) [K]
Definition: scale_const.F90:92
real(rp), parameter, public const_lhv0
latent heat of vaporizaion at 0C [J/kg]
Definition: scale_const.F90:77
module CONSTANT
Definition: scale_const.F90:14
integer, public ks
start point of inner domain: z, local

◆ atmos_hydrometeor_lhv_2d()

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

Definition at line 341 of file scale_atmos_sub_hydrometeor.F90.

References scale_const::const_lhv0, scale_const::const_tem00, scale_grid_index::ieb, scale_grid_index::isb, scale_grid_index::jeb, and scale_grid_index::jsb.

341  use scale_const, only: &
342  tem00 => const_tem00, &
343  lhv0 => const_lhv0
344  implicit none
345 
346  real(RP), intent(out) :: lhv (IA,JA)
347  real(RP), intent(in) :: temp(IA,JA)
348 
349  integer :: i, j
350  !---------------------------------------------------------------------------
351 
352  do j = jsb, jeb
353  do i = isb, ieb
354  lhv(i,j) = lhv0 + ( cp_vapor - cp_water ) * ( temp(i,j) - tem00 ) * thermodyn_emask
355  enddo
356  enddo
357 
358  return
integer, public jeb
real(rp), parameter, public const_tem00
temperature reference (0C) [K]
Definition: scale_const.F90:92
integer, public ieb
real(rp), parameter, public const_lhv0
latent heat of vaporizaion at 0C [J/kg]
Definition: scale_const.F90:77
module CONSTANT
Definition: scale_const.F90:14
integer, public isb
integer, public jsb

◆ atmos_hydrometeor_lhv_3d()

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

Definition at line 365 of file scale_atmos_sub_hydrometeor.F90.

References scale_const::const_lhv0, scale_const::const_tem00, scale_grid_index::ieb, scale_grid_index::isb, scale_grid_index::jeb, scale_grid_index::jsb, scale_grid_index::ke, and scale_grid_index::ks.

365  use scale_const, only: &
366  tem00 => const_tem00, &
367  lhv0 => const_lhv0
368  implicit none
369 
370  real(RP), intent(out) :: lhv (KA,IA,JA)
371  real(RP), intent(in) :: temp(KA,IA,JA)
372 
373  integer :: k, i, j
374  !---------------------------------------------------------------------------
375 
376  do j = jsb, jeb
377  do i = isb, ieb
378  do k = ks, ke
379  lhv(k,i,j) = lhv0 + ( cp_vapor - cp_water ) * ( temp(k,i,j) - tem00 ) * thermodyn_emask
380  enddo
381  enddo
382  enddo
383 
384  return
integer, public jeb
integer, public ke
end point of inner domain: z, local
real(rp), parameter, public const_tem00
temperature reference (0C) [K]
Definition: scale_const.F90:92
integer, public ieb
real(rp), parameter, public const_lhv0
latent heat of vaporizaion at 0C [J/kg]
Definition: scale_const.F90:77
module CONSTANT
Definition: scale_const.F90:14
integer, public ks
start point of inner domain: z, local
integer, public isb
integer, public jsb

◆ atmos_hydrometeor_lhs_0d()

subroutine scale_atmos_hydrometeor::atmos_hydrometeor_lhs_0d ( real(rp), intent(out)  lhs,
real(rp), intent(in)  temp 
)
Parameters
[out]lhslatent heat of sublimation [J/kg]
[in]temptemperature [K]

Definition at line 391 of file scale_atmos_sub_hydrometeor.F90.

References scale_const::const_lhs0, and scale_const::const_tem00.

391  use scale_const, only: &
392  tem00 => const_tem00, &
393  lhs0 => const_lhs0
394  implicit none
395 
396  real(RP), intent(out) :: lhs
397  real(RP), intent(in) :: temp
398  !---------------------------------------------------------------------------
399 
400  lhs = lhs0 + ( cp_vapor - cp_ice ) * ( temp - tem00 ) * thermodyn_emask
401 
402  return
real(rp), parameter, public const_tem00
temperature reference (0C) [K]
Definition: scale_const.F90:92
real(rp), parameter, public const_lhs0
latent heat of sublimation at 0C [J/kg]
Definition: scale_const.F90:79
module CONSTANT
Definition: scale_const.F90:14

◆ atmos_hydrometeor_lhs_1d()

subroutine scale_atmos_hydrometeor::atmos_hydrometeor_lhs_1d ( real(rp), dimension (ka), intent(out)  lhs,
real(rp), dimension(ka), intent(in)  temp 
)
Parameters
[out]lhslatent heat of sublimation [J/kg]
[in]temptemperature [K]

Definition at line 409 of file scale_atmos_sub_hydrometeor.F90.

References scale_const::const_lhs0, scale_const::const_tem00, scale_grid_index::ke, and scale_grid_index::ks.

409  use scale_const, only: &
410  tem00 => const_tem00, &
411  lhs0 => const_lhs0
412  implicit none
413 
414  real(RP), intent(out) :: lhs (KA)
415  real(RP), intent(in) :: temp(KA)
416 
417  integer :: k
418  !---------------------------------------------------------------------------
419 
420  do k = ks, ke
421  lhs(k) = lhs0 + ( cp_vapor - cp_ice ) * ( temp(k) - tem00 ) * thermodyn_emask
422  enddo
423 
424  return
integer, public ke
end point of inner domain: z, local
real(rp), parameter, public const_tem00
temperature reference (0C) [K]
Definition: scale_const.F90:92
real(rp), parameter, public const_lhs0
latent heat of sublimation at 0C [J/kg]
Definition: scale_const.F90:79
module CONSTANT
Definition: scale_const.F90:14
integer, public ks
start point of inner domain: z, local

◆ atmos_hydrometeor_lhs_2d()

subroutine scale_atmos_hydrometeor::atmos_hydrometeor_lhs_2d ( real(rp), dimension (ia,ja), intent(out)  lhs,
real(rp), dimension(ia,ja), intent(in)  temp 
)
Parameters
[out]lhslatent heat of sublimation [J/kg]
[in]temptemperature [K]

Definition at line 431 of file scale_atmos_sub_hydrometeor.F90.

References scale_const::const_lhs0, scale_const::const_tem00, scale_grid_index::ieb, scale_grid_index::isb, scale_grid_index::jeb, and scale_grid_index::jsb.

431  use scale_const, only: &
432  tem00 => const_tem00, &
433  lhs0 => const_lhs0
434  implicit none
435 
436  real(RP), intent(out) :: lhs (IA,JA)
437  real(RP), intent(in) :: temp(IA,JA)
438 
439  integer :: i, j
440  !---------------------------------------------------------------------------
441 
442  do j = jsb, jeb
443  do i = isb, ieb
444  lhs(i,j) = lhs0 + ( cp_vapor - cp_ice ) * ( temp(i,j) - tem00 ) * thermodyn_emask
445  enddo
446  enddo
447 
448  return
integer, public jeb
real(rp), parameter, public const_tem00
temperature reference (0C) [K]
Definition: scale_const.F90:92
integer, public ieb
real(rp), parameter, public const_lhs0
latent heat of sublimation at 0C [J/kg]
Definition: scale_const.F90:79
module CONSTANT
Definition: scale_const.F90:14
integer, public isb
integer, public jsb

◆ atmos_hydrometeor_lhs_3d()

subroutine scale_atmos_hydrometeor::atmos_hydrometeor_lhs_3d ( real(rp), dimension (ka,ia,ja), intent(out)  lhs,
real(rp), dimension(ka,ia,ja), intent(in)  temp 
)
Parameters
[out]lhslatent heat of sublimation [J/kg]
[in]temptemperature [K]

Definition at line 455 of file scale_atmos_sub_hydrometeor.F90.

References scale_const::const_lhs0, scale_const::const_tem00, scale_grid_index::ieb, scale_grid_index::isb, scale_grid_index::jeb, scale_grid_index::jsb, scale_grid_index::ke, and scale_grid_index::ks.

455  use scale_const, only: &
456  tem00 => const_tem00, &
457  lhs0 => const_lhs0
458  implicit none
459 
460  real(RP), intent(out) :: lhs (KA,IA,JA)
461  real(RP), intent(in) :: temp(KA,IA,JA)
462 
463  integer :: k, i, j
464  !---------------------------------------------------------------------------
465 
466  do j = jsb, jeb
467  do i = isb, ieb
468  do k = ks, ke
469  lhs(k,i,j) = lhs0 + ( cp_vapor - cp_ice ) * ( temp(k,i,j) - tem00 ) * thermodyn_emask
470  enddo
471  enddo
472  enddo
473 
474  return
integer, public jeb
integer, public ke
end point of inner domain: z, local
real(rp), parameter, public const_tem00
temperature reference (0C) [K]
Definition: scale_const.F90:92
integer, public ieb
real(rp), parameter, public const_lhs0
latent heat of sublimation at 0C [J/kg]
Definition: scale_const.F90:79
module CONSTANT
Definition: scale_const.F90:14
integer, public ks
start point of inner domain: z, local
integer, public isb
integer, public jsb

◆ atmos_hydrometeor_lhf_0d()

subroutine scale_atmos_hydrometeor::atmos_hydrometeor_lhf_0d ( real(rp), intent(out)  lhf,
real(rp), intent(in)  temp 
)
Parameters
[out]lhflatent heat of fusion [J/kg]
[in]temptemperature [K]

Definition at line 481 of file scale_atmos_sub_hydrometeor.F90.

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

481  use scale_const, only: &
482  tem00 => const_tem00, &
483  lhf0 => const_lhf0
484  implicit none
485 
486  real(RP), intent(out) :: lhf
487  real(RP), intent(in) :: temp
488  !---------------------------------------------------------------------------
489 
490  lhf = lhf0 + ( cp_water - cp_ice ) * ( temp - tem00 ) * thermodyn_emask
491 
492  return
real(rp), parameter, public const_tem00
temperature reference (0C) [K]
Definition: scale_const.F90:92
real(rp), public const_lhf0
latent heat of fusion at 0C [J/kg]
Definition: scale_const.F90:81
module CONSTANT
Definition: scale_const.F90:14

◆ atmos_hydrometeor_lhf_1d()

subroutine scale_atmos_hydrometeor::atmos_hydrometeor_lhf_1d ( real(rp), dimension (ka), intent(out)  lhf,
real(rp), dimension(ka), intent(in)  temp 
)
Parameters
[out]lhflatent heat of fusion [J/kg]
[in]temptemperature [K]

Definition at line 499 of file scale_atmos_sub_hydrometeor.F90.

References scale_const::const_lhf0, scale_const::const_tem00, scale_grid_index::ke, and scale_grid_index::ks.

499  use scale_const, only: &
500  tem00 => const_tem00, &
501  lhf0 => const_lhf0
502  implicit none
503 
504  real(RP), intent(out) :: lhf (KA)
505  real(RP), intent(in) :: temp(KA)
506 
507  integer :: k
508  !---------------------------------------------------------------------------
509 
510  do k = ks, ke
511  lhf(k) = lhf0 + ( cp_water - cp_ice ) * ( temp(k) - tem00 ) * thermodyn_emask
512  enddo
513 
514  return
integer, public ke
end point of inner domain: z, local
real(rp), parameter, public const_tem00
temperature reference (0C) [K]
Definition: scale_const.F90:92
real(rp), public const_lhf0
latent heat of fusion at 0C [J/kg]
Definition: scale_const.F90:81
module CONSTANT
Definition: scale_const.F90:14
integer, public ks
start point of inner domain: z, local

◆ atmos_hydrometeor_lhf_2d()

subroutine scale_atmos_hydrometeor::atmos_hydrometeor_lhf_2d ( real(rp), dimension (ia,ja), intent(out)  lhf,
real(rp), dimension(ia,ja), intent(in)  temp 
)
Parameters
[out]lhflatent heat of fusion [J/kg]
[in]temptemperature [K]

Definition at line 521 of file scale_atmos_sub_hydrometeor.F90.

References scale_const::const_lhf0, scale_const::const_tem00, scale_grid_index::ieb, scale_grid_index::isb, scale_grid_index::jeb, and scale_grid_index::jsb.

521  use scale_const, only: &
522  tem00 => const_tem00, &
523  lhf0 => const_lhf0
524  implicit none
525 
526  real(RP), intent(out) :: lhf (IA,JA)
527  real(RP), intent(in) :: temp(IA,JA)
528 
529  integer :: i, j
530  !---------------------------------------------------------------------------
531 
532  do j = jsb, jeb
533  do i = isb, ieb
534  lhf(i,j) = lhf0 + ( cp_water - cp_ice ) * ( temp(i,j) - tem00 ) * thermodyn_emask
535  enddo
536  enddo
537 
538  return
integer, public jeb
real(rp), parameter, public const_tem00
temperature reference (0C) [K]
Definition: scale_const.F90:92
real(rp), public const_lhf0
latent heat of fusion at 0C [J/kg]
Definition: scale_const.F90:81
integer, public ieb
module CONSTANT
Definition: scale_const.F90:14
integer, public isb
integer, public jsb

◆ atmos_hydrometeor_lhf_3d()

subroutine scale_atmos_hydrometeor::atmos_hydrometeor_lhf_3d ( real(rp), dimension (ka,ia,ja), intent(out)  lhf,
real(rp), dimension(ka,ia,ja), intent(in)  temp 
)
Parameters
[out]lhflatent heat of fusion [J/kg]
[in]temptemperature [K]

Definition at line 545 of file scale_atmos_sub_hydrometeor.F90.

References scale_const::const_lhf0, scale_const::const_tem00, scale_grid_index::ieb, scale_grid_index::isb, scale_grid_index::jeb, scale_grid_index::jsb, scale_grid_index::ke, and scale_grid_index::ks.

545  use scale_const, only: &
546  tem00 => const_tem00, &
547  lhf0 => const_lhf0
548  implicit none
549 
550  real(RP), intent(out) :: lhf (KA,IA,JA)
551  real(RP), intent(in) :: temp(KA,IA,JA)
552 
553  integer :: k, i, j
554  !---------------------------------------------------------------------------
555 
556  do j = jsb, jeb
557  do i = isb, ieb
558  do k = ks, ke
559  lhf(k,i,j) = lhf0 + ( cp_water - cp_ice ) * ( temp(k,i,j) - tem00 ) * thermodyn_emask
560  enddo
561  enddo
562  enddo
563 
564  return
integer, public jeb
integer, public ke
end point of inner domain: z, local
real(rp), parameter, public const_tem00
temperature reference (0C) [K]
Definition: scale_const.F90:92
real(rp), public const_lhf0
latent heat of fusion at 0C [J/kg]
Definition: scale_const.F90:81
integer, public ieb
module CONSTANT
Definition: scale_const.F90:14
integer, public ks
start point of inner domain: z, local
integer, public isb
integer, public jsb

◆ atmos_hydrometeor_entr_0d()

subroutine scale_atmos_hydrometeor::atmos_hydrometeor_entr_0d ( real(rp), intent(out)  entr,
real(rp), intent(in)  temp,
real(rp), intent(in)  pres,
real(rp), dimension(qa), intent(in)  q,
real(rp), dimension(qa), intent(in)  Rq 
)

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

Parameters
[out]entrentropy [J/K]
[in]temptemperature [K]
[in]prespressure [Pa]
[in]qmass concentration [kg/kg]
[in]rqgas constantt [J/kg/K]

Definition at line 575 of file scale_atmos_sub_hydrometeor.F90.

References scale_const::const_cpdry, scale_const::const_eps, scale_const::const_lhf0, scale_const::const_lhv0, scale_const::const_pre00, scale_const::const_psat0, scale_const::const_rdry, scale_const::const_rvap, scale_const::const_tem00, i_qv, scale_tracer::qa, qie, qis, qle, and qls.

575  use scale_const, only: &
576  eps => const_eps, &
577  pre00 => const_pre00, &
578  tem00 => const_tem00, &
579  rdry => const_rdry, &
580  cpdry => const_cpdry, &
581  rvap => const_rvap, &
582  lhv0 => const_lhv0, &
583  lhf0 => const_lhf0, &
584  psat0 => const_psat0
585  use scale_tracer, only: &
586  qa
587  implicit none
588 
589  real(RP), intent(out) :: entr
590  real(RP), intent(in) :: temp
591  real(RP), intent(in) :: pres
592  real(RP), intent(in) :: q(QA)
593  real(RP), intent(in) :: Rq(QA)
594 
595  real(RP) :: qdry, Rtot
596  real(RP) :: logT_T0, pres_dry, pres_vap
597 
598  integer :: iqw
599  !---------------------------------------------------------------------------
600 
601  logt_t0 = log( temp / tem00 )
602 
603  qdry = 1.0_rp
604  rtot = 0.0_rp
605  do iqw = 1, qa
606  qdry = qdry - q(iqw)
607  rtot = rtot + q(iqw) * rq(iqw)
608  enddo
609  rtot = rtot + rdry * qdry
610 
611  ! dry air + vapor
612  pres_dry = max( pres * qdry * rdry / rtot, eps )
613  entr = qdry * cpdry * logt_t0 &
614  - qdry * rdry * log( pres_dry / pre00 )
615 
616  if ( i_qv > 0 ) then
617  pres_vap = max( pres * q(i_qv) * rvap / rtot, eps )
618  entr = entr + q(i_qv) * cp_vapor * logt_t0 &
619  - q(i_qv) * rvap * log( pres_vap / psat0 ) &
620  + q(i_qv) * lhv0 / tem00
621  endif
622 
623  ! liquid water
624  if ( qls > 0 ) then
625  do iqw = qls, qle
626  entr = entr + q(iqw) * cp_water * logt_t0
627  enddo
628  endif
629 
630  ! ice
631  if ( qis > 0 ) then
632  do iqw = qis, qie
633  entr = entr + q(iqw) * cp_ice * logt_t0 &
634  - q(iqw) * lhf0 / tem00
635  enddo
636  endif
637 
638  return
real(rp), parameter, public const_psat0
saturate pressure of water vapor at 0C [Pa]
Definition: scale_const.F90:83
real(rp), public const_cpdry
specific heat (dry air,constant pressure) [J/kg/K]
Definition: scale_const.F90:58
real(rp), parameter, public const_tem00
temperature reference (0C) [K]
Definition: scale_const.F90:92
integer, public qa
real(rp), public const_lhf0
latent heat of fusion at 0C [J/kg]
Definition: scale_const.F90:81
real(rp), public const_rdry
specific gas constant (dry air) [J/kg/K]
Definition: scale_const.F90:57
module TRACER
real(rp), parameter, public const_lhv0
latent heat of vaporizaion at 0C [J/kg]
Definition: scale_const.F90:77
real(rp), public const_pre00
pressure reference [Pa]
Definition: scale_const.F90:90
real(rp), parameter, public const_rvap
specific gas constant (water vapor) [J/kg/K]
Definition: scale_const.F90:65
module CONSTANT
Definition: scale_const.F90:14
real(rp), public const_eps
small number
Definition: scale_const.F90:36

◆ atmos_hydrometeor_entr_2d()

subroutine scale_atmos_hydrometeor::atmos_hydrometeor_entr_2d ( real(rp), dimension(ia,ja), intent(out)  entr,
real(rp), dimension(ia,ja), intent(in)  temp,
real(rp), dimension(ia,ja), intent(in)  pres,
real(rp), dimension (ia,ja,qa), intent(in)  q,
real(rp), dimension (qa), intent(in)  Rq 
)

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

Parameters
[out]entrentropy [J/K]
[in]temptemperature [K]
[in]prespressure [Pa]
[in]qmass concentration [kg/kg]
[in]rqgas constant [J/kg/K]

Definition at line 649 of file scale_atmos_sub_hydrometeor.F90.

References scale_const::const_cpdry, scale_const::const_eps, scale_const::const_lhf0, scale_const::const_lhv0, scale_const::const_pre00, scale_const::const_psat0, scale_const::const_rdry, scale_const::const_rvap, scale_const::const_tem00, i_qv, scale_grid_index::ieb, scale_grid_index::isb, scale_grid_index::jeb, scale_grid_index::jsb, scale_tracer::qa, qie, qis, qle, and qls.

649  use scale_const, only: &
650  eps => const_eps, &
651  pre00 => const_pre00, &
652  tem00 => const_tem00, &
653  cpdry => const_cpdry, &
654  rdry => const_rdry, &
655  rvap => const_rvap, &
656  lhv0 => const_lhv0, &
657  lhf0 => const_lhf0, &
658  psat0 => const_psat0
659  use scale_tracer, only: &
660  qa
661  implicit none
662 
663  real(RP), intent(out) :: entr(IA,JA)
664  real(RP), intent(in) :: temp(IA,JA)
665  real(RP), intent(in) :: pres(IA,JA)
666  real(RP), intent(in) :: q (IA,JA,QA)
667  real(RP), intent(in) :: Rq (QA)
668 
669  real(RP) :: qdry, Rtot
670  real(RP) :: logT_T0, pres_dry, pres_vap
671 
672  integer :: i, j, iqw
673  !---------------------------------------------------------------------------
674 
675  ! dry air + vapor
676  do j = jsb, jeb
677  do i = isb, ieb
678 
679  logt_t0 = log( temp(i,j) / tem00 )
680 
681  qdry = 1.0_rp
682  rtot = 0.0_rp
683  do iqw = 1, qa
684  qdry = qdry - q(i,j,iqw)
685  rtot = rtot + q(i,j,iqw) * rq(iqw)
686  enddo
687  rtot = rtot + rdry * qdry
688 
689  ! dry air + vapor
690  pres_dry = max( pres(i,j) * qdry * rdry / rtot, eps )
691  entr(i,j) = qdry * cpdry * logt_t0 &
692  - qdry * rdry * log( pres_dry / pre00 )
693 
694  if ( i_qv > 0 ) then
695  pres_vap = max( pres(i,j) * q(i,j,i_qv) * rvap / rtot, eps )
696  entr(i,j) = entr(i,j) + q(i,j,i_qv) * cp_vapor * logt_t0 &
697  - q(i,j,i_qv) * rvap * log( pres_vap / psat0 ) &
698  + q(i,j,i_qv) * lhv0 / tem00
699  endif
700 
701  ! liquid water
702  if ( qls > 0 ) then
703  do iqw = qls, qle
704  entr(i,j) = entr(i,j) + q(i,j,iqw) * cp_water * logt_t0
705  enddo
706  endif
707 
708  ! ice
709  if ( qis > 0 ) then
710  do iqw = qis, qie
711  entr(i,j) = entr(i,j) + q(i,j,iqw) * cp_ice * logt_t0 &
712  - q(i,j,iqw) * lhf0 / tem00
713  enddo
714  endif
715 
716  enddo
717  enddo
718 
719  return
real(rp), parameter, public const_psat0
saturate pressure of water vapor at 0C [Pa]
Definition: scale_const.F90:83
real(rp), public const_cpdry
specific heat (dry air,constant pressure) [J/kg/K]
Definition: scale_const.F90:58
integer, public jeb
real(rp), parameter, public const_tem00
temperature reference (0C) [K]
Definition: scale_const.F90:92
real(rp), public const_lhf0
latent heat of fusion at 0C [J/kg]
Definition: scale_const.F90:81
real(rp), public const_rdry
specific gas constant (dry air) [J/kg/K]
Definition: scale_const.F90:57
integer, public ieb
module TRACER
real(rp), parameter, public const_lhv0
latent heat of vaporizaion at 0C [J/kg]
Definition: scale_const.F90:77
real(rp), public const_pre00
pressure reference [Pa]
Definition: scale_const.F90:90
real(rp), parameter, public const_rvap
specific gas constant (water vapor) [J/kg/K]
Definition: scale_const.F90:65
module CONSTANT
Definition: scale_const.F90:14
real(rp), public const_eps
small number
Definition: scale_const.F90:36
integer, public isb
integer, public jsb

◆ atmos_hydrometeor_entr_3d()

subroutine scale_atmos_hydrometeor::atmos_hydrometeor_entr_3d ( real(rp), dimension(ka,ia,ja), intent(out)  entr,
real(rp), dimension(ka,ia,ja), intent(in)  temp,
real(rp), dimension(ka,ia,ja), intent(in)  pres,
real(rp), dimension (ka,ia,ja,qa), intent(in)  q,
real(rp), dimension (qa), intent(in)  Rq 
)

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

Parameters
[out]entrentropy [J/K]
[in]temptemperature [K]
[in]prespressure [Pa]
[in]qmass concentration [kg/kg]
[in]rqgas constant [J/kg/K]

Definition at line 730 of file scale_atmos_sub_hydrometeor.F90.

References scale_const::const_cpdry, scale_const::const_eps, scale_const::const_lhf0, scale_const::const_lhv0, scale_const::const_pre00, scale_const::const_psat0, scale_const::const_rdry, scale_const::const_rvap, scale_const::const_tem00, i_qv, scale_grid_index::ieb, scale_grid_index::isb, scale_grid_index::jeb, scale_grid_index::jsb, scale_grid_index::ke, scale_grid_index::ks, scale_tracer::qa, qie, qis, qle, and qls.

730  use scale_const, only: &
731  eps => const_eps, &
732  pre00 => const_pre00, &
733  tem00 => const_tem00, &
734  cpdry => const_cpdry, &
735  rdry => const_rdry, &
736  rvap => const_rvap, &
737  lhv0 => const_lhv0, &
738  lhf0 => const_lhf0, &
739  psat0 => const_psat0
740  use scale_tracer, only: &
741  qa
742  implicit none
743 
744  real(RP), intent(out) :: entr(KA,IA,JA)
745  real(RP), intent(in) :: temp(KA,IA,JA)
746  real(RP), intent(in) :: pres(KA,IA,JA)
747  real(RP), intent(in) :: q (KA,IA,JA,QA)
748  real(RP), intent(in) :: Rq (QA)
749 
750  real(RP) :: qdry, Rtot
751  real(RP) :: logT_T0, pres_dry, pres_vap
752 
753  integer :: k, i, j, iqw
754  !---------------------------------------------------------------------------
755 
756  ! dry air + vapor
757  do j = jsb, jeb
758  do i = isb, ieb
759  do k = ks, ke
760 
761  logt_t0 = log( temp(k,i,j) / tem00 )
762 
763  qdry = 1.0_rp
764  rtot = 0.0_rp
765  do iqw = 1, qa
766  qdry = qdry - q(k,i,j,iqw)
767  rtot = rtot + q(k,i,j,iqw) * rq(iqw)
768  enddo
769  rtot = rtot + rdry * qdry
770 
771  ! dry air + vapor
772  pres_dry = max( pres(k,i,j) * qdry * rdry / rtot, eps )
773  entr(k,i,j) = qdry * cpdry * logt_t0 &
774  - qdry * rdry * log( pres_dry / pre00 )
775 
776  if ( i_qv > 0 ) then
777  pres_vap = max( pres(k,i,j) * q(k,i,j,i_qv) * rvap / rtot, eps )
778  entr(k,i,j) = entr(k,i,j) + q(k,i,j,i_qv) * cp_vapor * logt_t0 &
779  - q(k,i,j,i_qv) * rvap * log( pres_vap / psat0 ) &
780  + q(k,i,j,i_qv) * lhv0 / tem00
781  endif
782 
783  ! liquid water
784  if ( qls > 0 ) then
785  do iqw = qls, qle
786  entr(k,i,j) = entr(k,i,j) + q(k,i,j,iqw) * cp_water * logt_t0
787  enddo
788  endif
789 
790  ! ice
791  if ( qis > 0 ) then
792  do iqw = qis, qie
793  entr(k,i,j) = entr(k,i,j) + q(k,i,j,iqw) * cp_ice * logt_t0 &
794  - q(k,i,j,iqw) * lhf0 / tem00
795  enddo
796  endif
797 
798  enddo
799  enddo
800  enddo
801 
802  return
real(rp), parameter, public const_psat0
saturate pressure of water vapor at 0C [Pa]
Definition: scale_const.F90:83
real(rp), public const_cpdry
specific heat (dry air,constant pressure) [J/kg/K]
Definition: scale_const.F90:58
integer, public jeb
integer, public ke
end point of inner domain: z, local
real(rp), parameter, public const_tem00
temperature reference (0C) [K]
Definition: scale_const.F90:92
real(rp), public const_lhf0
latent heat of fusion at 0C [J/kg]
Definition: scale_const.F90:81
real(rp), public const_rdry
specific gas constant (dry air) [J/kg/K]
Definition: scale_const.F90:57
integer, public ieb
module TRACER
real(rp), parameter, public const_lhv0
latent heat of vaporizaion at 0C [J/kg]
Definition: scale_const.F90:77
real(rp), public const_pre00
pressure reference [Pa]
Definition: scale_const.F90:90
real(rp), parameter, public const_rvap
specific gas constant (water vapor) [J/kg/K]
Definition: scale_const.F90:65
module CONSTANT
Definition: scale_const.F90:14
integer, public ks
start point of inner domain: z, local
real(rp), public const_eps
small number
Definition: scale_const.F90:36
integer, public isb
integer, public jsb

◆ atmos_hydrometeor_diagnose_number_concentration()

subroutine, public scale_atmos_hydrometeor::atmos_hydrometeor_diagnose_number_concentration ( real(rp), dimension(:,:,:,:), intent(inout)  QTRC)

Definition at line 808 of file scale_atmos_sub_hydrometeor.F90.

References scale_const::const_pi, i_nc, i_ng, i_ni, i_nr, i_ns, i_qc, i_qg, i_qi, i_qr, and i_qs.

Referenced by scale_atmos_boundary::atmos_boundary_update(), mod_atmos_phy_cp_driver::atmos_phy_cp_driver(), and mod_realinput::parentatomsetup().

808  use scale_const, only: &
809  pi => const_pi
810  implicit none
811 
812  real(RP), intent(inout) :: QTRC(:,:,:,:)
813 
814  real(RP), parameter :: Dc = 20.e-6_rp ! typical particle diameter for cloud [m]
815  real(RP), parameter :: Dr = 200.e-6_rp ! typical particle diameter for rain [m]
816  real(RP), parameter :: Di = 80.e-6_rp ! typical particle diameter for ice [m]
817  real(RP), parameter :: Ds = 80.e-6_rp ! typical particle diameter for snow [m]
818  real(RP), parameter :: Dg = 200.e-6_rp ! typical particle diameter for grapel [m]
819  real(RP), parameter :: RHOw = 1000.0_rp ! typical density for warm particles [kg/m3]
820  real(RP), parameter :: RHOf = 100.0_rp ! typical density for frozen particles [kg/m3]
821  real(RP), parameter :: RHOg = 400.0_rp ! typical density for grapel particles [kg/m3]
822  real(RP), parameter :: b = 3.0_rp ! assume spherical form
823 
824  real(RP) :: piov6
825  !---------------------------------------------------------------------------
826 
827  piov6 = pi / 6.0_rp
828 
829  if ( i_nc > 0 ) qtrc(:,:,:,i_nc) = qtrc(:,:,:,i_qc) / ( (piov6*rhow) * dc**b )
830  if ( i_nr > 0 ) qtrc(:,:,:,i_nr) = qtrc(:,:,:,i_qr) / ( (piov6*rhow) * dr**b )
831  if ( i_ni > 0 ) qtrc(:,:,:,i_ni) = qtrc(:,:,:,i_qi) / ( (piov6*rhof) * di**b )
832  if ( i_ns > 0 ) qtrc(:,:,:,i_ns) = qtrc(:,:,:,i_qs) / ( (piov6*rhof) * ds**b )
833  if ( i_ng > 0 ) qtrc(:,:,:,i_ng) = qtrc(:,:,:,i_qg) / ( (piov6*rhog) * dg**b )
834 
835  return
module CONSTANT
Definition: scale_const.F90:14
real(rp), public const_pi
pi
Definition: scale_const.F90:34
Here is the caller graph for this function:

Variable Documentation

◆ i_qv

integer, public scale_atmos_hydrometeor::i_qv = -1

Definition at line 68 of file scale_atmos_sub_hydrometeor.F90.

Referenced by scale_atmos_adiabat::atmos_adiabat_cape(), scale_atmos_adiabat::atmos_adiabat_liftparcel(), scale_atmos_dyn_common::atmos_dyn_numfilter_coef_q(), atmos_hydrometeor_entr_0d(), atmos_hydrometeor_entr_2d(), atmos_hydrometeor_entr_3d(), atmos_hydrometeor_regist(), scale_atmos_hydrostatic::atmos_hydrostatic_setup(), scale_atmos_phy_ae_kajino13::atmos_phy_ae_kajino13(), scale_atmos_phy_ae_kajino13::atmos_phy_ae_kajino13_mkinit(), scale_atmos_phy_cp_kf::atmos_phy_cp_kf_setup(), scale_atmos_phy_mp_convert::atmos_phy_mp_bulk2bin(), scale_atmos_phy_mp_kessler::atmos_phy_mp_kessler(), scale_atmos_phy_mp_sn14::atmos_phy_mp_sn14(), scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08(), mod_atmos_phy_rd_driver::atmos_phy_rd_driver(), scale_atmos_phy_rd_mstrnx::atmos_phy_rd_mstrnx(), scale_atmos_phy_sf_bulk::atmos_phy_sf_bulk(), scale_atmos_phy_sf_const::atmos_phy_sf_const(), mod_atmos_phy_sf_driver::atmos_phy_sf_driver(), scale_atmos_phy_tb_mynn::atmos_phy_tb_mynn(), scale_atmos_refstate::atmos_refstate_update(), mod_atmos_vars::atmos_vars_history(), mod_atmos_vars::atmos_vars_monitor(), mod_atmos_vars::atmos_vars_total(), scale_atmos_phy_cp_kf::cp_kf_main(), mod_cpl_vars::cpl_getsfc_atm(), mod_cpl_vars::cpl_putatm(), mod_cpl_vars::cpl_vars_setup(), mod_mkinit::interporation_fact(), scale_atmos_phy_mp_sn14::mp_negativefilter(), scale_grid_nest::nest_setup(), mod_realinput_grads::parentatominputgrads(), mod_realinput_nicam::parentatominputnicam(), mod_realinput_wrfarw::parentatominputwrfarw(), mod_realinput::parentatomsetup(), and mod_mkinit::read_sounding().

68  integer, public :: I_QV = -1

◆ n_hyd

integer, parameter, public scale_atmos_hydrometeor::n_hyd = 6

◆ i_hc

integer, parameter, public scale_atmos_hydrometeor::i_hc = 1

◆ 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

◆ i_qc

integer, public scale_atmos_hydrometeor::i_qc = -1

◆ i_qr

integer, public scale_atmos_hydrometeor::i_qr = -1

◆ i_qi

integer, public scale_atmos_hydrometeor::i_qi = -1

◆ i_qs

integer, public scale_atmos_hydrometeor::i_qs = -1

◆ i_qg

integer, public scale_atmos_hydrometeor::i_qg = -1

◆ i_qh

integer, public scale_atmos_hydrometeor::i_qh = -1

Definition at line 84 of file scale_atmos_sub_hydrometeor.F90.

84  integer, public :: I_QH = -1

◆ i_nc

integer, public scale_atmos_hydrometeor::i_nc = -1

◆ i_nr

integer, public scale_atmos_hydrometeor::i_nr = -1

◆ i_ni

integer, public scale_atmos_hydrometeor::i_ni = -1

◆ i_ns

integer, public scale_atmos_hydrometeor::i_ns = -1

◆ i_ng

integer, public scale_atmos_hydrometeor::i_ng = -1

◆ i_nh

integer, public scale_atmos_hydrometeor::i_nh = -1

Definition at line 91 of file scale_atmos_sub_hydrometeor.F90.

91  integer, public :: I_NH = -1

◆ qhs

integer, public scale_atmos_hydrometeor::qhs = -1

◆ qhe

integer, public scale_atmos_hydrometeor::qhe = -2

◆ qls

integer, public scale_atmos_hydrometeor::qls = -1

◆ qle

integer, public scale_atmos_hydrometeor::qle = -2

◆ qis

integer, public scale_atmos_hydrometeor::qis = -1

◆ qie

integer, public scale_atmos_hydrometeor::qie = -2

◆ 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 104 of file scale_atmos_sub_hydrometeor.F90.

Referenced by atmos_hydrometeor_setup().

104  real(RP), public :: LHS

◆ lhf

real(rp), public scale_atmos_hydrometeor::lhf