SCALE-RM
scale_atmos_phy_mp.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
14 ! Warning: This file was generated from atmos-physics/microphysics/scale_atmos_phy_mp.F90.erb.
15 ! Do not edit this file.
16 !-------------------------------------------------------------------------------
18  !-----------------------------------------------------------------------------
19  !
20  !++ used modules
21  !
22  use scale_precision
23  use scale_stdio
24  use scale_prof
26  use scale_tracer
27  !-----------------------------------------------------------------------------
28  implicit none
29  private
30  !-----------------------------------------------------------------------------
31  !
32  !++ Public procedure
33  !
34  public :: atmos_phy_mp_config
35 
36  abstract interface
37  subroutine mp( &
38  DENS, &
39  MOMZ, &
40  MOMX, &
41  MOMY, &
42  RHOT, &
43  QTRC, &
44  CCN, &
45  EVAPORATE, &
46  SFLX_rain, &
47  SFLX_snow )
48  use scale_precision
50  use scale_tracer
51  real(RP), intent(inout) :: DENS(KA,IA,JA)
52  real(RP), intent(inout) :: MOMZ(KA,IA,JA)
53  real(RP), intent(inout) :: MOMX(KA,IA,JA)
54  real(RP), intent(inout) :: MOMY(KA,IA,JA)
55  real(RP), intent(inout) :: RHOT(KA,IA,JA)
56  real(RP), intent(inout) :: QTRC(KA,IA,JA,QA)
57  real(RP), intent(in) :: CCN(KA,IA,JA)
58  real(RP), intent(out) :: EVAPORATE(KA,IA,JA)
59  real(RP), intent(out) :: SFLX_rain(IA,JA)
60  real(RP), intent(out) :: SFLX_snow(IA,JA)
61  end subroutine mp
62  subroutine su
63  end subroutine su
64  subroutine cf( &
65  cldfrac, &
66  QTRC, &
67  mask_criterion )
68  use scale_precision
70  use scale_tracer
71  real(RP), intent(out) :: cldfrac(KA,IA,JA)
72  real(RP), intent(in) :: QTRC (KA,IA,JA,QA)
73  real(RP), intent(in) :: mask_criterion
74  end subroutine cf
75  subroutine er( &
76  Re, &
77  QTRC0, &
78  DENS0, &
79  TEMP0 )
80  use scale_precision
82  use scale_tracer
84  real(RP), intent(out) :: Re (KA,IA,JA,N_HYD) ! effective radius [cm]
85  real(RP), intent(in) :: QTRC0(KA,IA,JA,QA) ! tracer mass concentration [kg/kg]
86  real(RP), intent(in) :: DENS0(KA,IA,JA) ! density [kg/m3]
87  real(RP), intent(in) :: TEMP0(KA,IA,JA) ! temperature [K]
88  end subroutine er
89  subroutine mr( &
90  Qe, &
91  QTRC0 )
92  use scale_precision
94  use scale_tracer
96  real(RP), intent(out) :: Qe (KA,IA,JA,N_HYD) ! mixing ratio of each cateory [kg/kg]
97  real(RP), intent(in) :: QTRC0(KA,IA,JA,QA) ! tracer mass concentration [kg/kg]
98  end subroutine mr
99  end interface
100  procedure(mp), pointer :: atmos_phy_mp => null()
101  procedure(su), pointer :: atmos_phy_mp_setup => null()
102  procedure(cf), pointer :: atmos_phy_mp_cloudfraction => null()
103  procedure(er), pointer :: atmos_phy_mp_effectiveradius => null()
104  procedure(mr), pointer :: atmos_phy_mp_mixingratio => null()
105  public :: atmos_phy_mp
106  public :: atmos_phy_mp_setup
109  public :: atmos_phy_mp_mixingratio
110  !-----------------------------------------------------------------------------
111  !
112  !++ Public parameters & variables
113  !
114  integer, public :: qa_mp = 0 ! number of cloud microphysical tracers
115  integer, public :: qs_mp = -1 ! start index in QTRC
116  integer, public :: qe_mp = -1 ! end index in QTRC
117 
118  character(len=H_SHORT), pointer, public :: atmos_phy_mp_name(:)
119  character(len=H_MID), pointer, public :: atmos_phy_mp_desc(:)
120  character(len=H_SHORT), pointer, public :: atmos_phy_mp_unit(:)
121  real(RP), pointer, public :: atmos_phy_mp_dens(:)
122  !-----------------------------------------------------------------------------
123  !
124  !++ Private procedure
125  !
126  !-----------------------------------------------------------------------------
127  !
128  !++ Private parameters & variables
129  !
130  !-----------------------------------------------------------------------------
131 contains
132  !-----------------------------------------------------------------------------
134  !-----------------------------------------------------------------------------
135  subroutine atmos_phy_mp_config( MP_TYPE )
136  use scale_process, only: &
138  use scale_atmos_phy_mp_dry, only: &
149  use scale_atmos_phy_mp_kessler, only: &
160  use scale_atmos_phy_mp_tomita08, only: &
171  use scale_atmos_phy_mp_sn14, only: &
182  use scale_atmos_phy_mp_suzuki10, only: &
193  use scale_atmos_phy_mp_sdm, only: &
204  implicit none
205 
206  character(len=*), intent(in) :: mp_type
207  !---------------------------------------------------------------------------
208 
209  if( io_l ) write(io_fid_log,*) '*** => ', trim(mp_type), ' is selected.'
210 
211  select case( mp_type )
212  case( 'DRY' )
214  mp_type, & ! (in)
215  qa_mp, qs_mp ) ! (out)
225  case( 'KESSLER' )
227  mp_type, & ! (in)
228  qa_mp, qs_mp ) ! (out)
238  case( 'TOMITA08' )
240  mp_type, & ! (in)
241  qa_mp, qs_mp ) ! (out)
251  case( 'SN14' )
253  mp_type, & ! (in)
254  qa_mp, qs_mp ) ! (out)
264  case( 'SUZUKI10' )
266  mp_type, & ! (in)
267  qa_mp, qs_mp ) ! (out)
277  case( 'SDM' )
279  mp_type, & ! (in)
280  qa_mp, qs_mp ) ! (out)
290  end select
291 
292  qe_mp = qs_mp + qa_mp - 1
293 
294  return
295  end subroutine atmos_phy_mp_config
296 
297 end module scale_atmos_phy_mp
character(len=h_mid), dimension(qa_mp), target, public atmos_phy_mp_sn14_desc
procedure(cf), pointer, public atmos_phy_mp_cloudfraction
subroutine, public atmos_phy_mp_kessler(DENS, MOMZ, MOMX, MOMY, RHOT, QTRC, CCN, EVAPORATE, SFLX_rain, SFLX_snow)
Cloud Microphysics.
subroutine, public atmos_phy_mp_sn14_effectiveradius(Re, QTRC0, DENS0, TEMP0)
Calculate Effective Radius.
procedure(mp), pointer, public atmos_phy_mp
subroutine, public atmos_phy_mp_sdm_mixingratio(Qe, QTRC0)
Calculate mixing ratio of each category.
subroutine, public atmos_phy_mp_sdm_cloudfraction(cldfrac, QTRC, mask_criterion)
Calculate Cloud Fraction.
subroutine, public atmos_phy_mp_dry_cloudfraction(cldfrac, QTRC, mask_criterion)
Calculate Cloud Fraction.
subroutine, public atmos_phy_mp_suzuki10_config(MP_TYPE, QA, QS)
Config.
subroutine, public prc_mpistop
Abort MPI.
module ATMOSPHERE / Physics Cloud Microphysics
subroutine, public atmos_phy_mp_kessler_mixingratio(Qe, QTRC0)
Calculate mixing ratio of each category.
subroutine, public atmos_phy_mp_sn14_setup
Setup Cloud Microphysics.
subroutine, public atmos_phy_mp_tomita08_config(MP_TYPE, QA, QS)
Config.
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:61
module ATMOSPHERE / Physics Cloud Microphysics
subroutine, public atmos_phy_mp_kessler_cloudfraction(cldfrac, QTRC, mask_criterion)
Calculate Cloud Fraction.
module STDIO
Definition: scale_stdio.F90:12
character(len=h_short), dimension(qa_mp), target, public atmos_phy_mp_tomita08_unit
character(len=h_short), dimension(:), allocatable, target, public atmos_phy_mp_suzuki10_unit
character(len=h_short), dimension(qa_mp), target, public atmos_phy_mp_kessler_name
character(len=h_short), dimension(qa_mp), target, public atmos_phy_mp_dry_unit
subroutine, public atmos_phy_mp_dry(DENS, MOMZ, MOMX, MOMY, RHOT, QTRC, CCN, EVAPORATE, SFLX_rain, SFLX_snow)
Cloud Microphysics.
subroutine, public atmos_phy_mp_suzuki10_setup
Setup.
character(len=h_mid), dimension(qa_mp), target, public atmos_phy_mp_dry_desc
character(len=h_short), dimension(:), pointer, public atmos_phy_mp_name
real(rp), dimension(n_hyd), target, public atmos_phy_mp_sn14_dens
subroutine, public atmos_phy_mp_tomita08_mixingratio(Qe, QTRC0)
Calculate mixing ratio of each category.
subroutine, public atmos_phy_mp_kessler_effectiveradius(Re, QTRC0, DENS0, TEMP0)
Calculate Effective Radius.
character(len=h_short), dimension(qa_mp), target, public atmos_phy_mp_dry_name
subroutine, public atmos_phy_mp_kessler_config(MP_TYPE, QA, QS)
Setup.
character(len=h_short), dimension(qa_mp), target, public atmos_phy_mp_tomita08_name
subroutine, public atmos_phy_mp_sn14_cloudfraction(cldfrac, QTRC, mask_criterion)
Calculate Cloud Fraction.
subroutine, public atmos_phy_mp_suzuki10_mixingratio(Qe, QTRC0)
Calculate mixing ratio of each category.
module grid index
character(len=h_short), dimension(qa_mp), target, public atmos_phy_mp_kessler_unit
module TRACER
subroutine, public atmos_phy_mp_config(MP_TYPE)
Setup Cloud Microphysics.
character(len=h_short), dimension(qa_mp), target, public atmos_phy_mp_sn14_name
module ATMOSPHERE / Physics Cloud Microphysics
module ATMOSPHERE / Physics Cloud Microphysics
subroutine, public atmos_phy_mp_suzuki10_effectiveradius(Re, QTRC0, DENS0, TEMP0)
Calculate Effective Radius.
real(rp), dimension(n_hyd), target, public atmos_phy_mp_kessler_dens
subroutine, public atmos_phy_mp_dry_setup
Setup.
procedure(mr), pointer, public atmos_phy_mp_mixingratio
character(len=h_short), dimension(:), allocatable, target, public atmos_phy_mp_suzuki10_name
subroutine, public atmos_phy_mp_kessler_setup
Setup.
subroutine, public atmos_phy_mp_suzuki10_cloudfraction(cldfrac, QTRC, mask_criterion)
Calculate Cloud Fraction.
real(rp), dimension(n_hyd), target, public atmos_phy_mp_suzuki10_dens
module PROCESS
character(len=h_short), dimension(:), pointer, public atmos_phy_mp_unit
procedure(er), pointer, public atmos_phy_mp_effectiveradius
procedure(su), pointer, public atmos_phy_mp_setup
subroutine, public atmos_phy_mp_suzuki10(DENS, MOMZ, MOMX, MOMY, RHOT, QTRC, CCN, EVAPORATE, SFLX_rain, SFLX_snow)
Cloud Microphysics.
real(rp), dimension(n_hyd), target, public atmos_phy_mp_tomita08_dens
subroutine, public atmos_phy_mp_sdm_setup
Setup Cloud Microphysics.
subroutine, public atmos_phy_mp_sdm(DENS, MOMZ, MOMX, MOMY, RHOT, QTRC, CCN, EVAPORATE, SFLX_rain, SFLX_snow)
Cloud Microphysics.
subroutine, public atmos_phy_mp_sn14_config(MP_TYPE, QA, QS)
Configure.
real(rp), dimension(n_hyd), target, public atmos_phy_mp_dry_dens
subroutine, public atmos_phy_mp_sdm_config(MP_TYPE, QA, QS)
Confif.
subroutine, public atmos_phy_mp_tomita08(DENS, MOMZ, MOMX, MOMY, RHOT, QTRC, CCN, EVAPORATE, SFLX_rain, SFLX_snow)
Cloud Microphysics.
subroutine, public atmos_phy_mp_sn14_mixingratio(Qe, QTRC0)
Calculate mixing ratio of each category.
subroutine, public atmos_phy_mp_sn14(DENS, MOMZ, MOMX, MOMY, RHOT, QTRC, CCN, EVAPORATE, SFLX_rain, SFLX_snow)
Cloud Microphysics.
subroutine, public atmos_phy_mp_dry_config(MP_TYPE, QA, QS)
Configure.
character(len=h_mid), dimension(:), pointer, public atmos_phy_mp_desc
module profiler
Definition: scale_prof.F90:10
character(len=h_short), dimension(qa_mp), target, public atmos_phy_mp_sdm_unit
character(len=h_short), dimension(qa_mp), target, public atmos_phy_mp_sdm_name
character(len=h_short), dimension(qa_mp), target, public atmos_phy_mp_sn14_unit
module ATMOSPHERE / Physics Cloud Microphysics
character(len=h_mid), dimension(qa_mp), target, public atmos_phy_mp_kessler_desc
subroutine, public atmos_phy_mp_tomita08_cloudfraction(cldfrac, QTRC, mask_criterion)
Calculate Cloud Fraction.
module PRECISION
subroutine, public atmos_phy_mp_tomita08_setup
Setup.
character(len=h_mid), dimension(:), allocatable, target, public atmos_phy_mp_suzuki10_desc
character(len=h_mid), dimension(qa_mp), target, public atmos_phy_mp_sdm_desc
character(len=h_mid), dimension(qa_mp), target, public atmos_phy_mp_tomita08_desc
real(rp), dimension(n_hyd), target, public atmos_phy_mp_sdm_dens
integer, public io_fid_log
Log file ID.
Definition: scale_stdio.F90:56
integer, parameter, public n_hyd
module ATMOSPHERE / Physics Cloud Microphysics
module Spectran Bin Microphysics
real(rp), dimension(:), pointer, public atmos_phy_mp_dens
subroutine, public atmos_phy_mp_sdm_effectiveradius(Re, QTRC0, DENS0, TEMP0)
Calculate Effective Radius.
subroutine, public atmos_phy_mp_tomita08_effectiveradius(Re, QTRC0, DENS0, TEMP0)
Calculate Effective Radius.
subroutine, public atmos_phy_mp_dry_mixingratio(Qe, QTRC0)
Calculate mixing ratio of each category.
subroutine, public atmos_phy_mp_dry_effectiveradius(Re, QTRC0, DENS0, TEMP0)
Calculate Effective Radius.