SCALE-RM
Functions/Subroutines | Variables
scale_atmos_phy_mp_kessler Module Reference

module atmosphere / physics / microphysics / Kessler More...

Functions/Subroutines

subroutine, public atmos_phy_mp_kessler_setup
 ATMOS_PHY_MP_kessler_setup Setup. More...
 
subroutine, public atmos_phy_mp_kessler_adjustment (KA, KS, KE, IA, IS, IE, JA, JS, JE, DENS, PRES, dt, TEMP, QTRC, CPtot, CVtot, RHOE_t, EVAPORATE)
 ATMOS_PHY_MP_kessler_adjustment calculate state after saturation process. More...
 
subroutine, public atmos_phy_mp_kessler_terminal_velocity (KA, KS, KE, DENS0, RHOQ0, REFSTATE_dens_profile, vterm)
 Kessler-type warm rain microphysics (terminal velocity) More...
 
subroutine, public atmos_phy_mp_kessler_cloud_fraction (KA, KS, KE, IA, IS, IE, JA, JS, JE, QTRC, mask_criterion, cldfrac)
 Calculate Cloud Fraction. More...
 
subroutine, public atmos_phy_mp_kessler_effective_radius (KA, KS, KE, IA, IS, IE, JA, JS, JE, DENS0, TEMP0, QTRC0, Re)
 Calculate Effective Radius. More...
 
subroutine, public atmos_phy_mp_kessler_qtrc2qhyd (KA, KS, KE, IA, IS, IE, JA, JS, JE, QTRC, Qe)
 Calculate mass ratio of each category. More...
 
subroutine, public atmos_phy_mp_kessler_qhyd2qtrc (KA, KS, KE, IA, IS, IE, JA, JS, JE, Qe, QTRC)
 Calculate mass ratio of each category. More...
 

Variables

integer, parameter, public atmos_phy_mp_kessler_ntracers = QA_MP
 
integer, parameter, public atmos_phy_mp_kessler_nwaters = 2
 
integer, parameter, public atmos_phy_mp_kessler_nices = 0
 
character(len=h_short), dimension(qa_mp), parameter, public atmos_phy_mp_kessler_tracer_names = (/ 'QV', 'QC', 'QR' /)
 
character(len=h_mid), dimension(qa_mp), parameter, public atmos_phy_mp_kessler_tracer_descriptions = (/ 'Ratio of Water Vapor mass to total mass (Specific humidity)', 'Ratio of Cloud Water mass to total mass ', 'Ratio of Rain Water mass to total mass '/)
 
character(len=h_short), dimension(qa_mp), parameter, public atmos_phy_mp_kessler_tracer_units = (/ 'kg/kg', 'kg/kg', 'kg/kg' /)
 

Detailed Description

module atmosphere / physics / microphysics / Kessler

Description
Cloud Microphysics by Kessler-type parametarization Reference: Kessler(1969) Klemp and Wilhelmson(1978)
Author
Team SCALE
NAMELIST
  • No namelist group
History Output
namedescriptionunitvariable
Pcsat QC production term by satadjust kg/kg/s QC_t_sat

Function/Subroutine Documentation

◆ atmos_phy_mp_kessler_setup()

subroutine, public scale_atmos_phy_mp_kessler::atmos_phy_mp_kessler_setup ( )

ATMOS_PHY_MP_kessler_setup Setup.

Definition at line 90 of file scale_atmos_phy_mp_kessler.F90.

References scale_prc::prc_abort().

Referenced by mod_atmos_phy_mp_driver::atmos_phy_mp_driver_setup().

90  use scale_prc, only: &
91  prc_abort
92  implicit none
93 
94  !---------------------------------------------------------------------------
95 
96  log_newline
97  log_info("ATMOS_PHY_MP_kessler_setup",*) 'Setup'
98  log_info("ATMOS_PHY_MP_kessler_setup",*) 'KESSLER-type 1-moment bulk 3 category'
99 
100  if( couple_aerosol ) then
101  log_error("ATMOS_PHY_MP_kessler_setup",*) 'MP_aerosol_couple should be .false. for KESSLER type MP!'
102  call prc_abort
103  endif
104 
105  return
module PROCESS
Definition: scale_prc.F90:11
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:338
Here is the call graph for this function:
Here is the caller graph for this function:

◆ atmos_phy_mp_kessler_adjustment()

subroutine, public scale_atmos_phy_mp_kessler::atmos_phy_mp_kessler_adjustment ( integer, intent(in)  KA,
integer, intent(in)  KS,
integer, intent(in)  KE,
integer, intent(in)  IA,
integer, intent(in)  IS,
integer, intent(in)  IE,
integer, intent(in)  JA,
integer, intent(in)  JS,
integer, intent(in)  JE,
real(rp), dimension (ka,ia,ja), intent(in)  DENS,
real(rp), dimension (ka,ia,ja), intent(in)  PRES,
real(dp), intent(in)  dt,
real(rp), dimension (ka,ia,ja), intent(inout)  TEMP,
real(rp), dimension (ka,ia,ja,qa_mp), intent(inout)  QTRC,
real(rp), dimension (ka,ia,ja), intent(inout)  CPtot,
real(rp), dimension (ka,ia,ja), intent(inout)  CVtot,
real(rp), dimension (ka,ia,ja), intent(out)  RHOE_t,
real(rp), dimension(ka,ia,ja), intent(out)  EVAPORATE 
)

ATMOS_PHY_MP_kessler_adjustment calculate state after saturation process.

Definition at line 118 of file scale_atmos_phy_mp_kessler.F90.

References scale_const::const_dwatr, scale_const::const_eps, scale_const::const_pi, scale_atmos_hydrometeor::cp_vapor, scale_atmos_hydrometeor::cp_water, scale_atmos_hydrometeor::cv_vapor, scale_atmos_hydrometeor::cv_water, scale_atmos_hydrometeor::lhv, scale_prof::prof_rapend(), and scale_prof::prof_rapstart().

Referenced by mod_atmos_phy_mp_driver::atmos_phy_mp_driver_calc_tendency().

118  use scale_const, only: &
119  dwatr => const_dwatr, &
120  pi => const_pi
121  use scale_file_history, only: &
122  file_history_in
123  use scale_atmos_phy_mp_common, only: &
124  mp_saturation_adjustment => atmos_phy_mp_saturation_adjustment
125  implicit none
126 
127  integer, intent(in) :: ka, ks, ke
128  integer, intent(in) :: ia, is, ie
129  integer, intent(in) :: ja, js, je
130  real(RP), intent(in) :: dens (ka,ia,ja) ! density [kg/m3]
131  real(RP), intent(in) :: pres (ka,ia,ja) ! pressure [Pa]
132  real(DP), intent(in) :: dt ! time interval of microphysics [s]
133  real(RP), intent(inout) :: temp (ka,ia,ja) ! temperature [K]
134  real(RP), intent(inout) :: qtrc (ka,ia,ja,qa_mp) ! tracer mass concentration [kg/kg]
135  real(RP), intent(inout) :: cptot (ka,ia,ja) ! total specific heat capacity at constant pressure [J/kg/K]
136  real(RP), intent(inout) :: cvtot (ka,ia,ja) ! total specific heat capacity at constant volume [J/kg/K]
137  real(RP), intent(out) :: rhoe_t (ka,ia,ja) ! tendency of rhoe [J/m3/s]
138  real(RP), intent(out) :: evaporate(ka,ia,ja) ! number of evaporated cloud [#/m3/s]
139 
140  real(RP) :: qtrc_dummy(ka,ia,ja)
141  real(RP) :: rhoe_d_sat(ka,ia,ja)
142  real(RP) :: qc_t_sat (ka,ia,ja)
143 
144  integer :: k, i, j
145  !---------------------------------------------------------------------------
146 
147  log_progress(*) 'atmosphere / physics / microphysics / Kessler'
148 
149  !##### MP Main #####
150  call mp_kessler( &
151  ka, ks, ke, ia, is, ie, ja, js, je, &
152  dens(:,:,:), pres(:,:,:), & ! [IN]
153  dt, & ! [IN]
154  temp(:,:,:), qtrc(:,:,:,:), & ! [INOUT]
155  cptot(:,:,:), cvtot(:,:,:), & ! [INOUT]
156  rhoe_t(:,:,:) ) ! [OUT]
157 
158  ! save value before saturation adjustment
159  do j = js, je
160  do i = is, ie
161  do k = ks, ke
162  qc_t_sat(k,i,j) = qtrc(k,i,j,i_qc)
163  enddo
164  enddo
165  enddo
166 
167 !OCL XFILL
168  qtrc_dummy(:,:,:) = -1.0_rp
169 
170  call mp_saturation_adjustment( &
171  ka, ks, ke, ia, is, ie, ja, js, je, &
172  dens(:,:,:), & ! [IN]
173  flag_liquid, & ! [IN]
174  temp(:,:,:), & ! [INOUT]
175  qtrc(:,:,:,i_qv), & ! [INOUT]
176  qtrc(:,:,:,i_qc), qtrc_dummy(:,:,:), & ! [INOUT]
177  cptot(:,:,:), cvtot(:,:,:), & ! [INOUT]
178  rhoe_d_sat(:,:,:) ) ! [OUT]
179 
180  do j = js, je
181  do i = is, ie
182  do k = ks, ke
183  rhoe_t(k,i,j) = rhoe_t(k,i,j) + rhoe_d_sat(k,i,j) / dt
184  enddo
185  enddo
186  enddo
187  do j = js, je
188  do i = is, ie
189  do k = ks, ke
190  qc_t_sat(k,i,j) = ( qtrc(k,i,j,i_qc) - qc_t_sat(k,i,j) ) / dt
191  enddo
192  enddo
193  enddo
194 
195  call file_history_in( qc_t_sat(:,:,:), 'Pcsat', 'QC production term by satadjust', 'kg/kg/s' )
196 
197  do j = js, je
198  do i = is, ie
199  do k = ks, ke
200  evaporate(k,i,j) = max( -qc_t_sat(k,i,j), 0.0_rp ) & ! if negative, condensation
201  * dens(k,i,j) / (4.0_rp/3.0_rp*pi*dwatr*re_qc**3) ! mass -> number (assuming constant particle radius as re_qc)
202  enddo
203  enddo
204  enddo
205 
206  !##### END MP Main #####
207 
208  return
real(rp), parameter, public const_dwatr
density of water [kg/m3]
Definition: scale_const.F90:82
integer, public ia
of whole cells: x, local, with HALO
integer, public ja
of whole cells: y, local, with HALO
module ATMOSPHERE / Physics Cloud Microphysics - Common
integer, public is
start point of inner domain: x, local
integer, public ie
end point of inner domain: x, local
integer, public ke
end point of inner domain: z, local
integer, public je
end point of inner domain: y, local
integer, public ks
start point of inner domain: z, local
module CONSTANT
Definition: scale_const.F90:11
integer, public js
start point of inner domain: y, local
integer, public ka
of whole cells: z, local, with HALO
real(rp), public const_pi
pi
Definition: scale_const.F90:31
module file_history
Here is the call graph for this function:
Here is the caller graph for this function:

◆ atmos_phy_mp_kessler_terminal_velocity()

subroutine, public scale_atmos_phy_mp_kessler::atmos_phy_mp_kessler_terminal_velocity ( integer, intent(in)  KA,
integer, intent(in)  KS,
integer, intent(in)  KE,
real(rp), dimension (ka), intent(in)  DENS0,
real(rp), dimension (ka,qa_mp-1), intent(in)  RHOQ0,
real(rp), dimension(ka), intent(in)  REFSTATE_dens_profile,
real(rp), dimension (ka,qa_mp-1), intent(out)  vterm 
)

Kessler-type warm rain microphysics (terminal velocity)

Definition at line 356 of file scale_atmos_phy_mp_kessler.F90.

Referenced by mod_atmos_phy_mp_driver::atmos_phy_mp_driver_calc_tendency().

356  implicit none
357 
358  integer, intent(in) :: ka, ks, ke
359  real(RP), intent(in) :: dens0 (ka) ! density [kg/m3]
360  real(RP), intent(in) :: rhoq0 (ka,qa_mp-1) ! density of each hydrometeor [kg/m3]
361  real(RP), intent(in) :: refstate_dens_profile(ka) ! reference state density profile [kg/m3]
362  real(RP), intent(out) :: vterm (ka,qa_mp-1) ! terminal velocity of each hydrometeor [m/s]
363 
364  real(RP) :: qr
365  real(RP) :: zerosw
366 
367  integer :: k
368  !---------------------------------------------------------------------------
369 
370  do k = ks, ke
371  vterm(k,i_hyd_qc) = 0.0_rp
372  enddo
373 
374  do k = ks, ke
375  qr = rhoq0(k,i_hyd_qr) / dens0(k)
376  zerosw = 0.5_rp - sign(0.5_rp, qr - 1.e-12_rp )
377 
378  vterm(k,i_hyd_qr) = - 36.34_rp * ( dens0(k) * ( qr + zerosw ) )**0.1364_rp &
379  * refstate_dens_profile(ks) / refstate_dens_profile(k) * ( 1.0_rp - zerosw )
380  enddo
381 
382  return
integer, public ke
end point of inner domain: z, local
integer, public ks
start point of inner domain: z, local
integer, public ka
of whole cells: z, local, with HALO
Here is the caller graph for this function:

◆ atmos_phy_mp_kessler_cloud_fraction()

subroutine, public scale_atmos_phy_mp_kessler::atmos_phy_mp_kessler_cloud_fraction ( integer, intent(in)  KA,
integer, intent(in)  KS,
integer, intent(in)  KE,
integer, intent(in)  IA,
integer, intent(in)  IS,
integer, intent(in)  IE,
integer, intent(in)  JA,
integer, intent(in)  JS,
integer, intent(in)  JE,
real(rp), dimension(ka,ia,ja,qa_mp-1), intent(in)  QTRC,
real(rp), intent(in)  mask_criterion,
real(rp), dimension(ka,ia,ja), intent(out)  cldfrac 
)

Calculate Cloud Fraction.

Definition at line 392 of file scale_atmos_phy_mp_kessler.F90.

Referenced by mod_atmos_phy_mp_vars::atmos_phy_mp_vars_get_diagnostic().

392  implicit none
393  integer, intent(in) :: ka, ks, ke
394  integer, intent(in) :: ia, is, ie
395  integer, intent(in) :: ja, js, je
396 
397  real(RP), intent(in) :: qtrc(ka,ia,ja,qa_mp-1) ! hydrometeor mass concentration [kg/kg]
398  real(RP), intent(in) :: mask_criterion ! criterion of hydrometeor [kg/kg]
399 
400  real(RP), intent(out) :: cldfrac(ka,ia,ja) ! cloud fraction (0 or 1)
401 
402  real(RP) :: qhydro
403  integer :: k, i, j
404  !---------------------------------------------------------------------------
405 
406  do j = js, je
407  do i = is, ie
408  do k = ks, ke
409  qhydro = qtrc(k,i,j,i_hyd_qc) + qtrc(k,i,j,i_hyd_qr)
410  cldfrac(k,i,j) = 0.5_rp + sign(0.5_rp,qhydro-mask_criterion)
411  enddo
412  enddo
413  enddo
414 
415  return
integer, public ia
of whole cells: x, local, with HALO
integer, public ja
of whole cells: y, local, with HALO
integer, public is
start point of inner domain: x, local
integer, public ie
end point of inner domain: x, local
integer, public ke
end point of inner domain: z, local
integer, public je
end point of inner domain: y, local
integer, public ks
start point of inner domain: z, local
integer, public js
start point of inner domain: y, local
integer, public ka
of whole cells: z, local, with HALO
Here is the caller graph for this function:

◆ atmos_phy_mp_kessler_effective_radius()

subroutine, public scale_atmos_phy_mp_kessler::atmos_phy_mp_kessler_effective_radius ( integer, intent(in)  KA,
integer, intent(in)  KS,
integer, intent(in)  KE,
integer, intent(in)  IA,
integer, intent(in)  IS,
integer, intent(in)  IE,
integer, intent(in)  JA,
integer, intent(in)  JS,
integer, intent(in)  JE,
real(rp), dimension(ka,ia,ja), intent(in)  DENS0,
real(rp), dimension(ka,ia,ja), intent(in)  TEMP0,
real(rp), dimension(ka,ia,ja,qa_mp-1), intent(in)  QTRC0,
real(rp), dimension (ka,ia,ja,n_hyd), intent(out)  Re 
)

Calculate Effective Radius.

Definition at line 424 of file scale_atmos_phy_mp_kessler.F90.

References scale_atmos_hydrometeor::i_hc, scale_atmos_hydrometeor::i_hr, and scale_atmos_hydrometeor::n_hyd.

Referenced by mod_atmos_phy_mp_vars::atmos_phy_mp_vars_get_diagnostic().

424  use scale_atmos_hydrometeor, only: &
425  n_hyd, &
426  i_hc, &
427  i_hr
428  implicit none
429  integer, intent(in) :: ka, ks, ke
430  integer, intent(in) :: ia, is, ie
431  integer, intent(in) :: ja, js, je
432 
433  real(RP), intent(in) :: dens0(ka,ia,ja) ! density [kg/m3]
434  real(RP), intent(in) :: temp0(ka,ia,ja) ! temperature [K]
435  real(RP), intent(in) :: qtrc0(ka,ia,ja,qa_mp-1) ! hydrometeor mass concentration [kg/kg]
436 
437  real(RP), intent(out) :: re (ka,ia,ja,n_hyd) ! effective radius [cm]
438 
439  real(RP), parameter :: um2cm = 100.0_rp
440  !---------------------------------------------------------------------------
441 
442 !OCL XFILL
443  re(:,:,:,i_hc) = 8.e-6_rp * um2cm
444 !OCL XFILL
445  re(:,:,:,i_hr) = 100.e-6_rp * um2cm
446 !OCL XFILL
447  re(:,:,:,i_hr+1:) = 0.0_rp
448 
449  return
integer, public ia
of whole cells: x, local, with HALO
integer, parameter, public i_hr
liquid water rain
integer, public ja
of whole cells: y, local, with HALO
integer, public is
start point of inner domain: x, local
integer, public ie
end point of inner domain: x, local
module atmosphere / hydrometeor
integer, public ke
end point of inner domain: z, local
integer, public je
end point of inner domain: y, local
integer, public ks
start point of inner domain: z, local
integer, parameter, public i_hc
liquid water cloud
integer, public js
start point of inner domain: y, local
integer, public ka
of whole cells: z, local, with HALO
integer, parameter, public n_hyd
Here is the caller graph for this function:

◆ atmos_phy_mp_kessler_qtrc2qhyd()

subroutine, public scale_atmos_phy_mp_kessler::atmos_phy_mp_kessler_qtrc2qhyd ( integer, intent(in)  KA,
integer, intent(in)  KS,
integer, intent(in)  KE,
integer, intent(in)  IA,
integer, intent(in)  IS,
integer, intent(in)  IE,
integer, intent(in)  JA,
integer, intent(in)  JS,
integer, intent(in)  JE,
real(rp), dimension(ka,ia,ja,qa_mp-1), intent(in)  QTRC,
real(rp), dimension(ka,ia,ja,n_hyd), intent(out)  Qe 
)

Calculate mass ratio of each category.

Definition at line 458 of file scale_atmos_phy_mp_kessler.F90.

References scale_atmos_hydrometeor::i_hc, scale_atmos_hydrometeor::i_hr, and scale_atmos_hydrometeor::n_hyd.

Referenced by mod_atmos_phy_mp_vars::atmos_phy_mp_vars_get_diagnostic().

458  use scale_atmos_hydrometeor, only: &
459  n_hyd, &
460  i_hc, &
461  i_hr
462  implicit none
463  integer, intent(in) :: ka, ks, ke
464  integer, intent(in) :: ia, is, ie
465  integer, intent(in) :: ja, js, je
466 
467  real(RP), intent(in) :: qtrc(ka,ia,ja,qa_mp-1) ! hydrometeor mass concentration [kg/kg]
468 
469  real(RP), intent(out) :: qe(ka,ia,ja,n_hyd) ! mass ratio of each hydrometeor [kg/kg]
470  !---------------------------------------------------------------------------
471 
472 !OCL XFILL
473  qe(:,:,:,i_hc) = qtrc(:,:,:,i_hyd_qc)
474 !OCL XFILL
475  qe(:,:,:,i_hr) = qtrc(:,:,:,i_hyd_qr)
476 !OCL XFILL
477  qe(:,:,:,i_hr+1:) = 0.0_rp
478 
479  return
integer, public ia
of whole cells: x, local, with HALO
integer, parameter, public i_hr
liquid water rain
integer, public ja
of whole cells: y, local, with HALO
integer, public is
start point of inner domain: x, local
integer, public ie
end point of inner domain: x, local
module atmosphere / hydrometeor
integer, public ke
end point of inner domain: z, local
integer, public je
end point of inner domain: y, local
integer, public ks
start point of inner domain: z, local
integer, parameter, public i_hc
liquid water cloud
integer, public js
start point of inner domain: y, local
integer, public ka
of whole cells: z, local, with HALO
integer, parameter, public n_hyd
Here is the caller graph for this function:

◆ atmos_phy_mp_kessler_qhyd2qtrc()

subroutine, public scale_atmos_phy_mp_kessler::atmos_phy_mp_kessler_qhyd2qtrc ( integer, intent(in)  KA,
integer, intent(in)  KS,
integer, intent(in)  KE,
integer, intent(in)  IA,
integer, intent(in)  IS,
integer, intent(in)  IE,
integer, intent(in)  JA,
integer, intent(in)  JS,
integer, intent(in)  JE,
real(rp), dimension(ka,ia,ja,n_hyd), intent(in)  Qe,
real(rp), dimension(ka,ia,ja,qa_mp-1), intent(out)  QTRC 
)

Calculate mass ratio of each category.

Definition at line 488 of file scale_atmos_phy_mp_kessler.F90.

References scale_atmos_hydrometeor::i_hc, scale_atmos_hydrometeor::i_hr, and scale_atmos_hydrometeor::n_hyd.

Referenced by mod_atmos_phy_mp_driver::atmos_phy_mp_driver_qhyd2qtrc().

488  use scale_atmos_hydrometeor, only: &
489  n_hyd, &
490  i_hc, &
491  i_hr
492  implicit none
493  integer, intent(in) :: ka, ks, ke
494  integer, intent(in) :: ia, is, ie
495  integer, intent(in) :: ja, js, je
496 
497  real(RP), intent(in) :: qe(ka,ia,ja,n_hyd) ! mass ratio of each hydrometeor [kg/kg]
498 
499  real(RP), intent(out) :: qtrc(ka,ia,ja,qa_mp-1) ! hydrometeor mass concentration [kg/kg]
500 
501  !---------------------------------------------------------------------------
502 
503 !OCL XFILL
504  qtrc(:,:,:,i_hyd_qc) = qe(:,:,:,i_hc)
505 !OCL XFILL
506  qtrc(:,:,:,i_hyd_qr) = qe(:,:,:,i_hr)
507 
508  ! ignore ice water
509 
510  return
integer, public ia
of whole cells: x, local, with HALO
integer, parameter, public i_hr
liquid water rain
integer, public ja
of whole cells: y, local, with HALO
integer, public is
start point of inner domain: x, local
integer, public ie
end point of inner domain: x, local
module atmosphere / hydrometeor
integer, public ke
end point of inner domain: z, local
integer, public je
end point of inner domain: y, local
integer, public ks
start point of inner domain: z, local
integer, parameter, public i_hc
liquid water cloud
integer, public js
start point of inner domain: y, local
integer, public ka
of whole cells: z, local, with HALO
integer, parameter, public n_hyd
Here is the caller graph for this function:

Variable Documentation

◆ atmos_phy_mp_kessler_ntracers

integer, parameter, public scale_atmos_phy_mp_kessler::atmos_phy_mp_kessler_ntracers = QA_MP

Definition at line 44 of file scale_atmos_phy_mp_kessler.F90.

Referenced by mod_atmos_phy_mp_driver::atmos_phy_mp_driver_tracer_setup().

44  integer, parameter, public :: atmos_phy_mp_kessler_ntracers = qa_mp

◆ atmos_phy_mp_kessler_nwaters

integer, parameter, public scale_atmos_phy_mp_kessler::atmos_phy_mp_kessler_nwaters = 2

Definition at line 45 of file scale_atmos_phy_mp_kessler.F90.

Referenced by mod_atmos_phy_mp_driver::atmos_phy_mp_driver_tracer_setup().

45  integer, parameter, public :: atmos_phy_mp_kessler_nwaters = 2

◆ atmos_phy_mp_kessler_nices

integer, parameter, public scale_atmos_phy_mp_kessler::atmos_phy_mp_kessler_nices = 0

Definition at line 46 of file scale_atmos_phy_mp_kessler.F90.

Referenced by mod_atmos_phy_mp_driver::atmos_phy_mp_driver_tracer_setup().

46  integer, parameter, public :: atmos_phy_mp_kessler_nices = 0

◆ atmos_phy_mp_kessler_tracer_names

character(len=h_short), dimension(qa_mp), parameter, public scale_atmos_phy_mp_kessler::atmos_phy_mp_kessler_tracer_names = (/ 'QV', 'QC', 'QR' /)

Definition at line 47 of file scale_atmos_phy_mp_kessler.F90.

Referenced by mod_atmos_phy_mp_driver::atmos_phy_mp_driver_tracer_setup().

47  character(len=H_SHORT), parameter, public :: atmos_phy_mp_kessler_tracer_names(qa_mp) = (/ &
48  'QV', &
49  'QC', &
50  'QR' /)

◆ atmos_phy_mp_kessler_tracer_descriptions

character(len=h_mid), dimension(qa_mp), parameter, public scale_atmos_phy_mp_kessler::atmos_phy_mp_kessler_tracer_descriptions = (/ 'Ratio of Water Vapor mass to total mass (Specific humidity)', 'Ratio of Cloud Water mass to total mass ', 'Ratio of Rain Water mass to total mass '/)

Definition at line 51 of file scale_atmos_phy_mp_kessler.F90.

Referenced by mod_atmos_phy_mp_driver::atmos_phy_mp_driver_tracer_setup().

51  character(len=H_MID) , parameter, public :: atmos_phy_mp_kessler_tracer_descriptions(qa_mp) = (/ &
52  'Ratio of Water Vapor mass to total mass (Specific humidity)', &
53  'Ratio of Cloud Water mass to total mass ', &
54  'Ratio of Rain Water mass to total mass '/)

◆ atmos_phy_mp_kessler_tracer_units

character(len=h_short), dimension(qa_mp), parameter, public scale_atmos_phy_mp_kessler::atmos_phy_mp_kessler_tracer_units = (/ 'kg/kg', 'kg/kg', 'kg/kg' /)

Definition at line 55 of file scale_atmos_phy_mp_kessler.F90.

Referenced by mod_atmos_phy_mp_driver::atmos_phy_mp_driver_tracer_setup().

55  character(len=H_SHORT), parameter, public :: atmos_phy_mp_kessler_tracer_units(qa_mp) = (/ &
56  'kg/kg', &
57  'kg/kg', &
58  'kg/kg' /)