SCALE-RM
scale_atmos_phy_mp_kessler.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
18 !-------------------------------------------------------------------------------
19 #include "inc_openmp.h"
21  !-----------------------------------------------------------------------------
22  !
23  !++ used modules
24  !
25  use scale_precision
26  use scale_stdio
27  use scale_prof
29 
31  use scale_tracer, only: qa
32  !-----------------------------------------------------------------------------
33  implicit none
34  private
35  !-----------------------------------------------------------------------------
36  !
37  !++ Public procedure
38  !
40  public :: atmos_phy_mp_kessler
44 
45  !-----------------------------------------------------------------------------
46  !
47  !++ Public parameters & variables
48  !
49  real(RP), public, target :: atmos_phy_mp_dens(mp_qa) ! hydrometeor density [kg/m3]=[g/L]
50 
51  !-----------------------------------------------------------------------------
52  !
53  !++ Private procedure
54  !
55  private :: mp_kessler
56  private :: mp_kessler_vterm
57 
58  !-----------------------------------------------------------------------------
59  !
60  !++ Private parameters & variables
61  !
62  logical, private :: mp_donegative_fixer = .true. ! apply negative fixer?
63  logical, private :: mp_doprecipitation = .true. ! apply sedimentation (precipitation)?
64  logical, private :: mp_couple_aerosol = .false. ! Consider CCN effect ?
65 
66  real(RP), private, allocatable :: factor_vterm(:) ! collection factor for terminal velocity of QR
67 
68  logical, private :: first = .true.
69 
70  integer, private, save :: mp_ntmax_sedimentation = 1 ! number of time step for sedimentation
71  integer, private, save :: mp_nstep_sedimentation
72  real(RP), private, save :: mp_rnstep_sedimentation
73  real(DP), private, save :: mp_dtsec_sedimentation
74 
75  !-----------------------------------------------------------------------------
76 contains
77  !-----------------------------------------------------------------------------
79  subroutine atmos_phy_mp_kessler_setup( MP_TYPE )
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
156  end subroutine atmos_phy_mp_kessler_setup
157 
158  !-----------------------------------------------------------------------------
160  subroutine atmos_phy_mp_kessler( &
161  DENS, &
162  MOMZ, &
163  MOMX, &
164  MOMY, &
165  RHOT, &
166  QTRC, &
167  CCN, &
168  EVAPORATE, &
169  SFLX_rain, &
170  SFLX_snow )
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
325  end subroutine atmos_phy_mp_kessler
326 
327  !-----------------------------------------------------------------------------
329  subroutine mp_kessler( &
330  RHOE_t, &
331  QTRC_t, &
332  RHOE0, &
333  QTRC0, &
334  DENS0 )
335  use scale_const, only: &
336  lhv => const_lhv
337  use scale_time, only: &
339  use scale_atmos_thermodyn, only: &
340  thermodyn_temp_pres_e => atmos_thermodyn_temp_pres_e
341  use scale_atmos_saturation, only: &
342  saturation_dens2qsat_liq => atmos_saturation_dens2qsat_liq
343  implicit none
344 
345  real(RP), intent(inout) :: RHOE_t(ka,ia,ja) ! tendency rhoe [J/m3/s]
346  real(RP), intent(inout) :: QTRC_t(ka,ia,ja,qa) ! tendency tracer [kg/kg/s]
347  real(RP), intent(inout) :: RHOE0 (ka,ia,ja) ! density * internal energy [J/m3]
348  real(RP), intent(inout) :: QTRC0 (ka,ia,ja,qa) ! mass concentration [kg/kg]
349  real(RP), intent(in) :: DENS0 (ka,ia,ja) ! density [kg/m3]
350 
351  ! working
352  real(RP) :: QSAT (ka,ia,ja) ! saturated water vapor [kg/kg]
353  real(RP) :: TEMP0(ka,ia,ja) ! temperature [K]
354  real(RP) :: PRES0(ka,ia,ja) ! pressure [Pa]
355 
356  real(RP) :: dens (ka)
357  real(RP) :: rhoe (ka)
358  real(RP) :: temp (ka)
359  real(RP) :: pres (ka)
360  real(RP) :: qv (ka)
361  real(RP) :: qc (ka)
362  real(RP) :: qr (ka)
363  real(RP) :: qsatl(ka)
364  real(RP) :: Sliq (ka)
365 
366  ! tendency
367  real(RP) :: dq_evap(ka) ! tendency q (evaporation)
368  real(RP) :: dq_auto(ka) ! tendency q (autoconversion)
369  real(RP) :: dq_accr(ka) ! tendency q (accretion)
370 
371  real(RP) :: dqv, dqc, dqr
372  real(RP) :: vent_factor
373  real(RP) :: rdt
374 
375  integer :: k, i, j
376  !---------------------------------------------------------------------------
377 
378  call prof_rapstart('MP_kessler', 3)
379 
380  rdt = 1.0_rp / dt
381 
382  call thermodyn_temp_pres_e( temp0(:,:,:), & ! [OUT]
383  pres0(:,:,:), & ! [OUT]
384  dens0(:,:,:), & ! [IN]
385  rhoe0(:,:,:), & ! [IN]
386  qtrc0(:,:,:,:) ) ! [IN]
387 
388  call saturation_dens2qsat_liq( qsat(:,:,:), & ! [OUT]
389  temp0(:,:,:), & ! [IN]
390  dens0(:,:,:) ) ! [IN]
391 
392  !$omp parallel do private(i,j,k,dens,rhoe,temp,pres,qv,qc,qr,qsatl,Sliq,dq_evap,dq_auto,dq_accr,vent_factor,dqc,dqr,dqv) OMP_SCHEDULE_ collapse(2)
393  do j = js, je
394  do i = is, ie
395 
396  ! store to work
397  dens(ks:ke) = dens0(ks:ke,i,j)
398  rhoe(ks:ke) = rhoe0(ks:ke,i,j)
399  temp(ks:ke) = temp0(ks:ke,i,j)
400  pres(ks:ke) = pres0(ks:ke,i,j)
401  qv(ks:ke) = qtrc0(ks:ke,i,j,i_qv)
402  qc(ks:ke) = qtrc0(ks:ke,i,j,i_qc)
403  qr(ks:ke) = qtrc0(ks:ke,i,j,i_qr)
404  qsatl(ks:ke) = qsat(ks:ke,i,j)
405  sliq(ks:ke) = qv(ks:ke) / qsatl(ks:ke)
406 
407  ! Auto-conversion (QC->QR)
408  do k = ks, ke
409  dq_auto(k) = 1.e-3_rp * max( qc(k)-1.e-3_rp, 0.0_rp )
410  enddo
411 
412  ! Accretion (QC->QR)
413  do k = ks, ke
414  dq_accr(k) = 2.2_rp * qc(k) * qr(k)**0.875_rp
415  enddo
416 
417  ! Evaporation (QR->QV)
418  do k = ks, ke
419  vent_factor = 1.6_rp + 124.9_rp * ( dens(k)*qr(k) )**0.2046_rp
420 
421  dq_evap(k) = ( 1.0_rp-min(sliq(k),1.0_rp) ) / dens(k) * vent_factor &
422  * ( dens(k)*qr(k) )**0.525_rp / ( 5.4e5_rp + 2.55e8_rp / ( pres(k)*qsatl(k) ) )
423  enddo
424 
425  ! limiter
426  do k = ks, ke
427  dqc = (-dq_auto(k)-dq_accr(k) )
428 
429  qtrc_t(k,i,j,i_qc) = qtrc_t(k,i,j,i_qc) + max( dqc, -qc(k)*rdt )
430  enddo
431 
432  do k = ks, ke
433  dqr = ( dq_auto(k)+dq_accr(k)-dq_evap(k) )
434 
435  qtrc_t(k,i,j,i_qr) = qtrc_t(k,i,j,i_qr) + max( dqr, -qr(k)*rdt )
436  enddo
437 
438  do k = ks, ke
439  dqv = - ( qtrc_t(k,i,j,i_qc) &
440  + qtrc_t(k,i,j,i_qr) )
441 
442  qtrc_t(k,i,j,i_qv) = qtrc_t(k,i,j,i_qv) + max( dqv, -qv(k)*rdt )
443  enddo
444 
445  do k = ks, ke
446  rhoe_t(k,i,j) = rhoe_t(k,i,j) - dens(k) * ( lhv * qtrc_t(k,i,j,i_qv) )
447  enddo
448  enddo
449  enddo
450 
451  ! mass & energy update
452  !$omp parallel do private(i,j, k) OMP_SCHEDULE_ collapse(2)
453  do j = js, je
454  do i = is, ie
455  qtrc0(ks:ke,i,j,i_qv) = qtrc0(ks:ke,i,j,i_qv) + qtrc_t(ks:ke,i,j,i_qv) * dt
456  qtrc0(ks:ke,i,j,i_qc) = qtrc0(ks:ke,i,j,i_qc) + qtrc_t(ks:ke,i,j,i_qc) * dt
457  qtrc0(ks:ke,i,j,i_qr) = qtrc0(ks:ke,i,j,i_qr) + qtrc_t(ks:ke,i,j,i_qr) * dt
458 
459  rhoe0(ks:ke,i,j) = rhoe0(ks:ke,i,j) + rhoe_t(ks:ke,i,j) * dt
460  enddo
461  enddo
462 
463  call prof_rapend ('MP_kessler', 3)
464 
465  return
466  end subroutine mp_kessler
467 
468  !-----------------------------------------------------------------------------
470  subroutine mp_kessler_vterm( &
471  vterm, &
472  DENS0, &
473  QTRC0 )
474  use scale_atmos_refstate, only: &
475  refstate_dens => atmos_refstate_dens
476  implicit none
477 
478  real(RP), intent(out) :: vterm(ka,ia,ja,qa)
479  real(RP), intent(in) :: DENS0(ka,ia,ja)
480  real(RP), intent(in) :: QTRC0(ka,ia,ja,qa)
481 
482  real(RP) :: zerosw
483 
484  integer :: k, i, j
485  !---------------------------------------------------------------------------
486 
487  !$omp parallel do private(i,j,k,zerosw) OMP_SCHEDULE_ collapse(2)
488  do j = js, je
489  do i = is, ie
490  do k = ks, ke
491  vterm(k,i,j,i_qv) = 0.0_rp
492  vterm(k,i,j,i_qc) = 0.0_rp
493  enddo
494  enddo
495  enddo
496 
497  do j = js, je
498  do i = is, ie
499  do k = ks, ke
500  zerosw = 0.5_rp - sign(0.5_rp, qtrc0(k,i,j,i_qr) - 1.e-12_rp )
501  vterm(k,i,j,i_qr) = - 36.34_rp * ( dens0(k,i,j) * ( qtrc0(k,i,j,i_qr) + zerosw ) )**0.1364_rp &
502  * refstate_dens(ks,i,j) / refstate_dens(k,i,j) * ( 1.0_rp - zerosw )
503 ! * factor_vterm(k) * ( 1.0_RP - zerosw )
504  enddo
505  enddo
506  enddo
507 
508  return
509  end subroutine mp_kessler_vterm
510 
511  !-----------------------------------------------------------------------------
514  cldfrac, &
515  QTRC, &
516  mask_criterion )
518  use scale_tracer, only: &
519  qad => qa
520  implicit none
521 
522  real(RP), intent(out) :: cldfrac(ka,ia,ja)
523  real(RP), intent(in) :: QTRC (ka,ia,ja,qad)
524  real(RP), intent(in) :: mask_criterion
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-mask_criterion)
538  enddo
539  enddo
540  enddo
541 
542  return
544 
545  !-----------------------------------------------------------------------------
548  Re, &
549  QTRC0, &
550  DENS0, &
551  TEMP0 )
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
571 
572  !-----------------------------------------------------------------------------
575  Qe, &
576  QTRC0 )
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
591  end subroutine atmos_phy_mp_kessler_mixingratio
592 
integer, public is
start point of inner domain: x, local
subroutine, public atmos_phy_mp_kessler(DENS, MOMZ, MOMX, MOMY, RHOT, QTRC, CCN, EVAPORATE, SFLX_rain, SFLX_snow)
Cloud Microphysics.
integer, public je
end point of inner domain: y, local
module ATMOSPHERE / Saturation adjustment
subroutine, public prc_mpistop
Abort MPI.
subroutine, public atmos_phy_mp_precipitation(flux_rain, flux_snow, DENS, MOMZ, MOMX, MOMY, RHOE, QTRC, vterm, temp, dt)
precipitation transport
subroutine, public atmos_phy_mp_kessler_mixingratio(Qe, QTRC0)
Calculate mixing ratio of each category.
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
module ATMOSPHERE / Reference state
subroutine, public atmos_phy_mp_kessler_cloudfraction(cldfrac, QTRC, mask_criterion)
Calculate Cloud Fraction.
module STDIO
Definition: scale_stdio.F90:12
integer, public ke
end point of inner domain: z, local
integer, public qa
subroutine, public atmos_phy_mp_kessler_effectiveradius(Re, QTRC0, DENS0, TEMP0)
Calculate Effective Radius.
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
module ATMOSPHERE / Physics Cloud Microphysics
integer, public mp_qa
integer, public ka
of z whole cells (local, with HALO)
real(rp), dimension(:,:,:), allocatable, public atmos_refstate_dens
refernce density [kg/m3]
module TRACER / kessler
subroutine, public atmos_phy_mp_kessler_setup(MP_TYPE)
Setup.
module COMMUNICATION
Definition: scale_comm.F90:23
integer, public js
start point of inner domain: y, local
module TIME
Definition: scale_time.F90:15
module PROCESS
real(rp), public const_lhv
latent heat of vaporizaion for use
Definition: scale_const.F90:77
module CONSTANT
Definition: scale_const.F90:14
integer, public ks
start point of inner domain: z, local
subroutine, public atmos_phy_mp_negative_fixer(DENS, RHOT, QTRC)
Negative fixer.
module GRID (cartesian)
subroutine, public prof_rapstart(rapname_base, level)
Start raptime.
Definition: scale_prof.F90:132
module profiler
Definition: scale_prof.F90:10
integer, public ie
end point of inner domain: x, local
module ATMOSPHERE / Thermodynamics
logical, public io_lnml
output log or not? (for namelist, this process)
Definition: scale_stdio.F90:60
real(rp), dimension(mp_qa), target, public atmos_phy_mp_dens
module PRECISION
real(rp), dimension(:), allocatable, public grid_cdz
z-length of control volume [m]
module HISTORY
integer, public io_fid_conf
Config file ID.
Definition: scale_stdio.F90:55
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
subroutine, public prof_rapend(rapname_base, level)
Save raptime.
Definition: scale_prof.F90:178
integer, parameter, public rp
integer, public ja
of y whole cells (local, with HALO)