SCALE-RM
Functions/Subroutines | Variables
scale_atmos_phy_mp_tomita08 Module Reference

module ATMOSPHERE / Physics Cloud Microphysics More...

Functions/Subroutines

subroutine, public atmos_phy_mp_tomita08_config (MP_TYPE, QA, QS)
 Config. More...
 
subroutine, public atmos_phy_mp_tomita08_setup
 Setup. More...
 
subroutine, public atmos_phy_mp_tomita08 (DENS, MOMZ, MOMX, MOMY, RHOT, QTRC, CCN, EVAPORATE, SFLX_rain, SFLX_snow)
 Cloud Microphysics. More...
 
subroutine, public atmos_phy_mp_tomita08_cloudfraction (cldfrac, QTRC, mask_criterion)
 Calculate Cloud Fraction. More...
 
subroutine, public atmos_phy_mp_tomita08_effectiveradius (Re, QTRC0, DENS0, TEMP0)
 Calculate Effective Radius. More...
 
subroutine, public atmos_phy_mp_tomita08_mixingratio (Qe, QTRC0)
 Calculate mixing ratio of each category. More...
 

Variables

character(len=h_short), dimension(qa_mp), target, public atmos_phy_mp_tomita08_name
 
character(len=h_mid), dimension(qa_mp), target, public atmos_phy_mp_tomita08_desc
 
character(len=h_short), dimension(qa_mp), target, public atmos_phy_mp_tomita08_unit
 
real(rp), dimension(n_hyd), target, public atmos_phy_mp_tomita08_dens
 

Detailed Description

module ATMOSPHERE / Physics Cloud Microphysics

Description
Cloud Microphysics by Lin-type parametarization Reference: Tomita(2008)
Author
Team SCALE
History
  • 2013-03-25 (H.Yashiro) [new]
  • 2015-09-08 (Y.Sato) [add] Add evaporated cloud number concentration
NAMELIST
  • PARAM_ATMOS_PHY_MP
    nametypedefault valuecomment
    MP_DOPRECIPITATION logical .true. apply sedimentation (precipitation)?
    MP_DONEGATIVE_FIXER logical .true. apply negative fixer?
    MP_LIMIT_NEGATIVE real(RP) 1.0_RP Abort if abs(fixed negative vaue) > abs(MP_limit_negative)
    MP_DOEXPLICIT_ICEGEN logical .false. apply explicit ice generation?
    MP_NTMAX_SEDIMENTATION integer 1 number of time step for sedimentation
    MP_COUPLE_AEROSOL logical .false. apply CCN effect?

  • PARAM_ATMOS_PHY_MP_TOMITA08
    nametypedefault valuecomment
    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 mixing 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 mixing 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
Vterm_QG terminal velocity of QG m/s vterm
Vterm_QI terminal velocity of QI m/s vterm
Vterm_QR terminal velocity of QR m/s vterm
Vterm_QS terminal velocity of QS m/s vterm
pflux_QG precipitation flux of QG kg/m2/s FLX_hydro
pflux_QI precipitation flux of QI kg/m2/s FLX_hydro
pflux_QR precipitation flux of QR kg/m2/s FLX_hydro
pflux_QS precipitation flux of QS kg/m2/s FLX_hydro
Pgmlt individual tendency term in tomita08 kg/kg/s w3d

Function/Subroutine Documentation

◆ atmos_phy_mp_tomita08_config()

subroutine, public scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08_config ( character(len=*), intent(in)  MP_TYPE,
integer, intent(out)  QA,
integer, intent(out)  QS 
)

Config.

Definition at line 354 of file scale_atmos_phy_mp_tomita08.F90.

References scale_atmos_hydrometeor::atmos_hydrometeor_regist(), atmos_phy_mp_tomita08_desc, atmos_phy_mp_tomita08_name, atmos_phy_mp_tomita08_unit, scale_stdio::io_fid_log, scale_stdio::io_l, and scale_process::prc_mpistop().

Referenced by scale_atmos_phy_mp::atmos_phy_mp_config().

354  use scale_process, only: &
356  use scale_atmos_hydrometeor, only: &
358  implicit none
359 
360  character(len=*), intent(in) :: MP_TYPE
361  integer, intent(out) :: QA
362  integer, intent(out) :: QS
363  !---------------------------------------------------------------------------
364 
365  if( io_l ) write(io_fid_log,*)
366  if( io_l ) write(io_fid_log,*) '++++++ Module[Cloud Microphysics Tracer] / Categ[ATMOS PHYSICS] / Origin[SCALElib]'
367  if( io_l ) write(io_fid_log,*) '*** Tracers for Tomita (2008) 1-moment bulk 6 category'
368 
369  if ( mp_type /= 'TOMITA08' ) then
370  write(*,*) 'xxx ATMOS_PHY_MP_TYPE is not TOMITA08. Check!'
371  call prc_mpistop
372  endif
373 
374  call atmos_hydrometeor_regist( qs, & ! [OUT]
375  1, 2, 3, & ! [IN]
379 
380  qa = qa_mp
381  qs_mp = qs
382  qe_mp = qs + qa - 1
383 
384  i_qv = qs
385  i_qc = qs + i_mp_qc
386  i_qr = qs + i_mp_qr
387  i_qi = qs + i_mp_qi
388  i_qs = qs + i_mp_qs
389  i_qg = qs + i_mp_qg
390 
391  return
subroutine, public prc_mpistop
Abort MPI.
character(len=h_short), dimension(qa_mp), target, public atmos_phy_mp_tomita08_unit
character(len=h_short), dimension(qa_mp), target, public atmos_phy_mp_tomita08_name
module PROCESS
subroutine, public atmos_hydrometeor_regist(Q0, NV, NL, NI, NAME, DESC, UNIT, ADVC)
Regist tracer.
character(len=h_mid), dimension(qa_mp), target, public atmos_phy_mp_tomita08_desc
Here is the call graph for this function:
Here is the caller graph for this function:

◆ atmos_phy_mp_tomita08_setup()

subroutine, public scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08_setup ( )

Setup.

Definition at line 397 of file scale_atmos_phy_mp_tomita08.F90.

References atmos_phy_mp_tomita08_dens, scale_const::const_dice, scale_const::const_dwatr, scale_const::const_grav, scale_const::const_pi, scale_const::const_undef, scale_grid::grid_cdz, 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, scale_grid_index::ia, scale_grid_index::ie, scale_stdio::io_fid_conf, scale_stdio::io_fid_log, scale_stdio::io_fid_nml, scale_stdio::io_l, scale_stdio::io_nml, scale_grid_index::is, scale_grid_index::ja, scale_grid_index::je, scale_grid_index::js, scale_grid_index::ka, scale_process::prc_mpistop(), scale_specfunc::sf_gamma(), and scale_time::time_dtsec_atmos_phy_mp.

Referenced by scale_atmos_phy_mp::atmos_phy_mp_config().

397  use scale_process, only: &
399  use scale_const, only: &
400  undef => const_undef, &
401  pi => const_pi, &
402  grav => const_grav, &
403  dens_w => const_dwatr, &
404  dens_i => const_dice
405  use scale_specfunc, only: &
406  sf_gamma
407  use scale_time, only: &
409  use scale_grid, only: &
410  cdz => grid_cdz
411  implicit none
412 
413  real(RP) :: autoconv_nc = nc_ocn
414  logical :: enable_KK2000 = .false.
415  logical :: enable_RS2014 = .false.
416  logical :: enable_WDXZ2014 = .false.
417 
418  namelist / param_atmos_phy_mp / &
419  mp_doprecipitation, &
420  mp_donegative_fixer, &
421  mp_limit_negative, &
422  mp_doexplicit_icegen, &
423  mp_ntmax_sedimentation, &
424  mp_couple_aerosol
425 
426  namelist / param_atmos_phy_mp_tomita08 / &
427  autoconv_nc, &
428  enable_kk2000, &
429  enable_rs2014, &
430  enable_wdxz2014, &
431  n0r_def, &
432  n0s_def, &
433  n0g_def, &
434  dens_s, &
435  dens_g, &
436  re_qc, &
437  re_qi, &
438  drag_g, &
439  cr, &
440  cs, &
441  eiw, &
442  erw, &
443  esw, &
444  egw, &
445  eri, &
446  esi, &
447  egi, &
448  esr, &
449  egr, &
450  egs, &
451  gamma_sacr, &
452  gamma_gacs, &
453  mi, &
454  beta_saut, &
455  gamma_saut, &
456  qicrt_saut, &
457  beta_gaut, &
458  gamma_gaut, &
459  qscrt_gaut, &
460  fixed_re, &
461  const_rec, &
462  nofall_qr, &
463  nofall_qi, &
464  nofall_qs, &
465  nofall_qg
466 
467  real(RP), parameter :: max_term_vel = 10.0_rp !-- terminal velocity for calculate dt of sedimentation
468  integer :: nstep_max
469 
470  integer :: ierr
471  integer :: i, j
472  !---------------------------------------------------------------------------
473 
474  if( io_l ) write(io_fid_log,*)
475  if( io_l ) write(io_fid_log,*) '++++++ Module[Cloud Microphysics] / Categ[ATMOS PHYSICS] / Origin[SCALElib]'
476  if( io_l ) write(io_fid_log,*) '*** Tomita (2008) 1-moment bulk 6 category'
477 
478  allocate( w3d(ka,ia,ja,w_nmax) )
479  w3d(:,:,:,:) = 0.0_rp
480 
481  allocate( nc_def(ia,ja) )
482 
483  !--- read namelist
484  rewind(io_fid_conf)
485  read(io_fid_conf,nml=param_atmos_phy_mp,iostat=ierr)
486  if( ierr < 0 ) then !--- missing
487  if( io_l ) write(io_fid_log,*) '*** Not found namelist. Default used.'
488  elseif( ierr > 0 ) then !--- fatal error
489  write(*,*) 'xxx Not appropriate names in namelist PARAM_ATMOS_PHY_MP. Check!'
490  call prc_mpistop
491  endif
492  if( io_nml ) write(io_fid_nml,nml=param_atmos_phy_mp)
493 
494  if( io_l ) write(io_fid_log,*)
495  if( io_l ) write(io_fid_log,*) '*** Enable negative fixer? : ', mp_donegative_fixer
496  if( io_l ) write(io_fid_log,*) '*** Value limit of negative fixer (abs) : ', abs(mp_limit_negative)
497  if( io_l ) write(io_fid_log,*) '*** Enable sedimentation (precipitation)? : ', mp_doprecipitation
498  if( io_l ) write(io_fid_log,*) '*** Enable explicit ice generation? : ', mp_doexplicit_icegen
499 
500  nstep_max = ceiling( ( time_dtsec_atmos_phy_mp * max_term_vel ) / minval( cdz(:) ) )
501  mp_ntmax_sedimentation = max( mp_ntmax_sedimentation, nstep_max )
502 
503  mp_nstep_sedimentation = mp_ntmax_sedimentation
504  mp_rnstep_sedimentation = 1.0_rp / real(MP_ntmax_sedimentation,kind=rp)
505  mp_dtsec_sedimentation = time_dtsec_atmos_phy_mp * mp_rnstep_sedimentation
506 
507  if( io_l ) write(io_fid_log,*) '*** Timestep of sedimentation is divided into : ', mp_ntmax_sedimentation, 'step'
508  if( io_l ) write(io_fid_log,*) '*** DT of sedimentation : ', mp_dtsec_sedimentation, '[s]'
509 
510  !--- read namelist
511  rewind(io_fid_conf)
512  read(io_fid_conf,nml=param_atmos_phy_mp_tomita08,iostat=ierr)
513  if( ierr < 0 ) then !--- missing
514  if( io_l ) write(io_fid_log,*) '*** Not found namelist. Default used.'
515  elseif( ierr > 0 ) then !--- fatal error
516  write(*,*) 'xxx Not appropriate names in namelist PARAM_ATMOS_PHY_MP_TOMITA08. Check!'
517  call prc_mpistop
518  endif
519  if( io_nml ) write(io_fid_nml,nml=param_atmos_phy_mp_tomita08)
520 
521  if( io_l ) write(io_fid_log,*)
522  if( io_l ) write(io_fid_log,*) '*** density of the snow [kg/m3] : ', dens_s
523  if( io_l ) write(io_fid_log,*) '*** density of the graupel [kg/m3] : ', dens_g
524  if( io_l ) write(io_fid_log,*) '*** Nc for auto-conversion [num/m3]: ', autoconv_nc
525  if( io_l ) write(io_fid_log,*) '*** Use k-k scheme? : ', enable_kk2000
526  if( io_l ) write(io_fid_log,*) '*** Use Roh scheme? : ', enable_rs2014
527  if( io_l ) write(io_fid_log,*) '*** Use WDXZ scheme? : ', enable_wdxz2014
528  if( io_l ) write(io_fid_log,*)
529  if( io_l ) write(io_fid_log,*) '*** Use effective radius of ice for snow and graupel,'
530  if( io_l ) write(io_fid_log,*) ' and set rain transparent? : ', fixed_re
531  if( io_l ) write(io_fid_log,*) '*** Density of the ice is used for the calculation of '
532  if( io_l ) write(io_fid_log,*) ' optically effective volume of snow and graupel.'
533  if( io_l ) write(io_fid_log,*) '*** Surpress sedimentation of rain? : ', nofall_qr
534  if( io_l ) write(io_fid_log,*) '*** Surpress sedimentation of ice? : ', nofall_qi
535  if( io_l ) write(io_fid_log,*) '*** Surpress sedimentation of snow? : ', nofall_qs
536  if( io_l ) write(io_fid_log,*) '*** Surpress sedimentation of graupel? : ', nofall_qg
537  if( io_l ) write(io_fid_log,*)
538 
539  ! For the calculation of optically effective volume
540  atmos_phy_mp_tomita08_dens(:) = undef
541  atmos_phy_mp_tomita08_dens(i_hc) = dens_w
542  atmos_phy_mp_tomita08_dens(i_hr) = dens_w
543  atmos_phy_mp_tomita08_dens(i_hi) = dens_i
544  atmos_phy_mp_tomita08_dens(i_hs) = dens_i
545  atmos_phy_mp_tomita08_dens(i_hg) = dens_i
546 
547  do j = js, je
548  do i = is, ie
549  nc_def(i,j) = autoconv_nc
550  enddo
551  enddo
552 
553  !--- empirical coefficients A, B, C, D
554  ar = pi * dens_w / 6.0_rp
555  as = pi * dens_s / 6.0_rp
556  ag = pi * dens_g / 6.0_rp
557 
558  br = 3.0_rp
559  bs = 3.0_rp
560  bg = 3.0_rp
561 
562  cg = sqrt( ( 4.0_rp * dens_g * grav ) / ( 3.0_rp * dens00 * drag_g ) )
563 
564  dr = 0.50_rp
565  ds = 0.25_rp
566  dg = 0.50_rp
567 
568  if ( enable_rs2014 ) then ! overwrite parameters
569  mp_doexplicit_icegen = .true.
570 
571  sw_rs2014 = 1.0_rp
572  n0g_def = 4.e+8_rp
573  as = 0.069_rp
574  bs = 2.0_rp
575  esi = 0.25_rp
576  egi = 0.0_rp
577  egs = 0.0_rp
578  endif
579 
580  if ( mp_doexplicit_icegen ) then
581  only_liquid = .true.
582  sw_expice = 1.0_rp
583  else
584  only_liquid = .false.
585  sw_expice = 0.0_rp
586  endif
587 
588  if ( enable_kk2000 ) then
589  sw_kk2000 = 1.0_rp
590  else
591  sw_kk2000 = 0.0_rp
592  endif
593 
594  if ( enable_wdxz2014 ) then
595  sw_wdxz2014 = 1.0_rp
596  else
597  sw_wdxz2014 = 0.0_rp
598  endif
599 
600  gam = 1.0_rp ! =0!
601  gam_2 = 1.0_rp ! =1!
602  gam_3 = 2.0_rp ! =2!
603 
604  gam_1br = sf_gamma( 1.0_rp + br ) ! = 4!
605  gam_2br = sf_gamma( 2.0_rp + br ) ! = 5!
606  gam_3br = sf_gamma( 3.0_rp + br ) ! = 6!
607  gam_3dr = sf_gamma( 3.0_rp + dr )
608  gam_6dr = sf_gamma( 6.0_rp + dr )
609  gam_1brdr = sf_gamma( 1.0_rp + br + dr )
610  gam_5dr_h = sf_gamma( 0.5_rp * (5.0_rp+dr) )
611 
612  gam_1bs = sf_gamma( 1.0_rp + bs ) ! = 4!
613  gam_2bs = sf_gamma( 2.0_rp + bs ) ! = 5!
614  gam_3bs = sf_gamma( 3.0_rp + bs ) ! = 6!
615  gam_3ds = sf_gamma( 3.0_rp + ds )
616  gam_1bsds = sf_gamma( 1.0_rp + bs + ds )
617  gam_5ds_h = sf_gamma( 0.5_rp * (5.0_rp+ds) )
618 
619  gam_1bg = sf_gamma( 1.0_rp + bg ) ! = 4!
620  gam_3dg = sf_gamma( 3.0_rp + dg )
621  gam_1bgdg = sf_gamma( 1.0_rp + bg + dg)
622  gam_5dg_h = sf_gamma( 0.5_rp * (5.0_rp+dg) )
623 
624  ln10 = log(10.0_rp)
625 
626  return
integer, public is
start point of inner domain: x, local
integer, public je
end point of inner domain: y, local
subroutine, public prc_mpistop
Abort MPI.
real(rp), parameter, public const_dwatr
density of water [kg/m3]
Definition: scale_const.F90:84
real(dp), public time_dtsec_atmos_phy_mp
time interval of physics(microphysics) [sec]
Definition: scale_time.F90:41
real(rp) function, public sf_gamma(x)
Gamma function.
real(rp), parameter, public const_dice
density of ice [kg/m3]
Definition: scale_const.F90:85
real(rp), public const_undef
Definition: scale_const.F90:43
real(rp), public const_grav
standard acceleration of gravity [m/s2]
Definition: scale_const.F90:48
integer, public js
start point of inner domain: y, local
module TIME
Definition: scale_time.F90:15
module PROCESS
module SPECFUNC
module CONSTANT
Definition: scale_const.F90:14
real(rp), dimension(n_hyd), target, public atmos_phy_mp_tomita08_dens
module GRID (cartesian)
integer, public ie
end point of inner domain: x, local
real(rp), dimension(:), allocatable, public grid_cdz
z-length of control volume [m]
real(rp), public const_pi
pi
Definition: scale_const.F90:34
Here is the call graph for this function:
Here is the caller graph for this function:

◆ atmos_phy_mp_tomita08()

subroutine, public scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08 ( real(rp), dimension (ka,ia,ja), intent(inout)  DENS,
real(rp), dimension (ka,ia,ja), intent(inout)  MOMZ,
real(rp), dimension (ka,ia,ja), intent(inout)  MOMX,
real(rp), dimension (ka,ia,ja), intent(inout)  MOMY,
real(rp), dimension (ka,ia,ja), intent(inout)  RHOT,
real(rp), dimension (ka,ia,ja,qa), intent(inout)  QTRC,
real(rp), dimension (ka,ia,ja), intent(in)  CCN,
real(rp), dimension(ka,ia,ja), intent(out)  EVAPORATE,
real(rp), dimension(ia,ja), intent(out)  SFLX_rain,
real(rp), dimension(ia,ja), intent(out)  SFLX_snow 
)

Cloud Microphysics.

Definition at line 642 of file scale_atmos_phy_mp_tomita08.F90.

References scale_atmos_phy_mp_common::atmos_phy_mp_negative_fixer(), scale_atmos_phy_mp_common::atmos_phy_mp_precipitation(), scale_atmos_phy_mp_common::atmos_phy_mp_saturation_adjustment(), 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_atmos_hydrometeor::i_qc, scale_atmos_hydrometeor::i_qg, scale_atmos_hydrometeor::i_qi, scale_atmos_hydrometeor::i_qr, scale_atmos_hydrometeor::i_qs, scale_atmos_hydrometeor::i_qv, scale_grid_index::ie, scale_stdio::io_fid_log, scale_stdio::io_l, scale_grid_index::is, scale_grid_index::je, scale_grid_index::js, scale_grid_index::ka, scale_grid_index::ke, scale_grid_index::ks, scale_atmos_hydrometeor::lhf, scale_atmos_hydrometeor::lhv, scale_prof::prof_rapend(), scale_prof::prof_rapstart(), scale_tracer::qa, scale_time::time_dtsec_atmos_phy_mp, scale_tracer::tracer_cv, scale_tracer::tracer_mass, and scale_tracer::tracer_r.

Referenced by scale_atmos_phy_mp::atmos_phy_mp_config().

642  use scale_grid_index
643  use scale_const, only: &
644  dwatr => const_dwatr, &
645  pi => const_pi
646  use scale_time, only: &
647  dt_mp => time_dtsec_atmos_phy_mp
648  use scale_history, only: &
649  hist_in
650  use scale_tracer, only: &
651  qa, &
652  tracer_r, &
653  tracer_cv, &
655  use scale_atmos_thermodyn, only: &
656  thermodyn_rhoe => atmos_thermodyn_rhoe, &
657  thermodyn_rhot => atmos_thermodyn_rhot, &
658  thermodyn_temp_pres_e => atmos_thermodyn_temp_pres_e
659  use scale_atmos_phy_mp_common, only: &
660  mp_negative_fixer => atmos_phy_mp_negative_fixer, &
661  mp_precipitation => atmos_phy_mp_precipitation, &
662  mp_saturation_adjustment => atmos_phy_mp_saturation_adjustment
663  implicit none
664 
665  real(RP), intent(inout) :: DENS (KA,IA,JA)
666  real(RP), intent(inout) :: MOMZ (KA,IA,JA)
667  real(RP), intent(inout) :: MOMX (KA,IA,JA)
668  real(RP), intent(inout) :: MOMY (KA,IA,JA)
669  real(RP), intent(inout) :: RHOT (KA,IA,JA)
670  real(RP), intent(inout) :: QTRC (KA,IA,JA,QA)
671  real(RP), intent(in) :: CCN (KA,IA,JA)
672  real(RP), intent(out) :: EVAPORATE(KA,IA,JA) ! number of evaporated cloud [/m3/s]
673  real(RP), intent(out) :: SFLX_rain(IA,JA)
674  real(RP), intent(out) :: SFLX_snow(IA,JA)
675 
676  real(RP) :: RHOE_t (KA,IA,JA)
677  real(RP) :: QTRC_t (KA,IA,JA,QA)
678  real(RP) :: RHOE (KA,IA,JA)
679  real(RP) :: temp (KA,IA,JA)
680  real(RP) :: pres (KA,IA,JA)
681 
682  real(RP) :: vterm (KA,IA,JA,QA_MP-1) ! terminal velocity of each tracer [m/s]
683  real(RP) :: pflux (KA,IA,JA,QA_MP-1) ! precipitation flux of each tracer [kg/m2/s]
684  real(RP) :: FLX_hydro(KA,IA,JA,QA_MP-1)
685  real(RP) :: QC_t_sat (KA,IA,JA)
686  real(RP) :: QI_t_sat (KA,IA,JA)
687 
688  integer :: step
689  integer :: k, i, j
690  !---------------------------------------------------------------------------
691 
692  if( io_l ) write(io_fid_log,*) '*** Atmos physics step: Cloud microphysics(tomita08)'
693 
694  if ( mp_donegative_fixer ) then
695  call mp_negative_fixer( dens(:,:,:), & ! [INOUT]
696  rhot(:,:,:), & ! [INOUT]
697  qtrc(:,:,:,:), & ! [INOUT]
698  i_qv, & ! [IN]
699  mp_limit_negative ) ! [IN]
700  endif
701 
702  call thermodyn_rhoe( rhoe(:,:,:), & ! [OUT]
703  rhot(:,:,:), & ! [IN]
704  qtrc(:,:,:,:), & ! [IN]
705  tracer_cv(:), & ! [IN]
706  tracer_r(:), & ! [IN]
707  tracer_mass(:) ) ! [IN]
708 
709  !##### MP Main #####
710  call mp_tomita08( rhoe_t(:,:,:), & ! [OUT]
711  qtrc_t(:,:,:,:), & ! [OUT]
712  rhoe(:,:,:), & ! [INOUT]
713  qtrc(:,:,:,:), & ! [INOUT]
714  ccn(:,:,:), & ! [IN]
715  dens(:,:,:) ) ! [IN]
716 
717  vterm(:,:,:,:) = 0.0_rp
718  flx_hydro(:,:,:,:) = 0.0_rp
719 
720  if ( mp_doprecipitation ) then
721 
722  do step = 1, mp_nstep_sedimentation
723 
724  call thermodyn_temp_pres_e( temp(:,:,:), & ! [OUT]
725  pres(:,:,:), & ! [OUT]
726  dens(:,:,:), & ! [IN]
727  rhoe(:,:,:), & ! [IN]
728  qtrc(:,:,:,:), & ! [IN]
729  tracer_cv(:), & ! [IN]
730  tracer_r(:), & ! [IN]
731  tracer_mass(:) ) ! [IN]
732 
733  call mp_tomita08_vterm( vterm(:,:,:,:), & ! [OUT]
734  dens(:,:,:), & ! [IN]
735  temp(:,:,:), & ! [IN]
736  qtrc(:,:,:,:) ) ! [IN]
737 
738  call mp_precipitation( qa_mp, & ! [IN]
739  qs_mp, & ! [IN]
740  pflux(:,:,:,:), & ! [OUT]
741  vterm(:,:,:,:), & ! [INOUT]
742  dens(:,:,:), & ! [INOUT]
743  momz(:,:,:), & ! [INOUT]
744  momx(:,:,:), & ! [INOUT]
745  momy(:,:,:), & ! [INOUT]
746  rhoe(:,:,:), & ! [INOUT]
747  qtrc(:,:,:,:), & ! [INOUT]
748  temp(:,:,:), & ! [IN]
749  tracer_cv(:), & ! [IN]
750  mp_dtsec_sedimentation ) ! [IN]
751 
752  do j = js, je
753  do i = is, ie
754  do k = ks-1, ke-1
755  flx_hydro(k,i,j,i_mp_qr) = flx_hydro(k,i,j,i_mp_qr) + pflux(k,i,j,i_mp_qr) * mp_rnstep_sedimentation
756  flx_hydro(k,i,j,i_mp_qi) = flx_hydro(k,i,j,i_mp_qi) + pflux(k,i,j,i_mp_qi) * mp_rnstep_sedimentation
757  flx_hydro(k,i,j,i_mp_qs) = flx_hydro(k,i,j,i_mp_qs) + pflux(k,i,j,i_mp_qs) * mp_rnstep_sedimentation
758  flx_hydro(k,i,j,i_mp_qg) = flx_hydro(k,i,j,i_mp_qg) + pflux(k,i,j,i_mp_qg) * mp_rnstep_sedimentation
759  enddo
760  enddo
761  enddo
762 
763  enddo
764 
765  endif
766 
767  call hist_in( flx_hydro(:,:,:,i_mp_qr), 'pflux_QR', 'precipitation flux of QR', 'kg/m2/s', nohalo=.true. )
768  call hist_in( flx_hydro(:,:,:,i_mp_qi), 'pflux_QI', 'precipitation flux of QI', 'kg/m2/s', nohalo=.true. )
769  call hist_in( flx_hydro(:,:,:,i_mp_qs), 'pflux_QS', 'precipitation flux of QS', 'kg/m2/s', nohalo=.true. )
770  call hist_in( flx_hydro(:,:,:,i_mp_qg), 'pflux_QG', 'precipitation flux of QG', 'kg/m2/s', nohalo=.true. )
771 
772  call hist_in( vterm(:,:,:,i_mp_qr), 'Vterm_QR', 'terminal velocity of QR', 'm/s' )
773  call hist_in( vterm(:,:,:,i_mp_qi), 'Vterm_QI', 'terminal velocity of QI', 'm/s' )
774  call hist_in( vterm(:,:,:,i_mp_qs), 'Vterm_QS', 'terminal velocity of QS', 'm/s' )
775  call hist_in( vterm(:,:,:,i_mp_qg), 'Vterm_QG', 'terminal velocity of QG', 'm/s' )
776 
777  do j = js, je
778  do i = is, ie
779  do k = ks, ke
780  qc_t_sat(k,i,j) = qtrc(k,i,j,i_qc)
781  qi_t_sat(k,i,j) = qtrc(k,i,j,i_qi)
782  enddo
783  enddo
784  enddo
785 
786  call mp_saturation_adjustment( rhoe_t(:,:,:), & ! [INOUT]
787  qtrc_t(:,:,:,:), & ! [INOUT]
788  rhoe(:,:,:), & ! [INOUT]
789  qtrc(:,:,:,:), & ! [INOUT]
790  dens(:,:,:), & ! [IN]
791  i_qv, & ! [IN]
792  i_qc, & ! [IN]
793  i_qi, & ! [IN]
794  only_liquid ) ! [IN]
795 
796  do j = js, je
797  do i = is, ie
798  do k = ks, ke
799  qc_t_sat(k,i,j) = ( qtrc(k,i,j,i_qc) - qc_t_sat(k,i,j) ) / dt_mp
800  qi_t_sat(k,i,j) = ( qtrc(k,i,j,i_qi) - qi_t_sat(k,i,j) ) / dt_mp
801  enddo
802  enddo
803  enddo
804 
805  call hist_in( qc_t_sat(:,:,:), 'Pcsat', 'QC production term by satadjust', 'kg/kg/s' )
806  call hist_in( qi_t_sat(:,:,:), 'Pisat', 'QI production term by satadjust', 'kg/kg/s' )
807 
808  do j = js, je
809  do i = is, ie
810  do k = ks, ke
811  evaporate(k,i,j) = max( -qc_t_sat(k,i,j), 0.0_rp ) ! if negative, condensation
812  evaporate(k,i,j) = evaporate(k,i,j) * dens(k,i,j) / (4.0_rp/3.0_rp*pi*dwatr*re_qc**3) ! mass -> number (assuming constant particle radius as re_qc)
813  enddo
814  enddo
815  enddo
816 
817  !##### END MP Main #####
818 
819  call thermodyn_rhot( rhot(:,:,:), & ! [OUT]
820  rhoe(:,:,:), & ! [IN]
821  qtrc(:,:,:,:), & ! [IN]
822  tracer_cv(:), & ! [IN]
823  tracer_r(:), & ! [IN]
824  tracer_mass(:) ) ! [IN]
825 
826  if ( mp_donegative_fixer ) then
827  call mp_negative_fixer( dens(:,:,:), & ! [INOUT]
828  rhot(:,:,:), & ! [INOUT]
829  qtrc(:,:,:,:), & ! [INOUT]
830  i_qv, & ! [IN]
831  mp_limit_negative ) ! [IN]
832  endif
833 
834  do j = js, je
835  do i = is, ie
836  sflx_rain(i,j) = - ( flx_hydro(ks-1,i,j,i_mp_qr) )
837  sflx_snow(i,j) = - ( flx_hydro(ks-1,i,j,i_mp_qi) &
838  + flx_hydro(ks-1,i,j,i_mp_qs) &
839  + flx_hydro(ks-1,i,j,i_mp_qg) )
840  enddo
841  enddo
842 
843  return
integer, public is
start point of inner domain: x, local
integer, public je
end point of inner domain: y, local
real(rp), parameter, public const_dwatr
density of water [kg/m3]
Definition: scale_const.F90:84
subroutine, public atmos_phy_mp_negative_fixer(DENS, RHOT, QTRC, I_QV, limit_negative)
Negative fixer.
real(rp), dimension(qa_max), public tracer_r
real(dp), public time_dtsec_atmos_phy_mp
time interval of physics(microphysics) [sec]
Definition: scale_time.F90:41
integer, public ke
end point of inner domain: z, local
real(rp), dimension(qa_max), public tracer_cv
module ATMOSPHERE / Physics Cloud Microphysics - Common
module grid index
module TRACER
subroutine, public atmos_phy_mp_saturation_adjustment(RHOE_t, QTRC_t, RHOE0, QTRC0, DENS0, I_QV, I_QC, I_QI, flag_liquid)
Saturation adjustment.
integer, public js
start point of inner domain: y, local
module TIME
Definition: scale_time.F90:15
module CONSTANT
Definition: scale_const.F90:14
integer, public ks
start point of inner domain: z, local
integer, public ie
end point of inner domain: x, local
module ATMOSPHERE / Thermodynamics
module HISTORY
real(rp), public const_pi
pi
Definition: scale_const.F90:34
subroutine, public atmos_phy_mp_precipitation(QA_MP, QS_MP, qflx, vterm, DENS, MOMZ, MOMX, MOMY, RHOE, QTRC, temp, CVq, dt, vt_fixed)
precipitation transport
real(rp), dimension(qa_max), public tracer_mass
Here is the call graph for this function:
Here is the caller graph for this function:

◆ atmos_phy_mp_tomita08_cloudfraction()

subroutine, public scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08_cloudfraction ( real(rp), dimension(ka,ia,ja), intent(out)  cldfrac,
real(rp), dimension (ka,ia,ja,qa), intent(in)  QTRC,
real(rp), intent(in)  mask_criterion 
)

Calculate Cloud Fraction.

Definition at line 1964 of file scale_atmos_phy_mp_tomita08.F90.

References scale_grid_index::ie, scale_grid_index::is, scale_grid_index::je, scale_grid_index::js, scale_grid_index::ke, scale_grid_index::ks, and scale_tracer::qa.

Referenced by scale_atmos_phy_mp::atmos_phy_mp_config().

1964  use scale_grid_index
1965  use scale_tracer, only: &
1966  qa
1967  implicit none
1968 
1969  real(RP), intent(out) :: cldfrac(KA,IA,JA)
1970  real(RP), intent(in) :: QTRC (KA,IA,JA,QA)
1971  real(RP), intent(in) :: mask_criterion
1972 
1973  real(RP) :: qhydro
1974  integer :: k, i, j, iq
1975  !---------------------------------------------------------------------------
1976 
1977  do j = js, je
1978  do i = is, ie
1979  do k = ks, ke
1980  qhydro = 0.0_rp
1981  do iq = qs_mp+1, qe_mp
1982  qhydro = qhydro + qtrc(k,i,j,iq)
1983  enddo
1984  cldfrac(k,i,j) = 0.5_rp + sign(0.5_rp,qhydro-mask_criterion)
1985  enddo
1986  enddo
1987  enddo
1988 
1989  return
integer, public is
start point of inner domain: x, local
integer, public je
end point of inner domain: y, local
integer, public ke
end point of inner domain: z, local
module grid index
module TRACER
integer, public js
start point of inner domain: y, local
integer, public ks
start point of inner domain: z, local
integer, public ie
end point of inner domain: x, local
Here is the caller graph for this function:

◆ atmos_phy_mp_tomita08_effectiveradius()

subroutine, public scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08_effectiveradius ( real(rp), dimension (ka,ia,ja,n_hyd), intent(out)  Re,
real(rp), dimension(ka,ia,ja,qa), intent(in)  QTRC0,
real(rp), dimension(ka,ia,ja), intent(in)  DENS0,
real(rp), dimension(ka,ia,ja), intent(in)  TEMP0 
)

Calculate Effective Radius.

Definition at line 1999 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_grid_index::ie, scale_grid_index::is, scale_grid_index::je, scale_grid_index::js, scale_grid_index::ke, scale_grid_index::ks, scale_atmos_hydrometeor::n_hyd, and scale_tracer::qa.

Referenced by scale_atmos_phy_mp::atmos_phy_mp_config().

1999  use scale_grid_index
2000  use scale_const, only: &
2001  pi => const_pi, &
2002  dens_w => const_dwatr, &
2003  dens_i => const_dice, &
2004  tem00 => const_tem00
2005  use scale_tracer, only: &
2006  qa
2007  use scale_atmos_hydrometeor, only: &
2008  n_hyd
2009  implicit none
2010 
2011  real(RP), intent(out) :: Re (KA,IA,JA,N_HYD) ! effective radius [cm]
2012  real(RP), intent(in) :: QTRC0(KA,IA,JA,QA) ! tracer mass concentration [kg/kg]
2013  real(RP), intent(in) :: DENS0(KA,IA,JA) ! density [kg/m3]
2014  real(RP), intent(in) :: TEMP0(KA,IA,JA) ! temperature [K]
2015 
2016  real(RP) :: dens
2017  real(RP) :: temc
2018  real(RP) :: qc, qr, qs, qg
2019  real(RP) :: Nc(KA,IA,JA)
2020  real(RP) :: N0r, N0s, N0g
2021  real(RP) :: RLMDr, RLMDs, RLMDg
2022 
2023  real(RP), parameter :: um2cm = 100.0_rp
2024 
2025  !---< Roh and Satoh (2014) >---
2026  real(RP) :: tems, Xs2
2027  real(RP) :: coef_at(4), coef_bt(4)
2028  real(RP) :: loga_, b_, nm
2029 
2030  real(RP) :: zerosw
2031  integer :: k, i, j
2032  !---------------------------------------------------------------------------
2033 
2034  re(:,:,:,i_hi) = re_qi * um2cm
2035  re(:,:,:,i_hg+1:) = 0.0_rp
2036 
2037  if ( const_rec .or. fixed_re ) then
2038 
2039  re(:,:,:,i_hc) = re_qc * um2cm
2040 
2041  else
2042 
2043  if ( mp_couple_aerosol ) then
2044  do j = js, je
2045  do i = is, ie
2046  do k = ks, ke
2047  ! Nc(k,i,j) = max( CCN(k,i,j), Nc_def(i,j)*1.E+6_RP ) ! [#/m3] tentatively off the CCN effect
2048  nc(k,i,j) = nc_def(i,j) * 1.e+6_rp ! [#/m3]
2049  enddo
2050  enddo
2051  enddo
2052  else
2053  do j = js, je
2054  do i = is, ie
2055  do k = ks, ke
2056  nc(k,i,j) = nc_def(i,j) * 1.e+6_rp ! [#/m3]
2057  enddo
2058  enddo
2059  enddo
2060  endif
2061 
2062  do j = js, je
2063  do i = is, ie
2064  do k = ks, ke
2065  dens = dens0(k,i,j)
2066  qc = qtrc0(k,i,j,i_qc)
2067 
2068  re(k,i,j,i_hc) = 1.1_rp * ( dens * qc / nc(k,i,j) / ( 4.0_rp / 3.0_rp * pi * dens_w ) )**(1.0_rp/3.0_rp)
2069  re(k,i,j,i_hc) = min( 1.e-3_rp, max( 1.e-6_rp, re(k,i,j,i_hc) ) ) * um2cm
2070  enddo
2071  enddo
2072  enddo
2073 
2074  endif
2075 
2076  if ( fixed_re ) then
2077 
2078  re(:,:,:,i_hr) = 10000.e-6_rp * um2cm
2079  re(:,:,:,i_hs) = re_qi * um2cm
2080  re(:,:,:,i_hg) = re_qi * um2cm
2081 
2082  else
2083 
2084 #ifndef __GFORTRAN__
2085  !$omp parallel do default(none) &
2086  !$omp shared(JS,JE,IS,IE,KS,KE,DENS0,TEMP0,QTRC0,I_QR,I_QS,I_QG,sw_WDXZ2014,N0r_def) &
2087  !$omp shared(N0s_def,N0g_def,Ar,GAM_1br,Re,As,GAM_1bs) &
2088  !$omp shared(sw_RS2014,ln10,Ag,GAM_1bg) &
2089  !$omp private(i,j,k,dens,temc,qr,qs,qg,N0r,N0s,N0g,zerosw,RLMDr,RLMDs,Xs2,nm,tems,loga_) &
2090  !$omp private(b_,RLMDg,coef_at,coef_bt) &
2091  !$omp OMP_SCHEDULE_ collapse(2)
2092 #else
2093  !$omp parallel do default(shared) &
2094  !$omp private(i,j,k,dens,temc,qr,qs,qg,N0r,N0s,N0g,zerosw,RLMDr,RLMDs,Xs2,nm,tems,loga_) &
2095  !$omp private(b_,RLMDg,coef_at,coef_bt) &
2096  !$omp OMP_SCHEDULE_ collapse(2)
2097 #endif
2098  do j = js, je
2099  do i = is, ie
2100  do k = ks, ke
2101  dens = dens0(k,i,j)
2102  temc = temp0(k,i,j) - tem00
2103  qr = qtrc0(k,i,j,i_qr)
2104  qs = qtrc0(k,i,j,i_qs)
2105  qg = qtrc0(k,i,j,i_qg)
2106 
2107  ! intercept parameter N0
2108  n0r = ( 1.0_rp-sw_wdxz2014 ) * n0r_def & ! Marshall and Palmer (1948)
2109  + ( sw_wdxz2014 ) * 1.16e+5_rp * exp( log( max( dens*qr*1000.0_rp, 1.e-2_rp ) )*0.477_rp ) ! Wainwright et al. (2014)
2110  n0s = ( 1.0_rp-sw_wdxz2014 ) * n0s_def & ! Gunn and Marshall (1958)
2111  + ( sw_wdxz2014 ) * 4.58e+9_rp * exp( log( max( dens*qs*1000.0_rp, 1.e-2_rp ) )*0.788_rp ) ! Wainwright et al. (2014)
2112  n0g = ( 1.0_rp-sw_wdxz2014 ) * n0g_def & !
2113  + ( sw_wdxz2014 ) * 9.74e+8_rp * exp( log( max( dens*qg*1000.0_rp, 1.e-2_rp ) )*0.816_rp ) ! Wainwright et al. (2014)
2114 
2115  ! slope parameter lambda
2116  zerosw = 0.5_rp - sign(0.5_rp, qr - 1.e-12_rp )
2117  rlmdr = sqrt(sqrt( dens * qr / ( ar * n0r * gam_1br ) + zerosw )) * ( 1.0_rp-zerosw )
2118  ! Effective radius is defined by r3m/r2m = 1.5/lambda
2119  re(k,i,j,i_hr) = 1.5_rp * rlmdr * um2cm
2120 
2121 
2122  zerosw = 0.5_rp - sign(0.5_rp, qs - 1.e-12_rp )
2123  rlmds = sqrt(sqrt( dens * qs / ( as * n0s * gam_1bs ) + zerosw )) * ( 1.0_rp-zerosw )
2124  !---< modification by Roh and Satoh (2014) >---
2125  ! bimodal size distribution of snow
2126  zerosw = 0.5_rp - sign(0.5_rp, dens * qs - 1.e-12_rp )
2127  xs2 = dens * qs / as
2128 
2129  tems = min( -0.1_rp, temc )
2130  coef_at(1) = coef_a( 1) + tems * ( coef_a( 2) + tems * ( coef_a( 5) + tems * coef_a( 9) ) )
2131  coef_at(2) = coef_a( 3) + tems * ( coef_a( 4) + tems * coef_a( 7) )
2132  coef_at(3) = coef_a( 6) + tems * coef_a( 8)
2133  coef_at(4) = coef_a(10)
2134  coef_bt(1) = coef_b( 1) + tems * ( coef_b( 2) + tems * ( coef_b( 5) + tems * coef_b( 9) ) )
2135  coef_bt(2) = coef_b( 3) + tems * ( coef_b( 4) + tems * coef_b( 7) )
2136  coef_bt(3) = coef_b( 6) + tems * coef_b( 8)
2137  coef_bt(4) = coef_b(10)
2138 
2139  ! 1 + Bs(=2) moment
2140  nm = 3.0_rp
2141  loga_ = coef_at(1) + nm * ( coef_at(2) + nm * ( coef_at(3) + nm * coef_at(4) ) )
2142  b_ = coef_bt(1) + nm * ( coef_bt(2) + nm * ( coef_bt(3) + nm * coef_bt(4) ) )
2143 
2144  re(k,i,j,i_hs) = ( sw_rs2014 ) * 0.5_rp * exp(ln10*loga_) * exp(log(xs2+zerosw)*b_) * ( 1.0_rp-zerosw ) / ( xs2+zerosw ) * um2cm &
2145  + ( 1.0_rp-sw_rs2014 ) * 1.5_rp * rlmds * um2cm
2146 ! + ( 1.0_RP-sw_RS2014 ) * dens * qs / N0s / ( 2.0_RP / 3.0_RP * PI * dens_i ) / RLMDs**3 * um2cm
2147 
2148  zerosw = 0.5_rp - sign(0.5_rp, qg - 1.e-12_rp )
2149  rlmdg = sqrt(sqrt( dens * qg / ( ag * n0g * gam_1bg ) + zerosw )) * ( 1.0_rp-zerosw )
2150  re(k,i,j,i_hg) = 1.5_rp * rlmdg * um2cm
2151 ! Re(k,i,j,I_HG) = dens * qg / N0g / ( 2.0_RP / 3.0_RP * PI * dens_i ) / RLMDg**3 * um2cm
2152  enddo
2153  enddo
2154  enddo
2155 
2156  endif
2157 
2158  return
integer, public is
start point of inner domain: x, local
integer, public je
end point of inner domain: y, local
real(rp), parameter, public const_dwatr
density of water [kg/m3]
Definition: scale_const.F90:84
integer, public ke
end point of inner domain: z, local
real(rp), parameter, public const_tem00
temperature reference (0C) [K]
Definition: scale_const.F90:92
real(rp), parameter, public const_dice
density of ice [kg/m3]
Definition: scale_const.F90:85
module grid index
module TRACER
integer, public js
start point of inner domain: y, local
module CONSTANT
Definition: scale_const.F90:14
integer, public ks
start point of inner domain: z, local
integer, public ie
end point of inner domain: x, local
real(rp), public const_pi
pi
Definition: scale_const.F90:34
Here is the caller graph for this function:

◆ atmos_phy_mp_tomita08_mixingratio()

subroutine, public scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08_mixingratio ( real(rp), dimension (ka,ia,ja,n_hyd), intent(out)  Qe,
real(rp), dimension(ka,ia,ja,qa), intent(in)  QTRC0 
)

Calculate mixing ratio of each category.

Definition at line 2166 of file scale_atmos_phy_mp_tomita08.F90.

References scale_atmos_hydrometeor::n_hyd, and scale_tracer::qa.

Referenced by scale_atmos_phy_mp::atmos_phy_mp_config().

2166  use scale_grid_index
2167  use scale_tracer, only: &
2168  qa
2169  use scale_atmos_hydrometeor, only: &
2170  n_hyd
2171  implicit none
2172 
2173  real(RP), intent(out) :: Qe (KA,IA,JA,N_HYD) ! mixing ratio of each cateory [kg/kg]
2174  real(RP), intent(in) :: QTRC0(KA,IA,JA,QA) ! tracer mass concentration [kg/kg]
2175  !---------------------------------------------------------------------------
2176 
2177  qe(:,:,:,i_hc) = qtrc0(:,:,:,i_qc)
2178  qe(:,:,:,i_hr) = qtrc0(:,:,:,i_qr)
2179  qe(:,:,:,i_hi) = qtrc0(:,:,:,i_qi)
2180  qe(:,:,:,i_hs) = qtrc0(:,:,:,i_qs)
2181  qe(:,:,:,i_hg) = qtrc0(:,:,:,i_qg)
2182  qe(:,:,:,i_hg+1:) = 0.0_rp
2183 
2184  return
module grid index
module TRACER
Here is the caller graph for this function:

Variable Documentation

◆ atmos_phy_mp_tomita08_name

character(len=h_short), dimension(qa_mp), target, public scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08_name

Definition at line 60 of file scale_atmos_phy_mp_tomita08.F90.

Referenced by scale_atmos_phy_mp::atmos_phy_mp_config(), and atmos_phy_mp_tomita08_config().

60  character(len=H_SHORT), public, target :: ATMOS_PHY_MP_tomita08_NAME(QA_MP)

◆ atmos_phy_mp_tomita08_desc

character(len=h_mid), dimension(qa_mp), target, public scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08_desc

Definition at line 61 of file scale_atmos_phy_mp_tomita08.F90.

Referenced by scale_atmos_phy_mp::atmos_phy_mp_config(), and atmos_phy_mp_tomita08_config().

61  character(len=H_MID) , public, target :: ATMOS_PHY_MP_tomita08_DESC(QA_MP)

◆ atmos_phy_mp_tomita08_unit

character(len=h_short), dimension(qa_mp), target, public scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08_unit

Definition at line 62 of file scale_atmos_phy_mp_tomita08.F90.

Referenced by scale_atmos_phy_mp::atmos_phy_mp_config(), and atmos_phy_mp_tomita08_config().

62  character(len=H_SHORT), public, target :: ATMOS_PHY_MP_tomita08_UNIT(QA_MP)

◆ atmos_phy_mp_tomita08_dens

real(rp), dimension(n_hyd), target, public scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08_dens

Definition at line 64 of file scale_atmos_phy_mp_tomita08.F90.

Referenced by scale_atmos_phy_mp::atmos_phy_mp_config(), and atmos_phy_mp_tomita08_setup().

64  real(RP), public, target :: ATMOS_PHY_MP_tomita08_DENS(N_HYD) ! hydrometeor density [kg/m3]=[g/L]