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, flg_lt)
 ATMOS_PHY_MP_tomita08_setup Setup. More...
 
subroutine, public atmos_phy_mp_tomita08_finalize
 finalize More...
 
subroutine, public atmos_phy_mp_tomita08_putvar (in_drag_g, in_re_qc, in_re_qi, in_Cr, in_Cs)
 overwrite private variables More...
 
subroutine, public atmos_phy_mp_tomita08_getvar (out_drag_g, out_re_qc, out_re_qi, out_Cr, out_Cs)
 read private variables 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, flg_lt, d0_crg, v0_crg, dqcrg, beta_crg, QSPLT_in, Sarea, QTRC_crg)
 ATMOS_PHY_MP_tomita08_adjustment calculate state after saturation process. More...
 
subroutine, public atmos_phy_mp_tomita08_terminal_velocity (KA, KS, KE, DENS, TEMP, RHOQ, 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 .false. apply CCN effect?
    DO_EXPLICIT_ICEGEN logical .false. 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)
    ENABLE_HZDFHI2007 logical .false. use scheme by Heymsfield et al. (2007)
    CLOUD_TYPE character(len=12) "synoptic" ! "synoptic" or "crystal_face"
    ENABLE_TRAWLBG2017 logical .false. use scheme by Thornberry et al. (2017)
    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?
    ECOAL_GI real(RP)
    ECOAL_GS real(RP)

History Output
namedescriptionunitvariable
Pcsat QC production term by satulation adjustment kg/kg/s Pcsat
Pisat QI production term by satulation adjustment kg/kg/s Pisat
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,
logical, intent(in), optional  flg_lt 
)

ATMOS_PHY_MP_tomita08_setup Setup.

Definition at line 399 of file scale_atmos_phy_mp_tomita08.F90.

399  use scale_prc, only: &
400  prc_abort
401  use scale_const, only: &
402  pi => const_pi, &
403  grav => const_grav, &
404  dens_w => const_dwatr, &
405  dens_i => const_dice
406  use scale_specfunc, only: &
407  sf_gamma
408  use scale_file_history, only: &
410  implicit none
411 
412  integer, intent(in) :: KA, KS, KE
413  integer, intent(in) :: IA, IS, IE
414  integer, intent(in) :: JA, JS, JE
415 
416  logical, intent(in), optional :: flg_lt
417 
418  real(RP) :: autoconv_nc = nc_ocn
419 
420  namelist / param_atmos_phy_mp_tomita08 / &
421  do_couple_aerosol, &
422  do_explicit_icegen, &
423  autoconv_nc, &
424  enable_kk2000, &
425  enable_rs2014, &
426  enable_wdxz2014, &
427  enable_hzdfhi2007, &
428  cloud_type, &
429  enable_trawlbg2017, &
430  n0r_def, &
431  n0s_def, &
432  n0g_def, &
433  dens_s, &
434  dens_g, &
435  re_qc, &
436  re_qi, &
437  drag_g, &
438  cr, &
439  cs, &
440  eiw, &
441  erw, &
442  esw, &
443  egw, &
444  eri, &
445  esi, &
446  egi, &
447  esr, &
448  egr, &
449  egs, &
450  gamma_sacr, &
451  gamma_gacs, &
452  mi, &
453  beta_saut, &
454  gamma_saut, &
455  qicrt_saut, &
456  beta_gaut, &
457  gamma_gaut, &
458  qscrt_gaut, &
459  fixed_re, &
460  const_rec, &
461  nofall_qr, &
462  nofall_qi, &
463  nofall_qs, &
464  nofall_qg, &
465  ecoal_gi, &
466  ecoal_gs
467 
468  real(RP), parameter :: max_term_vel = 10.0_rp !-- terminal velocity for calculate dt of sedimentation
469 
470  logical :: flg_lt_
471  integer :: ierr
472  integer :: i, j, ip
473  !---------------------------------------------------------------------------
474 
475 
476  log_newline
477  log_info("ATMOS_PHY_MP_tomita08_setup",*) 'Setup'
478  log_info("ATMOS_PHY_MP_tomita08_setup",*) 'Tomita (2008) 1-moment bulk 6 category'
479 
480  allocate( w3d(ka,ia,ja,w_nmax) )
481  w3d(:,:,:,:) = 0.0_rp
482 
483  allocate( nc_def(ia,ja) )
484 
485  ecoal_gi = 0.0_rp
486  ecoal_gs = 0.0_rp
487 
488  !--- read namelist
489  rewind(io_fid_conf)
490  read(io_fid_conf,nml=param_atmos_phy_mp_tomita08,iostat=ierr)
491  if( ierr < 0 ) then !--- missing
492  log_info("ATMOS_PHY_MP_tomita08_setup",*) 'Not found namelist. Default used.'
493  elseif( ierr > 0 ) then !--- fatal error
494  log_error("ATMOS_PHY_MP_tomita08_setup",*) 'Not appropriate names in namelist PARAM_ATMOS_PHY_MP_TOMITA08. Check!'
495  call prc_abort
496  endif
497  log_nml(param_atmos_phy_mp_tomita08)
498 
499  log_newline
500  log_info("ATMOS_PHY_MP_tomita08_setup",*) 'density of the snow [kg/m3] : ', dens_s
501  log_info("ATMOS_PHY_MP_tomita08_setup",*) 'density of the graupel [kg/m3] : ', dens_g
502  log_info("ATMOS_PHY_MP_tomita08_setup",*) 'Nc for auto-conversion [num/m3]: ', autoconv_nc
503  log_info("ATMOS_PHY_MP_tomita08_setup",*) 'Use k-k scheme? : ', enable_kk2000
504  log_info("ATMOS_PHY_MP_tomita08_setup",*) 'Use Roh scheme? : ', enable_rs2014
505  log_info("ATMOS_PHY_MP_tomita08_setup",*) 'Use WDXZ scheme? : ', enable_wdxz2014
506  log_info("ATMOS_PHY_MP_tomita08_setup",*) 'Use HZDFHI scheme? : ', enable_hzdfhi2007
507  log_info("ATMOS_PHY_MP_tomita08_setup",*) 'Use TRAWLBG scheme? : ', enable_trawlbg2017
508  log_newline
509  log_info("ATMOS_PHY_MP_tomita08_setup",*) 'Use effective radius of ice for snow and graupel,'
510  log_info("ATMOS_PHY_MP_tomita08_setup",*) ' and set rain transparent? : ', fixed_re
511  log_info("ATMOS_PHY_MP_tomita08_setup",*) 'Density of the ice is used for the calculation of '
512  log_info("ATMOS_PHY_MP_tomita08_setup",*) ' optically effective volume of snow and graupel.'
513  log_info("ATMOS_PHY_MP_tomita08_setup",*) 'Surpress sedimentation of rain? : ', nofall_qr
514  log_info("ATMOS_PHY_MP_tomita08_setup",*) 'Surpress sedimentation of ice? : ', nofall_qi
515  log_info("ATMOS_PHY_MP_tomita08_setup",*) 'Surpress sedimentation of snow? : ', nofall_qs
516  log_info("ATMOS_PHY_MP_tomita08_setup",*) 'Surpress sedimentation of graupel? : ', nofall_qg
517  log_info("ATMOS_PHY_MP_tomita08_setup",*) 'Enable explicit ice generation? : ', do_explicit_icegen
518  log_newline
519 
520  do j = js, je
521  do i = is, ie
522  nc_def(i,j) = autoconv_nc
523  end do
524  end do
525 
526  !--- empirical coefficients A, B, C, D
527  ar = pi * dens_w / 6.0_rp
528  as = pi * dens_s / 6.0_rp
529  ag = pi * dens_g / 6.0_rp
530 
531  br = 3.0_rp
532  bs = 3.0_rp
533  bg = 3.0_rp
534 
535  cg = sqrt( ( 4.0_rp * dens_g * grav ) / ( 3.0_rp * dens00 * drag_g ) )
536 
537  dr = 0.50_rp
538  ds = 0.25_rp
539  dg = 0.50_rp
540 
541  if ( enable_rs2014 ) then ! overwrite parameters
542  do_explicit_icegen = .true.
543 
544  n0g_def = 4.e+8_rp
545  as = 0.069_rp
546  bs = 2.0_rp
547  esi = 0.25_rp
548  egi = 0.0_rp
549  egs = 0.0_rp
550  endif
551 
552  if ( enable_hzdfhi2007 ) then
553  select case( cloud_type )
554  case ( "synoptic" )
555  coef_a0 = 120.4_rp
556  coef_a1 = 2.21_rp
557  coef_b0 = 0.0487_rp
558  coef_b1 = 4.57e-4_rp
559  case ("crystal_face")
560  coef_a0 = 135.3_rp
561  coef_a1 = 1.93_rp
562  coef_b0 = 0.0058_rp
563  coef_b1 = -0.0024539_rp
564  case default
565  log_error("ATMOS_PHY_MP_tomita08_setup",*) 'cloud_type is invalid: ', trim(cloud_type)
566  call prc_abort
567  end select
568  end if
569 
570  if ( do_explicit_icegen ) then
571  only_liquid = .true.
572  sw_expice = 1.0_rp
573  else
574  only_liquid = .false.
575  sw_expice = 0.0_rp
576  endif
577 
578  gam = 1.0_rp ! =0!
579  gam_2 = 1.0_rp ! =1!
580  gam_3 = 2.0_rp ! =2!
581 
582  gam_1br = sf_gamma( 1.0_rp + br ) ! = 4!
583  gam_2br = sf_gamma( 2.0_rp + br ) ! = 5!
584  gam_3br = sf_gamma( 3.0_rp + br ) ! = 6!
585  gam_3dr = sf_gamma( 3.0_rp + dr )
586  gam_6dr = sf_gamma( 6.0_rp + dr )
587  gam_1brdr = sf_gamma( 1.0_rp + br + dr )
588  gam_5dr_h = sf_gamma( 0.5_rp * (5.0_rp+dr) )
589 
590  gam_1bs = sf_gamma( 1.0_rp + bs ) ! = 4!
591  gam_2bs = sf_gamma( 2.0_rp + bs ) ! = 5!
592  gam_3bs = sf_gamma( 3.0_rp + bs ) ! = 6!
593  gam_3ds = sf_gamma( 3.0_rp + ds )
594  gam_1bsds = sf_gamma( 1.0_rp + bs + ds )
595  gam_5ds_h = sf_gamma( 0.5_rp * (5.0_rp+ds) )
596 
597  gam_1bg = sf_gamma( 1.0_rp + bg ) ! = 4!
598  gam_3dg = sf_gamma( 3.0_rp + dg )
599  gam_1bgdg = sf_gamma( 1.0_rp + bg + dg)
600  gam_5dg_h = sf_gamma( 0.5_rp * (5.0_rp+dg) )
601 
602  ! history
603  do ip = 1, w_nmax
604  call file_history_reg( w_name(ip), 'individual tendency term in tomita08', 'kg/kg/s', & ! [IN]
605  hist_id(ip) ) ! [OUT]
606  end do
607  call file_history_reg( 'Pcsat', 'QC production term by satulation adjustment', 'kg/kg/s', hist_pcsat )
608  call file_history_reg( 'Pisat', 'QI production term by satulation adjustment', 'kg/kg/s', hist_pisat )
609 
610 
611  if( present(flg_lt) ) then
612  flg_lt_ = flg_lt
613  else
614  flg_lt_ = .false.
615  end if
616  if ( flg_lt_ ) then
617  if( enable_rs2014 ) then
618  write(*,*) 'xxx ROH and Satoh (2014) scheme is not supported for Lightning'
619  call prc_abort
620  endif
621 
622  if( ecoal_gi == 0.0_rp ) then
623  flg_ecoali = 0.0_rp
624  ecoal_gi = egi
625  else
626  flg_ecoali = 1.0_rp
627  endif
628  if( ecoal_gs == 0.0_rp ) then
629  flg_ecoals = 0.0_rp
630  ecoal_gs = egs
631  else
632  flg_ecoals = 1.0_rp
633  endif
634 
635  else
636  ecoal_gi = 1.0_rp
637  ecoal_gs = 1.0_rp
638  endif
639 
640  !$acc update device(Cr, Cs)
641  !$acc update device(Ar, As, Ag, Cg)
642  !$acc update device(nofall_qr, nofall_qi, nofall_qs, nofall_qg)
643  !$acc update device(N0r_def, N0s_def, N0g_def)
644  !$acc update device(Nc_def)
645  !$acc update device(GAM_1brdr, GAM_1bs, GAM_1bg, GAM_1br, GAM_1bsds, GAM_1bgdg)
646  !$acc update device(enable_RS2014)
647  !$acc update device(enable_WDXZ2014)
648  !$acc update device(enable_HZDFHI2007, coef_a0, coef_a1, coef_b0, coef_b1)
649 
650  return

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().

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

◆ atmos_phy_mp_tomita08_finalize()

subroutine, public scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08_finalize

finalize

Definition at line 656 of file scale_atmos_phy_mp_tomita08.F90.

656 
657  deallocate( w3d )
658  deallocate( nc_def )
659 
660  return

Referenced by mod_atmos_phy_mp_driver::atmos_phy_mp_driver_finalize().

Here is the caller graph for this function:

◆ atmos_phy_mp_tomita08_putvar()

subroutine, public scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08_putvar ( real(rp), intent(in), optional  in_drag_g,
real(rp), intent(in), optional  in_re_qc,
real(rp), intent(in), optional  in_re_qi,
real(rp), intent(in), optional  in_Cr,
real(rp), intent(in), optional  in_Cs 
)

overwrite private variables

Definition at line 671 of file scale_atmos_phy_mp_tomita08.F90.

671  use scale_const, only: &
672  grav => const_grav
673  implicit none
674 
675  real(RP), intent(in), optional :: in_drag_g
676  real(RP), intent(in), optional :: in_re_qc
677  real(RP), intent(in), optional :: in_re_qi
678  real(RP), intent(in), optional :: in_Cr
679  real(RP), intent(in), optional :: in_Cs
680  !---------------------------------------------------------------------------
681 
682  if( present(in_drag_g) ) drag_g = in_drag_g
683  if( present(in_re_qc ) ) re_qc = in_re_qc
684  if( present(in_re_qi ) ) re_qi = in_re_qi
685  if( present(in_cr ) ) cr = in_cr
686  if( present(in_cs ) ) cs = in_cs
687 
688  if( present(in_drag_g) ) then
689  cg = sqrt( ( 4.0_rp * dens_g * grav ) / ( 3.0_rp * dens00 * drag_g ) ) ! recalc. Cg
690  end if
691 
692  return

References scale_const::const_grav.

Referenced by mod_da_param_estimation::da_param_estimation_update().

Here is the caller graph for this function:

◆ atmos_phy_mp_tomita08_getvar()

subroutine, public scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08_getvar ( real(rp), intent(out), optional  out_drag_g,
real(rp), intent(out), optional  out_re_qc,
real(rp), intent(out), optional  out_re_qi,
real(rp), intent(out), optional  out_Cr,
real(rp), intent(out), optional  out_Cs 
)

read private variables

Definition at line 703 of file scale_atmos_phy_mp_tomita08.F90.

703  implicit none
704 
705  real(RP), intent(out), optional :: out_drag_g
706  real(RP), intent(out), optional :: out_re_qc
707  real(RP), intent(out), optional :: out_re_qi
708  real(RP), intent(out), optional :: out_Cr
709  real(RP), intent(out), optional :: out_Cs
710  !---------------------------------------------------------------------------
711 
712  if( present(out_drag_g) ) out_drag_g = drag_g
713  if( present(out_re_qc ) ) out_re_qc = re_qc
714  if( present(out_re_qi ) ) out_re_qi = re_qi
715  if( present(out_cr ) ) out_cr = cr
716  if( present(out_cs ) ) out_cs = cs
717 
718  return

Referenced by mod_da_param_estimation::da_param_estimation_update().

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,
logical, intent(in), optional  flg_lt,
real(rp), intent(in), optional  d0_crg,
real(rp), intent(in), optional  v0_crg,
real(rp), dimension(ka,ia,ja), intent(in), optional  dqcrg,
real(rp), dimension(ka,ia,ja), intent(in), optional  beta_crg,
real(rp), dimension(ka,ia,ja,3), intent(out), optional  QSPLT_in,
real(rp), dimension(ka,ia,ja,qa_mp-1), intent(out), optional  Sarea,
real(rp), dimension(ka,ia,ja,qa_mp-1), intent(inout), optional  QTRC_crg 
)

ATMOS_PHY_MP_tomita08_adjustment calculate state after saturation process.

Definition at line 735 of file scale_atmos_phy_mp_tomita08.F90.

735  use scale_const, only: &
736  dwatr => const_dwatr, &
737  pi => const_pi
738  use scale_file_history, only: &
739  file_history_query, &
740  file_history_put
741  use scale_atmos_phy_mp_common, only: &
742  mp_saturation_adjustment => atmos_phy_mp_saturation_adjustment
743  implicit none
744  integer, intent(in) :: KA, KS, KE
745  integer, intent(in) :: IA, IS, IE
746  integer, intent(in) :: JA, JS, JE
747 
748  real(RP), intent(in) :: DENS (KA,IA,JA)
749  real(RP), intent(in) :: PRES (KA,IA,JA)
750  real(RP), intent(in) :: CCN (KA,IA,JA)
751  real(DP), intent(in) :: dt
752 
753  real(RP), intent(inout) :: TEMP(KA,IA,JA)
754  real(RP), intent(inout) :: QTRC(KA,IA,JA,QA_MP)
755  real(RP), intent(inout) :: CPtot(KA,IA,JA)
756  real(RP), intent(inout) :: CVtot(KA,IA,JA)
757 
758  real(RP), intent(out) :: RHOE_t (KA,IA,JA)
759  real(RP), intent(out) :: EVAPORATE(KA,IA,JA) ! number of evaporated cloud [/m3/s]
760 
761  ! Optional for Lightning
762  logical, intent(in), optional :: flg_lt
763  real(RP), intent(in), optional :: d0_crg, v0_crg
764  real(RP), intent(in), optional :: dqcrg(KA,IA,JA)
765  real(RP), intent(in), optional :: beta_crg(KA,IA,JA)
766  real(RP), intent(out), optional :: Sarea(KA,IA,JA,QA_MP-1) ! Surface area of each hydrometeor
767  real(RP), intent(out), optional :: QSPLT_in(KA,IA,JA,3) ! Charge separation
768  real(RP), intent(inout), optional :: QTRC_crg(KA,IA,JA,QA_MP-1) ! Tracer for charge density
769 
770  real(RP) :: RHOE_d_sat(KA,IA,JA)
771 
772  real(RP) :: QC_t_sat(KA,IA,JA)
773  real(RP) :: QI_t_sat(KA,IA,JA)
774  logical :: hist_flag
775 
776  integer :: k, i, j
777  !---------------------------------------------------------------------------
778 
779  log_progress(*) 'atmosphere / physics / microphysics / Tomita08'
780 
781  !$acc data copy(TEMP, QTRC, CPtot, CVtot), &
782  !$acc copyin(KA, KS, KE, IA, IS, IE, JA, JS, JE, &
783  !$acc DENS, PRES, CCN, dt ) &
784  !$acc copyout(RHOE_t, EVAPORATE), &
785  !$acc create(RHOE_d_sat, QC_t_sat, QI_t_sat)
786 
787 
788  !##### MP Main #####
789  call mp_tomita08( &
790  ka, ks, ke, ia, is, ie, ja, js, je, & ! [IN]
791  dens(:,:,:), pres(:,:,:), ccn(:,:,:), & ! [IN]
792  dt, & ! [IN]
793  temp(:,:,:), qtrc(:,:,:,:), & ! [INOUT]
794  cptot(:,:,:), cvtot(:,:,:), & ! [INOUT]
795  rhoe_t(:,:,:), & ! [OUT]
796  flg_lt, d0_crg, v0_crg, & ! [IN:Optional]
797  dqcrg(:,:,:), beta_crg(:,:,:), & ! [IN:Optional]
798  qsplt_in(:,:,:,:), & ! [OUT:Optional]
799  sarea(:,:,:,:), & ! [OUT:Optional]
800  qtrc_crg(:,:,:,:) ) ! [INOUT:Optional]
801 
802  ! save value before saturation adjustment
803  !$omp parallel do collapse(2)
804  !$acc kernels
805  do j = js, je
806  do i = is, ie
807  do k = ks, ke
808  qc_t_sat(k,i,j) = qtrc(k,i,j,i_qc)
809  qi_t_sat(k,i,j) = qtrc(k,i,j,i_qi)
810  enddo
811  enddo
812  enddo
813  !$acc end kernels
814 
815  call mp_saturation_adjustment( &
816  ka, ks, ke, ia, is, ie, ja, js, je, &
817  dens(:,:,:), & ! [IN]
818  only_liquid, & ! [IN]
819  temp(:,:,:), & ! [INOUT]
820  qtrc(:,:,:,i_qv), & ! [INOUT]
821  qtrc(:,:,:,i_qc), qtrc(:,:,:,i_qi), & ! [INOUT]
822  cptot(:,:,:), cvtot(:,:,:), & ! [INOUT]
823  rhoe_d_sat(:,:,:) ) ! [OUT]
824 
825  !$omp parallel do collapse(2)
826  !$acc kernels
827  do j = js, je
828  do i = is, ie
829  do k = ks, ke
830  rhoe_t(k,i,j) = rhoe_t(k,i,j) + rhoe_d_sat(k,i,j) / dt
831  enddo
832  enddo
833  enddo
834  !$acc end kernels
835 
836  !$omp parallel do collapse(2)
837  !$acc kernels
838  do j = js, je
839  do i = is, ie
840  do k = ks, ke
841  qc_t_sat(k,i,j) = ( qtrc(k,i,j,i_qc) - qc_t_sat(k,i,j) ) / dt
842  evaporate(k,i,j) = max( -qc_t_sat(k,i,j), 0.0_rp ) & ! if negative, condensation
843  * dens(k,i,j) / (4.0_rp/3.0_rp*pi*dwatr*re_qc**3) ! mass -> number (assuming constant particle radius as re_qc)
844 
845  enddo
846  enddo
847  enddo
848  !$acc end kernels
849 
850  call file_history_query( hist_pcsat, hist_flag )
851  if ( hist_flag ) then
852  call file_history_put( hist_pcsat, qc_t_sat(:,:,:) )
853  end if
854 
855 
856  call file_history_query( hist_pisat, hist_flag )
857  if ( hist_flag ) then
858  !$omp parallel do collapse(2)
859  !$acc kernels
860  do j = js, je
861  do i = is, ie
862  do k = ks, ke
863  qi_t_sat(k,i,j) = ( qtrc(k,i,j,i_qi) - qi_t_sat(k,i,j) ) / dt
864  enddo
865  enddo
866  enddo
867  !$acc end kernels
868  call file_history_put( hist_pisat, qi_t_sat(:,:,:) )
869  end if
870 
871 
872  !$acc end data
873 
874  !##### END MP Main #####
875 
876  return

References atmos_phy_mp_tomita08_effective_radius(), scale_const::const_cl, scale_const::const_dice, 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::i_hs, scale_atmos_hydrometeor::lhf, scale_atmos_hydrometeor::lhv, scale_atmos_hydrometeor::n_hyd, scale_prof::prof_rapend(), and scale_prof::prof_rapstart().

Referenced by mod_atmos_phy_mp_driver::atmos_phy_mp_driver_calc_tendency().

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), value  KA,
integer, intent(in), value  KS,
integer, intent(in), value  KE,
real(rp), dimension(ka), intent(in)  DENS,
real(rp), dimension(ka), intent(in)  TEMP,
real(rp), dimension(ka,qa_mp-1), intent(in)  RHOQ,
real(rp), dimension(ka,qa_mp-1), intent(out)  vterm 
)

Lin-type cold rain microphysics (terminal velocity)

Definition at line 2257 of file scale_atmos_phy_mp_tomita08.F90.

2257  !$acc routine vector
2258  use scale_const, only: &
2259  tem00 => const_tem00
2260  implicit none
2261  integer, intent(in), value :: KA, KS, KE
2262 
2263  real(RP), intent(in) :: DENS(KA)
2264  real(RP), intent(in) :: TEMP(KA)
2265  real(RP), intent(in) :: RHOQ(KA,QA_MP-1)
2266 
2267  real(RP), intent(out) :: vterm(KA,QA_MP-1)
2268 
2269 #ifdef _OPENACC
2270  real(RP) :: temc_tv
2271  real(RP) :: qr_tv
2272  real(RP) :: qi_tv
2273  real(RP) :: qs_tv
2274  real(RP) :: qg_tv
2275  real(RP) :: rho_fact_tv
2276 #define temc_tv(k) temc_tv
2277 #define qr_tv(k) qr_tv
2278 #define qi_tv(k) qi_tv
2279 #define qs_tv(k) qs_tv
2280 #define qg_tv(k) qg_tv
2281 #define rho_fact_tv(k) rho_fact_tv
2282 #else
2283  real(RP) :: temc_tv(KA)
2284  real(RP) :: qr_tv(KA)
2285  real(RP) :: qi_tv(KA)
2286  real(RP) :: qs_tv(KA)
2287  real(RP) :: qg_tv(KA)
2288  real(RP) :: rho_fact_tv(KA) ! density factor
2289 #endif
2290 
2291 #ifdef _OPENACC
2292  real(RP) :: N0r_tv
2293  real(RP) :: N0s_tv
2294  real(RP) :: N0g_tv
2295 #define N0r_tv(k) N0r_tv
2296 #define N0s_tv(k) N0s_tv
2297 #define N0g_tv(k) N0g_tv
2298 #else
2299  real(RP) :: N0r_tv(KA)
2300  real(RP) :: N0s_tv(KA)
2301  real(RP) :: N0g_tv(KA)
2302 #endif
2303  real(RP) :: RLMDr, RLMDs, RLMDg
2304  real(RP) :: RLMDr_dr, RLMDs_ds, RLMDg_dg
2305 
2306  !---< Roh and Satoh (2014) >---
2307  real(RP) :: tems, Xs2
2308  real(RP) :: MOMs_0bs
2309 #ifdef _OPENACC
2310  real(RP) :: RMOMs_Vt_tv
2311 #define RMOMs_Vt_tv(k) RMOMs_Vt_tv
2312 #else
2313  real(RP) :: RMOMs_Vt_tv(KA)
2314 #endif
2315  real(RP) :: coef_at(4), coef_bt(4)
2316  real(RP) :: loga_, b_, nm
2317 
2318  real(RP) :: zerosw
2319  integer :: k
2320  !---------------------------------------------------------------------------
2321 
2322  !$acc loop private(coef_at, coef_bt)
2323  do k = ks, ke
2324  ! store to work
2325  temc_tv(k) = temp(k) - tem00
2326  qr_tv(k) = rhoq(k,i_hyd_qr) / dens(k)
2327  qi_tv(k) = rhoq(k,i_hyd_qi) / dens(k)
2328  qs_tv(k) = rhoq(k,i_hyd_qs) / dens(k)
2329  qg_tv(k) = rhoq(k,i_hyd_qg) / dens(k)
2330 
2331  rho_fact_tv(k) = sqrt( dens00 / dens(k) )
2332 #ifndef _OPENACC
2333  end do
2334 #endif
2335 
2336  if ( enable_wdxz2014 ) then
2337  ! Wainwright et al. (2014)
2338  ! intercept parameter N0
2339 #ifndef _OPENACC
2340  do k = ks, ke
2341 #endif
2342  n0r_tv(k) = 1.16e+5_rp * exp( log( max( dens(k)*qr_tv(k)*1000.0_rp, 1.e-2_rp ) )*0.477_rp )
2343  n0s_tv(k) = 4.58e+9_rp * exp( log( max( dens(k)*qs_tv(k)*1000.0_rp, 1.e-2_rp ) )*0.788_rp )
2344  n0g_tv(k) = 9.74e+8_rp * exp( log( max( dens(k)*qg_tv(k)*1000.0_rp, 1.e-2_rp ) )*0.816_rp )
2345 #ifndef _OPENACC
2346  end do
2347 #endif
2348  else
2349  ! intercept parameter N0
2350 #ifndef _OPENACC
2351  do k = ks, ke
2352 #endif
2353  n0r_tv(k) = n0r_def ! Marshall and Palmer (1948)
2354  n0s_tv(k) = n0s_def ! Gunn and Marshall (1958)
2355  n0g_tv(k) = n0g_def
2356 #ifndef _OPENACC
2357  end do
2358 #endif
2359  end if
2360 
2361  !---< terminal velocity >
2362  if ( enable_hzdfhi2007 ) then
2363 #ifndef _OPENACC
2364  do k = ks, ke
2365 #endif
2366  zerosw = 0.5_rp - sign(0.5_rp, qi_tv(k) - 1.e-8_rp )
2367  vterm(k,i_hyd_qi) = min( 0.0_rp, &
2368  - ( coef_a0 + coef_a1 * temc_tv(k) ) &
2369  * exp( log( dens(k)*qi_tv(k)*1000.0_rp+zerosw ) * ( coef_b0 + coef_b1 * temc_tv(k) ) ) &
2370  * 1e-2_rp * ( 1.0_rp-zerosw ) )
2371 #ifndef _OPENACC
2372  end do
2373 #endif
2374  else
2375 #ifndef _OPENACC
2376  do k = ks, ke
2377 #endif
2378  zerosw = 0.5_rp - sign(0.5_rp, qi_tv(k) - 1.e-8_rp )
2379  vterm(k,i_hyd_qi) = -3.29_rp * exp( log( dens(k)*qi_tv(k)+zerosw )*0.16_rp ) * ( 1.0_rp-zerosw )
2380 #ifndef _OPENACC
2381  end do
2382 #endif
2383  end if
2384 
2385 
2386  ! slope parameter lambda (Rain)
2387 #ifndef _OPENACC
2388  do k = ks, ke
2389 #endif
2390  zerosw = 0.5_rp - sign(0.5_rp, qr_tv(k) - 1.e-12_rp )
2391  rlmdr = sqrt(sqrt( dens(k) * qr_tv(k) / ( ar * n0r_tv(k) * gam_1br ) + zerosw )) * ( 1.0_rp-zerosw )
2392  rlmdr_dr = sqrt( rlmdr ) ! **Dr
2393  vterm(k,i_hyd_qr) = -cr * rho_fact_tv(k) * gam_1brdr / gam_1br * rlmdr_dr
2394 #ifndef _OPENACC
2395  end do
2396 #endif
2397 
2398 
2399  if ( enable_rs2014 ) then
2400  !---< modification by Roh and Satoh (2014) >---
2401  ! bimodal size distribution of snow
2402 #ifndef _OPENACC
2403  do k = ks, ke
2404 #endif
2405  zerosw = 0.5_rp - sign(0.5_rp, dens(k) * qs_tv(k) - 1.e-12_rp )
2406  xs2 = dens(k) * qs_tv(k) / as
2407 
2408  tems = min( -0.1_rp, temc_tv(k) )
2409  coef_at(1) = coef_a01 + tems * ( coef_a02 + tems * ( coef_a05 + tems * coef_a09 ) )
2410  coef_at(2) = coef_a03 + tems * ( coef_a04 + tems * coef_a07 )
2411  coef_at(3) = coef_a06 + tems * coef_a08
2412  coef_at(4) = coef_a10
2413  coef_bt(1) = coef_b01 + tems * ( coef_b02 + tems * ( coef_b05 + tems * coef_b09 ) )
2414  coef_bt(2) = coef_b03 + tems * ( coef_b04 + tems * coef_b07 )
2415  coef_bt(3) = coef_b06 + tems * coef_b08
2416  coef_bt(4) = coef_b10
2417  ! 0 + Bs(=2) moment
2418  nm = 2.0_rp
2419  loga_ = coef_at(1) + nm * ( coef_at(2) + nm * ( coef_at(3) + nm * coef_at(4) ) )
2420  b_ = coef_bt(1) + nm * ( coef_bt(2) + nm * ( coef_bt(3) + nm * coef_bt(4) ) )
2421  moms_0bs = exp( ln10 * loga_ + log(xs2+zerosw) * b_ ) * ( 1.0_rp-zerosw )
2422  ! Bs(=2) + Ds(=0.25) moment
2423  nm = 2.25_rp
2424  loga_ = coef_at(1) + nm * ( coef_at(2) + nm * ( coef_at(3) + nm * coef_at(4) ) )
2425  b_ = coef_bt(1) + nm * ( coef_bt(2) + nm * ( coef_bt(3) + nm * coef_bt(4) ) )
2426  rmoms_vt_tv(k) = exp( ln10 * loga_ + log(xs2+zerosw) * b_ ) * ( 1.0_rp-zerosw ) / ( moms_0bs + zerosw )
2427 #ifndef _OPENACC
2428  end do
2429 #endif
2430  else
2431  ! slope parameter lambda (Snow)
2432 #ifndef _OPENACC
2433  do k = ks, ke
2434 #endif
2435  zerosw = 0.5_rp - sign(0.5_rp, qs_tv(k) - 1.e-12_rp )
2436  rlmds = sqrt(sqrt( dens(k) * qs_tv(k) / ( as * n0s_tv(k) * gam_1bs ) + zerosw )) * ( 1.0_rp-zerosw )
2437  rlmds_ds = sqrt( sqrt(rlmds) ) ! **Ds
2438  rmoms_vt_tv(k) = gam_1bsds / gam_1bs * rlmds_ds
2439 #ifndef _OPENACC
2440  end do
2441 #endif
2442  end if
2443 
2444 #ifndef _OPENACC
2445  do k = ks, ke
2446 #endif
2447  vterm(k,i_hyd_qs) = -cs * rho_fact_tv(k) * rmoms_vt_tv(k)
2448 #ifndef _OPENACC
2449  end do
2450 
2451  do k = ks, ke
2452 #endif
2453  ! slope parameter lambda (Graupel)
2454  zerosw = 0.5_rp - sign(0.5_rp, qg_tv(k) - 1.e-12_rp )
2455  rlmdg = sqrt(sqrt( dens(k) * qg_tv(k) / ( ag * n0g_tv(k) * gam_1bg ) + zerosw )) * ( 1.0_rp-zerosw )
2456  rlmdg_dg = sqrt( rlmdg ) ! **Dg
2457  vterm(k,i_hyd_qg) = -cg * rho_fact_tv(k) * gam_1bgdg / gam_1bg * rlmdg_dg
2458 #ifndef _OPENACC
2459  enddo
2460 
2461 
2462 !OCL XFILL
2463  do k = ks, ke
2464 #endif
2465  vterm(k,i_hyd_qc) = 0.0_rp
2466 #ifndef _OPENACC
2467  end do
2468 #endif
2469 
2470  if ( nofall_qr ) then
2471 !OCL XFILL
2472 #ifndef _OPENACC
2473  do k = ks, ke
2474 #endif
2475  vterm(k,i_hyd_qr) = 0.0_rp
2476 #ifndef _OPENACC
2477  enddo
2478 #endif
2479  endif
2480 
2481  if ( nofall_qi ) then
2482 !OCL XFILL
2483 #ifndef _OPENACC
2484  do k = ks, ke
2485 #endif
2486  vterm(k,i_hyd_qi) = 0.0_rp
2487 #ifndef _OPENACC
2488  enddo
2489 #endif
2490  endif
2491 
2492  if ( nofall_qs ) then
2493 !OCL XFILL
2494 #ifndef _OPENACC
2495  do k = ks, ke
2496 #endif
2497  vterm(k,i_hyd_qs) = 0.0_rp
2498 #ifndef _OPENACC
2499  enddo
2500 #endif
2501  endif
2502 
2503  if ( nofall_qg ) then
2504 !OCL XFILL
2505 #ifndef _OPENACC
2506  do k = ks, ke
2507 #endif
2508  vterm(k,i_hyd_qg) = 0.0_rp
2509 #ifndef _OPENACC
2510  enddo
2511 #endif
2512  endif
2513 
2514 #ifdef _OPENACC
2515  end do
2516 #endif
2517 
2518  return

References scale_const::const_tem00.

Referenced by mod_atmos_phy_mp_driver::atmos_phy_mp_driver_calc_tendency().

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 2528 of file scale_atmos_phy_mp_tomita08.F90.

2528  implicit none
2529  integer, intent(in) :: KA, KS, KE
2530  integer, intent(in) :: IA, IS, IE
2531  integer, intent(in) :: JA, JS, JE
2532 
2533  real(RP), intent(in) :: QTRC(KA,IA,JA,QA_MP-1)
2534  real(RP), intent(in) :: mask_criterion
2535 
2536  real(RP), intent(out) :: cldfrac(KA,IA,JA)
2537 
2538  real(RP) :: qhydro
2539  integer :: k, i, j
2540  !---------------------------------------------------------------------------
2541 
2542  !$omp parallel do OMP_SCHEDULE_ &
2543  !$omp private(qhydro)
2544  !$acc kernels copyin(QTRC) copyout(cldfrac)
2545  do j = js, je
2546  do i = is, ie
2547  do k = ks, ke
2548  qhydro = qtrc(k,i,j,i_hyd_qc) + qtrc(k,i,j,i_hyd_qr) &
2549  + qtrc(k,i,j,i_hyd_qi) + qtrc(k,i,j,i_hyd_qs) + qtrc(k,i,j,i_hyd_qg)
2550  cldfrac(k,i,j) = 0.5_rp + sign(0.5_rp,qhydro-mask_criterion)
2551  enddo
2552  enddo
2553  enddo
2554  !$acc end kernels
2555 
2556  return

Referenced by mod_atmos_phy_mp_vars::atmos_phy_mp_vars_get_diagnostic().

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 2565 of file scale_atmos_phy_mp_tomita08.F90.

2565  use scale_const, only: &
2566  pi => const_pi, &
2567  dens_w => const_dwatr, &
2568  dens_i => const_dice, &
2569  tem00 => const_tem00
2570  use scale_atmos_hydrometeor, only: &
2571  n_hyd, &
2572  i_hc, &
2573  i_hr, &
2574  i_hi, &
2575  i_hs, &
2576  i_hg
2577  implicit none
2578  integer, intent(in) :: KA, KS, KE
2579  integer, intent(in) :: IA, IS, IE
2580  integer, intent(in) :: JA, JS, JE
2581 
2582  real(RP), intent(in) :: DENS0(KA,IA,JA)
2583  real(RP), intent(in) :: TEMP0(KA,IA,JA)
2584  real(RP), intent(in) :: QTRC0(KA,IA,JA,QA_MP-1)
2585 
2586  real(RP), intent(out) :: Re (KA,IA,JA,N_HYD) ! effective radius [cm]
2587  real(RP) :: dens(KA)
2588  real(RP) :: temc(KA)
2589  real(RP) :: qr(KA), qs(KA), qg(KA)
2590  real(RP) :: Nc(KA)
2591  real(RP) :: N0r(KA), N0s(KA), N0g(KA)
2592  real(RP) :: RLMDr, RLMDs, RLMDg
2593 
2594  real(RP), parameter :: m2cm = 100.0_rp
2595 
2596  !---< Roh and Satoh (2014) >---
2597  real(RP) :: tems, Xs2
2598  real(RP) :: coef_at(4), coef_bt(4)
2599  real(RP) :: loga_, b_, nm
2600 
2601  real(RP) :: zerosw, sw
2602  integer :: k, i, j, ih
2603  !---------------------------------------------------------------------------
2604 
2605  !$acc data copyin(DENS0,TEMP0,QTRC0) copyout(Re)
2606 
2607  if ( enable_trawlbg2017 ) then
2608  !$omp parallel do OMP_SCHEDULE_ &
2609  !$omp private(sw)
2610  !$acc kernels
2611  do j = js, je
2612  do i = is, ie
2613  do k = ks, ke
2614  sw = 0.5_rp + sign(0.5_rp, temp0(k,i,j) - 192.0_rp)
2615  re(k,i,j,i_hi) = ( ( 40.0_rp + 0.53_rp * ( temp0(k,i,j) - 192.0_rp ) ) * sw &
2616  + ( 12.0_rp + 28.0_rp * exp( 0.65_rp * ( temp0(k,i,j) - 192.0_rp ) ) ) * ( 1.0_rp-sw ) &
2617  ) * 1.0e-6_rp * m2cm
2618  end do
2619  end do
2620  end do
2621  !$acc end kernels
2622  else
2623  !$omp parallel do OMP_SCHEDULE_
2624  !$acc kernels
2625  do j = js, je
2626  do i = is, ie
2627  do k = ks, ke
2628  re(k,i,j,i_hi) = re_qi * m2cm
2629  end do
2630  end do
2631  end do
2632  !$acc end kernels
2633  end if
2634  !$omp parallel do OMP_SCHEDULE_
2635  !$acc kernels
2636  do ih = i_hg+1, n_hyd
2637  do j = js, je
2638  do i = is, ie
2639  do k = ks, ke
2640  re(k,i,j,ih) = 0.0_rp
2641  end do
2642  end do
2643  end do
2644  end do
2645  !$acc end kernels
2646 
2647  if ( const_rec .or. fixed_re ) then
2648 
2649  !$omp parallel do OMP_SCHEDULE_
2650  !$acc kernels
2651  do j = js, je
2652  do i = is, ie
2653  do k = ks, ke
2654  re(k,i,j,i_hc) = re_qc * m2cm
2655  end do
2656  end do
2657  end do
2658  !$acc end kernels
2659 
2660  else
2661 
2662  !$omp parallel do OMP_SCHEDULE_ &
2663  !$omp private(Nc)
2664  !$acc kernels
2665  do j = js, je
2666  !$acc loop private(Nc)
2667  do i = is, ie
2668  if ( do_couple_aerosol ) then
2669  do k = ks, ke
2670  ! Nc(k) = max( CCN(k,i,j), Nc_def(i,j)*1.E+6_RP ) ! [#/m3] tentatively off the CCN effect
2671  nc(k) = nc_def(i,j) * 1.e+6_rp ! [#/m3]
2672  enddo
2673  else
2674  do k = ks, ke
2675  nc(k) = nc_def(i,j) * 1.e+6_rp ! [#/m3]
2676  enddo
2677  endif
2678 
2679  do k = ks, ke
2680  re(k,i,j,i_hc) = 1.1_rp &
2681  * ( 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)
2682  re(k,i,j,i_hc) = min( 1.e-3_rp, max( 1.e-6_rp, re(k,i,j,i_hc) ) ) * m2cm
2683  enddo
2684  enddo
2685  enddo
2686  !$acc end kernels
2687 
2688  endif
2689 
2690  if ( fixed_re ) then
2691 
2692  !$omp parallel do OMP_SCHEDULE_
2693  !$acc kernels
2694  do j = js, je
2695  do i = is, ie
2696  do k = ks, ke
2697  re(k,i,j,i_hr) = 10000.e-6_rp * m2cm
2698  end do
2699  end do
2700  end do
2701  !$acc end kernels
2702 
2703  !$omp parallel do OMP_SCHEDULE_
2704  !$acc kernels
2705  do j = js, je
2706  do i = is, ie
2707  do k = ks, ke
2708  re(k,i,j,i_hs) = re_qi * m2cm
2709  end do
2710  end do
2711  end do
2712  !$acc end kernels
2713 
2714  !$omp parallel do OMP_SCHEDULE_
2715  !$acc kernels
2716  do j = js, je
2717  do i = is, ie
2718  do k = ks, ke
2719  re(k,i,j,i_hg) = re_qi * m2cm
2720  end do
2721  end do
2722  end do
2723  !$acc end kernels
2724 
2725  else
2726 
2727 #ifndef __GFORTRAN__
2728  !$omp parallel do default(none) &
2729  !$omp shared(JS,JE,IS,IE,KS,KE,DENS0,TEMP0,QTRC0,enable_WDXZ2014,N0r_def) &
2730  !$omp shared(N0s_def,N0g_def,Ar,GAM_1br,Re,As,GAM_1bs) &
2731  !$omp shared(enable_RS2014,Ag,GAM_1bg) &
2732  !$omp private(i,j,k,dens,temc,qr,qs,qg,N0r,N0s,N0g,zerosw,RLMDr,RLMDs,Xs2,nm,tems,loga_) &
2733  !$omp private(b_,RLMDg,coef_at,coef_bt) &
2734  !$omp OMP_SCHEDULE_ collapse(2)
2735 #else
2736  !$omp parallel do default(shared) &
2737  !$omp private(i,j,k,dens,temc,qr,qs,qg,N0r,N0s,N0g,zerosw,RLMDr,RLMDs,Xs2,nm,tems,loga_) &
2738  !$omp private(b_,RLMDg,coef_at,coef_bt) &
2739  !$omp OMP_SCHEDULE_ collapse(2)
2740 #endif
2741  !$acc kernels
2742  do j = js, je
2743  !$acc loop private(dens,temc,qr,qs,qg,N0r,N0s,N0g)
2744  do i = is, ie
2745 
2746  do k = ks, ke
2747  dens(k) = dens0(k,i,j)
2748  temc(k) = temp0(k,i,j) - tem00
2749  qr(k) = qtrc0(k,i,j,i_hyd_qr)
2750  qs(k) = qtrc0(k,i,j,i_hyd_qs)
2751  qg(k) = qtrc0(k,i,j,i_hyd_qg)
2752  end do
2753 
2754  ! intercept parameter N0
2755  if ( enable_wdxz2014 ) then
2756  ! Wainwright et al. (2014)
2757  do k = ks, ke
2758  n0r(k) = 1.16e+5_rp * exp( log( max( dens(k)*qr(k)*1000.0_rp, 1.e-2_rp ) )*0.477_rp )
2759  n0s(k) = 4.58e+9_rp * exp( log( max( dens(k)*qs(k)*1000.0_rp, 1.e-2_rp ) )*0.788_rp )
2760  n0g(k) = 9.74e+8_rp * exp( log( max( dens(k)*qg(k)*1000.0_rp, 1.e-2_rp ) )*0.816_rp )
2761  end do
2762  else
2763  do k = ks, ke
2764  n0r(k) = n0r_def ! Marshall and Palmer (1948)
2765  n0s(k) = n0s_def ! Gunn and Marshall (1958)
2766  n0g(k) = n0g_def
2767  end do
2768  end if
2769 
2770  ! slope parameter lambda
2771  do k = ks, ke
2772  zerosw = 0.5_rp - sign(0.5_rp, qr(k) - 1.e-12_rp )
2773  rlmdr = sqrt(sqrt( dens(k) * qr(k) / ( ar * n0r(k) * gam_1br ) + zerosw )) * ( 1.0_rp-zerosw )
2774  ! Effective radius is defined by r3m/r2m = 1.5/lambda
2775  re(k,i,j,i_hr) = 1.5_rp * rlmdr * m2cm
2776  end do
2777 
2778  if ( enable_rs2014 ) then
2779  !$acc loop private(coef_at, coef_bt)
2780  do k = ks, ke
2781  zerosw = 0.5_rp - sign(0.5_rp, qs(k) - 1.e-12_rp )
2782  !---< modification by Roh and Satoh (2014) >---
2783  ! bimodal size distribution of snow
2784  zerosw = 0.5_rp - sign(0.5_rp, dens(k) * qs(k) - 1.e-12_rp )
2785  xs2 = dens(k) * qs(k) / as
2786 
2787  tems = min( -0.1_rp, temc(k) )
2788  coef_at(1) = coef_a01 + tems * ( coef_a02 + tems * ( coef_a05 + tems * coef_a09 ) )
2789  coef_at(2) = coef_a03 + tems * ( coef_a04 + tems * coef_a07 )
2790  coef_at(3) = coef_a06 + tems * coef_a08
2791  coef_at(4) = coef_a10
2792  coef_bt(1) = coef_b01 + tems * ( coef_b02 + tems * ( coef_b05 + tems * coef_b09 ) )
2793  coef_bt(2) = coef_b03 + tems * ( coef_b04 + tems * coef_b07 )
2794  coef_bt(3) = coef_b06 + tems * coef_b08
2795  coef_bt(4) = coef_b10
2796  ! 1 + Bs(=2) moment
2797  nm = 3.0_rp
2798  loga_ = coef_at(1) + nm * ( coef_at(2) + nm * ( coef_at(3) + nm * coef_at(4) ) )
2799  b_ = coef_bt(1) + nm * ( coef_bt(2) + nm * ( coef_bt(3) + nm * coef_bt(4) ) )
2800 
2801  re(k,i,j,i_hs) = 0.5_rp * exp( ln10 * loga_ + log(xs2+zerosw) * b_ ) * ( 1.0_rp-zerosw ) / ( xs2+zerosw ) * m2cm
2802  end do
2803  else
2804  do k = ks, ke
2805  zerosw = 0.5_rp - sign(0.5_rp, qs(k) - 1.e-12_rp )
2806  rlmds = sqrt(sqrt( dens(k) * qs(k) / ( as * n0s(k) * gam_1bs ) + zerosw )) * ( 1.0_rp-zerosw )
2807  re(k,i,j,i_hs) = 1.5_rp * rlmds * m2cm
2808  ! Re(k,i,j,I_HS) = dens * qs / N0s / ( 2.0_RP / 3.0_RP * PI * dens_i ) / RLMDs**3 * m2cm
2809  end do
2810  end if
2811 
2812  do k = ks, ke
2813  zerosw = 0.5_rp - sign(0.5_rp, qg(k) - 1.e-12_rp )
2814  rlmdg = sqrt(sqrt( dens(k) * qg(k) / ( ag * n0g(k) * gam_1bg ) + zerosw )) * ( 1.0_rp-zerosw )
2815  re(k,i,j,i_hg) = 1.5_rp * rlmdg * m2cm
2816  ! Re(k,i,j,I_HG) = dens * qg / N0g / ( 2.0_RP / 3.0_RP * PI * dens_i ) / RLMDg**3 * m2cm
2817  enddo
2818 
2819  enddo
2820  enddo
2821  !$acc end kernels
2822 
2823  endif
2824 
2825  !$acc end data
2826 
2827  return

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 atmos_phy_mp_tomita08_adjustment(), and mod_atmos_phy_mp_vars::atmos_phy_mp_vars_get_diagnostic().

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 2836 of file scale_atmos_phy_mp_tomita08.F90.

2836  use scale_atmos_hydrometeor, only: &
2837  n_hyd, &
2838  i_hc, &
2839  i_hr, &
2840  i_hi, &
2841  i_hs, &
2842  i_hg
2843  implicit none
2844  integer, intent(in) :: KA, KS, KE
2845  integer, intent(in) :: IA, IS, IE
2846  integer, intent(in) :: JA, JS, JE
2847 
2848  real(RP), intent(in) :: QTRC(KA,IA,JA,QA_MP-1)
2849 
2850  real(RP), intent(out) :: Qe(KA,IA,JA,N_HYD) ! mass ratio of each cateory [kg/kg]
2851 
2852  integer :: k, i, j, ih
2853  !---------------------------------------------------------------------------
2854 
2855  !$acc data copyin(QTRC) copyout(Qe)
2856 
2857  !$omp parallel do OMP_SCHEDULE_
2858  !$acc kernels async
2859 !OCL XFILL
2860  do j = js, je
2861  do i = is, ie
2862  do k = ks, ke
2863  qe(k,i,j,i_hc) = qtrc(k,i,j,i_hyd_qc)
2864  end do
2865  end do
2866  end do
2867  !$acc end kernels
2868 
2869  !$omp parallel do OMP_SCHEDULE_
2870 !OCL XFILL
2871  !$acc kernels async
2872  do j = js, je
2873  do i = is, ie
2874  do k = ks, ke
2875  qe(k,i,j,i_hr) = qtrc(k,i,j,i_hyd_qr)
2876  end do
2877  end do
2878  end do
2879  !$acc end kernels
2880 
2881  !$omp parallel do OMP_SCHEDULE_
2882  !$acc kernels async
2883 !OCL XFILL
2884  do j = js, je
2885  do i = is, ie
2886  do k = ks, ke
2887  qe(k,i,j,i_hi) = qtrc(k,i,j,i_hyd_qi)
2888  end do
2889  end do
2890  end do
2891  !$acc end kernels
2892 
2893  !$omp parallel do OMP_SCHEDULE_
2894  !$acc kernels async
2895 !OCL XFILL
2896  do j = js, je
2897  do i = is, ie
2898  do k = ks, ke
2899  qe(k,i,j,i_hs) = qtrc(k,i,j,i_hyd_qs)
2900  end do
2901  end do
2902  end do
2903  !$acc end kernels
2904 
2905  !$omp parallel do OMP_SCHEDULE_
2906  !$acc kernels async
2907 !OCL XFILL
2908  do j = js, je
2909  do i = is, ie
2910  do k = ks, ke
2911  qe(k,i,j,i_hg) = qtrc(k,i,j,i_hyd_qg)
2912  end do
2913  end do
2914  end do
2915  !$acc end kernels
2916 
2917  !$omp parallel do OMP_SCHEDULE_
2918  !$acc kernels async
2919 !OCL XFILL
2920  do ih = i_hg+1, n_hyd
2921  do j = js, je
2922  do i = is, ie
2923  do k = ks, ke
2924  qe(k,i,j,ih) = 0.0_rp
2925  end do
2926  end do
2927  end do
2928  end do
2929  !$acc end kernels
2930 
2931  !$acc wait
2932 
2933  !$acc end data
2934 
2935  return

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().

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 2944 of file scale_atmos_phy_mp_tomita08.F90.

2944  use scale_atmos_hydrometeor, only: &
2945  n_hyd, &
2946  i_hc, &
2947  i_hr, &
2948  i_hi, &
2949  i_hs, &
2950  i_hg, &
2951  i_hh
2952  implicit none
2953 
2954  integer, intent(in) :: KA, KS, KE
2955  integer, intent(in) :: IA, IS, IE
2956  integer, intent(in) :: JA, JS, JE
2957  real(RP), intent(in) :: Qe (KA,IA,JA,N_HYD) ! mass ratio of each cateory [kg/kg]
2958  real(RP), intent(out) :: QTRC(KA,IA,JA,QA_MP-1)
2959 
2960  integer :: k, i, j
2961  !---------------------------------------------------------------------------
2962 
2963  !$acc data copyin(Qe) copyout(QTRC)
2964 
2965  !$omp parallel do OMP_SCHEDULE_
2966  !$acc kernels async
2967 !OCL XFILL
2968  do j = js, je
2969  do i = is, ie
2970  do k = ks, ke
2971  qtrc(k,i,j,i_hyd_qc) = qe(k,i,j,i_hc)
2972  end do
2973  end do
2974  end do
2975  !$acc end kernels
2976 
2977  !$omp parallel do OMP_SCHEDULE_
2978  !$acc kernels async
2979 !OCL XFILL
2980  do j = js, je
2981  do i = is, ie
2982  do k = ks, ke
2983  qtrc(k,i,j,i_hyd_qr) = qe(k,i,j,i_hr)
2984  end do
2985  end do
2986  end do
2987  !$acc end kernels
2988 
2989  !$omp parallel do OMP_SCHEDULE_
2990  !$acc kernels async
2991 !OCL XFILL
2992  do j = js, je
2993  do i = is, ie
2994  do k = ks, ke
2995  qtrc(k,i,j,i_hyd_qi) = qe(k,i,j,i_hi)
2996  end do
2997  end do
2998  end do
2999  !$acc end kernels
3000 
3001  !$omp parallel do OMP_SCHEDULE_
3002  !$acc kernels async
3003 !OCL XFILL
3004  do j = js, je
3005  do i = is, ie
3006  do k = ks, ke
3007  qtrc(k,i,j,i_hyd_qs) = qe(k,i,j,i_hs)
3008  end do
3009  end do
3010  end do
3011  !$acc end kernels
3012 
3013  !$omp parallel do OMP_SCHEDULE_
3014  !$acc kernels async
3015 !OCL XFILL
3016  do j = js, je
3017  do i = is, ie
3018  do k = ks, ke
3019  qtrc(k,i,j,i_hyd_qg) = qe(k,i,j,i_hg) + qe(k,i,j,i_hh)
3020  end do
3021  end do
3022  end do
3023  !$acc end kernels
3024 
3025  !$acc wait
3026 
3027  !$acc end data
3028 
3029  return

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().

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 46 of file scale_atmos_phy_mp_tomita08.F90.

46  integer, parameter, public :: ATMOS_PHY_MP_tomita08_ntracers = qa_mp

Referenced by mod_atmos_phy_mp_driver::atmos_phy_mp_driver_tracer_setup().

◆ atmos_phy_mp_tomita08_nwaters

integer, parameter, public scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08_nwaters = 2

Definition at line 47 of file scale_atmos_phy_mp_tomita08.F90.

47  integer, parameter, public :: ATMOS_PHY_MP_tomita08_nwaters = 2

Referenced by mod_atmos_phy_mp_driver::atmos_phy_mp_driver_tracer_setup().

◆ atmos_phy_mp_tomita08_nices

integer, parameter, public scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08_nices = 3

Definition at line 48 of file scale_atmos_phy_mp_tomita08.F90.

48  integer, parameter, public :: ATMOS_PHY_MP_tomita08_nices = 3

Referenced by mod_atmos_phy_mp_driver::atmos_phy_mp_driver_tracer_setup().

◆ 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 49 of file scale_atmos_phy_mp_tomita08.F90.

49  character(len=H_SHORT), parameter, public :: ATMOS_PHY_MP_tomita08_tracer_names(QA_MP) = (/ &
50  'QV', &
51  'QC', &
52  'QR', &
53  'QI', &
54  'QS', &
55  'QG' /)

Referenced by mod_atmos_phy_mp_driver::atmos_phy_mp_driver_tracer_setup().

◆ 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 56 of file scale_atmos_phy_mp_tomita08.F90.

56  character(len=H_MID) , parameter, public :: ATMOS_PHY_MP_tomita08_tracer_descriptions(QA_MP) = (/ &
57  'Ratio of Water Vapor mass to total mass (Specific humidity)', &
58  'Ratio of Cloud Water mass to total mass ', &
59  'Ratio of Rain Water mass to total mass ', &
60  'Ratio of Cloud Ice mass ratio to total mass ', &
61  'Ratio of Snow miass ratio to total mass ', &
62  'Ratio of Graupel mass ratio to total mass '/)

Referenced by mod_atmos_phy_mp_driver::atmos_phy_mp_driver_tracer_setup().

◆ 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 63 of file scale_atmos_phy_mp_tomita08.F90.

63  character(len=H_SHORT), parameter, public :: ATMOS_PHY_MP_tomita08_tracer_units(QA_MP) = (/ &
64  'kg/kg', &
65  'kg/kg', &
66  'kg/kg', &
67  'kg/kg', &
68  'kg/kg', &
69  'kg/kg' /)

Referenced by mod_atmos_phy_mp_driver::atmos_phy_mp_driver_tracer_setup().

scale_const::const_grav
real(rp), public const_grav
standard acceleration of gravity [m/s2]
Definition: scale_const.F90:49
scale_prc::prc_abort
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:350
scale_atmos_hydrometeor::i_hr
integer, parameter, public i_hr
liquid water rain
Definition: scale_atmos_hydrometeor.F90:98
scale_atmos_hydrometeor::i_hs
integer, parameter, public i_hs
snow
Definition: scale_atmos_hydrometeor.F90:100
scale_atmos_hydrometeor
module atmosphere / hydrometeor
Definition: scale_atmos_hydrometeor.F90:12
scale_atmos_hydrometeor::i_hh
integer, parameter, public i_hh
hail
Definition: scale_atmos_hydrometeor.F90:102
scale_file_history
module file_history
Definition: scale_file_history.F90:15
scale_prc
module PROCESS
Definition: scale_prc.F90:11
scale_atmos_hydrometeor::i_hi
integer, parameter, public i_hi
ice water cloud
Definition: scale_atmos_hydrometeor.F90:99
scale_const
module CONSTANT
Definition: scale_const.F90:11
scale_atmos_phy_mp_common
module ATMOSPHERE / Physics Cloud Microphysics - Common
Definition: scale_atmos_phy_mp_common.F90:13
scale_specfunc
module SPECFUNC
Definition: scale_specfunc.F90:14
scale_const::const_pi
real(rp), parameter, public const_pi
pi
Definition: scale_const.F90:32
scale_atmos_hydrometeor::i_hc
integer, parameter, public i_hc
liquid water cloud
Definition: scale_atmos_hydrometeor.F90:97
scale_const::const_dwatr
real(rp), parameter, public const_dwatr
density of water [kg/m3]
Definition: scale_const.F90:89
scale_const::const_tem00
real(rp), parameter, public const_tem00
temperature reference (0C) [K]
Definition: scale_const.F90:99
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:685
scale_const::const_dice
real(rp), parameter, public const_dice
density of ice [kg/m3]
Definition: scale_const.F90:90
scale_specfunc::sf_gamma
real(rp) function, public sf_gamma(x)
Gamma function.
Definition: scale_specfunc.F90:50
scale_atmos_hydrometeor::i_hg
integer, parameter, public i_hg
graupel
Definition: scale_atmos_hydrometeor.F90:101