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_setup
35 
36  !-----------------------------------------------------------------------------
37  !
38  !++ Public parameters & variables
39  !
40  !-----------------------------------------------------------------------------
41  !
42  !++ Private procedure
43  !
44  !-----------------------------------------------------------------------------
45  !
46  !++ Private parameters & variables
47  !
48  abstract interface
49  subroutine mp( &
50  DENS, &
51  MOMZ, &
52  MOMX, &
53  MOMY, &
54  RHOT, &
55  QTRC, &
56  CCN, &
57  EVAPORATE, &
58  SFLX_rain, &
59  SFLX_snow )
60  use scale_precision
62  use scale_tracer
63  real(RP), intent(inout) :: dens(ka,ia,ja)
64  real(RP), intent(inout) :: momz(ka,ia,ja)
65  real(RP), intent(inout) :: momx(ka,ia,ja)
66  real(RP), intent(inout) :: momy(ka,ia,ja)
67  real(RP), intent(inout) :: rhot(ka,ia,ja)
68  real(RP), intent(inout) :: qtrc(ka,ia,ja,qa)
69  real(RP), intent(in) :: ccn(ka,ia,ja)
70  real(RP), intent(out) :: evaporate(ka,ia,ja)
71  real(RP), intent(out) :: sflx_rain(ia,ja)
72  real(RP), intent(out) :: sflx_snow(ia,ja)
73  end subroutine mp
74  subroutine cf( &
75  cldfrac, &
76  QTRC, &
77  mask_criterion )
78  use scale_precision
80  use scale_tracer
81  real(RP), intent(out) :: cldfrac(ka,ia,ja)
82  real(RP), intent(in) :: qtrc (ka,ia,ja,qa)
83  real(RP), intent(in) :: mask_criterion
84  end subroutine cf
85  subroutine er( &
86  Re, &
87  QTRC0, &
88  DENS0, &
89  TEMP0 )
90  use scale_precision
92  use scale_tracer
93  real(RP), intent(out) :: re (ka,ia,ja,mp_qa) ! effective radius [cm]
94  real(RP), intent(in) :: qtrc0(ka,ia,ja,qa) ! tracer mass concentration [kg/kg]
95  real(RP), intent(in) :: dens0(ka,ia,ja) ! density [kg/m3]
96  real(RP), intent(in) :: temp0(ka,ia,ja) ! temperature [K]
97  end subroutine er
98  subroutine mr( &
99  Qe, &
100  QTRC0 )
101  use scale_precision
102  use scale_grid_index
103  use scale_tracer
104  real(RP), intent(out) :: qe (ka,ia,ja,mp_qa) ! mixing ratio of each cateory [kg/kg]
105  real(RP), intent(in) :: qtrc0(ka,ia,ja,qa) ! tracer mass concentration [kg/kg]
106  end subroutine mr
107  end interface
108  procedure(mp), pointer :: atmos_phy_mp => null()
109  procedure(cf), pointer :: atmos_phy_mp_cloudfraction => null()
110  procedure(er), pointer :: atmos_phy_mp_effectiveradius => null()
111  procedure(mr), pointer :: atmos_phy_mp_mixingratio => null()
112  public :: atmos_phy_mp
115  public :: atmos_phy_mp_mixingratio
116 
117  real(RP), pointer, public :: atmos_phy_mp_dens(:)
118 
119  !-----------------------------------------------------------------------------
120 contains
121  !-----------------------------------------------------------------------------
123  !-----------------------------------------------------------------------------
124  subroutine atmos_phy_mp_setup( MP_TYPE )
125  use scale_process, only: &
127 #define EXTM(pre, name, post) pre ## name ## post
128 #define NAME(pre, name, post) EXTM(pre, name, post)
129 #ifdef MP
130  use name(scale_atmos_phy_mp_, mp,), only: &
131  name(atmos_phy_mp_, mp, _setup), &
132  name(atmos_phy_mp_, mp,), &
133  name(atmos_phy_mp_, mp, _cloudfraction), &
134  name(atmos_phy_mp_, mp, _effectiveradius), &
135  name(atmos_phy_mp_, mp, _mixingratio), &
136  name(atmos_phy_mp_, mp, _dens) => atmos_phy_mp_dens
137 #else
138  use scale_atmos_phy_mp_dry, only: &
144  atmos_phy_mp_dry_dens => atmos_phy_mp_dens
145  use scale_atmos_phy_mp_kessler, only: &
151  atmos_phy_mp_kessler_dens => atmos_phy_mp_dens
152  use scale_atmos_phy_mp_tomita08, only: &
158  atmos_phy_mp_tomita08_dens => atmos_phy_mp_dens
159  use scale_atmos_phy_mp_sn14, only: &
165  atmos_phy_mp_sn14_dens => atmos_phy_mp_dens
166  use scale_atmos_phy_mp_suzuki10, only: &
172  atmos_phy_mp_suzuki10_dens => atmos_phy_mp_dens
173  use scale_atmos_phy_mp_sdm, only: &
179  atmos_phy_mp_sdm_dens => atmos_phy_mp_dens
180 #endif
181  implicit none
182 
183  character(len=*), intent(in) :: MP_TYPE
184  !---------------------------------------------------------------------------
185 
186 #ifdef MP
187  call name(atmos_phy_mp_, mp, _setup)( mp_type )
188  atmos_phy_mp => name(atmos_phy_mp_, mp,)
189  atmos_phy_mp_cloudfraction => name(atmos_phy_mp_, mp, _cloudfraction)
190  atmos_phy_mp_effectiveradius => name(atmos_phy_mp_, mp, _effectiveradius)
191  atmos_phy_mp_mixingratio => name(atmos_phy_mp_, mp, _mixingratio)
192  atmos_phy_mp_dens => name(atmos_phy_mp_, mp, _dens)
193 #else
194  select case ( mp_type )
195  case ( 'DRY' )
196  call atmos_phy_mp_dry_setup( mp_type )
201  atmos_phy_mp_dens => atmos_phy_mp_dry_dens
202  case ( 'KESSLER' )
203  call atmos_phy_mp_kessler_setup( mp_type )
208  atmos_phy_mp_dens => atmos_phy_mp_kessler_dens
209  case ( 'TOMITA08' )
210  call atmos_phy_mp_tomita08_setup( mp_type )
215  atmos_phy_mp_dens => atmos_phy_mp_tomita08_dens
216  case ( 'SN14' )
217  call atmos_phy_mp_sn14_setup( mp_type )
222  atmos_phy_mp_dens => atmos_phy_mp_sn14_dens
223  case ( 'SUZUKI10' )
224  call atmos_phy_mp_suzuki10_setup( mp_type )
229  atmos_phy_mp_dens => atmos_phy_mp_suzuki10_dens
230  case ( 'SDM' )
231  call atmos_phy_mp_sdm_setup( mp_type )
236  atmos_phy_mp_dens => atmos_phy_mp_sdm_dens
237  end select
238 #endif
239 
240  return
241  end subroutine atmos_phy_mp_setup
242 
243 end module scale_atmos_phy_mp
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 prc_mpistop
Abort MPI.
module ATMOSPHERE / Physics Cloud Microphysics
subroutine, public atmos_phy_mp_kessler_mixingratio(Qe, QTRC0)
Calculate mixing ratio of each category.
module ATMOSPHERE / Physics Cloud Microphysics
subroutine, public atmos_phy_mp_suzuki10_setup(MP_TYPE)
Setup.
subroutine, public atmos_phy_mp_kessler_cloudfraction(cldfrac, QTRC, mask_criterion)
Calculate Cloud Fraction.
module STDIO
Definition: scale_stdio.F90:12
integer, public qa
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_dry_setup(MP_TYPE)
Setup.
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.
subroutine, public atmos_phy_mp_sn14_cloudfraction(cldfrac, QTRC, mask_criterion)
Calculate Cloud Fraction.
subroutine, public atmos_phy_mp_setup(MP_TYPE)
Setup Cloud Microphysics.
subroutine, public atmos_phy_mp_suzuki10_mixingratio(Qe, QTRC0)
Calculate mixing ratio of each category.
module grid index
module TRACER
integer, public ia
of x whole cells (local, with HALO)
module ATMOSPHERE / Physics Cloud Microphysics
module ATMOSPHERE / Physics Cloud Microphysics
integer, public mp_qa
subroutine, public atmos_phy_mp_suzuki10_effectiveradius(Re, QTRC0, DENS0, TEMP0)
Calculate Effective Radius.
integer, public ka
of z whole cells (local, with HALO)
procedure(mr), pointer, public atmos_phy_mp_mixingratio
subroutine, public atmos_phy_mp_suzuki10_cloudfraction(cldfrac, QTRC, mask_criterion)
Calculate Cloud Fraction.
subroutine, public atmos_phy_mp_kessler_setup(MP_TYPE)
Setup.
module PROCESS
procedure(er), pointer, public atmos_phy_mp_effectiveradius
subroutine, public atmos_phy_mp_suzuki10(DENS, MOMZ, MOMX, MOMY, RHOT, QTRC, CCN, EVAPORATE, SFLX_rain, SFLX_snow)
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_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.
module profiler
Definition: scale_prof.F90:10
module ATMOSPHERE / Physics Cloud Microphysics
subroutine, public atmos_phy_mp_tomita08_cloudfraction(cldfrac, QTRC, mask_criterion)
Calculate Cloud Fraction.
module PRECISION
subroutine, public atmos_phy_mp_sdm_setup(MP_TYPE)
Setup Cloud Microphysics.
subroutine, public atmos_phy_mp_sn14_setup(MP_TYPE)
Setup Cloud Microphysics.
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_setup(MP_TYPE)
Setup.
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.
integer, public ja
of y whole cells (local, with HALO)
subroutine, public atmos_phy_mp_dry_effectiveradius(Re, QTRC0, DENS0, TEMP0)
Calculate Effective Radius.