SCALE-RM
Functions/Subroutines | Variables
scale_atmos_phy_mp_kessler Module Reference

module ATMOSPHERE / Physics Cloud Microphysics More...

Functions/Subroutines

subroutine, public atmos_phy_mp_kessler_setup (MP_TYPE)
 Setup. More...
 
subroutine, public atmos_phy_mp_kessler (DENS, MOMZ, MOMX, MOMY, RHOT, QTRC, CCN, EVAPORATE, SFLX_rain, SFLX_snow)
 Cloud Microphysics. More...
 
subroutine, public atmos_phy_mp_kessler_cloudfraction (cldfrac, QTRC)
 Calculate Cloud Fraction. More...
 
subroutine, public atmos_phy_mp_kessler_effectiveradius (Re, QTRC0, DENS0, TEMP0)
 Calculate Effective Radius. More...
 
subroutine, public atmos_phy_mp_kessler_mixingratio (Qe, QTRC0)
 Calculate mixing ratio of each category. More...
 

Variables

real(rp), dimension(mp_qa), target, public atmos_phy_mp_dens
 

Detailed Description

module ATMOSPHERE / Physics Cloud Microphysics

Description
Cloud Microphysics by Kessler-type parametarization Reference: Kessler(1969) Klemp and Wilhelmson(1978)
Author
Team SCALE
History
  • 2012-01-14 (Y.Miyamoto) [new]
  • 2012-03-23 (H.Yashiro) [mod] Explicit index parameter inclusion
  • 2012-12-23 (H.Yashiro) [mod] Reconstruction
  • 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_NTMAX_SEDIMENTATION integer 1 number of time step for sedimentation
    MP_COUPLE_AEROSOL logical .false. Consider CCN effect ?

History Output
namedescriptionunitvariable
Vterm_QR terminal velocity of QR m/s vterm

Function/Subroutine Documentation

◆ atmos_phy_mp_kessler_setup()

subroutine, public scale_atmos_phy_mp_kessler::atmos_phy_mp_kessler_setup ( character(len=*), intent(in)  MP_TYPE)

Setup.

Definition at line 80 of file scale_atmos_phy_mp_kessler.F90.

References atmos_phy_mp_dens, scale_const::const_dwatr, scale_grid::grid_cdz, scale_stdio::io_fid_conf, scale_stdio::io_fid_log, scale_stdio::io_l, scale_stdio::io_lnml, scale_grid_index::ka, scale_process::prc_mpistop(), and scale_time::time_dtsec_atmos_phy_mp.

Referenced by scale_atmos_phy_mp::atmos_phy_mp_setup().

80  use scale_process, only: &
82  use scale_const, only: &
84  use scale_time, only: &
86  use scale_grid, only: &
87  cdz => grid_cdz
88  implicit none
89 
90  character(len=*), intent(in) :: mp_type
91 
92  namelist / param_atmos_phy_mp / &
93  mp_doprecipitation, &
94  mp_donegative_fixer, &
95  mp_ntmax_sedimentation, &
96  mp_couple_aerosol
97 
98  real(RP), parameter :: max_term_vel = 10.0_rp !-- terminal velocity for calculate dt of sedimentation
99  integer :: nstep_max
100  integer :: ierr
101  !---------------------------------------------------------------------------
102 
103  if( io_l ) write(io_fid_log,*)
104  if( io_l ) write(io_fid_log,*) '++++++ Module[Cloud Microphysics] / Categ[ATMOS PHYSICS] / Origin[SCALElib]'
105  if( io_l ) write(io_fid_log,*) '*** KESSLER-type 1-moment bulk 3 category'
106 
107  if ( mp_type /= 'KESSLER' ) then
108  write(*,*) 'xxx ATMOS_PHY_MP_TYPE is not KESSLER. Check!'
109  call prc_mpistop
110  endif
111 
112  if ( i_qv <= 0 &
113  .OR. i_qc <= 0 &
114  .OR. i_qr <= 0 ) then
115  write(*,*) 'xxx KESSLER needs QV, QC, QR tracer. Check!'
116  call prc_mpistop
117  endif
118 
119  allocate( factor_vterm(ka) )
120 
121  !--- read namelist
122  rewind(io_fid_conf)
123  read(io_fid_conf,nml=param_atmos_phy_mp,iostat=ierr)
124  if( ierr < 0 ) then !--- missing
125  if( io_l ) write(io_fid_log,*) '*** Not found namelist. Default used.'
126  elseif( ierr > 0 ) then !--- fatal error
127  write(*,*) 'xxx Not appropriate names in namelist PARAM_ATMOS_PHY_MP. Check!'
128  call prc_mpistop
129  endif
130  if( io_lnml ) write(io_fid_log,nml=param_atmos_phy_mp)
131 
132  if( mp_couple_aerosol ) then
133  write(*,*) 'xxx MP_aerosol_couple should be .false. for KESSLER type MP!'
134  call prc_mpistop
135  endif
136 
137  if ( io_l ) write(io_fid_log,*)
138  if ( io_l ) write(io_fid_log,*) '*** Enable negative fixer? : ', mp_donegative_fixer
139  if ( io_l ) write(io_fid_log,*) '*** Enable sedimentation (precipitation)? : ', mp_doprecipitation
140 
141  atmos_phy_mp_dens(i_mp_qc) = const_dwatr
142  atmos_phy_mp_dens(i_mp_qr) = const_dwatr
143 
144  nstep_max = int( ( time_dtsec_atmos_phy_mp * max_term_vel ) / minval( cdz ) )
145  mp_ntmax_sedimentation = max( mp_ntmax_sedimentation, nstep_max )
146 
147  mp_nstep_sedimentation = mp_ntmax_sedimentation
148  mp_rnstep_sedimentation = 1.0_rp / real(mp_ntmax_sedimentation,kind=rp)
149  mp_dtsec_sedimentation = time_dtsec_atmos_phy_mp * mp_rnstep_sedimentation
150 
151  if ( io_l ) write(io_fid_log,*)
152  if ( io_l ) write(io_fid_log,*) '*** Timestep of sedimentation is divided into : ', mp_ntmax_sedimentation, ' step'
153  if ( io_l ) write(io_fid_log,*) '*** DT of sedimentation is : ', mp_dtsec_sedimentation, '[s]'
154 
155  return
subroutine, public prc_mpistop
Abort MPI.
real(rp), parameter, public const_dwatr
density of water [kg/m3]
Definition: scale_const.F90:87
real(dp), public time_dtsec_atmos_phy_mp
time interval of physics(microphysics) [sec]
Definition: scale_time.F90:41
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:59
integer, public ka
of z whole cells (local, with HALO)
real(rp), dimension(mp_qa), target, public atmos_phy_mp_dens
module TIME
Definition: scale_time.F90:15
module PROCESS
module CONSTANT
Definition: scale_const.F90:14
module GRID (cartesian)
logical, public io_lnml
output log or not? (for namelist, this process)
Definition: scale_stdio.F90:60
real(rp), dimension(:), allocatable, public grid_cdz
z-length of control volume [m]
integer, public io_fid_conf
Config file ID.
Definition: scale_stdio.F90:55
integer, public io_fid_log
Log file ID.
Definition: scale_stdio.F90:56
integer, parameter, public rp
Here is the call graph for this function:
Here is the caller graph for this function:

◆ atmos_phy_mp_kessler()

subroutine, public scale_atmos_phy_mp_kessler::atmos_phy_mp_kessler ( 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,qad), 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 171 of file scale_atmos_phy_mp_kessler.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_atmos_refstate::atmos_refstate_dens, scale_comm::comm_horizontal_mean(), scale_const::const_lhv, 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::ke, scale_grid_index::ks, scale_prof::prof_rapend(), scale_prof::prof_rapstart(), scale_tracer::qa, and scale_time::time_dtsec_atmos_phy_mp.

Referenced by scale_atmos_phy_mp::atmos_phy_mp_setup().

171  use scale_grid_index
172  use scale_comm, only: &
174  use scale_history, only: &
175  hist_in
176  use scale_tracer, only: &
177  qad => qa
178  use scale_atmos_thermodyn, only: &
179  thermodyn_rhoe => atmos_thermodyn_rhoe, &
180  thermodyn_rhot => atmos_thermodyn_rhot, &
181  thermodyn_temp_pres_e => atmos_thermodyn_temp_pres_e
182  use scale_atmos_phy_mp_common, only: &
183  mp_negative_fixer => atmos_phy_mp_negative_fixer, &
184  mp_precipitation => atmos_phy_mp_precipitation, &
185  mp_saturation_adjustment => atmos_phy_mp_saturation_adjustment
186  implicit none
187 
188  real(RP), intent(inout) :: dens(ka,ia,ja)
189  real(RP), intent(inout) :: momz(ka,ia,ja)
190  real(RP), intent(inout) :: momx(ka,ia,ja)
191  real(RP), intent(inout) :: momy(ka,ia,ja)
192  real(RP), intent(inout) :: rhot(ka,ia,ja)
193  real(RP), intent(inout) :: qtrc(ka,ia,ja,qad)
194  real(RP), intent(in) :: ccn(ka,ia,ja)
195  real(RP), intent(out) :: evaporate(ka,ia,ja) !--- number of evaporated cloud [/m3]
196  real(RP), intent(out) :: sflx_rain(ia,ja)
197  real(RP), intent(out) :: sflx_snow(ia,ja)
198 
199  real(RP) :: rhoe_t(ka,ia,ja)
200  real(RP) :: qtrc_t(ka,ia,ja,qad)
201  real(RP) :: rhoe (ka,ia,ja)
202  real(RP) :: temp (ka,ia,ja)
203  real(RP) :: pres (ka,ia,ja)
204 
205  real(RP) :: vterm (ka,ia,ja,qad) ! terminal velocity of each tracer [m/s]
206  real(RP) :: flx_rain(ka,ia,ja)
207  real(RP) :: flx_snow(ka,ia,ja)
208  real(RP) :: wflux_rain(ka,ia,ja)
209  real(RP) :: wflux_snow(ka,ia,ja)
210  integer :: step
211 
212  real(RP) :: rho_prof(ka) ! averaged profile of rho
213  integer :: k, i, j
214  !---------------------------------------------------------------------------
215 
216  if( io_l ) write(io_fid_log,*) '*** Physics step: Cloud microphysics(kessler)'
217 
218  if ( first ) then
219  ! Calculate collection factor for terminal velocity of QR
220  call comm_horizontal_mean( rho_prof(:), dens(:,:,:) )
221  rho_prof(:) = rho_prof(:) * 1.e-3_rp ! [kg/m3]->[g/cc]
222 
223  do k = ks, ke
224  factor_vterm(k) = sqrt( rho_prof(ks)/rho_prof(k) )
225  enddo
226  first = .false.
227  endif
228 
229  evaporate(:,:,:) = 0.0_rp
230 
231  if ( mp_donegative_fixer ) then
232  call mp_negative_fixer( dens(:,:,:), & ! [INOUT]
233  rhot(:,:,:), & ! [INOUT]
234  qtrc(:,:,:,:) ) ! [INOUT]
235  endif
236 
237  call thermodyn_rhoe( rhoe(:,:,:), & ! [OUT]
238  rhot(:,:,:), & ! [IN]
239  qtrc(:,:,:,:) ) ! [IN]
240 
241  !##### MP Main #####
242  rhoe_t(:,:,:) = 0.0_rp
243  qtrc_t(:,:,:,:) = 0.0_rp
244 
245  call mp_kessler( rhoe_t(:,:,:), & ! [INOUT]
246  qtrc_t(:,:,:,:), & ! [INOUT]
247  rhoe(:,:,:), & ! [INOUT]
248  qtrc(:,:,:,:), & ! [INOUT]
249  dens(:,:,:) ) ! [IN]
250 
251  if ( mp_doprecipitation ) then
252 
253  flx_rain(:,:,:) = 0.0_rp
254  flx_snow(:,:,:) = 0.0_rp
255 
256  vterm(:,:,:,:) = 0.0_rp
257 
258  do step = 1, mp_nstep_sedimentation
259 
260  call mp_kessler_vterm( vterm(:,:,:,:), & ! [OUT]
261  dens(:,:,:), & ! [IN]
262  qtrc(:,:,:,:) ) ! [IN]
263 
264  call thermodyn_temp_pres_e( temp(:,:,:), & ! [OUT]
265  pres(:,:,:), & ! [OUT]
266  dens(:,:,:), & ! [IN]
267  rhoe(:,:,:), & ! [IN]
268  qtrc(:,:,:,:) ) ! [IN]
269 
270  call mp_precipitation( wflux_rain(:,:,:), & ! [OUT]
271  wflux_snow(:,:,:), & ! [OUT]
272  dens(:,:,:), & ! [INOUT]
273  momz(:,:,:), & ! [INOUT]
274  momx(:,:,:), & ! [INOUT]
275  momy(:,:,:), & ! [INOUT]
276  rhoe(:,:,:), & ! [INOUT]
277  qtrc(:,:,:,:), & ! [INOUT]
278  vterm(:,:,:,:), & ! [IN]
279  temp(:,:,:), & ! [IN]
280  mp_dtsec_sedimentation ) ! [IN]
281 
282  do j = js, je
283  do i = is, ie
284  do k = ks-1, ke
285  flx_rain(k,i,j) = flx_rain(k,i,j) + wflux_rain(k,i,j) * mp_rnstep_sedimentation
286  flx_snow(k,i,j) = flx_snow(k,i,j) + wflux_snow(k,i,j) * mp_rnstep_sedimentation
287  enddo
288  enddo
289  enddo
290 
291  enddo
292 
293  else
294  vterm(:,:,:,:) = 0.0_rp
295 
296  flx_rain(:,:,:) = 0.0_rp
297  flx_snow(:,:,:) = 0.0_rp
298  endif
299 
300  call mp_saturation_adjustment( rhoe_t(:,:,:), & ! [INOUT]
301  qtrc_t(:,:,:,:), & ! [INOUT]
302  rhoe(:,:,:), & ! [INOUT]
303  qtrc(:,:,:,:), & ! [INOUT]
304  dens(:,:,:), & ! [IN]
305  flag_liquid = .true. ) ! [IN]
306 
307  call hist_in( vterm(:,:,:,i_qr), 'Vterm_QR', 'terminal velocity of QR', 'm/s' )
308 
309  !##### END MP Main #####
310 
311  call thermodyn_rhot( rhot(:,:,:), & ! [OUT]
312  rhoe(:,:,:), & ! [IN]
313  qtrc(:,:,:,:) ) ! [IN]
314 
315  if ( mp_donegative_fixer ) then
316  call mp_negative_fixer( dens(:,:,:), & ! [INOUT]
317  rhot(:,:,:), & ! [INOUT]
318  qtrc(:,:,:,:) ) ! [INOUT]
319  endif
320 
321  sflx_rain(:,:) = flx_rain(ks-1,:,:)
322  sflx_snow(:,:) = flx_snow(ks-1,:,:)
323 
324  return
integer, public is
start point of inner domain: x, local
integer, public je
end point of inner domain: y, local
subroutine, public atmos_phy_mp_precipitation(flux_rain, flux_snow, DENS, MOMZ, MOMX, MOMY, RHOE, QTRC, vterm, temp, dt)
precipitation transport
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:59
integer, public ke
end point of inner domain: z, local
integer, public qa
module ATMOSPHERE / Physics Cloud Microphysics - Common
module grid index
module TRACER
integer, public ia
of x whole cells (local, with HALO)
subroutine, public comm_horizontal_mean(varmean, var)
calculate horizontal mean (global total with communication)
Definition: scale_comm.F90:516
integer, public ka
of z whole cells (local, with HALO)
module COMMUNICATION
Definition: scale_comm.F90:23
integer, public js
start point of inner domain: y, local
integer, public ks
start point of inner domain: z, local
subroutine, public atmos_phy_mp_negative_fixer(DENS, RHOT, QTRC)
Negative fixer.
integer, public ie
end point of inner domain: x, local
module ATMOSPHERE / Thermodynamics
module HISTORY
subroutine, public atmos_phy_mp_saturation_adjustment(RHOE_t, QTRC_t, RHOE0, QTRC0, DENS0, flag_liquid)
Saturation adjustment.
integer, public io_fid_log
Log file ID.
Definition: scale_stdio.F90:56
integer, public ja
of y whole cells (local, with HALO)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ atmos_phy_mp_kessler_cloudfraction()

subroutine, public scale_atmos_phy_mp_kessler::atmos_phy_mp_kessler_cloudfraction ( real(rp), dimension(ka,ia,ja), intent(out)  cldfrac,
real(rp), dimension (ka,ia,ja,qad), intent(in)  QTRC 
)

Calculate Cloud Fraction.

Definition at line 516 of file scale_atmos_phy_mp_kessler.F90.

References scale_const::const_eps, 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_setup().

516  use scale_grid_index
517  use scale_tracer, only: &
518  qad => qa
519  use scale_const, only: &
520  eps => const_eps
521  implicit none
522 
523  real(RP), intent(out) :: cldfrac(ka,ia,ja)
524  real(RP), intent(in) :: qtrc (ka,ia,ja,qad)
525 
526  real(RP) :: qhydro
527  integer :: k, i, j, iq
528  !---------------------------------------------------------------------------
529 
530  do j = js, je
531  do i = is, ie
532  do k = ks, ke
533  qhydro = 0.0_rp
534  do iq = 1, mp_qa
535  qhydro = qhydro + qtrc(k,i,j,i_mp2all(iq))
536  enddo
537  cldfrac(k,i,j) = 0.5_rp + sign(0.5_rp,qhydro-eps)
538  enddo
539  enddo
540  enddo
541 
542  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
integer, public qa
module grid index
module TRACER
integer, public ia
of x whole cells (local, with HALO)
integer, public ka
of z whole cells (local, with HALO)
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_eps
small number
Definition: scale_const.F90:36
integer, public ja
of y whole cells (local, with HALO)
Here is the caller graph for this function:

◆ atmos_phy_mp_kessler_effectiveradius()

subroutine, public scale_atmos_phy_mp_kessler::atmos_phy_mp_kessler_effectiveradius ( real(rp), dimension (ka,ia,ja,mp_qad), intent(out)  Re,
real(rp), dimension(ka,ia,ja,qad), 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 552 of file scale_atmos_phy_mp_kessler.F90.

References scale_tracer::mp_qa, and scale_tracer::qa.

Referenced by scale_atmos_phy_mp::atmos_phy_mp_setup().

552  use scale_grid_index
553  use scale_tracer, only: &
554  qad => qa, &
555  mp_qad => mp_qa
556  implicit none
557 
558  real(RP), intent(out) :: re (ka,ia,ja,mp_qad) ! effective radius [cm]
559  real(RP), intent(in) :: qtrc0(ka,ia,ja,qad) ! tracer mass concentration [kg/kg]
560  real(RP), intent(in) :: dens0(ka,ia,ja) ! density [kg/m3]
561  real(RP), intent(in) :: temp0(ka,ia,ja) ! temperature [K]
562 
563  real(RP), parameter :: um2cm = 100.0_rp
564  !---------------------------------------------------------------------------
565 
566  re(:,:,:,i_mp_qc) = 8.e-6_rp * um2cm
567  re(:,:,:,i_mp_qr) = 100.e-6_rp * um2cm
568 
569  return
integer, public qa
module grid index
module TRACER
integer, public ia
of x whole cells (local, with HALO)
integer, public mp_qa
integer, public ka
of z whole cells (local, with HALO)
integer, public ja
of y whole cells (local, with HALO)
Here is the caller graph for this function:

◆ atmos_phy_mp_kessler_mixingratio()

subroutine, public scale_atmos_phy_mp_kessler::atmos_phy_mp_kessler_mixingratio ( real(rp), dimension (ka,ia,ja,mp_qad), intent(out)  Qe,
real(rp), dimension(ka,ia,ja,qad), intent(in)  QTRC0 
)

Calculate mixing ratio of each category.

Definition at line 577 of file scale_atmos_phy_mp_kessler.F90.

References scale_tracer::mp_qa, and scale_tracer::qa.

Referenced by scale_atmos_phy_mp::atmos_phy_mp_setup().

577  use scale_grid_index
578  use scale_tracer, only: &
579  qad => qa, &
580  mp_qad => mp_qa
581  implicit none
582 
583  real(RP), intent(out) :: qe (ka,ia,ja,mp_qad) ! mixing ratio of each cateory [kg/kg]
584  real(RP), intent(in) :: qtrc0(ka,ia,ja,qad) ! tracer mass concentration [kg/kg]
585  !---------------------------------------------------------------------------
586 
587  qe(:,:,:,i_mp_qc) = qtrc0(:,:,:,i_qc)
588  qe(:,:,:,i_mp_qr) = qtrc0(:,:,:,i_qr)
589 
590  return
integer, public qa
module grid index
module TRACER
integer, public ia
of x whole cells (local, with HALO)
integer, public mp_qa
integer, public ka
of z whole cells (local, with HALO)
integer, public ja
of y whole cells (local, with HALO)
Here is the caller graph for this function:

Variable Documentation

◆ atmos_phy_mp_dens

real(rp), dimension(mp_qa), target, public scale_atmos_phy_mp_kessler::atmos_phy_mp_dens

Definition at line 49 of file scale_atmos_phy_mp_kessler.F90.

Referenced by atmos_phy_mp_kessler_setup().

49  real(RP), public, target :: atmos_phy_mp_dens(mp_qa) ! hydrometeor density [kg/m3]=[g/L]
integer, public mp_qa
real(rp), dimension(mp_qa), target, public atmos_phy_mp_dens