SCALE-RM
Functions/Subroutines | Variables
scale_atmos_phy_mp_tomita08 Module Reference

module atmosphere / physics / microphysics / Tomita08 More...

Functions/Subroutines

subroutine, public atmos_phy_mp_tomita08_setup (KA, KS, KE, IA, IS, IE, JA, JS, JE)
 ATMOS_PHY_MP_tomita08_setup Setup. More...
 
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)
 ATMOS_PHY_MP_tomita08_adjustment calculate state after saturation process. More...
 
subroutine, public atmos_phy_mp_tomita08_terminal_velocity (KA, KS, KE, DENS0, TEMP0, RHOQ0, vterm)
 Lin-type cold rain microphysics (terminal velocity) More...
 
subroutine, public atmos_phy_mp_tomita08_cloud_fraction (KA, KS, KE, IA, IS, IE, JA, JS, JE, QTRC, mask_criterion, cldfrac)
 Calculate Cloud Fraction. More...
 
subroutine, public atmos_phy_mp_tomita08_effective_radius (KA, KS, KE, IA, IS, IE, JA, JS, JE, DENS0, TEMP0, QTRC0, Re)
 Calculate Effective Radius. More...
 
subroutine, public atmos_phy_mp_tomita08_qtrc2qhyd (KA, KS, KE, IA, IS, IE, JA, JS, JE, QTRC, Qe)
 Calculate mass ratio of each category. More...
 
subroutine, public atmos_phy_mp_tomita08_qhyd2qtrc (KA, KS, KE, IA, IS, IE, JA, JS, JE, Qe, QTRC)
 get mass ratio of each category More...
 

Variables

integer, parameter, public atmos_phy_mp_tomita08_ntracers = QA_MP
 
integer, parameter, public atmos_phy_mp_tomita08_nwaters = 2
 
integer, parameter, public atmos_phy_mp_tomita08_nices = 3
 
character(len=h_short), dimension(qa_mp), parameter, public atmos_phy_mp_tomita08_tracer_names = (/ 'QV', 'QC', 'QR', 'QI', 'QS', 'QG' /)
 
character(len=h_mid), dimension(qa_mp), parameter, public atmos_phy_mp_tomita08_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 ', 'Ratio of Cloud Ice mass ratio to total mass ', 'Ratio of Snow miass ratio to total mass ', 'Ratio of Graupel mass ratio to total mass '/)
 
character(len=h_short), dimension(qa_mp), parameter, public atmos_phy_mp_tomita08_tracer_units = (/ 'kg/kg', 'kg/kg', 'kg/kg', 'kg/kg', 'kg/kg', 'kg/kg' /)
 

Detailed Description

module atmosphere / physics / microphysics / Tomita08

Description
Cloud Microphysics by Lin-type parametarization Reference: Tomita(2008)
Author
Team SCALE
NAMELIST
  • PARAM_ATMOS_PHY_MP_TOMITA08
    nametypedefault valuecomment
    DO_COUPLE_AEROSOL logical apply CCN effect?
    DO_EXPLICIT_ICEGEN logical apply explicit ice generation?
    AUTOCONV_NC real(RP) Nc_ocn number concentration of cloud water [1/cc]
    ENABLE_KK2000 logical .false. use scheme by Khairoutdinov and Kogan (2000)
    ENABLE_RS2014 logical .false. use scheme by Roh and Satoh (2014)
    ENABLE_WDXZ2014 logical .false. use scheme by Wainwright et al. (2014)
    N0R_DEF real(RP) 8.E+6_RP intercept parameter for rain [1/m4]
    N0S_DEF real(RP) 3.E+6_RP intercept parameter for snow [1/m4]
    N0G_DEF real(RP) 4.E+6_RP intercept parameter for graupel [1/m4]
    DENS_S real(RP) 100.0_RP density of snow [kg/m3]
    DENS_G real(RP) 400.0_RP density of graupel [kg/m3]
    RE_QC real(RP) 8.E-6_RP effective radius for cloud water
    RE_QI real(RP) 40.E-6_RP effective radius for cloud ice
    DRAG_G real(RP) 0.6_RP drag coefficient for graupel
    CR real(RP) 130.0_RP
    CS real(RP) 4.84_RP
    EIW real(RP) 1.0_RP collection efficiency of cloud ice for cloud water
    ERW real(RP) 1.0_RP collection efficiency of rain for cloud water
    ESW real(RP) 1.0_RP collection efficiency of snow for cloud water
    EGW real(RP) 1.0_RP collection efficiency of graupel for cloud water
    ERI real(RP) 1.0_RP collection efficiency of rain for cloud ice
    ESI real(RP) 1.0_RP collection efficiency of snow for cloud ice
    EGI real(RP) 0.1_RP collection efficiency of graupel for cloud ice
    ESR real(RP) 1.0_RP collection efficiency of snow for rain
    EGR real(RP) 1.0_RP collection efficiency of graupel for rain
    EGS real(RP) 1.0_RP collection efficiency of graupel for snow
    GAMMA_SACR real(RP) 25.E-3_RP effect of low temperature for Esi
    GAMMA_GACS real(RP) 90.E-3_RP effect of low temperature for Egs
    MI real(RP) 4.19E-13_RP mass of one cloud ice crystal [kg]
    BETA_SAUT real(RP) 6.E-3_RP auto-conversion factor beta for ice
    GAMMA_SAUT real(RP) 60.E-3_RP auto-conversion factor gamma for ice
    QICRT_SAUT real(RP) 0.0_RP mass ratio threshold for Psaut [kg/kg]
    BETA_GAUT real(RP) 0.0_RP auto-conversion factor beta for snow
    GAMMA_GAUT real(RP) 90.E-3_RP auto-conversion factor gamma for snow
    QSCRT_GAUT real(RP) 6.E-4_RP mass ratio threshold for Pgaut [kg/kg]
    FIXED_RE logical .false. use ice's effective radius for snow and graupel, and set rain transparent?
    CONST_REC logical .true. use constant effective radius for cloud water?
    NOFALL_QR logical .false. surpress sedimentation of rain?
    NOFALL_QI logical .false. surpress sedimentation of ice?
    NOFALL_QS logical .false. surpress sedimentation of snow?
    NOFALL_QG logical .false. surpress sedimentation of graupel?

History Output
namedescriptionunitvariable
Pcsat QC production term by satadjust kg/kg/s QC_t_sat
Pisat QI production term by satadjust kg/kg/s QI_t_sat
Pgmlt individual tendency term in tomita08 kg/kg/s Pgmlt

Function/Subroutine Documentation

◆ atmos_phy_mp_tomita08_setup()

subroutine, public scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08_setup ( 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 
)

ATMOS_PHY_MP_tomita08_setup Setup.

Definition at line 348 of file scale_atmos_phy_mp_tomita08.F90.

References scale_const::const_dice, scale_const::const_dwatr, scale_const::const_grav, scale_const::const_pi, scale_file_history::file_history_reg(), scale_io::io_fid_conf, scale_prc::prc_abort(), and scale_specfunc::sf_gamma().

Referenced by mod_atmos_phy_mp_driver::atmos_phy_mp_driver_setup().

348  use scale_prc, only: &
349  prc_abort
350  use scale_const, only: &
351  pi => const_pi, &
352  grav => const_grav, &
353  dens_w => const_dwatr, &
354  dens_i => const_dice
355  use scale_specfunc, only: &
356  sf_gamma
357  use scale_file_history, only: &
359  implicit none
360 
361  integer, intent(in) :: ka, ks, ke
362  integer, intent(in) :: ia, is, ie
363  integer, intent(in) :: ja, js, je
364 
365  real(RP) :: autoconv_nc = nc_ocn
366 
367  namelist / param_atmos_phy_mp_tomita08 / &
368  do_couple_aerosol, &
369  do_explicit_icegen, &
370  autoconv_nc, &
371  enable_kk2000, &
372  enable_rs2014, &
373  enable_wdxz2014, &
374  n0r_def, &
375  n0s_def, &
376  n0g_def, &
377  dens_s, &
378  dens_g, &
379  re_qc, &
380  re_qi, &
381  drag_g, &
382  cr, &
383  cs, &
384  eiw, &
385  erw, &
386  esw, &
387  egw, &
388  eri, &
389  esi, &
390  egi, &
391  esr, &
392  egr, &
393  egs, &
394  gamma_sacr, &
395  gamma_gacs, &
396  mi, &
397  beta_saut, &
398  gamma_saut, &
399  qicrt_saut, &
400  beta_gaut, &
401  gamma_gaut, &
402  qscrt_gaut, &
403  fixed_re, &
404  const_rec, &
405  nofall_qr, &
406  nofall_qi, &
407  nofall_qs, &
408  nofall_qg
409 
410  real(RP), parameter :: max_term_vel = 10.0_rp !-- terminal velocity for calculate dt of sedimentation
411 
412  integer :: ierr
413  integer :: i, j, ip
414  !---------------------------------------------------------------------------
415 
416 
417  log_newline
418  log_info("ATMOS_PHY_MP_tomita08_setup",*) 'Setup'
419  log_info("ATMOS_PHY_MP_tomita08_setup",*) 'Tomita (2008) 1-moment bulk 6 category'
420 
421  allocate( w3d(ka,ia,ja,w_nmax) )
422  w3d(:,:,:,:) = 0.0_rp
423 
424  allocate( nc_def(ia,ja) )
425 
426 
427  !--- read namelist
428  rewind(io_fid_conf)
429  read(io_fid_conf,nml=param_atmos_phy_mp_tomita08,iostat=ierr)
430  if( ierr < 0 ) then !--- missing
431  log_info("ATMOS_PHY_MP_tomita08_setup",*) 'Not found namelist. Default used.'
432  elseif( ierr > 0 ) then !--- fatal error
433  log_error("ATMOS_PHY_MP_tomita08_setup",*) 'Not appropriate names in namelist PARAM_ATMOS_PHY_MP_TOMITA08. Check!'
434  call prc_abort
435  endif
436  log_nml(param_atmos_phy_mp_tomita08)
437 
438  log_newline
439  log_info("ATMOS_PHY_MP_tomita08_setup",*) 'density of the snow [kg/m3] : ', dens_s
440  log_info("ATMOS_PHY_MP_tomita08_setup",*) 'density of the graupel [kg/m3] : ', dens_g
441  log_info("ATMOS_PHY_MP_tomita08_setup",*) 'Nc for auto-conversion [num/m3]: ', autoconv_nc
442  log_info("ATMOS_PHY_MP_tomita08_setup",*) 'Use k-k scheme? : ', enable_kk2000
443  log_info("ATMOS_PHY_MP_tomita08_setup",*) 'Use Roh scheme? : ', enable_rs2014
444  log_info("ATMOS_PHY_MP_tomita08_setup",*) 'Use WDXZ scheme? : ', enable_wdxz2014
445  log_newline
446  log_info("ATMOS_PHY_MP_tomita08_setup",*) 'Use effective radius of ice for snow and graupel,'
447  log_info("ATMOS_PHY_MP_tomita08_setup",*) ' and set rain transparent? : ', fixed_re
448  log_info("ATMOS_PHY_MP_tomita08_setup",*) 'Density of the ice is used for the calculation of '
449  log_info("ATMOS_PHY_MP_tomita08_setup",*) ' optically effective volume of snow and graupel.'
450  log_info("ATMOS_PHY_MP_tomita08_setup",*) 'Surpress sedimentation of rain? : ', nofall_qr
451  log_info("ATMOS_PHY_MP_tomita08_setup",*) 'Surpress sedimentation of ice? : ', nofall_qi
452  log_info("ATMOS_PHY_MP_tomita08_setup",*) 'Surpress sedimentation of snow? : ', nofall_qs
453  log_info("ATMOS_PHY_MP_tomita08_setup",*) 'Surpress sedimentation of graupel? : ', nofall_qg
454  log_info("ATMOS_PHY_MP_tomita08_setup",*) 'Enable explicit ice generation? : ', do_explicit_icegen
455  log_newline
456 
457  do j = js, je
458  do i = is, ie
459  nc_def(i,j) = autoconv_nc
460  end do
461  end do
462 
463  !--- empirical coefficients A, B, C, D
464  ar = pi * dens_w / 6.0_rp
465  as = pi * dens_s / 6.0_rp
466  ag = pi * dens_g / 6.0_rp
467 
468  br = 3.0_rp
469  bs = 3.0_rp
470  bg = 3.0_rp
471 
472  cg = sqrt( ( 4.0_rp * dens_g * grav ) / ( 3.0_rp * dens00 * drag_g ) )
473 
474  dr = 0.50_rp
475  ds = 0.25_rp
476  dg = 0.50_rp
477 
478  if ( enable_rs2014 ) then ! overwrite parameters
479  do_explicit_icegen = .true.
480 
481  n0g_def = 4.e+8_rp
482  as = 0.069_rp
483  bs = 2.0_rp
484  esi = 0.25_rp
485  egi = 0.0_rp
486  egs = 0.0_rp
487  endif
488 
489  if ( do_explicit_icegen ) then
490  only_liquid = .true.
491  sw_expice = 1.0_rp
492  else
493  only_liquid = .false.
494  sw_expice = 0.0_rp
495  endif
496 
497  gam = 1.0_rp ! =0!
498  gam_2 = 1.0_rp ! =1!
499  gam_3 = 2.0_rp ! =2!
500 
501  gam_1br = sf_gamma( 1.0_rp + br ) ! = 4!
502  gam_2br = sf_gamma( 2.0_rp + br ) ! = 5!
503  gam_3br = sf_gamma( 3.0_rp + br ) ! = 6!
504  gam_3dr = sf_gamma( 3.0_rp + dr )
505  gam_6dr = sf_gamma( 6.0_rp + dr )
506  gam_1brdr = sf_gamma( 1.0_rp + br + dr )
507  gam_5dr_h = sf_gamma( 0.5_rp * (5.0_rp+dr) )
508 
509  gam_1bs = sf_gamma( 1.0_rp + bs ) ! = 4!
510  gam_2bs = sf_gamma( 2.0_rp + bs ) ! = 5!
511  gam_3bs = sf_gamma( 3.0_rp + bs ) ! = 6!
512  gam_3ds = sf_gamma( 3.0_rp + ds )
513  gam_1bsds = sf_gamma( 1.0_rp + bs + ds )
514  gam_5ds_h = sf_gamma( 0.5_rp * (5.0_rp+ds) )
515 
516  gam_1bg = sf_gamma( 1.0_rp + bg ) ! = 4!
517  gam_3dg = sf_gamma( 3.0_rp + dg )
518  gam_1bgdg = sf_gamma( 1.0_rp + bg + dg)
519  gam_5dg_h = sf_gamma( 0.5_rp * (5.0_rp+dg) )
520 
521  ln10 = log(10.0_rp)
522 
523  ! history
524  do ip = 1, w_nmax
525  call file_history_reg( w_name(ip), 'individual tendency term in tomita08', 'kg/kg/s', & ! [IN]
526  hist_id(ip) ) ! [OUT]
527  end do
528 
529  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
real(rp) function, public sf_gamma(x)
Gamma function.
integer, public ja
of whole cells: y, local, with HALO
integer, public io_fid_conf
Config file ID.
Definition: scale_io.F90:55
subroutine, public file_history_reg(name, desc, unit, itemid, standard_name, ndims, dim_type, cell_measures, fill_halo)
Register/Append variable to history file.
real(rp), parameter, public const_dice
density of ice [kg/m3]
Definition: scale_const.F90:83
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
module PROCESS
Definition: scale_prc.F90:11
integer, public je
end point of inner domain: y, local
real(rp), public const_grav
standard acceleration of gravity [m/s2]
Definition: scale_const.F90:46
integer, public ks
start point of inner domain: z, local
module SPECFUNC
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:338
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_tomita08_adjustment()

subroutine, public scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08_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(rp), dimension (ka,ia,ja), intent(in)  CCN,
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_tomita08_adjustment calculate state after saturation process.

Definition at line 543 of file scale_atmos_phy_mp_tomita08.F90.

References scale_const::const_cl, scale_const::const_dwatr, scale_const::const_eps, scale_const::const_lhf0, scale_const::const_lhs0, scale_const::const_lhv0, scale_const::const_pi, scale_const::const_pre00, scale_const::const_rvap, scale_const::const_tem00, scale_const::const_undef, scale_atmos_hydrometeor::cp_ice, scale_atmos_hydrometeor::cp_vapor, scale_atmos_hydrometeor::cp_water, scale_atmos_hydrometeor::cv_ice, scale_atmos_hydrometeor::cv_vapor, scale_atmos_hydrometeor::cv_water, scale_atmos_hydrometeor::lhf, scale_atmos_hydrometeor::lhv, scale_prof::prof_rapend(), and scale_prof::prof_rapstart().

Referenced by mod_atmos_phy_mp_driver::atmos_phy_mp_driver_calc_tendency().

543  use scale_const, only: &
544  dwatr => const_dwatr, &
545  pi => const_pi
546  use scale_file_history, only: &
547  file_history_in
548  use scale_atmos_phy_mp_common, only: &
549  mp_saturation_adjustment => atmos_phy_mp_saturation_adjustment
550  implicit none
551  integer, intent(in) :: ka, ks, ke
552  integer, intent(in) :: ia, is, ie
553  integer, intent(in) :: ja, js, je
554 
555  real(RP), intent(in) :: dens (ka,ia,ja)
556  real(RP), intent(in) :: pres (ka,ia,ja)
557  real(RP), intent(in) :: ccn (ka,ia,ja)
558  real(DP), intent(in) :: dt
559 
560  real(RP), intent(inout) :: temp(ka,ia,ja)
561  real(RP), intent(inout) :: qtrc(ka,ia,ja,qa_mp)
562  real(RP), intent(inout) :: cptot(ka,ia,ja)
563  real(RP), intent(inout) :: cvtot(ka,ia,ja)
564 
565  real(RP), intent(out) :: rhoe_t (ka,ia,ja)
566  real(RP), intent(out) :: evaporate(ka,ia,ja) ! number of evaporated cloud [/m3/s]
567 
568  real(RP) :: rhoe_d_sat(ka,ia,ja)
569 
570  real(RP) :: qc_t_sat(ka,ia,ja)
571  real(RP) :: qi_t_sat(ka,ia,ja)
572 
573  integer :: k, i, j
574  !---------------------------------------------------------------------------
575 
576  log_progress(*) 'atmosphere / physics / microphysics / Tomita08'
577 
578  !##### MP Main #####
579  call mp_tomita08( &
580  ka, ks, ke, ia, is, ie, ja, js, je, &
581  dens(:,:,:), pres(:,:,:), ccn(:,:,:), & ! [IN]
582  dt, & ! [IN]
583  temp(:,:,:), qtrc(:,:,:,:), & ! [INOUT]
584  cptot(:,:,:), cvtot(:,:,:), & ! [INOUT]
585  rhoe_t(:,:,:) ) ! [OUT]
586 
587  ! save value before saturation adjustment
588  do j = js, je
589  do i = is, ie
590  do k = ks, ke
591  qc_t_sat(k,i,j) = qtrc(k,i,j,i_qc)
592  qi_t_sat(k,i,j) = qtrc(k,i,j,i_qi)
593  enddo
594  enddo
595  enddo
596 
597  call mp_saturation_adjustment( &
598  ka, ks, ke, ia, is, ie, ja, js, je, &
599  dens(:,:,:), & ! [IN]
600  only_liquid, & ! [IN]
601  temp(:,:,:), & ! [INOUT]
602  qtrc(:,:,:,i_qv), & ! [INOUT]
603  qtrc(:,:,:,i_qc), qtrc(:,:,:,i_qi), & ! [INOUT]
604  cptot(:,:,:), cvtot(:,:,:), & ! [INOUT]
605  rhoe_d_sat(:,:,:) ) ! [OUT]
606 
607  do j = js, je
608  do i = is, ie
609  do k = ks, ke
610  rhoe_t(k,i,j) = rhoe_t(k,i,j) + rhoe_d_sat(k,i,j) / dt
611  enddo
612  enddo
613  enddo
614  do j = js, je
615  do i = is, ie
616  do k = ks, ke
617  qc_t_sat(k,i,j) = ( qtrc(k,i,j,i_qc) - qc_t_sat(k,i,j) ) / dt
618  enddo
619  enddo
620  enddo
621  do j = js, je
622  do i = is, ie
623  do k = ks, ke
624  qi_t_sat(k,i,j) = ( qtrc(k,i,j,i_qi) - qi_t_sat(k,i,j) ) / dt
625  enddo
626  enddo
627  enddo
628 
629  call file_history_in( qc_t_sat(:,:,:), 'Pcsat', 'QC production term by satadjust', 'kg/kg/s' )
630  call file_history_in( qi_t_sat(:,:,:), 'Pisat', 'QI production term by satadjust', 'kg/kg/s' )
631 
632  do j = js, je
633  do i = is, ie
634  do k = ks, ke
635  evaporate(k,i,j) = max( -qc_t_sat(k,i,j), 0.0_rp ) & ! if negative, condensation
636  * dens(k,i,j) / (4.0_rp/3.0_rp*pi*dwatr*re_qc**3) ! mass -> number (assuming constant particle radius as re_qc)
637  enddo
638  enddo
639  enddo
640 
641  !##### END MP Main #####
642 
643  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_tomita08_terminal_velocity()

subroutine, public scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08_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), intent(in)  TEMP0,
real(rp), dimension(ka,qa_mp-1), intent(in)  RHOQ0,
real(rp), dimension(ka,qa_mp-1), intent(out)  vterm 
)

Lin-type cold rain microphysics (terminal velocity)

Definition at line 1638 of file scale_atmos_phy_mp_tomita08.F90.

References scale_const::const_tem00.

Referenced by mod_atmos_phy_mp_driver::atmos_phy_mp_driver_calc_tendency().

1638  use scale_const, only: &
1639  tem00 => const_tem00
1640  implicit none
1641  integer, intent(in) :: ka, ks, ke
1642 
1643  real(RP), intent(in) :: dens0(ka)
1644  real(RP), intent(in) :: temp0(ka)
1645  real(RP), intent(in) :: rhoq0(ka,qa_mp-1)
1646 
1647  real(RP), intent(out) :: vterm(ka,qa_mp-1)
1648 
1649  real(RP) :: dens(ka)
1650  real(RP) :: temc(ka)
1651  real(RP) :: qr(ka), qi(ka), qs(ka), qg(ka)
1652 
1653  real(RP) :: rho_fact(ka) ! density factor
1654 
1655  real(RP) :: n0r(ka), n0s(ka), n0g(ka)
1656  real(RP) :: rlmdr, rlmds, rlmdg
1657  real(RP) :: rlmdr_dr, rlmds_ds, rlmdg_dg
1658 
1659  !---< Roh and Satoh (2014) >---
1660  real(RP) :: tems, xs2
1661  real(RP) :: moms_0bs, rmoms_vt(ka)
1662  real(RP) :: coef_at(4), coef_bt(4)
1663  real(RP) :: loga_, b_, nm
1664 
1665  real(RP) :: zerosw
1666  integer :: k
1667  !---------------------------------------------------------------------------
1668 
1669  do k = ks, ke
1670  ! store to work
1671  dens(k) = dens0(k)
1672  temc(k) = temp0(k) - tem00
1673  qr(k) = rhoq0(k,i_hyd_qr) / dens(k)
1674  qi(k) = rhoq0(k,i_hyd_qi) / dens(k)
1675  qs(k) = rhoq0(k,i_hyd_qs) / dens(k)
1676  qg(k) = rhoq0(k,i_hyd_qg) / dens(k)
1677 
1678  rho_fact(k) = sqrt( dens00 / dens(k) )
1679  end do
1680 
1681  if ( enable_wdxz2014 ) then
1682  ! Wainwright et al. (2014)
1683  do k = ks, ke
1684  ! intercept parameter N0
1685  n0r(k) = 1.16e+5_rp * exp( log( max( dens(k)*qr(k)*1000.0_rp, 1.e-2_rp ) )*0.477_rp )
1686  n0s(k) = 4.58e+9_rp * exp( log( max( dens(k)*qs(k)*1000.0_rp, 1.e-2_rp ) )*0.788_rp )
1687  n0g(k) = 9.74e+8_rp * exp( log( max( dens(k)*qg(k)*1000.0_rp, 1.e-2_rp ) )*0.816_rp )
1688  end do
1689  else
1690  do k = ks, ke
1691  ! intercept parameter N0
1692  n0r(k) = n0r_def ! Marshall and Palmer (1948)
1693  n0s(k) = n0s_def ! Gunn and Marshall (1958)
1694  n0g(k) = n0g_def
1695  end do
1696  end if
1697 
1698 
1699  do k = ks, ke
1700  !---< terminal velocity >
1701  zerosw = 0.5_rp - sign(0.5_rp, qi(k) - 1.e-8_rp )
1702  vterm(k,i_hyd_qi) = -3.29_rp * exp( log( dens(k)*qi(k)+zerosw )*0.16_rp ) * ( 1.0_rp-zerosw )
1703  end do
1704 
1705 
1706  do k = ks, ke
1707  ! slope parameter lambda (Rain)
1708  zerosw = 0.5_rp - sign(0.5_rp, qr(k) - 1.e-12_rp )
1709  rlmdr = sqrt(sqrt( dens(k) * qr(k) / ( ar * n0r(k) * gam_1br ) + zerosw )) * ( 1.0_rp-zerosw )
1710  rlmdr_dr = sqrt( rlmdr ) ! **Dr
1711  vterm(k,i_hyd_qr) = -cr * rho_fact(k) * gam_1brdr / gam_1br * rlmdr_dr
1712  end do
1713 
1714 
1715  if ( enable_rs2014 ) then
1716  !---< modification by Roh and Satoh (2014) >---
1717  ! bimodal size distribution of snow
1718  do k = ks, ke
1719  zerosw = 0.5_rp - sign(0.5_rp, dens(k) * qs(k) - 1.e-12_rp )
1720  xs2 = dens(k) * qs(k) / as
1721 
1722  tems = min( -0.1_rp, temc(k) )
1723  coef_at(1) = coef_a01 + tems * ( coef_a02 + tems * ( coef_a05 + tems * coef_a09 ) )
1724  coef_at(2) = coef_a03 + tems * ( coef_a04 + tems * coef_a07 )
1725  coef_at(3) = coef_a06 + tems * coef_a08
1726  coef_at(4) = coef_a10
1727  coef_bt(1) = coef_b01 + tems * ( coef_b02 + tems * ( coef_b05 + tems * coef_b09 ) )
1728  coef_bt(2) = coef_b03 + tems * ( coef_b04 + tems * coef_b07 )
1729  coef_bt(3) = coef_b06 + tems * coef_b08
1730  coef_bt(4) = coef_b10
1731  ! 0 + Bs(=2) moment
1732  nm = 2.0_rp
1733  loga_ = coef_at(1) + nm * ( coef_at(2) + nm * ( coef_at(3) + nm * coef_at(4) ) )
1734  b_ = coef_bt(1) + nm * ( coef_bt(2) + nm * ( coef_bt(3) + nm * coef_bt(4) ) )
1735  moms_0bs = exp(ln10*loga_) * exp(log(xs2+zerosw)*b_) * ( 1.0_rp-zerosw )
1736  ! Bs(=2) + Ds(=0.25) moment
1737  nm = 2.25_rp
1738  loga_ = coef_at(1) + nm * ( coef_at(2) + nm * ( coef_at(3) + nm * coef_at(4) ) )
1739  b_ = coef_bt(1) + nm * ( coef_bt(2) + nm * ( coef_bt(3) + nm * coef_bt(4) ) )
1740  rmoms_vt(k) = exp(ln10*loga_) * exp(log(xs2+zerosw)*b_) * ( 1.0_rp-zerosw ) / ( moms_0bs + zerosw )
1741  end do
1742  else
1743  ! slope parameter lambda (Snow)
1744  do k = ks, ke
1745  zerosw = 0.5_rp - sign(0.5_rp, qs(k) - 1.e-12_rp )
1746  rlmds = sqrt(sqrt( dens(k) * qs(k) / ( as * n0s(k) * gam_1bs ) + zerosw )) * ( 1.0_rp-zerosw )
1747  rlmds_ds = sqrt( sqrt(rlmds) ) ! **Ds
1748  rmoms_vt(k) = gam_1bsds / gam_1bs * rlmds_ds
1749  end do
1750  end if
1751 
1752  do k = ks, ke
1753  vterm(k,i_hyd_qs) = -cs * rho_fact(k) * rmoms_vt(k)
1754  end do
1755 
1756 
1757  do k = ks, ke
1758  ! slope parameter lambda (Graupel)
1759  zerosw = 0.5_rp - sign(0.5_rp, qg(k) - 1.e-12_rp )
1760  rlmdg = sqrt(sqrt( dens(k) * qg(k) / ( ag * n0g(k) * gam_1bg ) + zerosw )) * ( 1.0_rp-zerosw )
1761  rlmdg_dg = sqrt( rlmdg ) ! **Dg
1762  vterm(k,i_hyd_qg) = -cg * rho_fact(k) * gam_1bgdg / gam_1bg * rlmdg_dg
1763  enddo
1764 
1765 
1766 !OCL XFILL
1767  do k = ks, ke
1768  vterm(k,i_hyd_qc) = 0.0_rp
1769  end do
1770 
1771  if ( nofall_qr ) then
1772 !OCL XFILL
1773  do k = ks, ke
1774  vterm(k,i_hyd_qr) = 0.0_rp
1775  enddo
1776  endif
1777 
1778  if ( nofall_qi ) then
1779 !OCL XFILL
1780  do k = ks, ke
1781  vterm(k,i_hyd_qi) = 0.0_rp
1782  enddo
1783  endif
1784 
1785  if ( nofall_qs ) then
1786 !OCL XFILL
1787  do k = ks, ke
1788  vterm(k,i_hyd_qs) = 0.0_rp
1789  enddo
1790  endif
1791 
1792  if ( nofall_qg ) then
1793 !OCL XFILL
1794  do k = ks, ke
1795  vterm(k,i_hyd_qg) = 0.0_rp
1796  enddo
1797  endif
1798 
1799  vterm( 1:ks-1,:) = 0.0_rp
1800  vterm(ke+1:ka ,:) = 0.0_rp
1801 
1802  return
real(rp), parameter, public const_tem00
temperature reference (0C) [K]
Definition: scale_const.F90:90
integer, public ke
end point of inner domain: z, local
integer, public ks
start point of inner domain: z, local
module CONSTANT
Definition: scale_const.F90:11
integer, public ka
of whole cells: z, local, with HALO
Here is the caller graph for this function:

◆ atmos_phy_mp_tomita08_cloud_fraction()

subroutine, public scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08_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 1812 of file scale_atmos_phy_mp_tomita08.F90.

Referenced by mod_atmos_phy_mp_vars::atmos_phy_mp_vars_get_diagnostic().

1812  implicit none
1813  integer, intent(in) :: ka, ks, ke
1814  integer, intent(in) :: ia, is, ie
1815  integer, intent(in) :: ja, js, je
1816 
1817  real(RP), intent(in) :: qtrc(ka,ia,ja,qa_mp-1)
1818  real(RP), intent(in) :: mask_criterion
1819 
1820  real(RP), intent(out) :: cldfrac(ka,ia,ja)
1821 
1822  real(RP) :: qhydro
1823  integer :: k, i, j
1824  !---------------------------------------------------------------------------
1825 
1826  !$omp parallel do OMP_SCHEDULE_ &
1827  !$omp private(qhydro)
1828  do j = js, je
1829  do i = is, ie
1830  do k = ks, ke
1831  qhydro = qtrc(k,i,j,i_hyd_qc) + qtrc(k,i,j,i_hyd_qr) &
1832  + qtrc(k,i,j,i_hyd_qi) + qtrc(k,i,j,i_hyd_qs) + qtrc(k,i,j,i_hyd_qg)
1833  cldfrac(k,i,j) = 0.5_rp + sign(0.5_rp,qhydro-mask_criterion)
1834  enddo
1835  enddo
1836  enddo
1837 
1838  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_tomita08_effective_radius()

subroutine, public scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08_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 1847 of file scale_atmos_phy_mp_tomita08.F90.

References scale_const::const_dice, scale_const::const_dwatr, scale_const::const_pi, scale_const::const_tem00, scale_atmos_hydrometeor::i_hc, scale_atmos_hydrometeor::i_hg, scale_atmos_hydrometeor::i_hi, scale_atmos_hydrometeor::i_hr, scale_atmos_hydrometeor::i_hs, and scale_atmos_hydrometeor::n_hyd.

Referenced by mod_atmos_phy_mp_vars::atmos_phy_mp_vars_get_diagnostic().

1847  use scale_const, only: &
1848  pi => const_pi, &
1849  dens_w => const_dwatr, &
1850  dens_i => const_dice, &
1851  tem00 => const_tem00
1852  use scale_atmos_hydrometeor, only: &
1853  n_hyd, &
1854  i_hc, &
1855  i_hr, &
1856  i_hi, &
1857  i_hs, &
1858  i_hg
1859  implicit none
1860  integer, intent(in) :: ka, ks, ke
1861  integer, intent(in) :: ia, is, ie
1862  integer, intent(in) :: ja, js, je
1863 
1864  real(RP), intent(in) :: dens0(ka,ia,ja)
1865  real(RP), intent(in) :: temp0(ka,ia,ja)
1866  real(RP), intent(in) :: qtrc0(ka,ia,ja,qa_mp-1)
1867 
1868  real(RP), intent(out) :: re (ka,ia,ja,n_hyd) ! effective radius [cm]
1869  real(RP) :: dens(ka)
1870  real(RP) :: temc(ka)
1871  real(RP) :: qr(ka), qs(ka), qg(ka)
1872  real(RP) :: nc(ka)
1873  real(RP) :: n0r(ka), n0s(ka), n0g(ka)
1874  real(RP) :: rlmdr, rlmds, rlmdg
1875 
1876  real(RP), parameter :: um2cm = 100.0_rp
1877 
1878  !---< Roh and Satoh (2014) >---
1879  real(RP) :: tems, xs2
1880  real(RP) :: coef_at(4), coef_bt(4)
1881  real(RP) :: loga_, b_, nm
1882 
1883  real(RP) :: zerosw
1884  integer :: k, i, j, ih
1885  !---------------------------------------------------------------------------
1886 
1887  !$omp parallel do OMP_SCHEDULE_
1888  do j = js, je
1889  do i = is, ie
1890  do k = ks, ke
1891  re(k,i,j,i_hi) = re_qi * um2cm
1892  end do
1893  end do
1894  end do
1895  !$omp parallel do OMP_SCHEDULE_
1896  do ih = i_hg+1, n_hyd
1897  do j = js, je
1898  do i = is, ie
1899  do k = ks, ke
1900  re(k,i,j,ih) = 0.0_rp
1901  end do
1902  end do
1903  end do
1904  end do
1905 
1906  if ( const_rec .or. fixed_re ) then
1907 
1908  !$omp parallel do OMP_SCHEDULE_
1909  do j = js, je
1910  do i = is, ie
1911  do k = ks, ke
1912  re(k,i,j,i_hc) = re_qc * um2cm
1913  end do
1914  end do
1915  end do
1916 
1917  else
1918 
1919  !$omp parallel do OMP_SCHEDULE_ &
1920  !$omp private(Nc)
1921  do j = js, je
1922  do i = is, ie
1923  if ( do_couple_aerosol ) then
1924  do k = ks, ke
1925  ! Nc(k) = max( CCN(k,i,j), Nc_def(i,j)*1.E+6_RP ) ! [#/m3] tentatively off the CCN effect
1926  nc(k) = nc_def(i,j) * 1.e+6_rp ! [#/m3]
1927  enddo
1928  else
1929  do k = ks, ke
1930  nc(k) = nc_def(i,j) * 1.e+6_rp ! [#/m3]
1931  enddo
1932  endif
1933 
1934  do k = ks, ke
1935  re(k,i,j,i_hc) = 1.1_rp &
1936  * ( dens0(k,i,j) * qtrc0(k,i,j,i_hyd_qc) / nc(k) / ( 4.0_rp / 3.0_rp * pi * dens_w ) )**(1.0_rp/3.0_rp)
1937  re(k,i,j,i_hc) = min( 1.e-3_rp, max( 1.e-6_rp, re(k,i,j,i_hc) ) ) * um2cm
1938  enddo
1939  enddo
1940  enddo
1941 
1942  endif
1943 
1944  if ( fixed_re ) then
1945 
1946  !$omp parallel do OMP_SCHEDULE_
1947  do j = js, je
1948  do i = is, ie
1949  do k = ks, ke
1950  re(k,i,j,i_hr) = 10000.e-6_rp * um2cm
1951  end do
1952  end do
1953  end do
1954  !$omp parallel do OMP_SCHEDULE_
1955  do j = js, je
1956  do i = is, ie
1957  do k = ks, ke
1958  re(k,i,j,i_hs) = re_qi * um2cm
1959  end do
1960  end do
1961  end do
1962  !$omp parallel do OMP_SCHEDULE_
1963  do j = js, je
1964  do i = is, ie
1965  do k = ks, ke
1966  re(k,i,j,i_hg) = re_qi * um2cm
1967  end do
1968  end do
1969  end do
1970 
1971  else
1972 
1973 #ifndef __GFORTRAN__
1974  !$omp parallel do default(none) &
1975  !$omp shared(JS,JE,IS,IE,KS,KE,DENS0,TEMP0,QTRC0,enable_WDXZ2014,N0r_def) &
1976  !$omp shared(N0s_def,N0g_def,Ar,GAM_1br,Re,As,GAM_1bs) &
1977  !$omp shared(enable_RS2014,ln10,Ag,GAM_1bg) &
1978  !$omp private(i,j,k,dens,temc,qr,qs,qg,N0r,N0s,N0g,zerosw,RLMDr,RLMDs,Xs2,nm,tems,loga_) &
1979  !$omp private(b_,RLMDg,coef_at,coef_bt) &
1980  !$omp OMP_SCHEDULE_ collapse(2)
1981 #else
1982  !$omp parallel do default(shared) &
1983  !$omp private(i,j,k,dens,temc,qr,qs,qg,N0r,N0s,N0g,zerosw,RLMDr,RLMDs,Xs2,nm,tems,loga_) &
1984  !$omp private(b_,RLMDg,coef_at,coef_bt) &
1985  !$omp OMP_SCHEDULE_ collapse(2)
1986 #endif
1987  do j = js, je
1988  do i = is, ie
1989 
1990  do k = ks, ke
1991  dens(k) = dens0(k,i,j)
1992  temc(k) = temp0(k,i,j) - tem00
1993  qr(k) = qtrc0(k,i,j,i_hyd_qr)
1994  qs(k) = qtrc0(k,i,j,i_hyd_qs)
1995  qg(k) = qtrc0(k,i,j,i_hyd_qg)
1996  end do
1997 
1998  ! intercept parameter N0
1999  if ( enable_wdxz2014 ) then
2000  ! Wainwright et al. (2014)
2001  do k = ks, ke
2002  n0r(k) = 1.16e+5_rp * exp( log( max( dens(k)*qr(k)*1000.0_rp, 1.e-2_rp ) )*0.477_rp )
2003  n0s(k) = 4.58e+9_rp * exp( log( max( dens(k)*qs(k)*1000.0_rp, 1.e-2_rp ) )*0.788_rp )
2004  n0g(k) = 9.74e+8_rp * exp( log( max( dens(k)*qg(k)*1000.0_rp, 1.e-2_rp ) )*0.816_rp )
2005  end do
2006  else
2007  do k = ks, ke
2008  n0r(k) = n0r_def ! Marshall and Palmer (1948)
2009  n0s(k) = n0s_def ! Gunn and Marshall (1958)
2010  n0g(k) = n0g_def
2011  end do
2012  end if
2013 
2014  ! slope parameter lambda
2015  do k = ks, ke
2016  zerosw = 0.5_rp - sign(0.5_rp, qr(k) - 1.e-12_rp )
2017  rlmdr = sqrt(sqrt( dens(k) * qr(k) / ( ar * n0r(k) * gam_1br ) + zerosw )) * ( 1.0_rp-zerosw )
2018  ! Effective radius is defined by r3m/r2m = 1.5/lambda
2019  re(k,i,j,i_hr) = 1.5_rp * rlmdr * um2cm
2020  end do
2021 
2022  if ( enable_rs2014 ) then
2023  do k = ks, ke
2024  zerosw = 0.5_rp - sign(0.5_rp, qs(k) - 1.e-12_rp )
2025  !---< modification by Roh and Satoh (2014) >---
2026  ! bimodal size distribution of snow
2027  zerosw = 0.5_rp - sign(0.5_rp, dens(k) * qs(k) - 1.e-12_rp )
2028  xs2 = dens(k) * qs(k) / as
2029 
2030  tems = min( -0.1_rp, temc(k) )
2031  coef_at(1) = coef_a01 + tems * ( coef_a02 + tems * ( coef_a05 + tems * coef_a09 ) )
2032  coef_at(2) = coef_a03 + tems * ( coef_a04 + tems * coef_a07 )
2033  coef_at(3) = coef_a06 + tems * coef_a08
2034  coef_at(4) = coef_a10
2035  coef_bt(1) = coef_b01 + tems * ( coef_b02 + tems * ( coef_b05 + tems * coef_b09 ) )
2036  coef_bt(2) = coef_b03 + tems * ( coef_b04 + tems * coef_b07 )
2037  coef_bt(3) = coef_b06 + tems * coef_b08
2038  coef_bt(4) = coef_b10
2039  ! 1 + Bs(=2) moment
2040  nm = 3.0_rp
2041  loga_ = coef_at(1) + nm * ( coef_at(2) + nm * ( coef_at(3) + nm * coef_at(4) ) )
2042  b_ = coef_bt(1) + nm * ( coef_bt(2) + nm * ( coef_bt(3) + nm * coef_bt(4) ) )
2043 
2044  re(k,i,j,i_hs) = 0.5_rp * exp(ln10*loga_) * exp(log(xs2+zerosw)*b_) * ( 1.0_rp-zerosw ) / ( xs2+zerosw ) * um2cm
2045  end do
2046  else
2047  do k = ks, ke
2048  zerosw = 0.5_rp - sign(0.5_rp, qs(k) - 1.e-12_rp )
2049  rlmds = sqrt(sqrt( dens(k) * qs(k) / ( as * n0s(k) * gam_1bs ) + zerosw )) * ( 1.0_rp-zerosw )
2050  re(k,i,j,i_hs) = 1.5_rp * rlmds * um2cm
2051  ! Re(k,i,j,I_HS) = dens * qs / N0s / ( 2.0_RP / 3.0_RP * PI * dens_i ) / RLMDs**3 * um2cm
2052  end do
2053  end if
2054 
2055  do k = ks, ke
2056  zerosw = 0.5_rp - sign(0.5_rp, qg(k) - 1.e-12_rp )
2057  rlmdg = sqrt(sqrt( dens(k) * qg(k) / ( ag * n0g(k) * gam_1bg ) + zerosw )) * ( 1.0_rp-zerosw )
2058  re(k,i,j,i_hg) = 1.5_rp * rlmdg * um2cm
2059  ! Re(k,i,j,I_HG) = dens * qg / N0g / ( 2.0_RP / 3.0_RP * PI * dens_i ) / RLMDg**3 * um2cm
2060  enddo
2061 
2062  enddo
2063  enddo
2064 
2065  endif
2066 
2067  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, parameter, public i_hs
snow
integer, parameter, public i_hr
liquid water rain
integer, parameter, public i_hi
ice water cloud
real(rp), parameter, public const_tem00
temperature reference (0C) [K]
Definition: scale_const.F90:90
integer, public ja
of whole cells: y, local, with HALO
real(rp), parameter, public const_dice
density of ice [kg/m3]
Definition: scale_const.F90:83
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
module CONSTANT
Definition: scale_const.F90:11
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
real(rp), public const_pi
pi
Definition: scale_const.F90:31
integer, parameter, public n_hyd
integer, parameter, public i_hg
graupel
Here is the caller graph for this function:

◆ atmos_phy_mp_tomita08_qtrc2qhyd()

subroutine, public scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08_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 2136 of file scale_atmos_phy_mp_tomita08.F90.

References scale_atmos_hydrometeor::i_hc, scale_atmos_hydrometeor::i_hg, scale_atmos_hydrometeor::i_hi, scale_atmos_hydrometeor::i_hr, scale_atmos_hydrometeor::i_hs, and scale_atmos_hydrometeor::n_hyd.

Referenced by mod_atmos_phy_mp_vars::atmos_phy_mp_vars_get_diagnostic().

2136  use scale_atmos_hydrometeor, only: &
2137  n_hyd, &
2138  i_hc, &
2139  i_hr, &
2140  i_hi, &
2141  i_hs, &
2142  i_hg
2143  implicit none
2144  integer, intent(in) :: ka, ks, ke
2145  integer, intent(in) :: ia, is, ie
2146  integer, intent(in) :: ja, js, je
2147 
2148  real(RP), intent(in) :: qtrc(ka,ia,ja,qa_mp-1)
2149 
2150  real(RP), intent(out) :: qe(ka,ia,ja,n_hyd) ! mass ratio of each cateory [kg/kg]
2151 
2152  integer :: k, i, j, ih
2153  !---------------------------------------------------------------------------
2154 
2155  !$omp parallel do OMP_SCHEDULE_
2156 !OCL XFILL
2157  do j = js, je
2158  do i = is, ie
2159  do k = ks, ke
2160  qe(k,i,j,i_hc) = qtrc(k,i,j,i_hyd_qc)
2161  end do
2162  end do
2163  end do
2164  !$omp parallel do OMP_SCHEDULE_
2165 !OCL XFILL
2166  do j = js, je
2167  do i = is, ie
2168  do k = ks, ke
2169  qe(k,i,j,i_hr) = qtrc(k,i,j,i_hyd_qr)
2170  end do
2171  end do
2172  end do
2173  !$omp parallel do OMP_SCHEDULE_
2174 !OCL XFILL
2175  do j = js, je
2176  do i = is, ie
2177  do k = ks, ke
2178  qe(k,i,j,i_hi) = qtrc(k,i,j,i_hyd_qi)
2179  end do
2180  end do
2181  end do
2182  !$omp parallel do OMP_SCHEDULE_
2183 !OCL XFILL
2184  do j = js, je
2185  do i = is, ie
2186  do k = ks, ke
2187  qe(k,i,j,i_hs) = qtrc(k,i,j,i_hyd_qs)
2188  end do
2189  end do
2190  end do
2191  !$omp parallel do OMP_SCHEDULE_
2192 !OCL XFILL
2193  do j = js, je
2194  do i = is, ie
2195  do k = ks, ke
2196  qe(k,i,j,i_hg) = qtrc(k,i,j,i_hyd_qg)
2197  end do
2198  end do
2199  end do
2200  !$omp parallel do OMP_SCHEDULE_
2201 !OCL XFILL
2202  do ih = i_hg+1, n_hyd
2203  do j = js, je
2204  do i = is, ie
2205  do k = ks, ke
2206  qe(k,i,j,ih) = 0.0_rp
2207  end do
2208  end do
2209  end do
2210  end do
2211 
2212  return
integer, public ia
of whole cells: x, local, with HALO
integer, parameter, public i_hs
snow
integer, parameter, public i_hr
liquid water rain
integer, parameter, public i_hi
ice water cloud
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
integer, parameter, public i_hg
graupel
Here is the caller graph for this function:

◆ atmos_phy_mp_tomita08_qhyd2qtrc()

subroutine, public scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08_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 
)

get mass ratio of each category

Definition at line 2221 of file scale_atmos_phy_mp_tomita08.F90.

References scale_atmos_hydrometeor::i_hc, scale_atmos_hydrometeor::i_hg, scale_atmos_hydrometeor::i_hh, scale_atmos_hydrometeor::i_hi, scale_atmos_hydrometeor::i_hr, scale_atmos_hydrometeor::i_hs, and scale_atmos_hydrometeor::n_hyd.

Referenced by mod_atmos_phy_mp_driver::atmos_phy_mp_driver_qhyd2qtrc().

2221  use scale_atmos_hydrometeor, only: &
2222  n_hyd, &
2223  i_hc, &
2224  i_hr, &
2225  i_hi, &
2226  i_hs, &
2227  i_hg, &
2228  i_hh
2229  implicit none
2230 
2231  integer, intent(in) :: ka, ks, ke
2232  integer, intent(in) :: ia, is, ie
2233  integer, intent(in) :: ja, js, je
2234  real(RP), intent(in) :: qe (ka,ia,ja,n_hyd) ! mass ratio of each cateory [kg/kg]
2235  real(RP), intent(out) :: qtrc(ka,ia,ja,qa_mp-1)
2236 
2237  integer :: k, i, j
2238  !---------------------------------------------------------------------------
2239 
2240  !$omp parallel do OMP_SCHEDULE_
2241 !OCL XFILL
2242  do j = js, je
2243  do i = is, ie
2244  do k = ks, ke
2245  qtrc(k,i,j,i_hyd_qc) = qe(k,i,j,i_hc)
2246  end do
2247  end do
2248  end do
2249  !$omp parallel do OMP_SCHEDULE_
2250 !OCL XFILL
2251  do j = js, je
2252  do i = is, ie
2253  do k = ks, ke
2254  qtrc(k,i,j,i_hyd_qr) = qe(k,i,j,i_hr)
2255  end do
2256  end do
2257  end do
2258  !$omp parallel do OMP_SCHEDULE_
2259 !OCL XFILL
2260  do j = js, je
2261  do i = is, ie
2262  do k = ks, ke
2263  qtrc(k,i,j,i_hyd_qi) = qe(k,i,j,i_hi)
2264  end do
2265  end do
2266  end do
2267  !$omp parallel do OMP_SCHEDULE_
2268 !OCL XFILL
2269  do j = js, je
2270  do i = is, ie
2271  do k = ks, ke
2272  qtrc(k,i,j,i_hyd_qs) = qe(k,i,j,i_hs)
2273  end do
2274  end do
2275  end do
2276  !$omp parallel do OMP_SCHEDULE_
2277 !OCL XFILL
2278  do j = js, je
2279  do i = is, ie
2280  do k = ks, ke
2281  qtrc(k,i,j,i_hyd_qg) = qe(k,i,j,i_hg) + qe(k,i,j,i_hh)
2282  end do
2283  end do
2284  end do
2285 
2286  return
integer, public ia
of whole cells: x, local, with HALO
integer, parameter, public i_hs
snow
integer, parameter, public i_hr
liquid water rain
integer, parameter, public i_hi
ice water cloud
integer, public ja
of whole cells: y, local, with HALO
integer, parameter, public i_hh
hail
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
integer, parameter, public i_hg
graupel
Here is the caller graph for this function:

Variable Documentation

◆ atmos_phy_mp_tomita08_ntracers

integer, parameter, public scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08_ntracers = QA_MP

Definition at line 43 of file scale_atmos_phy_mp_tomita08.F90.

Referenced by mod_atmos_phy_mp_driver::atmos_phy_mp_driver_tracer_setup().

43  integer, parameter, public :: atmos_phy_mp_tomita08_ntracers = qa_mp

◆ atmos_phy_mp_tomita08_nwaters

integer, parameter, public scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08_nwaters = 2

Definition at line 44 of file scale_atmos_phy_mp_tomita08.F90.

Referenced by mod_atmos_phy_mp_driver::atmos_phy_mp_driver_tracer_setup().

44  integer, parameter, public :: atmos_phy_mp_tomita08_nwaters = 2

◆ atmos_phy_mp_tomita08_nices

integer, parameter, public scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08_nices = 3

Definition at line 45 of file scale_atmos_phy_mp_tomita08.F90.

Referenced by mod_atmos_phy_mp_driver::atmos_phy_mp_driver_tracer_setup().

45  integer, parameter, public :: atmos_phy_mp_tomita08_nices = 3

◆ atmos_phy_mp_tomita08_tracer_names

character(len=h_short), dimension(qa_mp), parameter, public scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08_tracer_names = (/ 'QV', 'QC', 'QR', 'QI', 'QS', 'QG' /)

Definition at line 46 of file scale_atmos_phy_mp_tomita08.F90.

Referenced by mod_atmos_phy_mp_driver::atmos_phy_mp_driver_tracer_setup().

46  character(len=H_SHORT), parameter, public :: atmos_phy_mp_tomita08_tracer_names(qa_mp) = (/ &
47  'QV', &
48  'QC', &
49  'QR', &
50  'QI', &
51  'QS', &
52  'QG' /)

◆ atmos_phy_mp_tomita08_tracer_descriptions

character(len=h_mid), dimension(qa_mp), parameter, public scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08_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 ', 'Ratio of Cloud Ice mass ratio to total mass ', 'Ratio of Snow miass ratio to total mass ', 'Ratio of Graupel mass ratio to total mass '/)

Definition at line 53 of file scale_atmos_phy_mp_tomita08.F90.

Referenced by mod_atmos_phy_mp_driver::atmos_phy_mp_driver_tracer_setup().

53  character(len=H_MID) , parameter, public :: atmos_phy_mp_tomita08_tracer_descriptions(qa_mp) = (/ &
54  'Ratio of Water Vapor mass to total mass (Specific humidity)', &
55  'Ratio of Cloud Water mass to total mass ', &
56  'Ratio of Rain Water mass to total mass ', &
57  'Ratio of Cloud Ice mass ratio to total mass ', &
58  'Ratio of Snow miass ratio to total mass ', &
59  'Ratio of Graupel mass ratio to total mass '/)

◆ atmos_phy_mp_tomita08_tracer_units

character(len=h_short), dimension(qa_mp), parameter, public scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08_tracer_units = (/ 'kg/kg', 'kg/kg', 'kg/kg', 'kg/kg', 'kg/kg', 'kg/kg' /)

Definition at line 60 of file scale_atmos_phy_mp_tomita08.F90.

Referenced by mod_atmos_phy_mp_driver::atmos_phy_mp_driver_tracer_setup().

60  character(len=H_SHORT), parameter, public :: atmos_phy_mp_tomita08_tracer_units(qa_mp) = (/ &
61  'kg/kg', &
62  'kg/kg', &
63  'kg/kg', &
64  'kg/kg', &
65  'kg/kg', &
66  'kg/kg' /)