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

module atmosphere / physics / cloud microphysics More...

Functions/Subroutines

subroutine, public atmos_phy_mp_driver_tracer_setup
 Config. More...
 
subroutine, public atmos_phy_mp_driver_setup
 Setup. More...
 
subroutine, public atmos_phy_mp_driver_adjustment
 adjustment More...
 
subroutine, public atmos_phy_mp_driver_calc_tendency (update_flag)
 calculate tendency More...
 
subroutine, public atmos_phy_mp_driver_qhyd2qtrc (KA, KS, KE, IA, IS, IE, JA, JS, JE, QV, QHYD, QTRC, QNUM)
 
subroutine, public atmos_phy_mp_driver_qhyd2qtrc_onlyqv (KA, KS, KE, IA, IS, IE, JA, JS, JE, QV, QHYD, QTRC, QNUM)
 

Variables

procedure(qhyd2qtrc), pointer, public atmos_phy_mp_user_qhyd2qtrc => NULL()
 

Detailed Description

module atmosphere / physics / cloud microphysics

Description
Cloud Microphysics driver
Author
Team SCALE
NAMELIST
  • PARAM_ATMOS_PHY_MP
    nametypedefault valuecomment
    MP_DO_PRECIPITATION logical .true. > apply sedimentation (precipitation)?
    MP_DO_NEGATIVE_FIXER logical .true. > apply negative fixer?
    MP_LIMIT_NEGATIVE real(RP) 0.1_RP > Abort if abs(fixed negative vaue) > abs(MP_limit_negative)
    MP_NTMAX_SEDIMENTATION integer 1 > number of time step for sedimentation
    MP_MAX_TERM_VEL real(RP) 10.0_RP > terminal velocity for calculate dt of sedimentation
    MP_CLDFRAC_THLESHOLD real(RP) > thleshold for cloud fraction

History Output
namedescriptionunitvariable
"DENS_t_MP_NF" "vapor supply by the negative fixer" "kg/m3/s" "DENS_t_MP_NF"
"ENGI_t_MP_NF" "internal energy supply by the negative fixer" "J/m3/s" "ENGI_t_MP_NF"
"RHOH_MP_NF" "sensible heat by the negative fixer" "J/m3/s" "RHOH_MP_NF"
Vterm_{TRACER_NAME} terminal velocity of {TRACER_NAME};
{TRACER_NAME} depends on the physics schemes, e.g., QV, QC, QR.
m/s {'Vterm_'//trim}
DENS_t_MP tendency DENS in MP kg/m3/s DENS_t_MP
EVAPORATE evaporated cloud number num/m3/s EVAPORATE
MOMZ_t_MP tendency MOMZ in MP kg/m2/s2 MOMZ_t_MP
PREC_MP surface precipitation rate by MP kg/m2/s precip
QSPLT_G Charge split of QG by Non-inductive process fC/m3/s QSPLT_in
QSPLT_I Charge split of QI by Non-inductive process fC/m3/s QSPLT_in
QSPLT_S Charge split of QS by Non-inductive process fC/m3/s QSPLT_in
RAIN_MP surface rain rate by MP kg/m2/s SFLX_rain
RHOH_MP diabatic heating rate in MP J/m3/s RHOH_MP
RHOU_t_MP tendency RHOU in MP kg/m2/s2 RHOU_t_MP
RHOV_t_MP tendency RHOV in MP kg/m2/s2 RHOV_t_MP
SNOW_MP surface snow rate by MP kg/m2/s SFLX_snow
{TRACER_NAME}_t_MP tendency rho*{TRACER_NAME} in MP;
{TRACER_NAME} depends on the physics schemes, e.g., QV, QC, QR.
fC/m3/s RHOC_t_MP
QCRG_H Ratio of charge density of Hail fC/kg QCRG_H
QH Ratio of Hail Water mass to total mass kg/kg QH

Function/Subroutine Documentation

◆ atmos_phy_mp_driver_tracer_setup()

subroutine, public mod_atmos_phy_mp_driver::atmos_phy_mp_driver_tracer_setup

Config.

Definition at line 93 of file mod_atmos_phy_mp_driver.F90.

93  use scale_prc, only: &
94  prc_abort
95  use scale_tracer, only: &
97  use scale_atmos_hydrometeor, only: &
99  use scale_atmos_phy_mp_kessler, only: &
106  use scale_atmos_phy_mp_tomita08, only: &
113  use scale_atmos_phy_mp_sn14, only: &
120  use scale_atmos_phy_mp_suzuki10, only: &
129  use mod_atmos_admin, only: &
131  atmos_sw_phy_mp, &
133  use mod_atmos_phy_mp_vars, only: &
134  qa_mp, &
135  qs_mp, &
136  qe_mp
137  implicit none
138 
139  integer :: QS2
140  !---------------------------------------------------------------------------
141 
142  log_newline
143  log_info("ATMOS_PHY_MP_driver_tracer_setup",*) 'Setup'
144 
145  if ( atmos_sw_phy_mp ) then
146  select case ( atmos_phy_mp_type )
147  case ( 'KESSLER' )
154  qs_mp ) ! [OUT]
156  if ( atmos_sw_phy_lt ) then
157  log_error("ATMOS_PHY_MP_driver_tracer_setup",*) 'ATMOS_PHY_MP_TYPE is invalud for Lightning: ', trim(atmos_phy_mp_type)
158  log_error("ATMOS_PHY_MP_driver_tracer_setup",*) 'ATMOS_PHY_MP_TYPE is should be TOMITA08 or SN14 or SUZUKI10 for Lightning: '
159  call prc_abort
160  endif
161  case ( 'TOMITA08' )
168  qs_mp ) ! [OUT]
170  case( 'SN14' )
172  atmos_phy_mp_sn14_nwaters, & ! [IN]
173  atmos_phy_mp_sn14_nices, & ! [IN]
174  atmos_phy_mp_sn14_tracer_names(1:6), & ! [IN]
176  atmos_phy_mp_sn14_tracer_units(1:6), & ! [IN]
177  qs_mp ) ! [OUT]
178  call tracer_regist( qs2, & ! [OUT]
179  5, & ! [IN]
180  atmos_phy_mp_sn14_tracer_names(7:11), & ! [IN]
182  atmos_phy_mp_sn14_tracer_units(7:11) ) ! [IN]
184  case ( 'SUZUKI10' )
192  qs_mp ) ! [OUT]
193  if( atmos_phy_mp_suzuki10_nccn > 0 ) then
194  call tracer_regist( qs2, & ! [OUT]
198  + 2 : ), & ! [IN]
201  + 2 : ), & ! [IN]
204  + 2 : ) ) ! [IN]
205  end if
207  case default
208  log_error("ATMOS_PHY_MP_driver_tracer_setup",*) 'ATMOS_PHY_MP_TYPE is invalud: ', trim(atmos_phy_mp_type)
209  call prc_abort
210  end select
211 
212  qe_mp = qs_mp + qa_mp - 1
213 
214  else
215  qa_mp = 0
216  qs_mp = -1
217  qe_mp = -2
218  end if
219 
220  return

References scale_atmos_hydrometeor::atmos_hydrometeor_regist(), scale_atmos_phy_mp_kessler::atmos_phy_mp_kessler_nices, scale_atmos_phy_mp_kessler::atmos_phy_mp_kessler_ntracers, scale_atmos_phy_mp_kessler::atmos_phy_mp_kessler_nwaters, scale_atmos_phy_mp_kessler::atmos_phy_mp_kessler_tracer_descriptions, scale_atmos_phy_mp_kessler::atmos_phy_mp_kessler_tracer_names, scale_atmos_phy_mp_kessler::atmos_phy_mp_kessler_tracer_units, scale_atmos_phy_mp_sn14::atmos_phy_mp_sn14_nices, scale_atmos_phy_mp_sn14::atmos_phy_mp_sn14_ntracers, scale_atmos_phy_mp_sn14::atmos_phy_mp_sn14_nwaters, scale_atmos_phy_mp_sn14::atmos_phy_mp_sn14_tracer_descriptions, scale_atmos_phy_mp_sn14::atmos_phy_mp_sn14_tracer_names, scale_atmos_phy_mp_sn14::atmos_phy_mp_sn14_tracer_units, scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_nccn, scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_nices, scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_ntracers, scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_nwaters, scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_tracer_descriptions, scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_tracer_names, scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_tracer_setup(), scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_tracer_units, scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08_nices, scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08_ntracers, scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08_nwaters, scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08_tracer_descriptions, scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08_tracer_names, scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08_tracer_units, mod_atmos_admin::atmos_phy_mp_type, mod_atmos_admin::atmos_sw_phy_lt, mod_atmos_admin::atmos_sw_phy_mp, scale_prc::prc_abort(), mod_atmos_phy_mp_vars::qa_mp, mod_atmos_phy_mp_vars::qe_mp, mod_atmos_phy_mp_vars::qs_mp, and scale_tracer::tracer_regist().

Referenced by mod_atmos_driver::atmos_driver_tracer_setup().

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

◆ atmos_phy_mp_driver_setup()

subroutine, public mod_atmos_phy_mp_driver::atmos_phy_mp_driver_setup

Setup.

Definition at line 226 of file mod_atmos_phy_mp_driver.F90.

226  use scale_prc, only: &
227  prc_abort
228  use scale_atmos_grid_cartesc, only: &
230  use scale_const, only: &
231  eps => const_eps
232  use scale_time, only: &
234  use scale_atmos_hydrometeor, only: &
235  n_hyd
236  use scale_atmos_phy_mp_kessler, only: &
238  use scale_atmos_phy_mp_tomita08, only: &
240  use scale_atmos_phy_mp_sn14, only: &
242  use scale_atmos_phy_mp_suzuki10, only: &
244  use scale_file_history, only: &
246  use scale_monitor, only: &
247  monitor_reg, &
248  monitor_put
249  use mod_atmos_admin, only: &
252  use mod_atmos_phy_mp_vars, only: &
254  sflx_rain => atmos_phy_mp_sflx_rain, &
255  sflx_snow => atmos_phy_mp_sflx_snow, &
256  sflx_engi => atmos_phy_mp_sflx_engi
257  use mod_atmos_phy_mp_vars, only: &
258  qs_mp, &
259  qe_mp
260  use mod_atmos_phy_lt_vars, only: &
261  flg_lt
262  implicit none
263 
264  namelist / param_atmos_phy_mp / &
265  mp_do_precipitation, &
266  mp_do_negative_fixer, &
267  mp_limit_negative, &
268  mp_ntmax_sedimentation, &
269  mp_max_term_vel, &
270  mp_cldfrac_thleshold
271 
272  character(len=H_SHORT) :: w_name(N_HYD)
273  character(len=H_MID) :: w_longname(N_HYD)
274  character(len=H_SHORT) :: w_unit(N_HYD)
275  data w_name / 'QC', &
276  'QR', &
277  'QI', &
278  'QS', &
279  'QG', &
280  'QH' /
281  data w_longname / &
282  'Ratio of Cloud Water mass to total mass', &
283  'Ratio of Rain Water mass to total mass', &
284  'Ratio of Ice Water mass to total mass', &
285  'Ratio of Snow Water mass to total mass', &
286  'Ratio of Graupel Water mass to total mass', &
287  'Ratio of Hail Water mass to total mass' /
288  data w_unit / &
289  'kg/kg', &
290  'kg/kg', &
291  'kg/kg', &
292  'kg/kg', &
293  'kg/kg', &
294  'kg/kg' /
295 
296  character(len=H_SHORT) :: w_crg_name(N_HYD)
297  character(len=H_MID) :: w_crg_longname(N_HYD)
298  character(len=H_SHORT) :: w_crg_unit(N_HYD)
299  data w_crg_name / 'QCRG_C', &
300  'QCRG_R', &
301  'QCRG_I', &
302  'QCRG_S', &
303  'QCRG_G', &
304  'QCRG_H' /
305  data w_crg_longname / &
306  'Ratio of charge density of Cloud Water', &
307  'Ratio of charge density of Rain Water', &
308  'Ratio of charge density of Ice', &
309  'Ratio of charge density of Snow', &
310  'Ratio of charge density of Graupel', &
311  'Ratio of charge density of Hail' /
312  data w_crg_unit / &
313  'fC/kg', &
314  'fC/kg', &
315  'fC/kg', &
316  'fC/kg', &
317  'fC/kg', &
318  'fC/kg' /
319 
320  real(RP) :: ZERO(KA,IA,JA)
321 
322  integer :: nstep_max
323 
324  integer :: iq
325  integer :: ierr
326  !---------------------------------------------------------------------------
327 
328  log_newline
329  log_info("ATMOS_PHY_MP_driver_setup",*) 'Setup'
330 
331  if ( atmos_sw_phy_mp ) then
332 
333  mp_cldfrac_thleshold = eps
334 
335  !--- read namelist
336  rewind(io_fid_conf)
337  read(io_fid_conf,nml=param_atmos_phy_mp,iostat=ierr)
338  if( ierr < 0 ) then !--- missing
339  log_info("ATMOS_PHY_MP_driver_setup",*) 'Not found namelist. Default used.'
340  elseif( ierr > 0 ) then !--- fatal error
341  log_error("ATMOS_PHY_MP_driver_setup",*) 'Not appropriate names in namelist PARAM_ATMOS_PHY_MP. Check!'
342  call prc_abort
343  endif
344  log_nml(param_atmos_phy_mp)
345 
346  nstep_max = ceiling( ( time_dtsec_atmos_phy_mp * mp_max_term_vel ) / minval( cdz(:) ) )
347  mp_ntmax_sedimentation = max( mp_ntmax_sedimentation, nstep_max )
348 
349  mp_nstep_sedimentation = mp_ntmax_sedimentation
350  mp_rnstep_sedimentation = 1.0_rp / real(mp_ntmax_sedimentation,kind=rp)
351  mp_dtsec_sedimentation = time_dtsec_atmos_phy_mp * mp_rnstep_sedimentation
352 
353  atmos_phy_mp_cldfrac_thleshold = mp_cldfrac_thleshold
354 
355  log_newline
356  log_info("ATMOS_PHY_MP_driver_setup",*) 'Enable negative fixer? : ', mp_do_negative_fixer
357  log_info("ATMOS_PHY_MP_driver_setup",*) 'Value limit of negative fixer (abs) : ', abs(mp_limit_negative)
358  log_info("ATMOS_PHY_MP_driver_setup",*) 'Enable sedimentation (precipitation)? : ', mp_do_precipitation
359  log_info("ATMOS_PHY_MP_driver_setup",*) 'Timestep of sedimentation is divided into : ', mp_ntmax_sedimentation, 'step'
360  log_info("ATMOS_PHY_MP_driver_setup",*) 'DT of sedimentation : ', mp_dtsec_sedimentation, '[s]'
361 
362  select case ( atmos_phy_mp_type )
363  case ( 'KESSLER' )
365  case ( 'TOMITA08' )
367  ka, ks, ke, ia, is, ie, ja, js, je, &
368  flg_lt )
369  case ( 'SN14' )
371  ka, ia, ja )
372  case ( 'SUZUKI10' )
374  ka, ia, ja, &
375  flg_lt )
376  allocate( hist_hyd_id(n_hyd) )
377  do iq = 1, n_hyd
378  call file_history_reg( w_name(iq), w_longname(iq), w_unit(iq), & ! [IN]
379  hist_hyd_id(iq) ) ! [OUT]
380  enddo
381  if( flg_lt ) then
382  allocate( hist_crg_id(n_hyd) )
383  do iq = 1, n_hyd
384  call file_history_reg( w_crg_name(iq), w_crg_longname(iq), w_crg_unit(iq), & ! [IN]
385  hist_crg_id(iq) ) ! [OUT]
386  enddo
387  endif
388  end select
389 
390  ! history putput
391  if ( mp_do_precipitation ) then
392 
393  allocate( hist_vterm_id(qs_mp+1:qe_mp) )
394  do iq = qs_mp+1, qe_mp
395  call file_history_reg( 'Vterm_'//trim(tracer_name(iq)), 'terminal velocity of '//trim(tracer_name(iq)), 'm/s', hist_vterm_id(iq) )
396  end do
397 
398  else
399 
400  sflx_rain(:,:) = 0.0_rp
401  sflx_snow(:,:) = 0.0_rp
402  sflx_engi(:,:) = 0.0_rp
403 
404  end if
405 
406  ! monitor
407  if ( mp_do_negative_fixer ) then
408  call file_history_reg( "RHOH_MP_NF", "sensible heat by the negative fixer", "J/m3/s", & ! [IN]
409  hist_nf_rhoh_id ) ! [OUT]
410  call file_history_reg( "DENS_t_MP_NF", "vapor supply by the negative fixer", "kg/m3/s", & ! [IN]
411  hist_nf_dens_id ) ! [OUT]
412  call file_history_reg( "ENGI_t_MP_NF", "internal energy supply by the negative fixer", "J/m3/s", & ! [IN]
413  hist_nf_engi_id ) ! [OUT]
414  call monitor_reg( "QTOTTND_NF", "vapor supply by the negative fixer", "kg", & ! [IN]
415  monit_nf_mass_id, & ! [OUT]
416  is_tendency=.true. ) ! [IN]
417  call monitor_reg( "ENGITND_NF", "internal energy supply by the negative fixer", "J", & ! [IN]
418  monit_nf_engi_id, & ! [OUT]
419  is_tendency=.true. ) ! [IN]
420  zero(:,:,:) = 0.0_rp
421  call monitor_put( monit_nf_mass_id, zero(:,:,:) )
422  call monitor_put( monit_nf_engi_id, zero(:,:,:) )
423  end if
424 
425  else
426 
427  log_info("ATMOS_PHY_MP_driver_setup",*) 'this component is never called.'
428  log_info("ATMOS_PHY_MP_driver_setup",*) 'SFLX_rain and SFLX_snow is set to zero.'
429  sflx_rain(:,:) = 0.0_rp
430  sflx_snow(:,:) = 0.0_rp
431  sflx_engi(:,:) = 0.0_rp
432 
433  endif
434 
435  return

References scale_atmos_grid_cartesc::atmos_grid_cartesc_cdz, mod_atmos_phy_mp_vars::atmos_phy_mp_cldfrac_thleshold, scale_atmos_phy_mp_kessler::atmos_phy_mp_kessler_setup(), mod_atmos_phy_mp_vars::atmos_phy_mp_sflx_engi, mod_atmos_phy_mp_vars::atmos_phy_mp_sflx_rain, mod_atmos_phy_mp_vars::atmos_phy_mp_sflx_snow, scale_atmos_phy_mp_sn14::atmos_phy_mp_sn14_setup(), scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_setup(), scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08_setup(), mod_atmos_admin::atmos_phy_mp_type, mod_atmos_admin::atmos_sw_phy_mp, scale_const::const_eps, scale_file_history::file_history_reg(), mod_atmos_phy_lt_vars::flg_lt, scale_atmos_grid_cartesc_index::ia, scale_atmos_grid_cartesc_index::ie, scale_io::io_fid_conf, scale_atmos_grid_cartesc_index::is, scale_atmos_grid_cartesc_index::ja, scale_atmos_grid_cartesc_index::je, scale_atmos_grid_cartesc_index::js, scale_atmos_grid_cartesc_index::ka, scale_atmos_grid_cartesc_index::ke, scale_atmos_grid_cartesc_index::ks, scale_monitor::monitor_reg(), scale_atmos_hydrometeor::n_hyd, scale_prc::prc_abort(), mod_atmos_phy_mp_vars::qe_mp, mod_atmos_phy_mp_vars::qs_mp, scale_precision::rp, scale_time::time_dtsec_atmos_phy_mp, and scale_tracer::tracer_name.

Referenced by mod_atmos_driver::atmos_driver_setup().

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

◆ atmos_phy_mp_driver_adjustment()

subroutine, public mod_atmos_phy_mp_driver::atmos_phy_mp_driver_adjustment

adjustment

Definition at line 441 of file mod_atmos_phy_mp_driver.F90.

441  use scale_const, only: &
442  pre00 => const_pre00
443  use scale_atmos_hydrometeor, only: &
445  i_qv, &
446  qla, &
447  qia, &
448  qhs, &
449  qhe
450  use scale_atmos_phy_mp_common, only: &
452  use mod_atmos_vars, only: &
453  temp, &
454  cvtot, &
455  cptot, &
456  dens, &
457  rhot, &
458  qtrc
459  use mod_atmos_phy_mp_vars, only: &
460  qa_mp
461  use scale_time, only: &
462  dt => time_dtsec
463  use scale_file_history, only: &
464  file_history_query, &
465  file_history_put
466  use scale_monitor, only: &
467  monitor_put
468 
469  real(RP) :: RHOH (KA,IA,JA)
470  real(RP) :: DENS_d(KA,IA,JA)
471  real(RP) :: ENGI_d(KA,IA,JA)
472  real(RP) :: Rtot
473 
474  logical :: do_put_rhoh
475  logical :: do_put_dens
476  logical :: do_put_engi
477 
478  integer :: k, i, j, iq
479 
480  if ( mp_do_negative_fixer .and. (.not. atmos_hydrometeor_dry) ) then
481 
482  call file_history_query( hist_nf_rhoh_id, do_put_rhoh )
483  call file_history_query( hist_nf_dens_id, do_put_dens )
484  call file_history_query( hist_nf_engi_id, do_put_engi )
485 
486  if ( monit_nf_mass_id > 0 .or. monit_nf_engi_id > 0 .or. &
487  do_put_rhoh .or. do_put_dens .or. do_put_engi ) then
489  ka, ks, ke, ia, 1, ia, ja, 1, ja, qla, qia, &
490  mp_limit_negative, & ! [IN]
491  dens(:,:,:), temp(:,:,:), & ! [INOUT]
492  cvtot(:,:,:), cptot(:,:,:), & ! [INOUT]
493  qtrc(:,:,:,i_qv), qtrc(:,:,:,qhs:qhe), & ! [INOUT]
494  rhoh = rhoh, & ! [OUT, optional]
495  dens_diff = dens_d, engi_diff = engi_d ) ! [OUT, optional]
496  else
498  ka, ks, ke, ia, 1, ia, ja, 1, ja, qla, qia, &
499  mp_limit_negative, & ! [IN]
500  dens(:,:,:), temp(:,:,:), & ! [INOUT]
501  cvtot(:,:,:), cptot(:,:,:), & ! [INOUT]
502  qtrc(:,:,:,i_qv), qtrc(:,:,:,qhs:qhe) ) ! [INOUT]
503  end if
504 
505  !$omp parallel private(Rtot)
506 
507  !$omp do
508  do j = js, je
509  do i = is, ie
510  do k = ks, ke
511  rtot = cptot(k,i,j) - cvtot(k,i,j)
512  rhot(k,i,j) = pre00 / rtot * ( dens(k,i,j) * temp(k,i,j) * rtot / pre00 )**( cvtot(k,i,j) / cptot(k,i,j) )
513  end do
514  end do
515  end do
516  !$omp end do nowait
517 
518  ! for non-mass tracers, such as number density
519  do iq = qhe+1, qa_mp
520  !$omp do
521  do j = js, je
522  do i = is, ie
523  do k = ks, ke
524  qtrc(k,i,j,iq) = max( qtrc(k,i,j,iq), 0.0_rp )
525  end do
526  end do
527  end do
528  end do
529 
530  !$omp end parallel
531 
532  if ( do_put_rhoh ) then
533  !$omp parallel do
534  do j = js, je
535  do i = is, ie
536  do k = ks, ke
537  rhoh(k,i,j) = rhoh(k,i,j) / dt
538  end do
539  end do
540  end do
541  call file_history_put( hist_nf_rhoh_id, rhoh(:,:,:) )
542  end if
543  if ( monit_nf_mass_id > 0 .or. do_put_dens ) then
544  !$omp parallel do
545  do j = js, je
546  do i = is, ie
547  do k = ks, ke
548  dens_d(k,i,j) = dens_d(k,i,j) / dt
549  end do
550  end do
551  end do
552  call file_history_put( hist_nf_dens_id, dens_d(:,:,:) )
553  call monitor_put( monit_nf_mass_id, dens_d(:,:,:) )
554  end if
555  if ( monit_nf_engi_id > 0 .or. do_put_engi ) then
556  !$omp parallel do
557  do j = js, je
558  do i = is, ie
559  do k = ks, ke
560  engi_d(k,i,j) = engi_d(k,i,j) / dt
561  end do
562  end do
563  end do
564  call file_history_put( hist_nf_engi_id, engi_d(:,:,:) )
565  call monitor_put( monit_nf_engi_id, engi_d(:,:,:) )
566  end if
567 
568  end if
569 
570  return

References scale_atmos_hydrometeor::atmos_hydrometeor_dry, scale_atmos_phy_mp_common::atmos_phy_mp_negative_fixer(), scale_const::const_pre00, mod_atmos_vars::cptot, mod_atmos_vars::cvtot, mod_atmos_vars::dens, scale_atmos_hydrometeor::i_qv, scale_atmos_grid_cartesc_index::ia, scale_atmos_grid_cartesc_index::ie, scale_atmos_grid_cartesc_index::is, scale_atmos_grid_cartesc_index::ja, scale_atmos_grid_cartesc_index::je, scale_atmos_grid_cartesc_index::js, scale_tracer::k, scale_atmos_grid_cartesc_index::ka, scale_atmos_grid_cartesc_index::ke, scale_atmos_grid_cartesc_index::ks, mod_atmos_phy_mp_vars::qa_mp, scale_atmos_hydrometeor::qhe, scale_atmos_hydrometeor::qhs, scale_atmos_hydrometeor::qia, scale_atmos_hydrometeor::qla, mod_atmos_vars::qtrc, mod_atmos_vars::rhot, mod_atmos_vars::temp, and scale_time::time_dtsec.

Referenced by mod_atmos_driver::atmos_driver_update().

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

◆ atmos_phy_mp_driver_calc_tendency()

subroutine, public mod_atmos_phy_mp_driver::atmos_phy_mp_driver_calc_tendency ( logical, intent(in)  update_flag)

calculate tendency

Definition at line 576 of file mod_atmos_phy_mp_driver.F90.

576  use scale_const, only: &
577  pre00 => const_pre00
578  use scale_time, only: &
579  dt_mp => time_dtsec_atmos_phy_mp, &
580  dt_lt => time_dtsec_atmos_phy_lt
581  use scale_atmos_grid_cartesc_real, only: &
582  real_cz => atmos_grid_cartesc_real_cz, &
583  real_fz => atmos_grid_cartesc_real_fz, &
592  use scale_statistics, only: &
594  statistics_total
595  use scale_file_history, only: &
596  file_history_in, &
597  file_history_query, &
598  file_history_put
599  use scale_atmos_hydrometeor, only: &
600  lhf, &
601  i_qv, &
602  n_hyd, &
603  qha, &
604  qhs, &
605  qhe, &
606  qla, &
607  qls, &
608  qle, &
609  qia
610  use scale_atmos_phy_mp_common, only: &
614  use scale_atmos_refstate, only: &
615  refstate_dens => atmos_refstate_dens
616  use scale_atmos_phy_mp_kessler, only: &
619  use scale_atmos_phy_mp_tomita08, only: &
622  use scale_atmos_phy_mp_sn14, only: &
625  use scale_atmos_phy_mp_suzuki10, only: &
630  use scale_atmos_phy_lt_sato2019, only: &
632  use scale_file_history, only: &
633  file_history_query, &
634  file_history_put
635  use mod_atmos_admin, only: &
639  use mod_atmos_vars, only: &
640  dens => dens_av, &
641  momz => momz_av, &
642  u, &
643  v, &
644  w, &
645  pott, &
646  qtrc => qtrc_av, &
647  dens_t => dens_tp, &
648  momz_t => momz_tp, &
649  rhou_t => rhou_tp, &
650  rhov_t => rhov_tp, &
651 ! RHOT_t => RHOT_tp, &
652  rhoq_t => rhoq_tp, &
653  rhoh => rhoh_p, &
654  temp, &
655  pres, &
656  qdry, &
657  cvtot, &
658  cptot, &
659  exner, &
660  rhot => rhot_av
661  use mod_atmos_phy_mp_vars, only: &
662  qa_mp, &
663  qs_mp, &
664  qe_mp, &
665  dens_t_mp => atmos_phy_mp_dens_t, &
666  momz_t_mp => atmos_phy_mp_momz_t, &
667 ! RHOT_t_MP => ATMOS_PHY_MP_RHOT_t, &
668  rhou_t_mp => atmos_phy_mp_rhou_t, &
669  rhov_t_mp => atmos_phy_mp_rhov_t, &
670  rhoq_t_mp => atmos_phy_mp_rhoq_t, &
671  rhoc_t_mp => atmos_phy_mp_rhoc_t, &
672  rhoh_mp => atmos_phy_mp_rhoh, &
673  evaporate => atmos_phy_mp_evaporate, &
674  sflx_rain => atmos_phy_mp_sflx_rain, &
675  sflx_snow => atmos_phy_mp_sflx_snow, &
676  sflx_engi => atmos_phy_mp_sflx_engi
677  use mod_atmos_phy_ae_vars, only: &
678  ccn_t => atmos_phy_ae_ccn_t
679  use mod_atmos_phy_lt_vars, only: &
680  qa_lt, &
681  qs_lt, &
682  qe_lt, &
683  d0_crg, &
684  v0_crg, &
685  flg_lt, &
686  sarea => atmos_phy_lt_sarea
687  implicit none
688 
689  logical, intent(in) :: update_flag
690 
691  real(RP) :: RHOE_t(KA,IA,JA)
692  real(RP) :: TEMP1 (KA,IA,JA)
693  real(RP) :: CPtot1(KA,IA,JA)
694  real(RP) :: CVtot1(KA,IA,JA)
695  real(RP) :: CCN (KA,IA,JA)
696  real(RP) :: vterm (KA,QS_MP+1:QE_MP)
697  real(RP), target :: QTRC1(KA,IA,JA,QS_MP:QE_MP)
698 
699  real(RP) :: FLX_hydro(KA)
700  real(RP) :: DENS2 (KA)
701  real(RP) :: TEMP2 (KA)
702  real(RP) :: PRES2 (KA)
703  real(RP) :: CPtot2 (KA)
704  real(RP) :: CVtot2 (KA)
705  real(RP) :: RHOE (KA)
706  real(RP) :: RHOE2 (KA)
707  real(RP) :: RHOQ (KA,QS_MP+1:QE_MP)
708  real(RP) :: RHOQ2 (KA,QS_MP+1:QE_MP)
709  real(RP) :: mflux (KA)
710  real(RP) :: sflux (2)
711  real(RP) :: eflux
712 
713  real(RP) :: FZ (KA)
714  real(RP) :: FDZ (KA)
715  real(RP) :: RFDZ(KA)
716  real(RP) :: RCDZ(KA)
717 
718  real(RP) :: CPtot_t(KA,IA,JA), CVtot_t(KA,IA,JA)
719  real(RP) :: CP_t, CV_t
720 
721  real(RP) :: precip (IA,JA)
722 
723  real(RP) :: Qe(KA,IA,JA,N_HYD)
724  logical :: HIST_sw(N_HYD)
725  real(RP) :: Qecrg(KA,IA,JA,N_HYD)
726  logical :: HIST_crg_sw(N_HYD)
727 
728  ! for history output
729  real(RP), allocatable :: vterm_hist(:,:,:,:)
730  integer :: hist_vterm_idx(QS_MP+1:QE_MP)
731  logical :: flag
732  integer :: ih
733 
734  integer :: k, i, j, iq
735  integer :: step
736 
737  real(RP) :: QTRC1_crg(KA,IA,JA,QS_LT:QE_LT)
738  real(RP) :: RHOQ2_crg(KA,QS_LT:QE_LT)
739  real(RP) :: mflux_crg(KA), sflux_crg(2), eflux_crg
740  real(RP) :: QSPLT_in(KA,IA,JA,3)
741  real(RP) :: dqcrg(KA,IA,JA), beta_crg(KA,IA,JA)
742  real(RP) :: dummy(KA)
743  !---------------------------------------------------------------------------
744 
745  if ( update_flag ) then
746 
747  ccn(:,:,:) = ccn_t(:,:,:) * dt_mp
748 
749  select case ( atmos_phy_mp_type )
750  case ( 'KESSLER' )
751 !OCL XFILL
752  temp1(:,:,:) = temp(:,:,:)
753 !OCL XFILL
754  qtrc1(:,:,:,qs_mp:qe_mp) = qtrc(:,:,:,qs_mp:qe_mp)
755 !OCL XFILL
756  cvtot1(:,:,:) = cvtot(:,:,:)
757 !OCL XFILL
758  cptot1(:,:,:) = cptot(:,:,:)
759 
761  ka, ks, ke, ia, is, ie, ja, js, je, &
762  dens(:,:,:), pres(:,:,:), dt_mp, & ! [IN]
763  temp1(:,:,:), qtrc1(:,:,:,qs_mp:qe_mp), cptot1(:,:,:), cvtot1(:,:,:), & ! [INOUT]
764  rhoe_t(:,:,:), evaporate(:,:,:) ) ! [OUT]
765 
766  do iq = qs_mp, qe_mp
767  do j = js, je
768  do i = is, ie
769  do k = ks, ke
770  rhoq_t_mp(k,i,j,iq) = ( qtrc1(k,i,j,iq) - qtrc(k,i,j,iq) ) * dens(k,i,j) / dt_mp
771  enddo
772  enddo
773  enddo
774  enddo
775 
776  do j = js, je
777  do i = is, ie
778  do k = ks, ke
779  cptot_t(k,i,j) = ( cptot1(k,i,j) - cptot(k,i,j) ) / dt_mp
780  cvtot_t(k,i,j) = ( cvtot1(k,i,j) - cvtot(k,i,j) ) / dt_mp
781  end do
782  end do
783  end do
784 
785  case ( 'TOMITA08' )
786 !OCL XFILL
787  do j = js, je
788  do i = is, ie
789  do k = ks, ke
790  temp1(k,i,j) = temp(k,i,j)
791  end do
792  end do
793  end do
794 !OCL XFILL
795  do iq = qs_mp, qe_mp
796  do j = js, je
797  do i = is, ie
798  do k = ks, ke
799  qtrc1(k,i,j,iq) = qtrc(k,i,j,iq)
800  end do
801  end do
802  end do
803  end do
804 !OCL XFILL
805  do j = js, je
806  do i = is, ie
807  do k = ks, ke
808  cvtot1(k,i,j) = cvtot(k,i,j)
809  end do
810  end do
811  end do
812 !OCL XFILL
813  do j = js, je
814  do i = is, ie
815  do k = ks, ke
816  cptot1(k,i,j) = cptot(k,i,j)
817  end do
818  end do
819  end do
820 
821  if( flg_lt ) then
822 !OCL XFILL
823  do iq = qs_lt, qe_lt
824  do j = js, je
825  do i = is, ie
826  do k = ks, ke
827  qtrc1_crg(k,i,j,iq) = qtrc(k,i,j,iq)
828  end do
829  end do
830  end do
831  end do
832 
834  ka, ks, ke, ia, is, ie, ja, js, je, & ! [IN]
835  qla, & ! [IN]
836  temp1(:,:,:), dens(:,:,:), & ! [IN]
837  qtrc1(:,:,:,qls:qle), & ! [IN]
838  dqcrg(:,:,:), beta_crg(:,:,:) ) ! [OUT]
840  ka, ks, ke, ia, is, ie, ja, js, je, &
841  dens(:,:,:), pres(:,:,:), ccn(:,:,:), dt_mp, & ! [IN]
842  temp1(:,:,:), qtrc1(:,:,:,qs_mp:qe_mp), cptot1(:,:,:), cvtot1(:,:,:), & ! [INOUT]
843  rhoe_t(:,:,:), evaporate(:,:,:), & ! [OUT]
844  flg_lt, d0_crg, v0_crg, dqcrg(:,:,:), beta_crg(:,:,:), & ! [IN:optional]
845  qsplt_in(:,:,:,:), sarea(:,:,:,:), & ! [OUT:optional]
846  qtrc1_crg(:,:,:,qs_lt:qe_lt) ) ! [INOUT:optional]
847  else
849  ka, ks, ke, ia, is, ie, ja, js, je, &
850  dens(:,:,:), pres(:,:,:), ccn(:,:,:), dt_mp, & ! [IN]
851  temp1(:,:,:), qtrc1(:,:,:,qs_mp:qe_mp), cptot1(:,:,:), cvtot1(:,:,:), & ! [INOUT]
852  rhoe_t(:,:,:), evaporate(:,:,:) ) ! [OUT]
853  endif
854 
855  do iq = qs_mp, qe_mp
856  do j = js, je
857  do i = is, ie
858  do k = ks, ke
859  rhoq_t_mp(k,i,j,iq) = ( qtrc1(k,i,j,iq) - qtrc(k,i,j,iq) ) * dens(k,i,j) / dt_mp
860  enddo
861  enddo
862  enddo
863  enddo
864 
865  do j = js, je
866  do i = is, ie
867  do k = ks, ke
868  cptot_t(k,i,j) = ( cptot1(k,i,j) - cptot(k,i,j) ) / dt_mp
869  cvtot_t(k,i,j) = ( cvtot1(k,i,j) - cvtot(k,i,j) ) / dt_mp
870  end do
871  end do
872  end do
873 
874  if( flg_lt ) then
875  do iq = qs_lt, qe_lt
876  do j = js, je
877  do i = is, ie
878  do k = ks, ke
879  rhoc_t_mp(k,i,j,iq) = ( qtrc1_crg(k,i,j,iq) - qtrc(k,i,j,iq) ) * dens(k,i,j) / dt_mp
880  enddo
881  enddo
882  enddo
883  enddo
884  endif
885 
886  case ( 'SN14' )
887 
888  if( flg_lt ) then
890  ka, ks, ke, ia, is, ie, ja, js, je, & ! [IN]
891  qla, & ! [IN]
892  temp(:,:,:), dens(:,:,:), & ! [IN]
893  qtrc(:,:,:,qls:qle), & ! [IN]
894  dqcrg(:,:,:), beta_crg(:,:,:) ) ! [OUT]
896  ka, ks, ke, ia, is, ie, ja, js, je, &
897  dens(:,:,:), w(:,:,:), qtrc(:,:,:,qs_mp:qe_mp), pres(:,:,:), temp(:,:,:), & ! [IN]
898  qdry(:,:,:), cptot(:,:,:), cvtot(:,:,:), ccn(:,:,:), dt_mp, real_cz(:,:,:), real_fz(:,:,:), & ! [IN]
899  rhoq_t_mp(:,:,:,qs_mp:qe_mp), rhoe_t(:,:,:), cptot_t(:,:,:), cvtot_t(:,:,:), evaporate(:,:,:),& ! [OUT]
900  flg_lt, d0_crg, v0_crg, dqcrg(:,:,:), beta_crg(:,:,:), & ! [IN:optional]
901  qtrc(:,:,:,qs_lt:qe_lt), & ! [IN:optional]
902  qsplt_in(:,:,:,:), sarea(:,:,:,:), & ! [OUT:optional]
903  rhoc_t_mp(:,:,:,qs_lt:qe_lt) ) ! [OUT:optional]
904  else
906  ka, ks, ke, ia, is, ie, ja, js, je, &
907  dens(:,:,:), w(:,:,:), qtrc(:,:,:,qs_mp:qe_mp), pres(:,:,:), temp(:,:,:), & ! [IN]
908  qdry(:,:,:), cptot(:,:,:), cvtot(:,:,:), ccn(:,:,:), dt_mp, real_cz(:,:,:), real_fz(:,:,:), & ! [IN]
909  rhoq_t_mp(:,:,:,qs_mp:qe_mp), rhoe_t(:,:,:), cptot_t(:,:,:), cvtot_t(:,:,:), evaporate(:,:,:) ) ! [OUT]
910  endif
911 
912  case ( 'SUZUKI10' )
913 
914  if( flg_lt ) then
916  ka, ks, ke, ia, is, ie, ja, js, je, & ! [IN]
917  qla, & ! [IN]
918  temp(:,:,:), dens(:,:,:), & ! [IN]
919  qtrc(:,:,:,qls:qle), & ! [IN]
920  dqcrg(:,:,:), beta_crg(:,:,:) ) ! [OUT]
921  call atmos_phy_mp_suzuki10_tendency( ka, ks, ke, ia, is, ie, ja, js, je, kijmax, &
922  dt_mp, & ! [IN]
923  dens(:,:,:), pres(:,:,:), temp(:,:,:), & ! [IN]
924  qtrc(:,:,:,qs_mp:qe_mp), qdry(:,:,:), & ! [IN]
925  cptot(:,:,:), cvtot(:,:,:), & ! [IN]
926  ccn(:,:,:), & ! [IN]
927  rhoq_t_mp(:,:,:,qs_mp:qe_mp), & ! [OUT]
928  rhoe_t(:,:,:), & ! [OUT]
929  cptot_t(:,:,:), cvtot_t(:,:,:), & ! [OUT]
930  evaporate(:,:,:), & ! [OUT]
931  flg_lt, d0_crg, v0_crg, & ! [IN:optional]
932  dqcrg(:,:,:), beta_crg(:,:,:), & ! [IN:optional]
933  qtrc(:,:,:,qs_lt:qe_lt), & ! [IN:optional]
934  qsplt_in(:,:,:,:), sarea(:,:,:,:), & ! [OUT:optional]
935  rhoc_t_mp(:,:,:,qs_lt:qe_lt) ) ! [OUT:optional]
936  else
937  call atmos_phy_mp_suzuki10_tendency( ka, ks, ke, ia, is, ie, ja, js, je, kijmax, &
938  dt_mp, & ! [IN]
939  dens(:,:,:), pres(:,:,:), temp(:,:,:), & ! [IN]
940  qtrc(:,:,:,qs_mp:qe_mp), qdry(:,:,:), & ! [IN]
941  cptot(:,:,:), cvtot(:,:,:), & ! [IN]
942  ccn(:,:,:), & ! [IN]
943  rhoq_t_mp(:,:,:,qs_mp:qe_mp), & ! [OUT]
944  rhoe_t(:,:,:), & ! [OUT]
945  cptot_t(:,:,:), cvtot_t(:,:,:), & ! [OUT]
946  evaporate(:,:,:) ) ! [OUT]
947  endif
948 
949  call atmos_phy_mp_suzuki10_qtrc2qhyd( ka, ks, ke, ia, is, ie, ja, js, je, & ! [IN]
950  qtrc(:,:,:,qs_mp+1:qe_mp), & ! [IN]
951  qe(:,:,:,:) ) ! [OUT]
952 
953  do iq = 1, n_hyd
954  call file_history_query( hist_hyd_id(iq), hist_sw(iq) )
955  enddo
956  do iq = 1, n_hyd
957  if( hist_sw(iq) ) then
958  call file_history_put( hist_hyd_id(iq), qe(:,:,:,iq) )
959  endif
960  enddo
961 
962  if( flg_lt ) then
963 
964  call atmos_phy_mp_suzuki10_crg_qtrc2qhyd( ka, ks, ke, ia, is, ie, ja, js, je, & ! [IN]
965  qtrc(:,:,:,qs_lt:qe_lt), & ! [IN]
966  qecrg(:,:,:,:) ) ! [OUT]
967 
968  do iq = 1, n_hyd
969  call file_history_query( hist_crg_id(iq), hist_crg_sw(iq) )
970  enddo
971  do iq = 1, n_hyd
972  if( hist_crg_sw(iq) ) then
973  call file_history_put( hist_crg_id(iq), qecrg(:,:,:,iq) )
974  endif
975  enddo
976 
977  endif
978 
979  end select
980 
981 
982  do j = js, je
983  do i = is, ie
984  do k = ks, ke
985  rhoh_mp(k,i,j) = rhoe_t(k,i,j) &
986  - ( cptot_t(k,i,j) + log( pres(k,i,j) / pre00 ) * ( cvtot(k,i,j) / cptot(k,i,j) * cptot_t(k,i,j) - cvtot_t(k,i,j) ) ) &
987  * dens(k,i,j) * temp(k,i,j)
988 ! RHOT_t_MP(k,i,j) = RHOE_t(k,i,j) / ( EXNER(k,i,j) * CPtot(k,i,j) ) &
989 ! - RHOT(k,i,j) * CPtot_t(k,i,j) / CPtot(k,i,j) &
990 ! - RHOT(k,i,j) * CVtot(k,i,j) / ( CPtot(k,i,j) - CVtot(k,i,j) ) &
991 ! * log( EXNER(k,i,j) ) * ( CPtot_t(k,i,j) / CPtot(k,i,j) - CVtot_t(k,i,j) / CVtot(k,i,j) )
992  end do
993  end do
994  end do
995 
996  if ( mp_do_precipitation ) then
997 
998  call prof_rapstart('MP_Precipitation', 2)
999 
1000  ! prepare for history output
1001  hist_vterm_idx(:) = -1
1002  ih = 0
1003  do iq = qs_mp+1, qe_mp
1004  call file_history_query( hist_vterm_id(iq), flag )
1005  if ( flag ) then
1006  ih = ih + 1
1007  hist_vterm_idx(iq) = ih
1008  end if
1009  end do
1010  if ( ih > 0 ) then
1011  allocate( vterm_hist(ka,ia,ja,ih) )
1012  vterm_hist(:,:,:,:) = 0.0_rp
1013  end if
1014 
1015  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
1016  !$omp shared (KA,KS,KE,IS,IE,JS,JE,QS_MP,QE_MP,QHA,QHS,QHE,QLA,QIA, &
1017  !$omp QA_LT,QS_LT,QE_LT, &
1018  !$omp PRE00,LHF, &
1019  !$omp ATMOS_PHY_MP_TYPE, ATMOS_PHY_PRECIP_TYPE, &
1020  !$omp dt_MP,MP_NSTEP_SEDIMENTATION,MP_DTSEC_SEDIMENTATION,MP_RNSTEP_SEDIMENTATION, &
1021  !$omp REAL_CZ,REAL_FZ, &
1022  !$omp DENS,MOMZ,U,V,RHOT,TEMP,PRES,QTRC,CPtot,CVtot,EXNER, &
1023  !$omp DENS_t_MP,MOMZ_t_MP,RHOU_t_MP,RHOV_t_MP,RHOQ_t_MP,RHOH_MP, &
1024  !$omp SFLX_rain,SFLX_snow,SFLX_ENGI, &
1025  !$omp REFSTATE_dens, &
1026  !$omp flg_lt,RHOC_t_MP, &
1027  !$omp vterm_hist,hist_vterm_idx) &
1028  !$omp private(i,j,k,iq,step, &
1029  !$omp FZ,FDZ,RFDZ,RCDZ, &
1030  !$omp DENS2,TEMP2,PRES2,CPtot2,CVtot2,RHOE,RHOE2,RHOQ,RHOQ2, &
1031  !$omp RHOQ2_crg,mflux_crg,sflux_crg,eflux_crg, &
1032  !$omp vterm,mflux,sflux,eflux,FLX_hydro,CP_t,CV_t)
1033  do j = js, je
1034  do i = is, ie
1035 
1036  fz(1:ka) = real_fz(1:ka,i,j)
1037 
1038  fdz(ks-1) = real_cz(ks,i,j) - real_fz(ks-1,i,j)
1039  rfdz(ks-1) = 1.0_rp / fdz(ks-1)
1040  do k = ks, ke
1041  fdz(k) = real_cz(k+1,i,j) - real_cz(k ,i,j)
1042  rfdz(k) = 1.0_rp / fdz(k)
1043  rcdz(k) = 1.0_rp / ( real_fz(k ,i,j) - real_fz(k-1,i,j) )
1044  enddo
1045 
1046  do k = ks, ke
1047  dens2(k) = dens(k,i,j)
1048  temp2(k) = temp(k,i,j)
1049  pres2(k) = pres(k,i,j)
1050  cptot2(k) = cptot(k,i,j)
1051  cvtot2(k) = cvtot(k,i,j)
1052  rhoe(k) = temp(k,i,j) * cvtot(k,i,j) * dens2(k)
1053  rhoe2(k) = rhoe(k)
1054  end do
1055  do iq = qs_mp+1, qe_mp
1056  do k = ks, ke
1057  rhoq(k,iq) = dens2(k) * qtrc(k,i,j,iq) + rhoq_t_mp(k,i,j,iq) * dt_mp
1058  rhoq2(k,iq) = rhoq(k,iq)
1059  end do
1060  end do
1061 
1062  if( flg_lt ) then
1063  do iq = qs_lt, qe_lt
1064  do k = ks, ke
1065  rhoq2_crg(k,iq) = dens2(k) * qtrc(k,i,j,iq)
1066  end do
1067  end do
1068  endif
1069 
1070  sflx_rain(i,j) = 0.0_rp
1071  sflx_snow(i,j) = 0.0_rp
1072  sflx_engi(i,j) = 0.0_rp
1073  flx_hydro(:) = 0.0_rp
1074  do step = 1, mp_nstep_sedimentation
1075 
1076  select case ( atmos_phy_mp_type )
1077  case ( 'KESSLER' )
1079  ka, ks, ke, &
1080  dens2(:), rhoq2(:,:), & ! [IN]
1081  refstate_dens(:,i,j), & ! [IN]
1082  vterm(:,:) ) ! [OUT]
1083  case ( 'TOMITA08' )
1085  ka, ks, ke, &
1086  dens2(:), temp2(:), rhoq2(:,:), & ! [IN]
1087  vterm(:,:) ) ! [OUT]
1088  case ( 'SN14' )
1090  ka, ks, ke, &
1091  dens2(:), temp2(:), rhoq2(:,:), pres2(:), & ! [IN]
1092  vterm(:,:) ) ! [OUT]
1093  case ( 'SUZUKI10' )
1095  ka, & ! [IN]
1096  vterm(:,:) ) ! [OUT]
1097  case default
1098  vterm(:,:) = 0.0_rp ! tentative
1099  end select
1100 
1101  ! store to history output
1102  do iq = qs_mp+1, qe_mp
1103  if ( hist_vterm_idx(iq) > 0 ) then
1104  do k = ks, ke
1105  vterm_hist(k,i,j,hist_vterm_idx(iq)) = vterm_hist(k,i,j,hist_vterm_idx(iq)) &
1106  + vterm(k,iq) * mp_rnstep_sedimentation
1107  end do
1108  end if
1109  end do
1110 
1111  select case ( atmos_phy_precip_type )
1112  case ( 'Upwind-Euler' )
1114  ka, ks, ke, qe_mp-qs_mp, qla, qia, &
1115  temp2(:), vterm(:,:), & ! [IN]
1116  fdz(:), rcdz(:), & ! [IN]
1117  mp_dtsec_sedimentation, & ! [IN]
1118  i, j, & ! [IN]
1119  dens2(:), rhoq2(:,:), & ! [INOUT]
1120  cptot2(:), cvtot2(:), & ! [INOUT]
1121  rhoe2(:), & ! [INOUT]
1122  mflux(:), sflux(:), & ! [OUT]
1123  eflux ) ! [OUT]
1124  case ( 'Semilag' )
1126  ka, ks, ke, qe_mp-qs_mp, qla, qia, &
1127  temp2(:), vterm(:,:), & ! [IN]
1128  fz(:), fdz(:), rcdz(:), & ! [IN]
1129  mp_dtsec_sedimentation, & ! [IN]
1130  i, j, & ! [IN]
1131  dens2(:), rhoq2(:,:), & ! [INOUT]
1132  cptot2(:), cvtot2(:), & ! [INOUT]
1133  rhoe2(:), & ! [INOUT]
1134  mflux(:), sflux(:), & ! [OUT]
1135  eflux ) ! [OUT]
1136  case default
1137  mflux(:) = 0.0_rp
1138  sflux(:) = 0.0_rp
1139  eflux = 0.0_rp
1140  end select
1141 
1142  if( flg_lt ) then
1143 
1144  select case ( atmos_phy_precip_type )
1145  case ( 'Upwind-Euler' )
1147  ka, ks, ke, qa_lt, 0, 0, & ! no mass tracer for charge density
1148  temp2(:), vterm(:,qhs:qhe), & ! [IN]
1149  fdz(:), rcdz(:), & ! [IN]
1150  mp_dtsec_sedimentation, & ! [IN]
1151  i, j, & ! [IN]
1152  dens2(:), rhoq2_crg(:,:), & ! [INOUT]
1153  cptot2(:), cvtot2(:), & ! [INOUT]
1154  rhoe2(:), & ! [INOUT]
1155  mflux_crg(:), sflux_crg(:), & ! [OUT] dummy
1156  eflux_crg ) ! [OUT] dummy
1157  case ( 'Semilag' )
1159  ka, ks, ke, qa_lt, 0, 0, & ! no mass tracer for charge density
1160  temp2(:), vterm(:,qhs:qhe), & ! [IN]
1161  fz(:), fdz(:), rcdz(:), & ! [IN]
1162  mp_dtsec_sedimentation, & ! [IN]
1163  i, j, & ! [IN]
1164  dens2(:), rhoq2_crg(:,:), & ! [INOUT]
1165  cptot2(:), cvtot2(:), & ! [INOUT]
1166  rhoe2(:), & ! [INOUT]
1167  mflux_crg(:), sflux_crg(:), & ! [OUT] dummy
1168  eflux_crg ) ! [OUT] dummy
1169  case default
1170  mflux_crg(:) = 0.0_rp
1171  sflux_crg(:) = 0.0_rp
1172  eflux_crg = 0.0_rp
1173  end select
1174 
1175  endif
1176 
1177  do k = ks, ke
1178  temp2(k) = rhoe2(k) / ( dens2(k) * cvtot2(k) )
1179  end do
1180 
1181  do k = ks-1, ke-1
1182  flx_hydro(k) = flx_hydro(k) + mflux(k) * mp_rnstep_sedimentation
1183  enddo
1184 
1185  sflx_rain(i,j) = sflx_rain(i,j) - sflux(1) * mp_rnstep_sedimentation
1186  sflx_snow(i,j) = sflx_snow(i,j) - sflux(2) * mp_rnstep_sedimentation
1187  sflx_engi(i,j) = sflx_engi(i,j) - eflux * mp_rnstep_sedimentation
1188 
1189  enddo
1190  sflx_engi(i,j) = sflx_engi(i,j) - sflx_snow(i,j) * lhf ! moist internal energy
1191 
1192 !OCL XFILL
1193  do k = ks, ke
1194  dens_t_mp(k,i,j) = ( dens2(k) - dens(k,i,j) ) / dt_mp
1195  end do
1196 
1197  do k = ks, ke
1198  cp_t = ( cptot2(k) - cptot(k,i,j) ) / dt_mp
1199  cv_t = ( cvtot2(k) - cvtot(k,i,j) ) / dt_mp
1200  rhoh_mp(k,i,j) = rhoh_mp(k,i,j) &
1201  + ( rhoe2(k) - rhoe(k) ) / dt_mp &
1202  - ( cp_t + log( pres(k,i,j) / pre00 ) * ( cvtot(k,i,j) / cptot(k,i,j) * cp_t - cv_t ) ) &
1203  * dens(k,i,j) * temp(k,i,j)
1204 ! RHOT_t_MP(k,i,j) = RHOT_t_MP(k,i,j) &
1205 ! + ( RHOE2(k) - RHOE(k) ) / ( dt_MP * EXNER(k,i,j) * CPtot(k,i,j) ) &
1206 ! - RHOT(k,i,j) * CP_t / CPtot(k,i,j) &
1207 ! - RHOT(k,i,j) * CVtot(k,i,j) / ( CPtot(k,i,j) - CVtot(k,i,j) ) &
1208 ! * log( EXNER(k,i,j) ) * ( CP_t / CPtot(k,i,j) - CV_t / CVtot(k,i,j) )
1209  end do
1210 
1211  do iq = qs_mp+1, qe_mp
1212  do k = ks, ke
1213  rhoq_t_mp(k,i,j,iq) = rhoq_t_mp(k,i,j,iq) &
1214  + ( rhoq2(k,iq) - rhoq(k,iq) ) / dt_mp
1215  enddo
1216  enddo
1217 
1218  if( flg_lt ) then
1219  do iq = qs_lt, qe_lt
1220  do k = ks, ke
1221  rhoc_t_mp(k,i,j,iq) = rhoc_t_mp(k,i,j,iq) &
1222  + ( rhoq2_crg(k,iq) - dens(k,i,j) * qtrc(k,i,j,iq) ) / dt_mp
1223  enddo
1224  enddo
1225  endif
1226 
1228  ka, ks, ke, &
1229  dens(:,i,j), momz(:,i,j), u(:,i,j), v(:,i,j), & ! [IN]
1230  flx_hydro(:), & ! [IN]
1231  rcdz(:), rfdz(:), & ! [IN]
1232  momz_t_mp(:,i,j), rhou_t_mp(:,i,j), rhov_t_mp(:,i,j) ) ! [OUT]
1233 
1234  enddo
1235  enddo
1236 
1237  ! history output
1238  do iq = qs_mp+1, qe_mp
1239  if ( hist_vterm_idx(iq) > 0 ) &
1240  call file_history_put( hist_vterm_id(iq), vterm_hist(:,:,:,hist_vterm_idx(iq)) )
1241  end do
1242  if ( allocated( vterm_hist ) ) deallocate( vterm_hist )
1243 
1244  call prof_rapend ('MP_Precipitation', 2)
1245 
1246  end if
1247 
1248 !OCL XFILL
1249  do j = js, je
1250  do i = is, ie
1251  precip(i,j) = sflx_rain(i,j) + sflx_snow(i,j)
1252  end do
1253  end do
1254 
1255  call file_history_in( sflx_rain(:,:), 'RAIN_MP', 'surface rain rate by MP', 'kg/m2/s', fill_halo=.true. )
1256  call file_history_in( sflx_snow(:,:), 'SNOW_MP', 'surface snow rate by MP', 'kg/m2/s', fill_halo=.true. )
1257  call file_history_in( precip(:,:), 'PREC_MP', 'surface precipitation rate by MP', 'kg/m2/s', fill_halo=.true. )
1258  call file_history_in( evaporate(:,:,:), 'EVAPORATE', 'evaporated cloud number', 'num/m3/s', fill_halo=.true. )
1259 
1260  call file_history_in( dens_t_mp(:,:,:), 'DENS_t_MP', 'tendency DENS in MP', 'kg/m3/s' , fill_halo=.true. )
1261  call file_history_in( momz_t_mp(:,:,:), 'MOMZ_t_MP', 'tendency MOMZ in MP', 'kg/m2/s2' , fill_halo=.true. )
1262  call file_history_in( rhou_t_mp(:,:,:), 'RHOU_t_MP', 'tendency RHOU in MP', 'kg/m2/s2' , fill_halo=.true. )
1263  call file_history_in( rhov_t_mp(:,:,:), 'RHOV_t_MP', 'tendency RHOV in MP', 'kg/m2/s2' , fill_halo=.true. )
1264 ! call FILE_HISTORY_in( RHOT_t_MP(:,:,:), 'RHOT_t_MP', 'tendency RHOT in MP', 'K*kg/m3/s', fill_halo=.true. )
1265  call file_history_in( rhoh_mp(:,:,:), 'RHOH_MP', 'diabatic heating rate in MP', 'J/m3/s', fill_halo=.true. )
1266  do iq = qs_mp, qe_mp
1267  call file_history_in( rhoq_t_mp(:,:,:,iq), trim(tracer_name(iq))//'_t_MP', &
1268  'tendency rho*'//trim(tracer_name(iq))//' in MP', 'kg/m3/s', fill_halo=.true. )
1269  enddo
1270 
1271  if ( flg_lt ) then
1272  call file_history_in( qsplt_in(:,:,:,1), 'QSPLT_G', 'Charge split of QG by Non-inductive process', 'fC/m3/s', fill_halo=.true. )
1273  call file_history_in( qsplt_in(:,:,:,2), 'QSPLT_I', 'Charge split of QI by Non-inductive process', 'fC/m3/s', fill_halo=.true. )
1274  call file_history_in( qsplt_in(:,:,:,3), 'QSPLT_S', 'Charge split of QS by Non-inductive process', 'fC/m3/s', fill_halo=.true. )
1275 
1276  do iq = qs_lt, qe_lt
1277  call file_history_in( rhoc_t_mp(:,:,:,iq), trim(tracer_name(iq))//'_t_MP', &
1278  'tendency rho*'//trim(tracer_name(iq))//' in MP', 'fC/m3/s', fill_halo=.true. )
1279  enddo
1280 
1281  end if
1282 
1283  endif
1284 
1285  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
1286  !$omp shared(KS,KE,IS,IE,JS,JE, &
1287  !$omp DENS_t_MP,MOMZ_t_MP,RHOU_t_MP,RHOV_t_MP,RHOH_MP, &
1288  !$omp DENS_t,MOMZ_t,RHOU_t,RHOV_t,RHOH)
1289  do j = js, je
1290  do i = is, ie
1291  do k = ks, ke
1292  dens_t(k,i,j) = dens_t(k,i,j) + dens_t_mp(k,i,j)
1293  momz_t(k,i,j) = momz_t(k,i,j) + momz_t_mp(k,i,j)
1294  rhou_t(k,i,j) = rhou_t(k,i,j) + rhou_t_mp(k,i,j)
1295  rhov_t(k,i,j) = rhov_t(k,i,j) + rhov_t_mp(k,i,j)
1296 ! RHOT_t(k,i,j) = RHOT_t(k,i,j) + RHOT_t_MP(k,i,j)
1297  rhoh(k,i,j) = rhoh(k,i,j) + rhoh_mp(k,i,j)
1298  enddo
1299  enddo
1300  enddo
1301 
1302  do iq = qs_mp, qe_mp
1303  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ &
1304  !$omp shared(JS,JE,IS,IE,KS,KE,RHOQ_t,iq,RHOQ_t_MP)
1305  do j = js, je
1306  do i = is, ie
1307  do k = ks, ke
1308  rhoq_t(k,i,j,iq) = rhoq_t(k,i,j,iq) + rhoq_t_mp(k,i,j,iq)
1309  enddo
1310  enddo
1311  enddo
1312  enddo
1313 
1314  if( flg_lt ) then
1315  do iq = qs_lt, qe_lt
1316  do j = js, je
1317  do i = is, ie
1318  do k = ks, ke
1319  rhoq_t(k,i,j,iq) = rhoq_t(k,i,j,iq) &
1320  + rhoc_t_mp(k,i,j,iq)
1321  enddo
1322  enddo
1323  enddo
1324  enddo
1325  endif
1326 
1327  if ( statistics_checktotal ) then
1328  call statistics_total( ka, ks, ke, ia, is, ie, ja, js, je, &
1329  dens_t_mp(:,:,:), 'DENS_t_MP', &
1330  atmos_grid_cartesc_real_vol(:,:,:), &
1332  call statistics_total( ka, ks, ke, ia, is, ie, ja, js, je, &
1333  momz_t_mp(:,:,:), 'MOMZ_t_MP', &
1336  call statistics_total( ka, ks, ke, ia, is, ie, ja, js, je, &
1337  rhoh_mp(:,:,:), 'RHOH_MP', &
1338  atmos_grid_cartesc_real_vol(:,:,:), &
1340 ! call STATISTICS_total( KA, KS, KE, IA, IS, IE, JA, JS, JE, &
1341 ! RHOT_t_MP(:,:,:), 'RHOT_t_MP', &
1342 ! ATMOS_GRID_CARTESC_REAL_VOL(:,:,:), &
1343 ! ATMOS_GRID_CARTESC_REAL_TOTVOL )
1344 
1345  do iq = qs_mp, qe_mp
1346  call statistics_total( ka, ks, ke, ia, is, ie, ja, js, je, &
1347  rhoq_t_mp(:,:,:,iq), trim(tracer_name(iq))//'_t_MP', &
1348  atmos_grid_cartesc_real_vol(:,:,:), &
1350  enddo
1351  endif
1352 
1353  return

References scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_cz, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_fz, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_totvol, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_totvolwxy, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_totvolzuy, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_totvolzxv, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_vol, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_volwxy, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_volzuy, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_volzxv, mod_atmos_phy_ae_vars::atmos_phy_ae_ccn_t, mod_atmos_phy_lt_vars::atmos_phy_lt_sarea, scale_atmos_phy_lt_sato2019::atmos_phy_lt_sato2019_select_dqcrg_from_lut(), mod_atmos_phy_mp_vars::atmos_phy_mp_dens_t, scale_atmos_phy_mp_kessler::atmos_phy_mp_kessler_adjustment(), scale_atmos_phy_mp_kessler::atmos_phy_mp_kessler_terminal_velocity(), mod_atmos_phy_mp_vars::atmos_phy_mp_momz_t, scale_atmos_phy_mp_common::atmos_phy_mp_precipitation_momentum(), scale_atmos_phy_mp_common::atmos_phy_mp_precipitation_semilag(), scale_atmos_phy_mp_common::atmos_phy_mp_precipitation_upwind(), scale_atmos_phy_mp_sn14::atmos_phy_mp_sn14_tendency(), 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_suzuki10::atmos_phy_mp_suzuki10_qtrc2qhyd(), scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_tendency(), scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_terminal_velocity(), scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08_adjustment(), scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08_terminal_velocity(), mod_atmos_admin::atmos_phy_mp_type, mod_atmos_admin::atmos_phy_precip_type, scale_atmos_refstate::atmos_refstate_dens, mod_atmos_admin::atmos_sw_phy_lt, scale_const::const_pre00, mod_atmos_phy_lt_vars::d0_crg, mod_atmos_vars::dens, mod_atmos_vars::dens_av, mod_atmos_vars::dens_tp, mod_atmos_phy_lt_vars::flg_lt, scale_atmos_hydrometeor::i_qv, scale_atmos_grid_cartesc_index::ia, scale_atmos_grid_cartesc_index::ie, scale_atmos_grid_cartesc_index::is, scale_atmos_grid_cartesc_index::ja, scale_atmos_grid_cartesc_index::je, scale_atmos_grid_cartesc_index::js, scale_tracer::k, scale_atmos_grid_cartesc_index::ka, scale_atmos_grid_cartesc_index::ke, scale_atmos_grid_cartesc_index::kijmax, scale_atmos_grid_cartesc_index::ks, scale_atmos_hydrometeor::lhf, mod_atmos_vars::momz, mod_atmos_vars::momz_av, mod_atmos_vars::momz_tp, scale_atmos_hydrometeor::n_hyd, mod_atmos_vars::pott, scale_prof::prof_rapend(), scale_prof::prof_rapstart(), mod_atmos_phy_lt_vars::qa_lt, mod_atmos_phy_mp_vars::qa_mp, mod_atmos_phy_lt_vars::qe_lt, mod_atmos_phy_mp_vars::qe_mp, scale_atmos_hydrometeor::qha, scale_atmos_hydrometeor::qhe, scale_atmos_hydrometeor::qhs, scale_atmos_hydrometeor::qia, scale_atmos_hydrometeor::qla, scale_atmos_hydrometeor::qle, scale_atmos_hydrometeor::qls, mod_atmos_phy_lt_vars::qs_lt, mod_atmos_phy_mp_vars::qs_mp, mod_atmos_vars::qtrc, mod_atmos_vars::qtrc_av, mod_atmos_vars::rhou_tp, mod_atmos_vars::rhov_tp, scale_statistics::statistics_checktotal, scale_time::time_dtsec_atmos_phy_lt, scale_time::time_dtsec_atmos_phy_mp, scale_tracer::tracer_name, mod_atmos_vars::u, mod_atmos_vars::v, mod_atmos_phy_lt_vars::v0_crg, and mod_atmos_vars::w.

Referenced by mod_atmos_driver::atmos_driver_calc_tendency().

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

◆ atmos_phy_mp_driver_qhyd2qtrc()

subroutine, public mod_atmos_phy_mp_driver::atmos_phy_mp_driver_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), intent(in)  QV,
real(rp), dimension(ka,ia,ja,n_hyd), intent(in)  QHYD,
real(rp), dimension(ka,ia,ja,qa_mp), intent(out)  QTRC,
real(rp), dimension(ka,ia,ja,n_hyd), intent(in), optional  QNUM 
)

Definition at line 1361 of file mod_atmos_phy_mp_driver.F90.

1361  use scale_prc, only: &
1362  prc_abort
1363  use scale_atmos_hydrometeor, only: &
1364  n_hyd
1365  use mod_atmos_phy_mp_vars, only: &
1366  qa_mp
1367  use mod_atmos_admin, only: &
1369  use scale_atmos_phy_mp_kessler, only: &
1371  use scale_atmos_phy_mp_tomita08, only: &
1373  use scale_atmos_phy_mp_sn14, only: &
1375  use scale_atmos_phy_mp_suzuki10, only: &
1377  integer, intent(in) :: KA, KS, KE
1378  integer, intent(in) :: IA, IS, IE
1379  integer, intent(in) :: JA, JS, JE
1380 
1381  real(RP), intent(in) :: QV (KA,IA,JA)
1382  real(RP), intent(in) :: QHYD(KA,IA,JA,N_HYD)
1383 
1384  real(RP), intent(out) :: QTRC(KA,IA,JA,QA_MP)
1385 
1386  real(RP), intent(in), optional :: QNUM(KA,IA,JA,N_HYD)
1387 
1388  integer :: k, i, j
1389 
1390  select case( atmos_phy_mp_type )
1391  case ( "NONE" )
1392  if ( associated( atmos_phy_mp_user_qhyd2qtrc ) ) then
1393  call atmos_phy_mp_user_qhyd2qtrc( ka, ks, ke, ia, is, ie, ja, js, je, &
1394  qv(:,:,:), qhyd(:,:,:,:), & ! [IN]
1395  qtrc(:,:,:,:), & ! [OUT]
1396  qnum=qnum ) ! [IN]
1397  end if
1398  case ( "KESSLER" )
1399  !$omp parallel do OMP_SCHEDULE_
1400  do j = js, je
1401  do i = is, ie
1402  do k = ks, ke
1403  qtrc(k,i,j,1) = qv(k,i,j)
1404  end do
1405  end do
1406  end do
1407  call atmos_phy_mp_kessler_qhyd2qtrc( ka, ks, ke, ia, is, ie, ja, js, je, &
1408  qhyd(:,:,:,:), & ! [IN]
1409  qtrc(:,:,:,2:) ) ! [OUT]
1410  case ( "TOMITA08" )
1411  !$omp parallel do OMP_SCHEDULE_
1412  do j = js, je
1413  do i = is, ie
1414  do k = ks, ke
1415  qtrc(k,i,j,1) = qv(k,i,j)
1416  end do
1417  end do
1418  end do
1419  call atmos_phy_mp_tomita08_qhyd2qtrc( ka, ks, ke, ia, is, ie, ja, js, je, &
1420  qhyd(:,:,:,:), & ! [IN]
1421  qtrc(:,:,:,2:) ) ! [OUT]
1422  case ( "SN14" )
1423  !$omp parallel do OMP_SCHEDULE_
1424  do j = js, je
1425  do i = is, ie
1426  do k = ks, ke
1427  qtrc(k,i,j,1) = qv(k,i,j)
1428  end do
1429  end do
1430  end do
1431  call atmos_phy_mp_sn14_qhyd2qtrc( ka, ks, ke, ia, is, ie, ja, js, je, &
1432  qhyd(:,:,:,:), & ! [IN]
1433  qtrc(:,:,:,2:), & ! [OUT]
1434  qnum=qnum ) ! [IN]
1435  case ( "SUZUKI10" )
1436  !$omp parallel do OMP_SCHEDULE_
1437  do j = js, je
1438  do i = is, ie
1439  do k = ks, ke
1440  qtrc(k,i,j,1) = qv(k,i,j)
1441  end do
1442  end do
1443  end do
1444  call atmos_phy_mp_suzuki10_qhyd2qtrc( ka, ks, ke, ia, is, ie, ja, js, je, &
1445  qhyd(:,:,:,:), & ! [IN]
1446  qtrc(:,:,:,2:), & ! [OUT]
1447  qnum=qnum ) ! [IN]
1448  case default
1449  log_error("ATMOS_PHY_MP_driver_qhyd2qtrc",*) 'ATMOS_PHY_MP_TYPE (', trim(atmos_phy_mp_type), ') is not supported'
1450  call prc_abort
1451  end select
1452 
1453  return

References scale_atmos_phy_mp_kessler::atmos_phy_mp_kessler_qhyd2qtrc(), scale_atmos_phy_mp_sn14::atmos_phy_mp_sn14_qhyd2qtrc(), scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_qhyd2qtrc(), scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08_qhyd2qtrc(), mod_atmos_admin::atmos_phy_mp_type, atmos_phy_mp_user_qhyd2qtrc, scale_atmos_grid_cartesc_index::ia, scale_atmos_grid_cartesc_index::ie, scale_atmos_grid_cartesc_index::is, scale_atmos_grid_cartesc_index::ja, scale_atmos_grid_cartesc_index::je, scale_atmos_grid_cartesc_index::js, scale_tracer::k, scale_atmos_grid_cartesc_index::ka, scale_atmos_grid_cartesc_index::ke, scale_atmos_grid_cartesc_index::ks, scale_atmos_hydrometeor::n_hyd, scale_prc::prc_abort(), and mod_atmos_phy_mp_vars::qa_mp.

Referenced by mod_atmos_bnd_driver::atmos_boundary_driver_send(), mod_atmos_phy_cp_driver::atmos_phy_cp_driver_calc_tendency(), mod_mkinit::mkinit(), mod_realinput_scale::parentatmosinputscale(), and mod_realinput::realinput_surface().

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

◆ atmos_phy_mp_driver_qhyd2qtrc_onlyqv()

subroutine, public mod_atmos_phy_mp_driver::atmos_phy_mp_driver_qhyd2qtrc_onlyqv ( 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)  QV,
real(rp), dimension(ka,ia,ja,n_hyd), intent(in)  QHYD,
real(rp), dimension(ka,ia,ja,qa_mp), intent(out)  QTRC,
real(rp), dimension(ka,ia,ja,n_hyd), intent(in), optional  QNUM 
)

Definition at line 1461 of file mod_atmos_phy_mp_driver.F90.

1461  use scale_atmos_hydrometeor, only: &
1462  n_hyd
1463  use mod_atmos_phy_mp_vars, only: &
1464  qa_mp
1465  integer, intent(in) :: KA, KS, KE
1466  integer, intent(in) :: IA, IS, IE
1467  integer, intent(in) :: JA, JS, JE
1468 
1469  real(RP), intent(in) :: QV (KA,IA,JA)
1470  real(RP), intent(in) :: QHYD(KA,IA,JA,N_HYD)
1471 
1472  real(RP), intent(out) :: QTRC(KA,IA,JA,QA_MP)
1473 
1474  real(RP), intent(in), optional :: QNUM(KA,IA,JA,N_HYD)
1475 
1476  integer :: k, i, j
1477 
1478  !$omp parallel do
1479  do j = js, je
1480  do i = is, ie
1481  do k = ks, ke
1482  qtrc(k,i,j,1) = qv(k,i,j)
1483  end do
1484  end do
1485  end do
1486 
1487  return

References scale_atmos_grid_cartesc_index::ie, scale_atmos_grid_cartesc_index::is, scale_atmos_grid_cartesc_index::je, scale_atmos_grid_cartesc_index::js, scale_tracer::k, scale_atmos_grid_cartesc_index::ke, scale_atmos_grid_cartesc_index::ks, scale_atmos_hydrometeor::n_hyd, and mod_atmos_phy_mp_vars::qa_mp.

Referenced by mod_atmos_driver::atmos_driver_tracer_setup().

Here is the caller graph for this function:

Variable Documentation

◆ atmos_phy_mp_user_qhyd2qtrc

procedure(qhyd2qtrc), pointer, public mod_atmos_phy_mp_driver::atmos_phy_mp_user_qhyd2qtrc => NULL()

Definition at line 54 of file mod_atmos_phy_mp_driver.F90.

54  procedure(qhyd2qtrc), pointer :: ATMOS_PHY_MP_USER_qhyd2qtrc => null()

Referenced by mod_atmos_driver::atmos_driver_tracer_setup(), and atmos_phy_mp_driver_qhyd2qtrc().

mod_atmos_phy_mp_vars::atmos_phy_mp_sflx_engi
real(rp), dimension(:,:), allocatable, public atmos_phy_mp_sflx_engi
Definition: mod_atmos_phy_mp_vars.F90:75
mod_atmos_vars::momz_av
real(rp), dimension(:,:,:), pointer, public momz_av
Definition: mod_atmos_vars.F90:90
scale_statistics
module Statistics
Definition: scale_statistics.F90:11
scale_atmos_phy_mp_kessler::atmos_phy_mp_kessler_terminal_velocity
subroutine, public atmos_phy_mp_kessler_terminal_velocity(KA, KS, KE, DENS0, RHOQ0, REFSTATE_dens_profile, vterm)
Kessler-type warm rain microphysics (terminal velocity)
Definition: scale_atmos_phy_mp_kessler.F90:356
mod_atmos_phy_mp_vars
module Atmosphere / Physics Cloud Microphysics
Definition: mod_atmos_phy_mp_vars.F90:12
scale_prc::prc_abort
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:342
scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08_terminal_velocity
subroutine, public atmos_phy_mp_tomita08_terminal_velocity(KA, KS, KE, DENS0, TEMP0, RHOQ0, vterm)
Lin-type cold rain microphysics (terminal velocity)
Definition: scale_atmos_phy_mp_tomita08.F90:2036
scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08_qhyd2qtrc
subroutine, public atmos_phy_mp_tomita08_qhyd2qtrc(KA, KS, KE, IA, IS, IE, JA, JS, JE, Qe, QTRC)
get mass ratio of each category
Definition: scale_atmos_phy_mp_tomita08.F90:2619
scale_atmos_hydrometeor::qia
integer, public qia
Definition: scale_atmos_hydrometeor.F90:122
mod_atmos_phy_mp_vars::qa_mp
integer, public qa_mp
Definition: mod_atmos_phy_mp_vars.F90:77
scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08_ntracers
integer, parameter, public atmos_phy_mp_tomita08_ntracers
Definition: scale_atmos_phy_mp_tomita08.F90:43
scale_atmos_phy_mp_sn14::atmos_phy_mp_sn14_tracer_names
character(len=h_short), dimension(qa_mp), parameter, public atmos_phy_mp_sn14_tracer_names
Definition: scale_atmos_phy_mp_sn14.F90:114
scale_atmos_grid_cartesc::atmos_grid_cartesc_cdz
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_cdz
z-length of control volume [m]
Definition: scale_atmos_grid_cartesC.F90:42
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_totvolzxv
real(rp), public atmos_grid_cartesc_real_totvolzxv
total volume (zxv, local) [m3]
Definition: scale_atmos_grid_cartesC_real.F90:90
mod_atmos_admin::atmos_phy_precip_type
character(len=h_short), public atmos_phy_precip_type
Definition: mod_atmos_admin.F90:46
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_cz
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_cz
geopotential height [m] (zxy)
Definition: scale_atmos_grid_cartesC_real.F90:38
mod_atmos_vars::pott
real(rp), dimension(:,:,:), allocatable, target, public pott
Definition: mod_atmos_vars.F90:132
scale_atmos_phy_mp_tomita08
module atmosphere / physics / microphysics / Tomita08
Definition: scale_atmos_phy_mp_tomita08.F90:13
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_volwxy
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_volwxy
control volume (wxy) [m3]
Definition: scale_atmos_grid_cartesC_real.F90:84
scale_atmos_phy_mp_kessler::atmos_phy_mp_kessler_adjustment
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.
Definition: scale_atmos_phy_mp_kessler.F90:118
mod_atmos_vars::qtrc_av
real(rp), dimension(:,:,:,:), pointer, public qtrc_av
Definition: mod_atmos_vars.F90:94
scale_atmos_hydrometeor::qhs
integer, public qhs
Definition: scale_atmos_hydrometeor.F90:115
scale_atmos_refstate::atmos_refstate_dens
real(rp), dimension(:,:,:), allocatable, public atmos_refstate_dens
refernce density [kg/m3]
Definition: scale_atmos_refstate.F90:40
mod_atmos_phy_mp_vars::atmos_phy_mp_momz_t
real(rp), dimension(:,:,:), allocatable, public atmos_phy_mp_momz_t
Definition: mod_atmos_phy_mp_vars.F90:64
scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08_nwaters
integer, parameter, public atmos_phy_mp_tomita08_nwaters
Definition: scale_atmos_phy_mp_tomita08.F90:44
scale_atmos_phy_mp_sn14::atmos_phy_mp_sn14_ntracers
integer, parameter, public atmos_phy_mp_sn14_ntracers
Definition: scale_atmos_phy_mp_sn14.F90:111
mod_atmos_vars::rhov_tp
real(rp), dimension(:,:,:), allocatable, public rhov_tp
Definition: mod_atmos_vars.F90:117
mod_atmos_admin
module ATMOS admin
Definition: mod_atmos_admin.F90:11
scale_time::time_dtsec_atmos_phy_mp
real(dp), public time_dtsec_atmos_phy_mp
time interval of physics(microphysics) [sec]
Definition: scale_time.F90:38
scale_atmos_phy_mp_sn14::atmos_phy_mp_sn14_qhyd2qtrc
subroutine, public atmos_phy_mp_sn14_qhyd2qtrc(KA, KS, KE, IA, IS, IE, JA, JS, JE, Qe, QTRC, QNUM)
Definition: scale_atmos_phy_mp_sn14.F90:995
mod_atmos_phy_lt_vars
module Atmosphere / Physics Chemistry
Definition: mod_atmos_phy_lt_vars.F90:12
scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08_tracer_units
character(len=h_short), dimension(qa_mp), parameter, public atmos_phy_mp_tomita08_tracer_units
Definition: scale_atmos_phy_mp_tomita08.F90:60
mod_atmos_phy_mp_vars::qs_mp
integer, public qs_mp
Definition: mod_atmos_phy_mp_vars.F90:78
scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_crg_qtrc2qhyd
subroutine, public atmos_phy_mp_suzuki10_crg_qtrc2qhyd(KA, KS, KE, IA, IS, IE, JA, JS, JE, QTRC0, Qecrg)
get charge density ratio of each category
Definition: scale_atmos_phy_mp_suzuki10.F90:1901
scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_qtrc2qhyd
subroutine, public atmos_phy_mp_suzuki10_qtrc2qhyd(KA, KS, KE, IA, IS, IE, JA, JS, JE, QTRC0, Qe)
Calculate mass ratio of each category.
Definition: scale_atmos_phy_mp_suzuki10.F90:1567
scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_tendency
subroutine, public atmos_phy_mp_suzuki10_tendency(KA, KS, KE, IA, IS, IE, JA, JS, JE, KIJMAX, dt, DENS, PRES, TEMP, QTRC, QDRY, CPtot, CVtot, CCN, RHOQ_t, RHOE_t, CPtot_t, CVtot_t, EVAPORATE, flg_lt, d0_crg, v0_crg, dqcrg, beta_crg, QTRC_crg, QSPLT_in, Sarea, RHOC_t_mp)
Cloud Microphysics.
Definition: scale_atmos_phy_mp_suzuki10.F90:893
mod_atmos_vars::cptot
real(rp), dimension(:,:,:), allocatable, target, public cptot
Definition: mod_atmos_vars.F90:142
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_totvolwxy
real(rp), public atmos_grid_cartesc_real_totvolwxy
total volume (wxy, local) [m3]
Definition: scale_atmos_grid_cartesC_real.F90:88
scale_const::const_eps
real(rp), public const_eps
small number
Definition: scale_const.F90:33
mod_atmos_admin::atmos_sw_phy_mp
logical, public atmos_sw_phy_mp
Definition: mod_atmos_admin.F90:52
scale_atmos_hydrometeor
module atmosphere / hydrometeor
Definition: scale_atmos_hydrometeor.F90:12
scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08_setup
subroutine, public atmos_phy_mp_tomita08_setup(KA, KS, KE, IA, IS, IE, JA, JS, JE, flg_lt)
ATMOS_PHY_MP_tomita08_setup Setup.
Definition: scale_atmos_phy_mp_tomita08.F90:355
scale_atmos_refstate
module atmosphere / reference state
Definition: scale_atmos_refstate.F90:12
mod_atmos_phy_mp_vars::atmos_phy_mp_sflx_rain
real(rp), dimension(:,:), allocatable, public atmos_phy_mp_sflx_rain
Definition: mod_atmos_phy_mp_vars.F90:73
scale_atmos_phy_mp_kessler::atmos_phy_mp_kessler_tracer_names
character(len=h_short), dimension(qa_mp), parameter, public atmos_phy_mp_kessler_tracer_names
Definition: scale_atmos_phy_mp_kessler.F90:47
scale_atmos_hydrometeor::atmos_hydrometeor_dry
logical, public atmos_hydrometeor_dry
Definition: scale_atmos_hydrometeor.F90:97
mod_atmos_phy_lt_vars::qe_lt
integer, public qe_lt
Definition: mod_atmos_phy_lt_vars.F90:60
mod_atmos_phy_lt_vars::atmos_phy_lt_sarea
real(rp), dimension(:,:,:,:), allocatable, public atmos_phy_lt_sarea
Definition: mod_atmos_phy_lt_vars.F90:56
mod_atmos_vars::rhot
real(rp), dimension(:,:,:), allocatable, target, public rhot
Definition: mod_atmos_vars.F90:79
scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08_adjustment
subroutine, public atmos_phy_mp_tomita08_adjustment(KA, KS, KE, IA, IS, IE, JA, JS, JE, DENS, PRES, CCN, dt, TEMP, QTRC, CPtot, CVtot, RHOE_t, EVAPORATE, flg_lt, d0_crg, v0_crg, dqcrg, beta_crg, QSPLT_in, Sarea, QTRC_crg)
ATMOS_PHY_MP_tomita08_adjustment calculate state after saturation process.
Definition: scale_atmos_phy_mp_tomita08.F90:588
scale_atmos_grid_cartesc_real
module Atmosphere GRID CartesC Real(real space)
Definition: scale_atmos_grid_cartesC_real.F90:11
mod_atmos_vars::qtrc
real(rp), dimension(:,:,:,:), allocatable, target, public qtrc
Definition: mod_atmos_vars.F90:80
scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_terminal_velocity
subroutine, public atmos_phy_mp_suzuki10_terminal_velocity(KA, vterm_o)
get terminal velocity
Definition: scale_atmos_phy_mp_suzuki10.F90:1349
scale_file_history
module file_history
Definition: scale_file_history.F90:15
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_volzuy
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_volzuy
control volume (zuy) [m3]
Definition: scale_atmos_grid_cartesC_real.F90:85
scale_atmos_phy_mp_sn14
module ATMOSPHERE / Physics Cloud Microphysics
Definition: scale_atmos_phy_mp_sn14.F90:51
scale_atmos_phy_mp_sn14::atmos_phy_mp_sn14_tracer_descriptions
character(len=h_mid), dimension(qa_mp), parameter, public atmos_phy_mp_sn14_tracer_descriptions
Definition: scale_atmos_phy_mp_sn14.F90:126
mod_atmos_phy_mp_vars::atmos_phy_mp_dens_t
real(rp), dimension(:,:,:), allocatable, public atmos_phy_mp_dens_t
Definition: mod_atmos_phy_mp_vars.F90:63
scale_prc
module PROCESS
Definition: scale_prc.F90:11
scale_atmos_hydrometeor::qle
integer, public qle
Definition: scale_atmos_hydrometeor.F90:120
mod_atmos_vars::rhou_tp
real(rp), dimension(:,:,:), allocatable, public rhou_tp
Definition: mod_atmos_vars.F90:116
scale_atmos_phy_mp_common::atmos_phy_mp_negative_fixer
subroutine, public atmos_phy_mp_negative_fixer(KA, KS, KE, IA, IS, IE, JA, JS, JE, QLA, QIA, limit_negative, DENS, TEMP, CVtot, CPtot, QV, QTRC, RHOH, DENS_diff, ENGI_diff)
ATMOS_PHY_MP_negative_fixer negative fixer.
Definition: scale_atmos_phy_mp_common.F90:69
scale_atmos_phy_mp_kessler::atmos_phy_mp_kessler_nwaters
integer, parameter, public atmos_phy_mp_kessler_nwaters
Definition: scale_atmos_phy_mp_kessler.F90:45
scale_atmos_phy_mp_common::atmos_phy_mp_precipitation_semilag
subroutine, public atmos_phy_mp_precipitation_semilag(KA, KS, KE, QHA, QLA, QIA, TEMP, vterm, FZ, FDZ, RCDZ, dt, i, j, DENS, RHOQ, CPtot, CVtot, RHOE, mflx, sflx, esflx)
Definition: scale_atmos_phy_mp_common.F90:548
mod_atmos_vars::dens
real(rp), dimension(:,:,:), allocatable, target, public dens
Definition: mod_atmos_vars.F90:75
scale_atmos_phy_mp_kessler::atmos_phy_mp_kessler_ntracers
integer, parameter, public atmos_phy_mp_kessler_ntracers
Definition: scale_atmos_phy_mp_kessler.F90:44
scale_atmos_hydrometeor::qhe
integer, public qhe
Definition: scale_atmos_hydrometeor.F90:116
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_volzxv
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_volzxv
control volume (zxv) [m3]
Definition: scale_atmos_grid_cartesC_real.F90:86
scale_const
module CONSTANT
Definition: scale_const.F90:11
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_vol
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_vol
control volume (zxy) [m3]
Definition: scale_atmos_grid_cartesC_real.F90:83
scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_nwaters
integer, public atmos_phy_mp_suzuki10_nwaters
Definition: scale_atmos_phy_mp_suzuki10.F90:56
mod_atmos_phy_lt_vars::v0_crg
real(rp), parameter, public v0_crg
Definition: mod_atmos_phy_lt_vars.F90:63
mod_atmos_vars::momz
real(rp), dimension(:,:,:), allocatable, target, public momz
Definition: mod_atmos_vars.F90:76
mod_atmos_admin::atmos_sw_phy_lt
logical, public atmos_sw_phy_lt
Definition: mod_atmos_admin.F90:60
scale_atmos_phy_mp_common::atmos_phy_mp_precipitation_momentum
subroutine, public atmos_phy_mp_precipitation_momentum(KA, KS, KE, DENS, MOMZ, U, V, mflx, RCDZ, RFDZ, MOMZ_t, RHOU_t, RHOV_t)
Definition: scale_atmos_phy_mp_common.F90:756
scale_atmos_phy_mp_sn14::atmos_phy_mp_sn14_terminal_velocity
subroutine, public atmos_phy_mp_sn14_terminal_velocity(KA, KS, KE, DENS, TEMP, RHOQ, PRES, vterm)
ATMOS_PHY_MP_sn14_terminal_velocity Calculate terminal velocity.
Definition: scale_atmos_phy_mp_sn14.F90:1231
scale_atmos_phy_mp_common
module ATMOSPHERE / Physics Cloud Microphysics - Common
Definition: scale_atmos_phy_mp_common.F90:13
scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_qhyd2qtrc
subroutine, public atmos_phy_mp_suzuki10_qhyd2qtrc(KA, KS, KE, IA, IS, IE, JA, JS, JE, Qe, QTRC, QNUM)
get mass ratio of each category
Definition: scale_atmos_phy_mp_suzuki10.F90:1703
scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_tracer_units
character(len=h_short), dimension(:), allocatable, public atmos_phy_mp_suzuki10_tracer_units
Definition: scale_atmos_phy_mp_suzuki10.F90:64
mod_atmos_vars::v
real(rp), dimension(:,:,:), allocatable, target, public v
Definition: mod_atmos_vars.F90:130
mod_atmos_vars::w
real(rp), dimension(:,:,:), allocatable, target, public w
Definition: mod_atmos_vars.F90:128
mod_atmos_vars::momz_tp
real(rp), dimension(:,:,:), allocatable, public momz_tp
Definition: mod_atmos_vars.F90:115
scale_atmos_phy_mp_common::atmos_phy_mp_precipitation_upwind
subroutine, public atmos_phy_mp_precipitation_upwind(KA, KS, KE, QHA, QLA, QIA, TEMP, vterm, FDZ, RCDZ, dt, i, j, DENS, RHOQ, CPtot, CVtot, RHOE, mflx, sflx, esflx)
Definition: scale_atmos_phy_mp_common.F90:423
scale_atmos_phy_mp_suzuki10
module Spectran Bin Microphysics
Definition: scale_atmos_phy_mp_suzuki10.F90:23
mod_atmos_vars::temp
real(rp), dimension(:,:,:), allocatable, target, public temp
Definition: mod_atmos_vars.F90:133
scale_monitor::monitor_reg
subroutine, public monitor_reg(name, desc, unit, itemid, ndims, dim_type, is_tendency)
Search existing item, or matching check between requested and registered item.
Definition: scale_monitor.F90:241
mod_atmos_vars::dens_tp
real(rp), dimension(:,:,:), allocatable, public dens_tp
Definition: mod_atmos_vars.F90:114
scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08_tracer_names
character(len=h_short), dimension(qa_mp), parameter, public atmos_phy_mp_tomita08_tracer_names
Definition: scale_atmos_phy_mp_tomita08.F90:46
scale_atmos_hydrometeor::qha
integer, public qha
Definition: scale_atmos_hydrometeor.F90:114
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_totvolzuy
real(rp), public atmos_grid_cartesc_real_totvolzuy
total volume (zuy, local) [m3]
Definition: scale_atmos_grid_cartesC_real.F90:89
scale_atmos_phy_mp_kessler
module atmosphere / physics / microphysics / Kessler
Definition: scale_atmos_phy_mp_kessler.F90:14
scale_time
module TIME
Definition: scale_time.F90:11
scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_tracer_descriptions
character(len=h_mid), dimension(:), allocatable, public atmos_phy_mp_suzuki10_tracer_descriptions
Definition: scale_atmos_phy_mp_suzuki10.F90:63
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
mod_atmos_phy_bl_vars::qe
integer, public qe
Definition: mod_atmos_phy_bl_vars.F90:45
scale_atmos_phy_mp_kessler::atmos_phy_mp_kessler_tracer_descriptions
character(len=h_mid), dimension(qa_mp), parameter, public atmos_phy_mp_kessler_tracer_descriptions
Definition: scale_atmos_phy_mp_kessler.F90:51
scale_tracer
module TRACER
Definition: scale_tracer.F90:12
scale_atmos_phy_mp_sn14::atmos_phy_mp_sn14_tendency
subroutine, public atmos_phy_mp_sn14_tendency(KA, KS, KE, IA, IS, IE, JA, JS, JE, DENS, W, QTRC, PRES, TEMP, Qdry, CPtot, CVtot, CCN, dt, cz, fz, RHOQ_t, RHOE_t, CPtot_t, CVtot_t, EVAPORATE, flg_lt, d0_crg, v0_crg, dqcrg, beta_crg, QTRC_crg, QSPLT_in, Sarea, RHOQcrg_t)
ATMOS_PHY_MP_sn14_tendency calculate tendency.
Definition: scale_atmos_phy_mp_sn14.F90:672
scale_atmos_hydrometeor::i_qv
integer, public i_qv
Definition: scale_atmos_hydrometeor.F90:77
mod_atmos_phy_lt_vars::qa_lt
integer, public qa_lt
Definition: mod_atmos_phy_lt_vars.F90:58
mod_atmos_phy_lt_vars::d0_crg
real(rp), parameter, public d0_crg
Definition: mod_atmos_phy_lt_vars.F90:62
mod_atmos_vars::dens_av
real(rp), dimension(:,:,:), pointer, public dens_av
Definition: mod_atmos_vars.F90:89
scale_atmos_phy_mp_kessler::atmos_phy_mp_kessler_qhyd2qtrc
subroutine, public atmos_phy_mp_kessler_qhyd2qtrc(KA, KS, KE, IA, IS, IE, JA, JS, JE, Qe, QTRC)
Calculate mass ratio of each category.
Definition: scale_atmos_phy_mp_kessler.F90:488
mod_atmos_vars::u
real(rp), dimension(:,:,:), allocatable, target, public u
Definition: mod_atmos_vars.F90:129
scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_tracer_setup
subroutine, public atmos_phy_mp_suzuki10_tracer_setup
Config.
Definition: scale_atmos_phy_mp_suzuki10.F90:299
scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_ntracers
integer, public atmos_phy_mp_suzuki10_ntracers
Definition: scale_atmos_phy_mp_suzuki10.F90:55
scale_statistics::statistics_checktotal
logical, public statistics_checktotal
calc&report variable totals to logfile?
Definition: scale_statistics.F90:64
scale_time::time_dtsec
real(dp), public time_dtsec
time interval of model [sec]
Definition: scale_time.F90:33
mod_atmos_vars
module ATMOSPHERIC Variables
Definition: mod_atmos_vars.F90:12
scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_tracer_names
character(len=h_short), dimension(:), allocatable, public atmos_phy_mp_suzuki10_tracer_names
Definition: scale_atmos_phy_mp_suzuki10.F90:62
scale_atmos_phy_mp_kessler::atmos_phy_mp_kessler_nices
integer, parameter, public atmos_phy_mp_kessler_nices
Definition: scale_atmos_phy_mp_kessler.F90:46
mod_atmos_phy_mp_vars::atmos_phy_mp_sflx_snow
real(rp), dimension(:,:), allocatable, public atmos_phy_mp_sflx_snow
Definition: mod_atmos_phy_mp_vars.F90:74
mod_atmos_phy_lt_vars::flg_lt
logical, public flg_lt
Definition: mod_atmos_phy_lt_vars.F90:64
scale_atmos_phy_mp_sn14::atmos_phy_mp_sn14_nices
integer, parameter, public atmos_phy_mp_sn14_nices
Definition: scale_atmos_phy_mp_sn14.F90:113
scale_atmos_hydrometeor::lhf
real(rp), public lhf
latent heat of fusion for use [J/kg]
Definition: scale_atmos_hydrometeor.F90:128
scale_atmos_hydrometeor::atmos_hydrometeor_regist
subroutine, public atmos_hydrometeor_regist(NL, NI, NAME, DESC, UNIT, Q0, ADVC)
ATMOS_HYDROMETEOR_regist Regist tracer.
Definition: scale_atmos_hydrometeor.F90:237
mod_atmos_phy_ae_vars
module ATMOSPHERE / Physics Aerosol Microphysics
Definition: mod_atmos_phy_ae_vars.F90:12
scale_atmos_phy_lt_sato2019
module atmosphere / physics / lightninh / SATO2019
Definition: scale_atmos_phy_lt_sato2019.F90:12
mod_atmos_admin::atmos_phy_mp_type
character(len=h_short), public atmos_phy_mp_type
Definition: mod_atmos_admin.F90:36
scale_file_history::file_history_reg
subroutine, public file_history_reg(name, desc, unit, itemid, standard_name, ndims, dim_type, cell_measures, fill_halo)
Register/Append variable to history file.
Definition: scale_file_history.F90:650
scale_atmos_phy_mp_sn14::atmos_phy_mp_sn14_nwaters
integer, parameter, public atmos_phy_mp_sn14_nwaters
Definition: scale_atmos_phy_mp_sn14.F90:112
mod_atmos_phy_mp_vars::qe_mp
integer, public qe_mp
Definition: mod_atmos_phy_mp_vars.F90:79
mod_atmos_phy_mp_vars::atmos_phy_mp_cldfrac_thleshold
real(rp), public atmos_phy_mp_cldfrac_thleshold
Definition: mod_atmos_phy_mp_vars.F90:61
scale_time::time_dtsec_atmos_phy_lt
real(dp), public time_dtsec_atmos_phy_lt
time interval of physics(lightning ) [sec]
Definition: scale_time.F90:45
scale_atmos_phy_mp_sn14::atmos_phy_mp_sn14_tracer_units
character(len=h_short), dimension(qa_mp), parameter, public atmos_phy_mp_sn14_tracer_units
Definition: scale_atmos_phy_mp_sn14.F90:138
mod_atmos_phy_ae_vars::atmos_phy_ae_ccn_t
real(rp), dimension(:,:,:), allocatable, public atmos_phy_ae_ccn_t
Definition: mod_atmos_phy_ae_vars.F90:63
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_fz
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_fz
geopotential height [m] (wxy)
Definition: scale_atmos_grid_cartesC_real.F90:42
scale_atmos_phy_lt_sato2019::atmos_phy_lt_sato2019_select_dqcrg_from_lut
subroutine, public atmos_phy_lt_sato2019_select_dqcrg_from_lut(KA, KS, KE, IA, IS, IE, JA, JS, JE, NLIQ, TEMP, DENS, QLIQ, dqcrg, beta_crg)
Select cwc-temp point on LUT.
Definition: scale_atmos_phy_lt_sato2019.F90:2621
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_totvol
real(rp), public atmos_grid_cartesc_real_totvol
total volume (zxy, local) [m3]
Definition: scale_atmos_grid_cartesC_real.F90:87
scale_atmos_grid_cartesc
module atmosphere / grid / cartesC
Definition: scale_atmos_grid_cartesC.F90:12
scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08_nices
integer, parameter, public atmos_phy_mp_tomita08_nices
Definition: scale_atmos_phy_mp_tomita08.F90:45
scale_atmos_hydrometeor::qla
integer, public qla
Definition: scale_atmos_hydrometeor.F90:118
scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_setup
subroutine, public atmos_phy_mp_suzuki10_setup(KA, IA, JA, flg_lt)
Setup.
Definition: scale_atmos_phy_mp_suzuki10.F90:391
scale_atmos_hydrometeor::qls
integer, public qls
Definition: scale_atmos_hydrometeor.F90:119
scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_nices
integer, public atmos_phy_mp_suzuki10_nices
Definition: scale_atmos_phy_mp_suzuki10.F90:57
scale_atmos_phy_mp_sn14::atmos_phy_mp_sn14_setup
subroutine, public atmos_phy_mp_sn14_setup(KA, IA, JA)
ATMOS_PHY_MP_sn14_setup setup.
Definition: scale_atmos_phy_mp_sn14.F90:597
scale_atmos_phy_mp_kessler::atmos_phy_mp_kessler_setup
subroutine, public atmos_phy_mp_kessler_setup
ATMOS_PHY_MP_kessler_setup Setup.
Definition: scale_atmos_phy_mp_kessler.F90:90
scale_const::const_pre00
real(rp), public const_pre00
pressure reference [Pa]
Definition: scale_const.F90:88
scale_atmos_hydrometeor::n_hyd
integer, parameter, public n_hyd
Definition: scale_atmos_hydrometeor.F90:79
scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_nccn
integer, public atmos_phy_mp_suzuki10_nccn
Definition: scale_atmos_phy_mp_suzuki10.F90:58
mod_atmos_vars::cvtot
real(rp), dimension(:,:,:), allocatable, target, public cvtot
Definition: mod_atmos_vars.F90:141
scale_atmos_phy_mp_kessler::atmos_phy_mp_kessler_tracer_units
character(len=h_short), dimension(qa_mp), parameter, public atmos_phy_mp_kessler_tracer_units
Definition: scale_atmos_phy_mp_kessler.F90:55
scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08_tracer_descriptions
character(len=h_mid), dimension(qa_mp), parameter, public atmos_phy_mp_tomita08_tracer_descriptions
Definition: scale_atmos_phy_mp_tomita08.F90:53
mod_atmos_phy_lt_vars::qs_lt
integer, public qs_lt
Definition: mod_atmos_phy_lt_vars.F90:59
scale_monitor
module MONITOR
Definition: scale_monitor.F90:12