SCALE-RM
scale_atmos_phy_mp_sn14.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
48 !-------------------------------------------------------------------------------
49 
50 #ifdef PROFILE_FAPP
51 #define PROFILE_START(name) call fapp_start(name, 1, 1)
52 #define PROFILE_STOP(name) call fapp_stop (name, 1, 1)
53 #elif defined(PROFILE_FINEPA)
54 #define PROFILE_START(name) call start_collection(name)
55 #define PROFILE_STOP(name) call stop_collection (name)
56 #else
57 #define PROFILE_START(name)
58 #define PROFILE_STOP(name)
59 #endif
60 
61 #include "scalelib.h"
63  !-----------------------------------------------------------------------------
64  !
65  !++ used modules
66  !
67  use scale_precision
68  use scale_io
69  use scale_prof
70 
71  use scale_const, only: &
72  grav => const_grav, &
73  pi => const_pi, &
74  undef8 => const_undef8, &
75  rdry => const_rdry, &
76  cpdry => const_cpdry, &
77  cvdry => const_cvdry, &
78  p00 => const_pre00, &
79  t00 => const_tem00, &
80  rvap => const_rvap, &
81  cpvap => const_cpvap, &
82  cvvap => const_cvvap, &
83  cl => const_cl, &
84  ci => const_ci, &
85  lhv => const_lhv00, &
86  lhf => const_lhf00, &
87  lhv0 => const_lhv0, &
88  lhf0 => const_lhf0, &
89  lhs0 => const_lhs0, &
90  lhv00 => const_lhv00, &
91  lhf00 => const_lhf00, &
92  psat0 => const_psat0, &
93  emelt => const_emelt, &
94  dwatr => const_dwatr
95 
96  use scale_atmos_hydrometeor, only: &
97  n_hyd
98 
99  !-----------------------------------------------------------------------------
100  implicit none
101  private
102  !-----------------------------------------------------------------------------
103  !
104  !++ Public procedure
105  !
106  public :: atmos_phy_mp_sn14_setup
114 
115  !-----------------------------------------------------------------------------
116  !
117  !++ Public parameters & variables
118  !
119  integer, public, parameter :: qa_mp = 11
120 
121  integer, parameter, public :: atmos_phy_mp_sn14_ntracers = qa_mp
122  integer, parameter, public :: atmos_phy_mp_sn14_nwaters = 2
123  integer, parameter, public :: atmos_phy_mp_sn14_nices = 3
124  character(len=H_SHORT), parameter, public :: atmos_phy_mp_sn14_tracer_names(qa_mp) = (/ &
125  'QV', &
126  'QC', &
127  'QR', &
128  'QI', &
129  'QS', &
130  'QG', &
131  'NC', &
132  'NR', &
133  'NI', &
134  'NS', &
135  'NG' /)
136  character(len=H_MID) , parameter, public :: atmos_phy_mp_sn14_tracer_descriptions(qa_mp) = (/ &
137  'Ratio of Water Vapor mass to total mass (Specific humidity)', &
138  'Ratio of Cloud Water mass to total mass ', &
139  'Ratio of Rain Water mass to total mass ', &
140  'Ratio of Cloud Ice mass ratio to total mass ', &
141  'Ratio of Snow miass ratio to total mass ', &
142  'Ratio of Graupel mass ratio to total mass ', &
143  'Cloud Water Number Density ', &
144  'Rain Water Number Density ', &
145  'Cloud Ice Number Density ', &
146  'Snow Number Density ', &
147  'Graupel Number Density '/)
148  character(len=H_SHORT), parameter, public :: atmos_phy_mp_sn14_tracer_units(qa_mp) = (/ &
149  'kg/kg ', &
150  'kg/kg ', &
151  'kg/kg ', &
152  'kg/kg ', &
153  'kg/kg ', &
154  'kg/kg ', &
155  'num/kg', &
156  'num/kg', &
157  'num/kg', &
158  'num/kg', &
159  'num/kg' /)
160 
161  !-----------------------------------------------------------------------------
162  !
163  !++ Private procedure
164  !
165  private :: mp_sn14_init
166  private :: mp_sn14
167 
168  !-----------------------------------------------------------------------------
169  !
170  !++ Private parameters
171  !
172  integer, private, parameter :: i_qv = 1
173  integer, private, parameter :: i_qc = 2
174  integer, private, parameter :: i_qr = 3
175  integer, private, parameter :: i_qi = 4
176  integer, private, parameter :: i_qs = 5
177  integer, private, parameter :: i_qg = 6
178  integer, private, parameter :: i_nc = 7
179  integer, private, parameter :: i_nr = 8
180  integer, private, parameter :: i_ni = 9
181  integer, private, parameter :: i_ns = 10
182  integer, private, parameter :: i_ng = 11
183 
184  integer, private, parameter :: hydro_max = 5
185 
186  integer, private, parameter :: i_mp_qc = 1
187  integer, private, parameter :: i_mp_qr = 2
188  integer, private, parameter :: i_mp_qi = 3
189  integer, private, parameter :: i_mp_qs = 4
190  integer, private, parameter :: i_mp_qg = 5
191  integer, private, parameter :: i_mp_nc = 6
192  integer, private, parameter :: i_mp_nr = 7
193  integer, private, parameter :: i_mp_ni = 8
194  integer, private, parameter :: i_mp_ns = 9
195  integer, private, parameter :: i_mp_ng = 10
196 
197  ! production rate
198  ! nucleation
199  integer, parameter :: i_lcccn = 1
200  integer, parameter :: i_ncccn = 2
201  integer, parameter :: i_liccn = 3
202  integer, parameter :: i_niccn = 4
203  ! freezing
204  integer, parameter :: i_lchom = 5
205  integer, parameter :: i_nchom = 6
206  integer, parameter :: i_lchet = 7
207  integer, parameter :: i_nchet = 8
208  integer, parameter :: i_lrhet = 9
209  integer, parameter :: i_nrhet = 10
210  ! melting
211  integer, parameter :: i_limlt = 11
212  integer, parameter :: i_nimlt = 12
213  integer, parameter :: i_lsmlt = 13
214  integer, parameter :: i_nsmlt = 14
215  integer, parameter :: i_lgmlt = 15
216  integer, parameter :: i_ngmlt = 16
217  ! vapor deposition
218  integer, parameter :: i_lrdep = 17
219  integer, parameter :: i_nrdep = 18
220  integer, parameter :: i_lidep = 19
221  integer, parameter :: i_nidep = 20
222  integer, parameter :: i_lsdep = 21
223  integer, parameter :: i_nsdep = 22
224  integer, parameter :: i_lgdep = 23
225  integer, parameter :: i_ngdep = 24
226  integer, parameter :: i_lcdep = 25
227 ! integer, parameter :: I_NCdep = 26
228  ! warm collection process
229  ! auto-conversion
230  integer, parameter :: i_lcaut = 26
231  integer, parameter :: i_ncaut = 27
232  integer, parameter :: i_nraut = 28
233  ! accretion
234  integer, parameter :: i_lcacc = 29
235  integer, parameter :: i_ncacc = 30
236  ! self-colletion, break-up
237  integer, parameter :: i_nrslc = 31
238  integer, parameter :: i_nrbrk = 32
239 
240  ! partial conversion(ice, snow => graupel)
241  integer, parameter :: i_licon = 33
242  integer, parameter :: i_nicon = 34
243  integer, parameter :: i_lscon = 35
244  integer, parameter :: i_nscon = 36
245  ! enhanced melting due to
246  integer, parameter :: i_liacm = 37 ! ice-cloud
247  integer, parameter :: i_niacm = 38
248  integer, parameter :: i_liarm = 39 ! ice-rain
249  integer, parameter :: i_niarm = 40
250  integer, parameter :: i_lsacm = 41 ! snow-cloud
251  integer, parameter :: i_nsacm = 42
252  integer, parameter :: i_lsarm = 43 ! snow-rain
253  integer, parameter :: i_nsarm = 44
254  integer, parameter :: i_lgacm = 45 ! graupel-cloud
255  integer, parameter :: i_ngacm = 46
256  integer, parameter :: i_lgarm = 47 ! graupel-rain
257  integer, parameter :: i_ngarm = 48
258  ! ice multiplication by splintering
259  integer, parameter :: i_lgspl = 49
260  integer, parameter :: i_lsspl = 50
261  integer, parameter :: i_nispl = 51
262 
263  integer, parameter :: pq_max = 51
264 
265  ! production rate of mixed-phase collection process
266  ! PXXacYY2ZZ means XX collect YY produce ZZ
267  integer, parameter :: i_liaclc2li = 1 ! cloud-ice
268  integer, parameter :: i_niacnc2ni = 2
269  integer, parameter :: i_lsaclc2ls = 3 ! cloud-snow(cloud change)
270  integer, parameter :: i_nsacnc2ns = 4
271  integer, parameter :: i_lgaclc2lg = 5 ! cloud-graupel
272  integer, parameter :: i_ngacnc2ng = 6
273  integer, parameter :: i_lracli2lg_i = 7 ! rain-ice(ice change)
274  integer, parameter :: i_nracni2ng_i = 8
275  integer, parameter :: i_lracli2lg_r = 9 ! rain-ice(rain change)
276  integer, parameter :: i_nracni2ng_r = 10
277  integer, parameter :: i_lracls2lg_s = 11 ! rain-snow(snow change)
278  integer, parameter :: i_nracns2ng_s = 12
279  integer, parameter :: i_lracls2lg_r = 13 ! rain-snow(rain change)
280  integer, parameter :: i_nracns2ng_r = 14
281  integer, parameter :: i_lraclg2lg = 15 ! rain-graupel(rain change)
282  integer, parameter :: i_nracng2ng = 16
283  integer, parameter :: i_liacli2ls = 17 ! ice-ice
284  integer, parameter :: i_niacni2ns = 18
285  integer, parameter :: i_liacls2ls = 19 ! ice-snow(ice change)
286  integer, parameter :: i_niacns2ns = 20
287  integer, parameter :: i_nsacns2ns = 21 ! snow-snow
288  integer, parameter :: i_ngacng2ng = 22 ! graupel-graupel
289  integer, parameter :: i_lgacls2lg = 23 ! snow-graupel
290  integer, parameter :: i_ngacns2ng = 24
291 
292  integer, parameter :: pac_max = 24
293 
294 
295 
296  character(len=H_SHORT), save :: wlabel(hydro_max)
297 
298  ! empirical value from Meyers etal.(1991), 1[/liter] = 1.d3[/m3]
299  real(RP), private, parameter :: nqmin(hydro_max) = (/ 1.e+4_rp, 1.0_rp, 1.0_rp, 1.e-4_rp, 1.e-4_rp /) ! [1/m3]
300  ! refer to Seifert(2002) (Dr. Thesis, Table.5.1)
301  ! max mass, for D_min=79um, 2mm, 5mm, 1cm, 1cm
302  real(RP), private, parameter :: xqmax(hydro_max) = (/ 2.6e-10_rp, 5.0e-6_rp, 1.377e-6_rp, 7.519e-6_rp, 4.90e-5_rp /)! [kg]
303  ! SB06, Table 1.
304  ! min mass, for D_min=2um, 79um, 10um, 20um, 100um
305  real(RP), private, parameter :: xqmin(hydro_max) = (/ 4.20e-15_rp, 2.60e-10_rp, 3.382e-13_rp, 1.847e-12_rp, 1.230e-10_rp /)! [kg]
306 
307 
308 
309  ! for all processes
310  ! SB06, Table 1.
311  real(RP), private, parameter :: xc_min = 4.20e-15_rp ! [kg] : min mass, D_min=2um
312  real(RP), private, parameter :: xr_min = 2.60e-10_rp ! [kg] : min mass, D_min=79um
313  real(RP), private, parameter :: xi_min = 3.382e-13_rp ! [kg] : min mass, D_min=10um
314  real(RP), private, parameter :: xs_min = 1.847e-12_rp ! [kg] : min mass, D_min=20um
315  real(RP), private, parameter :: xg_min = 1.230e-10_rp ! [kg] : min mass, D_min=100um
316  ! refer to Seifert(2002) (Dr. Thesis, Table.5.1)
317  real(RP), private, parameter :: xc_max = 2.6e-10_rp ! [kg] : max, D_max=79um
318  real(RP), private, parameter :: xr_max = 5.00e-6_rp ! [kg] : max, D_max=2mm
319  real(RP), private, parameter :: xi_max = 1.377e-6_rp ! [kg] : max, D_max=5mm
320  real(RP), private, parameter :: xs_max = 7.519e-6_rp ! [kg] : max, D_max=1cm
321  real(RP), private, parameter :: xg_max = 4.900e-5_rp ! [kg] : max, D_max=1cm
322  ! filter similar to Ikawa et al.(1991) sec.3.5
323  real(RP), private, parameter :: xmin_filter= xc_min
324  ! filter of effective radius(1 micron)
325  real(RP), private, parameter :: rmin_re= 1.e-6_rp
326  !
327  ! SB06(95),(96)
328  real(RP), private, parameter :: n0r_min= 2.5e+5_rp ! [m-4]: min intercept parameter of rain
329  real(RP), private, parameter :: n0r_max= 2.0e+7_rp ! [m-4]: max
330  real(RP), private, parameter :: lambdar_min= 1.e+3_rp ! [m-1]: min slope parameter of rain
331  real(RP), private, parameter :: lambdar_max= 1.e+4_rp ! [m-1]: max
332  ! empirical value from Meyers etal.(1991), 1[/liter] = 1.d3[/m3]
333  real(RP), private, parameter :: nc_min = 1.e+4_rp ! [m-3] empirical T.Mitsui
334  real(RP), private, parameter :: nr_min = 1.0_rp ! [m-3] 1/1000 [/liter]
335  real(RP), private, parameter :: ni_min = 1.0_rp ! [m-3]
336  real(RP), private, parameter :: ns_min = 1.e-4_rp ! [m-3]
337  real(RP), private, parameter :: ng_min = 1.e-4_rp ! [m-3]
338  ! empirical filter
339  real(RP), private, parameter :: lc_min = xc_min*nc_min
340  real(RP), private, parameter :: lr_min = xr_min*nr_min
341  real(RP), private, parameter :: li_min = xi_min*ni_min
342  real(RP), private, parameter :: ls_min = xs_min*ns_min
343  real(RP), private, parameter :: lg_min = xg_min*ng_min
344  !
345  real(RP), private, parameter :: x_sep = 2.6e-10_rp ! boundary mass between cloud and rain
346  !
347  real(RP), private, parameter :: tem_min=100.0_rp
348  real(RP), private, parameter :: rho_min=1.e-5_rp ! 3.e-3 is lower limit recognized in many experiments.
349  real(RP), private, parameter :: rhoi = 916.70_rp
350  !
351  integer, private, save :: ntmax_phase_change = 1
352  integer, private, save :: ntmax_collection = 1
353  !
354  !--- standard density
355  real(RP), private, parameter :: rho_0 = 1.280_rp
356  !--- max number of Nc( activatable aerosol number concentration )
357  real(RP), allocatable, private, save :: nc_uplim_d(:,:,:)
358  !
359  !--- thermal conductivity of air
360  real(RP), private, parameter :: ka0 = 2.428e-2_rp
361 
362  real(RP), private, parameter :: dka_dt = 7.47e-5_rp
363 
364  !====== Ka = Ka0 + temc*dKa_dT
365  !
366  !--- Dynamic viscosity
367  real(RP), private, parameter :: mua0 = 1.718e-5_rp
368 
369  real(RP), private, parameter :: dmua_dt = 5.28e-8_rp
370 
371  !====== mua = mua0 + temc*dmua_dT
372  !
373  real(RP), private, save :: xc_ccn = 1.e-12_rp ! [kg]
374  real(RP), private, save :: xi_ccn = 1.e-12_rp ! [kg] ! [move] 11/08/30 T.Mitsui
375  !
376  ! capacity of diffusional growth
377  ! ( dependent of their geometries )
378  real(RP), private, save :: cap(hydro_max)
379  !
380  ! constants for Diameter-Mass relation
381  ! D = a * x^b
382  real(RP), private, save :: a_m(hydro_max)
383  real(RP), private, save :: b_m(hydro_max)
384  ! constants for Terminal velocity-Mass relation
385  ! vt = alpha * x^beta * f
386  real(RP), private, save :: alpha_v(hydro_max,2)
387  real(RP), private, save :: beta_v(hydro_max,2)
388  real(RP), private, save :: alpha_vn(hydro_max,2) !
389  real(RP), private, save :: beta_vn(hydro_max,2) !
390  real(RP), private, save :: gamma_v(hydro_max)
391  ! Aerodynamical factor for correction of terminal velocity.(Heymsfield and Iaquinta, 2000)
392  ! vt(tem,pre) = vt0 * (pre/pre0)**a_pre0 * (tem/tem0)**a_tem0
393  real(RP), private, parameter :: pre0_vt = 300.e+2_rp ! 300hPa
394  real(RP), private, parameter :: tem0_vt = 233.0_rp ! -40degC
395  real(RP), private, parameter :: a_pre0_vt = -0.1780_rp
396  real(RP), private, parameter :: a_tem0_vt = -0.3940_rp
397  ! Parameters to determine Droplet Size Distribution
398  ! as a General Gamma Distribution
399  ! f(x) = A x^nu exp(-lambda x^mu )
400  ! for Marshall Palmer Distribution ( popular for rain )
401  ! mu=1/3, nu=-2/3
402  ! for Gamma Distribution ( popular for cloud )
403  ! mu=1
404  real(RP), private, save :: nu(hydro_max)
405  real(RP), private, save :: mu(hydro_max)
406  ! Mitchell(1996), JAS, vol.53, No.12, pp.1710-1723
407  ! area = a_area*D^b_area
408  ! area = ax_area*x^bx_area
409  ! Auer and Veal(1970), JAS, vol.27, pp.919-pp.926
410  ! height = a_h*x^b_h( based on h=a_ar*D^b_ar, ar:aspect ratio)
411  real(RP), private, save :: a_area(hydro_max) !
412  real(RP), private, save :: b_area(hydro_max) !
413  real(RP), private, save :: ax_area(hydro_max) !
414  real(RP), private, save :: bx_area(hydro_max) !
415  ! parameters for radius of equivalent area
416  ! r_ea = a_rea*x**b_rea
417  real(RP), private, save :: a_rea(hydro_max) !
418  real(RP), private, save :: b_rea(hydro_max) !
419  real(RP), private, save :: a_rea2(hydro_max) !
420  real(RP), private, save :: b_rea2(hydro_max) !
421  real(RP), private, save :: a_rea3(hydro_max) !
422  real(RP), private, save :: b_rea3(hydro_max) !
423  !
424  real(RP), private, save :: a_d2vt(hydro_max) !
425  real(RP), private, save :: b_d2vt(hydro_max) !
426  ! coefficient of x^2 moment of DSD
427  ! Z = integral x*x*f(x) dx
428  ! = coef_m2*N*(L/N)^2
429  real(RP), private, save :: coef_m2(hydro_max)
430  ! radar reflectivity coefficient defined by diameter
431  real(RP), private, save :: coef_d6(hydro_max) !
432  ! volume coefficient defined by diameter
433  real(RP), private, save :: coef_d3(hydro_max) !
434  ! coefficient of weighted mean diameter
435  real(RP), private, save :: coef_d(hydro_max)
436  ! coefficient of weighted mean d*d*v
437  real(RP), private, save :: coef_d2v(hydro_max) !
438  ! coefficient of moment of d*d*v
439  real(RP), private, save :: coef_md2v(hydro_max) !
440  !
441  ! for effective radius(spherical particle)
442  real(RP), private, save :: coef_r2(hydro_max)
443  real(RP), private, save :: coef_r3(hydro_max)
444  real(RP), private, save :: coef_re(hydro_max)
445  ! for effective radius(hexagonal plate)
446  real(RP), private, save :: coef_rea2(hydro_max) !
447  real(RP), private, save :: coef_rea3(hydro_max) !
448  logical, private, save :: opt_m96_ice=.true. !
449  logical, private, save :: opt_m96_column_ice=.false. !
450  !
451  ! coefficeint of weighted mean terminal velocity
452  ! vt0 is number weighted and
453  ! vt1 is mass weighted.
454  real(RP), private, save :: coef_vt0(hydro_max,2)
455  real(RP), private, save :: coef_vt1(hydro_max,2)
456  real(RP), private, save :: coef_deplc
457  real(RP), private, save :: coef_dave_n(hydro_max) !
458  real(RP), private, save :: coef_dave_l(hydro_max) !
459  ! diameter of terminal velocity branch
460  !
461  real(RP), private, save :: d0_ni=261.76e-6_rp !
462  real(RP), private, save :: d0_li=398.54e-6_rp
463  real(RP), private, parameter :: d0_ns=270.03e-6_rp !
464  real(RP), private, parameter :: d0_ls=397.47e-6_rp !
465  real(RP), private, parameter :: d0_ng=269.08e-6_rp !
466  real(RP), private, parameter :: d0_lg=376.36e-6_rp !
467  !
468  real(RP), private, parameter :: coef_vtr_ar1=9.65_rp ! coef. for large branch
469  ! original parameter of Rogers etal.(1993)
470  real(RP), private, parameter :: coef_vtr_br1=10.43_rp ! ...
471  real(RP), private, parameter :: coef_vtr_cr1=600.0_rp ! ...
472  real(RP), private, parameter :: coef_vtr_ar2=4.e+3_rp ! coef. for small branch
473  real(RP), private, parameter :: coef_vtr_br2=12.e+3_rp ! ...
474  real(RP), private, parameter :: d_vtr_branch=0.745e-3_rp ! 0.745 mm (diameter dividing 2-branches)
475  ! equilibrium diameter of rain break-up
476  real(RP), private, parameter :: dr_eq = 1.10e-3_rp ! eqilibrium diameter, Seifert 2008(36)
477  ! coefficient of General Gamma.
478  ! f(x) = A x^nu exp(-lambda x^mu )
479  ! lambda = coef_lambda * (L/N)^{-mu}
480  ! A = coef_A*N*lambda^slope_A
481  real(RP), private, save :: coef_a(hydro_max)
482  real(RP), private, save :: coef_lambda(hydro_max)
483 ! real(RP), private, save :: slope_A(HYDRO_MAX)
484  ! coefficeint of weighted ventilation effect.
485  ! large, and small branch is by PK97(13-60),(13-61),(13-88),(13-89)
486  real(RP), private, save :: ah_vent (hydro_max,2) !
487  real(RP), private, save :: bh_vent (hydro_max,2) !
488  real(RP), private, save :: ah_vent0 (hydro_max,2) !
489  real(RP), private, save :: bh_vent0 (hydro_max,2) !
490  real(RP), private, save :: ah_vent1 (hydro_max,2) !
491  real(RP), private, save :: bh_vent1 (hydro_max,2) !
492  ! coefficient of collision growth
493  real(RP), private, save :: delta_b0 (hydro_max)
494  real(RP), private, save :: delta_b1 (hydro_max)
495  real(RP), private, save :: delta_ab0(hydro_max,hydro_max)
496  real(RP), private, save :: delta_ab1(hydro_max,hydro_max)
497  !
498  real(RP), private, save :: theta_b0 (hydro_max)
499  real(RP), private, save :: theta_b1 (hydro_max)
500  real(RP), private, save :: theta_ab0(hydro_max,hydro_max)
501  real(RP), private, save :: theta_ab1(hydro_max,hydro_max)
502  !
503  logical, private, save :: opt_debug=.false.
504  !
505  logical, private, save :: opt_debug_tem=.false.
506  logical, private, save :: opt_debug_inc=.true.
507  logical, private, save :: opt_debug_act=.true.
508  logical, private, save :: opt_debug_ree=.true.
509  logical, private, save :: opt_debug_bcs=.true.
510 
511  logical, private, save :: mp_doautoconversion = .true.
512  logical, private, save :: mp_couple_aerosol = .false. ! apply CCN effect?
513  real(RP), private, save :: mp_ssw_lim = 1.e+1_rp
514 
515  !-----------------------------------------------------------------------------
516 contains
517  !-----------------------------------------------------------------------------
521  subroutine atmos_phy_mp_sn14_setup( &
522  & KA, IA, JA )
523  use scale_prc, only: &
524  prc_abort
525  implicit none
526 
527  integer, intent(in) :: KA
528  integer, intent(in) :: IA
529  integer, intent(in) :: JA
530 
531  namelist / param_atmos_phy_mp_sn14 / &
532  mp_doautoconversion, &
533  mp_ssw_lim, &
534  mp_couple_aerosol
535 
536  integer :: ierr
537  !---------------------------------------------------------------------------
538 
539  log_newline
540  log_info("ATMOS_PHY_MP_sn14_setup",*) 'Setup'
541  log_info("ATMOS_PHY_MP_sn14_setup",*) 'Seiki and Nakajima (2014) 2-moment bulk 6 category'
542 
543  !--- read namelist
544  rewind(io_fid_conf)
545  read(io_fid_conf,nml=param_atmos_phy_mp_sn14,iostat=ierr)
546  if( ierr < 0 ) then !--- missing
547  log_info("ATMOS_PHY_MP_sn14_setup",*) 'Not found namelist. Default used.'
548  elseif( ierr > 0 ) then !--- fatal error
549  log_error("ATMOS_PHY_MP_sn14_setup",*) 'Not appropriate names in namelist PARAM_ATMOS_PHY_MP_SN14. Check!'
550  call prc_abort
551  endif
552  log_nml(param_atmos_phy_mp_sn14)
553 
554  wlabel(1) = "CLOUD"
555  wlabel(2) = "RAIN"
556  wlabel(3) = "ICE"
557  wlabel(4) = "SNOW"
558  wlabel(5) = "GRAUPEL"
559 
560  call mp_sn14_init
561 
562  allocate(nc_uplim_d(1,ia,ja))
563  nc_uplim_d(:,:,:) = 150.e6_rp
564 
565  return
566  end subroutine atmos_phy_mp_sn14_setup
567 
568  !-----------------------------------------------------------------------------
572  subroutine atmos_phy_mp_sn14_tendency( &
573  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
574  DENS, &
575  W, &
576  QTRC, &
577  PRES, &
578  TEMP, &
579  Qdry, &
580  CPtot, &
581  CVtot, &
582  CCN, &
583  dt, &
584  cz, &
585  fz, &
586  RHOQ_t, &
587  RHOE_t, &
588  CPtot_t, &
589  CVtot_t, &
590  EVAPORATE )
591  implicit none
592 
593  integer, intent(in) :: KA, KS, KE
594  integer, intent(in) :: IA, IS, IE
595  integer, intent(in) :: JA, JS, JE
596 
597  real(RP), intent(in) :: DENS (ka,ia,ja)
598  real(RP), intent(in) :: W (ka,ia,ja)
599  real(RP), intent(in) :: QTRC (ka,ia,ja,qa_mp)
600  real(RP), intent(in) :: PRES(ka,ia,ja)
601  real(RP), intent(in) :: TEMP(ka,ia,ja)
602  real(RP), intent(in) :: Qdry(ka,ia,ja)
603  real(RP), intent(in) :: CPtot(ka,ia,ja)
604  real(RP), intent(in) :: CVtot(ka,ia,ja)
605  real(RP), intent(in) :: CCN (ka,ia,ja)
606  real(DP), intent(in) :: dt
607  real(RP), intent(in) :: cz( ka,ia,ja)
608  real(RP), intent(in) :: fz(0:ka,ia,ja)
609 
610  real(RP), intent(out) :: RHOQ_t (ka,ia,ja,qa_mp)
611  real(RP), intent(out) :: RHOE_t (ka,ia,ja)
612  real(RP), intent(out) :: CPtot_t(ka,ia,ja)
613  real(RP), intent(out) :: CVtot_t(ka,ia,ja)
614  real(RP), intent(out) :: EVAPORATE(ka,ia,ja) !--- number of evaporated cloud [/m3]
615 
616  !---------------------------------------------------------------------------
617 
618  log_progress(*) 'atmosphere / physics / microphysics / SN14'
619 
620 #ifdef PROFILE_FIPP
621  call fipp_start()
622 #endif
623 
624  !##### MP Main #####
625  call mp_sn14 ( &
626  ka, ks, ke, ia, is, ie, ja, js, je, &
627  dens(:,:,:), w(:,:,:), qtrc(:,:,:,:), pres(:,:,:), temp(:,:,:), & ! (in)
628  qdry(:,:,:), cptot(:,:,:), cvtot(:,:,:), ccn(:,:,:), & ! (in)
629  real(dt,RP), cz(:,:,:), fz(:,:,:), & ! (in)
630  RHOQ_t(:,:,:,:), RHOE_t(:,:,:), CPtot_t(:,:,:), CVtot_t(:,:,:), & ! (out)
631  EVAPORATE(:,:,:) ) ! (out)
632 
633 #ifdef PROFILE_FIPP
634  call fipp_stop()
635 #endif
636 
637  return
638  end subroutine atmos_phy_mp_sn14_tendency
639 
640  !-----------------------------------------------------------------------------
645  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
646  QTRC, &
647  mask_criterion, &
648  cldfrac )
649  implicit none
650  integer, intent(in) :: KA, KS, KE
651  integer, intent(in) :: IA, IS, IE
652  integer, intent(in) :: JA, JS, JE
653 
654  real(RP), intent(in) :: QTRC (ka,ia,ja,qa_mp-1)
655  real(RP), intent(in) :: mask_criterion
656 
657  real(RP), intent(out) :: cldfrac(ka,ia,ja)
658 
659  real(RP) :: qhydro
660  integer :: k, i, j, iq
661  !---------------------------------------------------------------------------
662 
663  do j = js, je
664  do i = is, ie
665  do k = ks, ke
666  qhydro = 0.0_rp
667  do iq = 1, qa_mp-1
668  qhydro = qhydro + qtrc(k,i,j,iq)
669  enddo
670  cldfrac(k,i,j) = 0.5_rp + sign(0.5_rp,qhydro-mask_criterion)
671  enddo
672  enddo
673  enddo
674 
675  return
676  end subroutine atmos_phy_mp_sn14_cloud_fraction
677  !-----------------------------------------------------------------------------
682  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
683  DENS0, TEMP0, QTRC0, &
684  Re )
686  n_hyd, &
687  i_hc, &
688  i_hr, &
689  i_hi, &
690  i_hs, &
691  i_hg
692  implicit none
693  integer, intent(in) :: KA, KS, KE
694  integer, intent(in) :: IA, IS, IE
695  integer, intent(in) :: JA, JS, JE
696 
697  real(RP), intent(in) :: DENS0(ka,ia,ja) ! density [kg/m3]
698  real(RP), intent(in) :: TEMP0(ka,ia,ja) ! temperature [K]
699  real(RP), intent(in) :: QTRC0(ka,ia,ja,i_qc:i_ng) ! tracer mass concentration [kg/kg]
700 
701  real(RP), intent(out) :: Re (ka,ia,ja,n_hyd) ! effective radius [cm]
702 
703  ! mass concentration[kg/m3] and mean particle mass[kg]
704  real(RP) :: xc(ka,ia,ja)
705  real(RP) :: xr(ka,ia,ja)
706  real(RP) :: xi(ka,ia,ja)
707  real(RP) :: xs(ka,ia,ja)
708  real(RP) :: xg(ka,ia,ja)
709  ! diameter of average mass[kg/m3]
710  real(RP) :: dc_ave(ka,ia,ja)
711  real(RP) :: dr_ave(ka,ia,ja)
712  ! radius of average mass
713  real(RP) :: rc, rr
714  ! 2nd. and 3rd. order moment of DSD
715  real(RP) :: ri2m(ka,ia,ja), ri3m(ka,ia,ja)
716  real(RP) :: rs2m(ka,ia,ja), rs3m(ka,ia,ja)
717  real(RP) :: rg2m(ka,ia,ja), rg3m(ka,ia,ja)
718 
719  real(RP) :: coef_Fuetal1998
720  ! r2m_min is minimum value(moment of 1 particle with 1 micron)
721  real(RP), parameter :: r2m_min=1.e-12_rp
722  real(RP), parameter :: um2cm = 100.0_rp
723 
724  real(RP) :: limitsw, zerosw
725  integer :: k, i, j
726  !---------------------------------------------------------------------------
727 
728  ! mean particle mass[kg]
729  do j = js, je
730  do i = is, ie
731  do k = ks, ke
732  xc(k,i,j) = min( xc_max, max( xc_min, dens0(k,i,j)*qtrc0(k,i,j,i_qc)/(qtrc0(k,i,j,i_nc)+nc_min) ) )
733  xr(k,i,j) = min( xr_max, max( xr_min, dens0(k,i,j)*qtrc0(k,i,j,i_qr)/(qtrc0(k,i,j,i_nr)+nr_min) ) )
734  xi(k,i,j) = min( xi_max, max( xi_min, dens0(k,i,j)*qtrc0(k,i,j,i_qi)/(qtrc0(k,i,j,i_ni)+ni_min) ) )
735  xs(k,i,j) = min( xs_max, max( xs_min, dens0(k,i,j)*qtrc0(k,i,j,i_qs)/(qtrc0(k,i,j,i_ns)+ns_min) ) )
736  xg(k,i,j) = min( xg_max, max( xg_min, dens0(k,i,j)*qtrc0(k,i,j,i_qg)/(qtrc0(k,i,j,i_ng)+ng_min) ) )
737  enddo
738  enddo
739  enddo
740 
741  ! diameter of average mass : SB06 eq.(32)
742  do j = js, je
743  do i = is, ie
744  do k = ks, ke
745  dc_ave(k,i,j) = a_m(i_mp_qc) * xc(k,i,j)**b_m(i_mp_qc)
746  dr_ave(k,i,j) = a_m(i_mp_qr) * xr(k,i,j)**b_m(i_mp_qr)
747  enddo
748  enddo
749  enddo
750 
751  ! cloud effective radius
752  do j = js, je
753  do i = is, ie
754  do k = ks, ke
755  rc = 0.5_rp * dc_ave(k,i,j)
756  limitsw = 0.5_rp + sign(0.5_rp, rc-rmin_re )
757  re(k,i,j,i_hc) = coef_re(i_mp_qc) * rc * limitsw * um2cm
758  enddo
759  enddo
760  enddo
761 
762  ! rain effective radius
763  do j = js, je
764  do i = is, ie
765  do k = ks, ke
766  rr = 0.5_rp * dr_ave(k,i,j)
767  limitsw = 0.5_rp + sign(0.5_rp, rr-rmin_re )
768  re(k,i,j,i_hr) = coef_re(i_mp_qr) * rr * limitsw * um2cm
769  enddo
770  enddo
771  enddo
772 
773  do j = js, je
774  do i = is, ie
775  do k = ks, ke
776  ri2m(k,i,j) = pi * coef_rea2(i_mp_qi) * qtrc0(k,i,j,i_ni) * a_rea2(i_mp_qi) * xi(k,i,j)**b_rea2(i_mp_qi)
777  rs2m(k,i,j) = pi * coef_rea2(i_mp_qs) * qtrc0(k,i,j,i_ns) * a_rea2(i_mp_qs) * xs(k,i,j)**b_rea2(i_mp_qs)
778  rg2m(k,i,j) = pi * coef_rea2(i_mp_qg) * qtrc0(k,i,j,i_ng) * a_rea2(i_mp_qg) * xg(k,i,j)**b_rea2(i_mp_qg)
779  enddo
780  enddo
781  enddo
782 
783  ! Fu(1996), eq.(3.11) or Fu et al.(1998), eq.(2.5)
784  coef_fuetal1998 = 3.0_rp / (4.0_rp*rhoi)
785  do j = js, je
786  do i = is, ie
787  do k = ks, ke
788  ri3m(k,i,j) = coef_fuetal1998 * qtrc0(k,i,j,i_ni) * xi(k,i,j)
789  rs3m(k,i,j) = coef_fuetal1998 * qtrc0(k,i,j,i_ns) * xs(k,i,j)
790  rg3m(k,i,j) = coef_fuetal1998 * qtrc0(k,i,j,i_ng) * xg(k,i,j)
791  enddo
792  enddo
793  enddo
794 
795  ! ice effective radius
796  do j = js, je
797  do i = is, ie
798  do k = ks, ke
799  zerosw = 0.5_rp - sign(0.5_rp, ri2m(k,i,j) - r2m_min )
800  re(k,i,j,i_hi) = ri3m(k,i,j) / ( ri2m(k,i,j) + zerosw ) * ( 1.0_rp - zerosw ) * um2cm
801  enddo
802  enddo
803  enddo
804 
805  ! snow effective radius
806  do j = js, je
807  do i = is, ie
808  do k = ks, ke
809  zerosw = 0.5_rp - sign(0.5_rp, rs2m(k,i,j) - r2m_min )
810  re(k,i,j,i_hs) = rs3m(k,i,j) / ( rs2m(k,i,j) + zerosw ) * ( 1.0_rp - zerosw ) * um2cm
811  enddo
812  enddo
813  enddo
814 
815  ! graupel effective radius
816  do j = js, je
817  do i = is, ie
818  do k = ks, ke
819  zerosw = 0.5_rp - sign(0.5_rp, rg2m(k,i,j) - r2m_min )
820  re(k,i,j,i_hg) = rg3m(k,i,j) / ( rg2m(k,i,j) + zerosw ) * ( 1.0_rp - zerosw ) * um2cm
821  enddo
822  enddo
823  enddo
824 
825  re(:,:,:,i_hg+1:) = 0.0_rp
826 
827  return
829  !-----------------------------------------------------------------------------
833  subroutine atmos_phy_mp_sn14_qtrc2qhyd( &
834  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
835  QTRC0, &
836  Qe )
838  n_hyd, &
839  i_hc, &
840  i_hr, &
841  i_hi, &
842  i_hs, &
843  i_hg
844  implicit none
845  integer, intent(in) :: KA, KS, KE
846  integer, intent(in) :: IA, IS, IE
847  integer, intent(in) :: JA, JS, JE
848 
849  real(RP), intent(in) :: QTRC0(ka,ia,ja,qa_mp-1) ! tracer mass concentration [kg/kg]
850 
851  real(RP), intent(out) :: Qe (ka,ia,ja,n_hyd) ! mixing ratio of each cateory [kg/kg]
852  !---------------------------------------------------------------------------
853 
854 !OCL XFILL
855  qe(:,:,:,i_hc) = qtrc0(:,:,:,i_mp_qc)
856 !OCL XFILL
857  qe(:,:,:,i_hr) = qtrc0(:,:,:,i_mp_qr)
858 !OCL XFILL
859  qe(:,:,:,i_hi) = qtrc0(:,:,:,i_mp_qi)
860 !OCL XFILL
861  qe(:,:,:,i_hs) = qtrc0(:,:,:,i_mp_qs)
862 !OCL XFILL
863  qe(:,:,:,i_hg) = qtrc0(:,:,:,i_mp_qg)
864 !OCL XFILL
865  qe(:,:,:,i_hg+1:) = 0.0_rp
866 
867  return
868  end subroutine atmos_phy_mp_sn14_qtrc2qhyd
869 
871  subroutine atmos_phy_mp_sn14_qtrc2nhyd( &
872  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
873  QTRC0, &
874  Ne )
876  n_hyd, &
877  i_hc, &
878  i_hr, &
879  i_hi, &
880  i_hs, &
881  i_hg
882  implicit none
883  integer, intent(in) :: KA, KS, KE
884  integer, intent(in) :: IA, IS, IE
885  integer, intent(in) :: JA, JS, JE
886 
887  real(RP), intent(in) :: QTRC0(ka,ia,ja,qa_mp-1) ! tracer mass concentration [kg/kg]
888 
889  real(RP), intent(out) :: Ne (ka,ia,ja,n_hyd) ! number density of each cateory [1/m3]
890  !---------------------------------------------------------------------------
891 
892 !OCL XFILL
893  ne(:,:,:,i_hc) = qtrc0(:,:,:,i_mp_nc)
894 !OCL XFILL
895  ne(:,:,:,i_hr) = qtrc0(:,:,:,i_mp_nr)
896 !OCL XFILL
897  ne(:,:,:,i_hi) = qtrc0(:,:,:,i_mp_ni)
898 !OCL XFILL
899  ne(:,:,:,i_hs) = qtrc0(:,:,:,i_mp_ns)
900 !OCL XFILL
901  ne(:,:,:,i_hg) = qtrc0(:,:,:,i_mp_ng)
902 !OCL XFILL
903  ne(:,:,:,i_hg+1:) = 0.0_rp
904 
905  return
906  end subroutine atmos_phy_mp_sn14_qtrc2nhyd
907 
908  subroutine atmos_phy_mp_sn14_qhyd2qtrc( &
909  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
910  Qe, &
911  QTRC, &
912  QNUM )
913  use scale_const, only: &
914  pi => const_pi
915  use scale_atmos_hydrometeor, only: &
916  n_hyd, &
917  i_hc, &
918  i_hr, &
919  i_hi, &
920  i_hs, &
921  i_hg, &
922  i_hh
923  implicit none
924  integer, intent(in) :: KA, KS, KE
925  integer, intent(in) :: IA, IS, IE
926  integer, intent(in) :: JA, JS, JE
927 
928  real(RP), intent(in) :: Qe(ka,ia,ja,n_hyd) ! mass ratio of each cateory [kg/kg]
929 
930  real(RP), intent(out) :: QTRC(ka,ia,ja,qa_mp-1) ! tracer mass concentration [kg/kg]
931 
932  real(RP), intent(in), optional :: QNUM(ka,ia,ja,n_hyd)
933 
934  real(RP), parameter :: Dc = 20.e-6_rp ! typical particle diameter for cloud [m]
935  real(RP), parameter :: Dr = 200.e-6_rp ! typical particle diameter for rain [m]
936  real(RP), parameter :: Di = 80.e-6_rp ! typical particle diameter for ice [m]
937  real(RP), parameter :: Ds = 80.e-6_rp ! typical particle diameter for snow [m]
938  real(RP), parameter :: Dg = 200.e-6_rp ! typical particle diameter for grapel [m]
939  real(RP), parameter :: RHOw = 1000.0_rp ! typical density for warm particles [kg/m3]
940  real(RP), parameter :: RHOf = 100.0_rp ! typical density for frozen particles [kg/m3]
941  real(RP), parameter :: RHOg = 400.0_rp ! typical density for grapel particles [kg/m3]
942  real(RP), parameter :: b = 3.0_rp ! assume spherical form
943 
944  real(RP) :: piov6
945 
946  integer :: k, i, j
947  !---------------------------------------------------------------------------
948 
949 
950 !OCL XFILL
951  !$omp parallel do
952  do j = js, je
953  do i = is, ie
954  do k = ks, ke
955  qtrc(k,i,j,i_mp_qc) = qe(k,i,j,i_hc)
956  end do
957  end do
958  end do
959 
960 !OCL XFILL
961  !$omp parallel do
962  do j = js, je
963  do i = is, ie
964  do k = ks, ke
965  qtrc(k,i,j,i_mp_qr) = qe(k,i,j,i_hr)
966  end do
967  end do
968  end do
969 
970 !OCL XFILL
971  !$omp parallel do
972  do j = js, je
973  do i = is, ie
974  do k = ks, ke
975  qtrc(k,i,j,i_mp_qi) = qe(k,i,j,i_hi)
976  end do
977  end do
978  end do
979 
980 !OCL XFILL
981  !$omp parallel do
982  do j = js, je
983  do i = is, ie
984  do k = ks, ke
985  qtrc(k,i,j,i_mp_qs) = qe(k,i,j,i_hs)
986  end do
987  end do
988  end do
989 
990 !OCL XFILL
991  !$omp parallel do
992  do j = js, je
993  do i = is, ie
994  do k = ks, ke
995  qtrc(k,i,j,i_mp_qg) = qe(k,i,j,i_hg) + qe(k,i,j,i_hh)
996  end do
997  end do
998  end do
999 
1000  if ( present(qnum) ) then
1001 
1002 !OCL XFILL
1003  !$omp parallel do
1004  do j = js, je
1005  do i = is, ie
1006  do k = ks, ke
1007  qtrc(k,i,j,i_mp_nc) = qnum(k,i,j,i_hc)
1008  end do
1009  end do
1010  end do
1011 
1012 !OCL XFILL
1013  !$omp parallel do
1014  do j = js, je
1015  do i = is, ie
1016  do k = ks, ke
1017  qtrc(k,i,j,i_mp_nr) = qnum(k,i,j,i_hr)
1018  end do
1019  end do
1020  end do
1021 
1022 !OCL XFILL
1023  !$omp parallel do
1024  do j = js, je
1025  do i = is, ie
1026  do k = ks, ke
1027  qtrc(k,i,j,i_mp_ni) = qnum(k,i,j,i_hi)
1028  end do
1029  end do
1030  end do
1031 
1032 !OCL XFILL
1033  !$omp parallel do
1034  do j = js, je
1035  do i = is, ie
1036  do k = ks, ke
1037  qtrc(k,i,j,i_mp_ns) = qnum(k,i,j,i_hs)
1038  end do
1039  end do
1040  end do
1041 
1042 !OCL XFILL
1043  !$omp parallel do
1044  do j = js, je
1045  do i = is, ie
1046  do k = ks, ke
1047  qtrc(k,i,j,i_mp_ng) = qnum(k,i,j,i_hg) + qnum(k,i,j,i_hh)
1048  end do
1049  end do
1050  end do
1051 
1052  else
1053  piov6 = pi / 6.0_rp
1054 
1055 !OCL XFILL
1056  !$omp parallel do
1057  do j = js, je
1058  do i = is, ie
1059  do k = ks, ke
1060  qtrc(k,i,j,i_mp_nc) = qtrc(k,i,j,i_mp_qc) / ( (piov6*rhow) * dc**b )
1061  end do
1062  end do
1063  end do
1064 
1065 !OCL XFILL
1066  !$omp parallel do
1067  do j = js, je
1068  do i = is, ie
1069  do k = ks, ke
1070  qtrc(k,i,j,i_mp_nr) = qtrc(k,i,j,i_mp_qr) / ( (piov6*rhow) * dr**b )
1071  end do
1072  end do
1073  end do
1074 
1075 !OCL XFILL
1076  !$omp parallel do
1077  do j = js, je
1078  do i = is, ie
1079  do k = ks, ke
1080  qtrc(k,i,j,i_mp_ni) = qtrc(k,i,j,i_mp_qi) / ( (piov6*rhof) * di**b )
1081  end do
1082  end do
1083  end do
1084 
1085 !OCL XFILL
1086  !$omp parallel do
1087  do j = js, je
1088  do i = is, ie
1089  do k = ks, ke
1090  qtrc(k,i,j,i_mp_ns) = qtrc(k,i,j,i_mp_qs) / ( (piov6*rhof) * ds**b )
1091  end do
1092  end do
1093  end do
1094 
1095 !OCL XFILL
1096  !$omp parallel do
1097  do j = js, je
1098  do i = is, ie
1099  do k = ks, ke
1100  qtrc(k,i,j,i_mp_ng) = qtrc(k,i,j,i_mp_qg) / ( (piov6*rhog) * dg**b )
1101  end do
1102  end do
1103  end do
1104 
1105  end if
1106 
1107  return
1108  end subroutine atmos_phy_mp_sn14_qhyd2qtrc
1109 
1110 
1111  ! private
1112  !-----------------------------------------------------------------------------
1113  subroutine mp_sn14_init
1114  use scale_prc, only: &
1115  prc_abort
1116  use scale_specfunc, only: &
1117  gammafunc => sf_gamma
1118  implicit none
1119 
1120  real(RP) :: w1(hydro_max)
1121  real(RP) :: w2(hydro_max)
1122  real(RP) :: w3(hydro_max)
1123  real(RP) :: w4(hydro_max)
1124  real(RP) :: w5(hydro_max)
1125  real(RP) :: w6(hydro_max)
1126  real(RP) :: w7(hydro_max)
1127  real(RP) :: w8(hydro_max)
1128 
1129  ! work for calculation of capacity, Mitchell and Arnott (1994) , eq.(9)
1130  real(RP) :: ar_ice_fix = 0.7_rp
1131  real(RP) :: wcap1, wcap2
1132  ! work for ventilation coefficient
1133  logical :: flag_vent0(hydro_max), flag_vent1(hydro_max)
1134  integer :: ierr
1135  integer :: iw, ia, ib
1136  integer :: n
1137  !
1138  namelist / param_atmos_phy_mp_sn14_init / &
1139  opt_debug, &
1140  opt_debug_tem, &
1141  opt_debug_inc, &
1142  opt_debug_act, &
1143  opt_debug_ree, &
1144  opt_debug_bcs, &
1145  ntmax_phase_change, &
1146  ntmax_collection
1147  !
1148  namelist / param_atmos_phy_mp_sn14_particles / &
1149  a_m, b_m, alpha_v, beta_v, gamma_v, &
1150  alpha_vn, beta_vn, &
1151  a_area, b_area, cap, &
1152  nu, mu, &
1153  opt_m96_column_ice, &
1154  opt_m96_ice, &
1155  ar_ice_fix
1156  real(RP), parameter :: eps_gamma=1.e-30_rp
1157 
1158  a_m(:) = undef8
1159  b_m(:) = undef8
1160  alpha_v(:,:) = undef8
1161  beta_v(:,:) = undef8
1162  alpha_vn(:,:) = undef8
1163  beta_vn(:,:) = undef8
1164  gamma_v(:) = undef8
1165  a_d2vt(:) = undef8
1166  b_d2vt(:) = undef8
1167  a_area(:) = undef8
1168  b_area(:) = undef8
1169  ax_area(:) = undef8
1170  bx_area(:) = undef8
1171  a_rea(:) = undef8
1172  b_rea(:) = undef8
1173  a_rea2(:) = undef8
1174  b_rea2(:) = undef8
1175  a_rea3(:) = undef8
1176  b_rea3(:) = undef8
1177  nu(:) = undef8
1178  mu(:) = undef8
1179  cap(:) = undef8
1180  coef_m2(:) = undef8
1181  coef_dave_n(:) = undef8
1182  coef_dave_l(:) = undef8
1183  coef_d(:) = undef8
1184  coef_d3(:) = undef8
1185  coef_d6(:) = undef8
1186  coef_d2v(:) = undef8
1187  coef_md2v(:) = undef8
1188  coef_r2(:) = undef8
1189  coef_r3(:) = undef8
1190  coef_re(:) = undef8
1191  coef_rea2(:) = undef8
1192  coef_rea3(:) = undef8
1193  coef_a(:) = undef8
1194 ! slope_A(:) = UNDEF8
1195  coef_lambda(:) = undef8
1196  coef_vt0(:,:) = undef8
1197  coef_vt1(:,:) = undef8
1198  delta_b0(:) = undef8
1199  delta_b1(:) = undef8
1200  delta_ab0(:,:) = undef8
1201  delta_ab1(:,:) = undef8
1202  theta_b0(:) = undef8
1203  theta_b1(:) = undef8
1204  theta_ab0(:,:) = undef8
1205  theta_ab1(:,:) = undef8
1206  !
1207  ah_vent(:,:) = undef8
1208  ah_vent0(:,:) = undef8
1209  ah_vent1(:,:) = undef8
1210  bh_vent(:,:) = undef8
1211  bh_vent0(:,:) = undef8
1212  bh_vent1(:,:) = undef8
1213 
1214  !--- read namelist
1215  rewind(io_fid_conf)
1216  read(io_fid_conf,nml=param_atmos_phy_mp_sn14_init,iostat=ierr)
1217 
1218  if( ierr < 0 ) then !--- missing
1219  log_info("ATMOS_PHY_MP_sn14_init",*) 'Not found namelist. Default used.'
1220  elseif( ierr > 0 ) then !--- fatal error
1221  log_error("ATMOS_PHY_MP_sn14_init",*) 'Not appropriate names in namelist PARAM_ATMOS_PHY_MP_SN14_init. Check!'
1222  call prc_abort
1223  endif
1224  log_nml(param_atmos_phy_mp_sn14_init)
1225 
1226  !
1227  ! default setting
1228  !
1229  ! Area parameters with mks unit originated by Mitchell(1996)
1230  a_area(i_mp_qc) = pi/4.0_rp ! sphere
1231  a_area(i_mp_qr) = pi/4.0_rp ! sphere
1232  a_area(i_mp_qi) = 0.65_rp*1.e-4_rp*100.0_rp**(2.00_rp) ! Mitchell(1996), Hexagonal Plate
1233  a_area(i_mp_qs) = 0.2285_rp*1.e-4_rp*100.0_rp**(1.88_rp) ! Mitchell(1996), Aggregates
1234  a_area(i_mp_qg) = 0.50_rp*1.e-4_rp*100.0_rp**(2.0_rp) ! Mitchell(1996), Lump Graupel
1235  b_area(i_mp_qc) = 2.0_rp
1236  b_area(i_mp_qr) = 2.0_rp
1237  b_area(i_mp_qi) = 2.0_rp
1238  b_area(i_mp_qs) = 1.88_rp
1239  b_area(i_mp_qg) = 2.0_rp
1240  !
1241  ! Seifert and Beheng(2006), Table. 1 or List of symbols
1242  !----------------------------------------------------------
1243  ! Diameter-Mass relationship
1244  ! D = a * x^b
1245  a_m(i_mp_qc) = 0.124_rp
1246  a_m(i_mp_qr) = 0.124_rp
1247  a_m(i_mp_qi) = 0.217_rp
1248  a_m(i_mp_qs) = 8.156_rp
1249  a_m(i_mp_qg) = 0.190_rp
1250  b_m(i_mp_qc) = 1.0_rp/3.0_rp
1251  b_m(i_mp_qr) = 1.0_rp/3.0_rp
1252  b_m(i_mp_qi) = 0.302_rp
1253  b_m(i_mp_qs) = 0.526_rp
1254  b_m(i_mp_qg) = 0.323_rp
1255  !----------------------------------------------------------
1256  ! Terminal velocity-Mass relationship
1257  ! vt = alpha * x^beta * (rho0/rho)^gamma
1258  alpha_v(i_mp_qc,:)= 3.75e+5_rp
1259  alpha_v(i_mp_qr,:)= 159.0_rp ! not for sedimantation
1260  alpha_v(i_mp_qi,:)= 317.0_rp
1261  alpha_v(i_mp_qs,:)= 27.70_rp
1262  alpha_v(i_mp_qg,:)= 40.0_rp
1263  beta_v(i_mp_qc,:) = 2.0_rp/3.0_rp
1264  beta_v(i_mp_qr,:) = 0.266_rp ! not for sedimantation
1265  beta_v(i_mp_qi,:) = 0.363_rp
1266  beta_v(i_mp_qs,:) = 0.216_rp
1267  beta_v(i_mp_qg,:) = 0.230_rp
1268  gamma_v(i_mp_qc) = 1.0_rp
1269  ! This is high Reynolds number limit(Beard 1980)
1270  gamma_v(i_mp_qr) = 1.0_rp/2.0_rp
1271  gamma_v(i_mp_qi) = 1.0_rp/2.0_rp
1272  gamma_v(i_mp_qs) = 1.0_rp/2.0_rp
1273  gamma_v(i_mp_qg) = 1.0_rp/2.0_rp
1274  !----------------------------------------------------------
1275  ! DSD parameters
1276  ! f(x) = A x^nu exp( -lambda x^mu )
1277  ! Gamma Disribution : mu=1 , nu:arbitrary
1278  ! Marshall-Palmer Distribution: mu=1/3, nu:-2/3
1279  ! In the case of MP, f(D) dD = f(x)dx
1280  ! f(x) = c * f(D)/D^2 (c:coefficient)
1281  nu(i_mp_qc) = 1.0_rp ! arbitrary for Gamma
1282  nu(i_mp_qr) = -1.0_rp/3.0_rp ! nu(diameter)=1, equilibrium condition.
1283  nu(i_mp_qi) = 1.0_rp !
1284  nu(i_mp_qs) = 1.0_rp !
1285  nu(i_mp_qg) = 1.0_rp !
1286  !
1287  mu(i_mp_qc) = 1.0_rp ! Gamma
1288  mu(i_mp_qr) = 1.0_rp/3.0_rp ! Marshall Palmer
1289  mu(i_mp_qi) = 1.0_rp/3.0_rp !
1290  mu(i_mp_qs) = 1.0_rp/3.0_rp !
1291  mu(i_mp_qg) = 1.0_rp/3.0_rp !
1292  !----------------------------------------------------------
1293  ! Geomeries for diffusion growth
1294  ! Pruppacher and Klett(1997), (13-77)-(13-80) and
1295  ! originally derived by McDonald(1963b)
1296  ! sphere: cap=2
1297  ! plate : cap=pi
1298  ! needle with aspect ratio a/b
1299  ! : cap=log(2*a/b)
1300  cap(i_mp_qc) = 2.0_rp ! sphere
1301  cap(i_mp_qr) = 2.0_rp ! sphere
1302  cap(i_mp_qi) = pi ! hexagonal plate
1303  cap(i_mp_qs) = 2.0_rp ! mix aggregates
1304  cap(i_mp_qg) = 2.0_rp ! lump
1305  !
1306  alpha_vn(:,:) = alpha_v(:,:)
1307  beta_vn(:,:) = beta_v(:,:)
1308  !------------------------------------------------------------------------
1309  !
1310  ! additional setting
1311  !
1312 
1313  !--- read namelist
1314  rewind(io_fid_conf)
1315  read(io_fid_conf,nml=param_atmos_phy_mp_sn14_particles,iostat=ierr)
1316 
1317  if( ierr < 0 ) then !--- missing
1318  log_info("ATMOS_PHY_MP_sn14_init",*) 'Not found namelist. Default used.'
1319  elseif( ierr > 0 ) then !--- fatal error
1320  log_error("ATMOS_PHY_MP_sn14_init",*) 'Not appropriate names in namelist PARAM_ATMOS_PHY_MP_SN14_particles. Check!'
1321  call prc_abort
1322  endif
1323  log_nml(param_atmos_phy_mp_sn14_particles)
1324 
1325  ! [Add] 10/08/03 T.Mitsui
1326  ! particles shapes are
1327  if( opt_m96_ice ) then
1328  ! ice is randomly oriented Hexagonal plate (Auer and Veal 1970, Takano and Liou 1995, Mitchell 1996)
1329  ! snow is assemblages of planar polycrystals(Mitchell 1996)
1330  ! graupel is Lump graupel(R4b) is assumed(Mitchell 1996)
1331  a_area(i_mp_qi) = 0.120284936_rp
1332  a_area(i_mp_qs) = 0.131488_rp
1333  a_area(i_mp_qg) = 0.5_rp
1334  b_area(i_mp_qi) = 1.850000_rp
1335  b_area(i_mp_qs) = 1.880000_rp
1336  b_area(i_mp_qg) = 2.0_rp
1337  a_m(i_mp_qi) = 1.23655360084766_rp
1338  a_m(i_mp_qs) = a_m(i_mp_qi)
1339  a_m(i_mp_qg) = 0.346111225718402_rp
1340  b_m(i_mp_qi) = 0.408329930583912_rp
1341  b_m(i_mp_qs) = b_m(i_mp_qi)
1342  b_m(i_mp_qg) = 0.357142857142857_rp
1343  !
1344  if( opt_m96_column_ice )then
1345  d0_ni=240.49e-6_rp ! this is column
1346  d0_li=330.09e-6_rp ! this is column
1347  a_area(i_mp_qi)= (0.684_rp*1.e-4_rp)*10.0_rp**(2.0_rp*2.00_rp)
1348  b_area(i_mp_qi)= 2.0_rp
1349  a_m(i_mp_qi) = 0.19834046116844_rp
1350  b_m(i_mp_qi) = 0.343642611683849_rp
1351  ! [Add] 11/08/30 T.Mitsui
1352  ! approximated by the capacity for prolate spheroid with constant aspect ratio
1353  wcap1 = sqrt(1.0_rp-ar_ice_fix**2)
1354  wcap2 = log( (1.0_rp+wcap1)/ar_ice_fix )
1355  cap(i_mp_qi) = 2.0_rp*wcap2/wcap1
1356  !
1357  end if
1358  !
1359  ! These value are derived by least-square fitting in the range
1360  ! qi [100um:1000um] in diameter
1361  ! qs [100um:1000um] in diameter
1362  ! qg [200um:2000um] in diameter
1363  ! small branch , large branch
1364  alpha_v(i_mp_qi,:) =(/ 5798.60107421875_rp, 167.347076416016_rp/)
1365  alpha_vn(i_mp_qi,:) =(/ 12408.177734375_rp, 421.799865722656_rp/)
1366  if( opt_m96_column_ice )then
1367  alpha_v(i_mp_qi,:) = (/2901.0_rp, 32.20_rp/)
1368  alpha_vn(i_mp_qi,:) = (/9675.2_rp, 64.16_rp/)
1369  end if
1370  alpha_v(i_mp_qs,:) =(/ 15173.3916015625_rp, 305.678619384766_rp/)
1371  alpha_vn(i_mp_qs,:) =(/ 29257.1601562500_rp, 817.985717773438_rp/)
1372  alpha_v(i_mp_qg,:) =(/ 15481.6904296875_rp, 311.642242431641_rp/)
1373  alpha_vn(i_mp_qg,:) =(/ 27574.6562500000_rp, 697.536132812500_rp/)
1374  !
1375  beta_v(i_mp_qi,:) =(/ 0.504873454570770_rp, 0.324817866086960_rp/)
1376  beta_vn(i_mp_qi,:) =(/ 0.548495233058929_rp, 0.385287821292877_rp/)
1377  if( opt_m96_column_ice )then
1378  beta_v(i_mp_qi,:) =(/ 0.465552181005478_rp, 0.223826110363007_rp/)
1379  beta_vn(i_mp_qi,:) =(/ 0.530453503131866_rp, 0.273761242628098_rp/)
1380  end if
1381  beta_v(i_mp_qs,:) =(/ 0.528109610080719_rp, 0.329863965511322_rp/)
1382  beta_vn(i_mp_qs,:) =(/ 0.567154467105865_rp, 0.393876969814301_rp/)
1383  beta_v(i_mp_qg,:) =(/ 0.534656763076782_rp, 0.330253750085831_rp/)
1384  beta_vn(i_mp_qg,:) =(/ 0.570551633834839_rp, 0.387124240398407_rp/)
1385  end if
1386  !
1387  ! area-diameter relation => area-mass relation
1388  ax_area(:) = a_area(:)*a_m(:)**b_area(:)
1389  bx_area(:) = b_area(:)*b_m(:)
1390  !
1391  ! radius of equivalent area - m ass relation
1392  ! pi*rea**2 = ax_area*x**bx_area
1393  a_rea(:) = sqrt(ax_area(:)/pi)
1394  b_rea(:) = bx_area(:)/2.0_rp
1395  a_rea2(:) = a_rea(:)**2
1396  b_rea2(:) = b_rea(:)*2.0_rp
1397  a_rea3(:) = a_rea(:)**3
1398  b_rea3(:) = b_rea(:)*3.0_rp
1399  !
1400  a_d2vt(:)=alpha_v(:,2)*(1.0_rp/alpha_v(:,2))**(beta_v(:,2)/b_m(:))
1401  b_d2vt(:)=(beta_v(:,2)/b_m(:))
1402  !
1403  ! Calculation of Moment Coefficient
1404  !
1405  w1(:) = 0.0_rp
1406  w2(:) = 0.0_rp
1407  w3(:) = 0.0_rp
1408  w4(:) = 0.0_rp
1409  w5(:) = 0.0_rp
1410  w6(:) = 0.0_rp
1411  w7(:) = 0.0_rp
1412  w8(:) = 0.0_rp
1413  !-------------------------------------------------------
1414  ! moment coefficient
1415  ! SB06 (82)
1416  ! M^n /= coef_mn * N * (L/N)**n
1417  ! M^2 = Z = coef_m2 * N *(L/N)**2
1418  ! a*M^b = a*integral x^b f(x) dx = ave D
1419  do iw=1, hydro_max
1420  n = 2
1421  w1(iw) = gammafunc( (n+nu(iw)+1.0_rp)/mu(iw) )
1422  w2(iw) = gammafunc( (nu(iw)+1.0_rp)/mu(iw) )
1423  w3(iw) = gammafunc( (nu(iw)+2.0_rp)/mu(iw) )
1424  coef_m2(iw) = w1(iw)/w2(iw)*( w2(iw)/w3(iw) )**n
1425  !
1426  w4(iw) = gammafunc( (b_m(iw)+nu(iw)+1.0_rp)/mu(iw) )
1427  coef_d(iw) = a_m(iw) * w4(iw)/w2(iw)*( w2(iw)/w3(iw) )**b_m(iw)
1428  w5(iw) = gammafunc( (2.0_rp*b_m(iw)+beta_v(iw,2)+nu(iw)+1.0_rp)/mu(iw) )
1429  w6(iw) = gammafunc( (3.0_rp*b_m(iw)+beta_v(iw,2)+nu(iw)+1.0_rp)/mu(iw) )
1430  coef_d2v(iw) = a_m(iw) * w6(iw)/w5(iw)* ( w2(iw)/w3(iw) )**b_m(iw)
1431  coef_md2v(iw)= w5(iw)/w2(iw)* ( w2(iw)/w3(iw) )**(2.0_rp*b_m(iw)+beta_v(iw,2))
1432  ! 09/04/14 [Add] T.Mitsui, volume and radar reflectivity
1433  w7(iw) = gammafunc( (3.0_rp*b_m(iw)+nu(iw)+1.0_rp)/mu(iw) )
1434  coef_d3(iw) = a_m(iw)**3 * w7(iw)/w2(iw)*( w2(iw)/w3(iw) )**(3.0_rp*b_m(iw))
1435  w8(iw) = gammafunc( (6.0_rp*b_m(iw)+nu(iw)+1.0_rp)/mu(iw) )
1436  coef_d6(iw) = a_m(iw)**6 * w8(iw)/w2(iw)*( w2(iw)/w3(iw) )**(6.0_rp*b_m(iw))
1437  end do
1438  !
1439  coef_deplc = coef_d(i_mp_qc)/a_m(i_mp_qc)
1440  !-------------------------------------------------------
1441  ! coefficient of 2nd and 3rd moments for effective radius
1442  ! for spherical particle
1443  do iw=1, hydro_max
1444  ! integ r^2 f(x)dx
1445  w1(iw) = gammafunc( (2.0_rp*b_m(iw)+nu(iw)+1.0_rp)/mu(iw) )
1446  w2(iw) = gammafunc( (nu(iw)+1.0_rp)/mu(iw) )
1447  w3(iw) = gammafunc( (nu(iw)+2.0_rp)/mu(iw) )
1448  ! integ r^3 f(x)dx
1449  w4(iw) = gammafunc( (3.0_rp*b_m(iw)+nu(iw)+1.0_rp)/mu(iw) )
1450  !
1451  coef_r2(iw)=w1(iw)/w2(iw)*( w2(iw)/w3(iw) )**(2.0_rp*b_m(iw))
1452  coef_r3(iw)=w4(iw)/w2(iw)*( w2(iw)/w3(iw) )**(3.0_rp*b_m(iw))
1453  coef_re(iw)=coef_r3(iw)/coef_r2(iw)
1454  !
1455  end do
1456  !-------------------------------------------------------
1457  ! coefficient for effective radius of equivalent area and
1458  ! coefficient for volume of equivalent area
1459  do iw=1, hydro_max
1460  w1(iw) = gammafunc( (nu(iw)+1.0_rp)/mu(iw) )
1461  w2(iw) = gammafunc( (nu(iw)+2.0_rp)/mu(iw) )
1462  w3(iw) = gammafunc( (b_rea2(iw)+nu(iw)+1.0_rp)/mu(iw) )
1463  w4(iw) = gammafunc( (b_rea3(iw)+nu(iw)+1.0_rp)/mu(iw) )
1464  !
1465  coef_rea2(iw) = w3(iw)/w1(iw)*( w1(iw)/w2(iw) )**b_rea2(iw)
1466  coef_rea3(iw) = w4(iw)/w1(iw)*( w1(iw)/w2(iw) )**b_rea3(iw)
1467  end do
1468  !-------------------------------------------------------
1469  ! coefficient of gamma-distribution
1470  ! SB06(80)
1471  do iw=1, hydro_max
1472  w1(iw) = gammafunc( (nu(iw) + 1.0_rp)/mu(iw) )
1473  w2(iw) = gammafunc( (nu(iw) + 2.0_rp)/mu(iw) )
1474  coef_a(iw) = mu(iw)/w1(iw)
1475 ! slope_A(iw) = w1(iw)
1476  coef_lambda(iw) = (w1(iw)/w2(iw))**(-mu(iw))
1477  end do
1478  !-------------------------------------------------------
1479  ! coefficient for terminal velocity in sedimentation
1480  ! SB06(78)
1481  do ia=1,2
1482  do iw=1, hydro_max
1483  n = 0
1484  w1(iw) = gammafunc( (beta_vn(iw,ia) + nu(iw) + 1.0_rp + n)/mu(iw) )
1485  w2(iw) = gammafunc( ( nu(iw) + 1.0_rp + n)/mu(iw) )
1486  w3(iw) = gammafunc( (nu(iw) + 1.0_rp)/mu(iw) )
1487  w4(iw) = gammafunc( (nu(iw) + 2.0_rp)/mu(iw) )
1488  ! coefficient of terminal velocity for number
1489  coef_vt0(iw,ia)=alpha_vn(iw,ia)*w1(iw)/w2(iw)*(w3(iw)/w4(iw))**beta_vn(iw,ia)
1490  n = 1
1491  w1(iw) = gammafunc( (beta_v(iw,ia) + nu(iw) + 1.0_rp + n)/mu(iw) )
1492  w2(iw) = gammafunc( ( nu(iw) + 1.0_rp + n)/mu(iw) )
1493  ! coefficient of terminal velocity for mass
1494  coef_vt1(iw,ia)=alpha_v(iw,ia)*w1(iw)/w2(iw)*(w3(iw)/w4(iw))**beta_v(iw,ia)
1495  end do
1496  end do
1497  ! coefficient for weighted diameter used in calculation of terminal velocity
1498  do iw=1, hydro_max
1499  w1(iw) = gammafunc( ( b_m(iw) + nu(iw) + 1.0_rp)/mu(iw) )
1500  w2(iw) = gammafunc( (1.0_rp + b_m(iw) + nu(iw) + 1.0_rp)/mu(iw) )
1501  w3(iw) = gammafunc( (nu(iw) + 1.0_rp)/mu(iw) )
1502  w4(iw) = gammafunc( (nu(iw) + 2.0_rp)/mu(iw) )
1503  coef_dave_n(iw) = (w1(iw)/w3(iw))*(w3(iw)/w4(iw))** b_m(iw)
1504  coef_dave_l(iw) = (w2(iw)/w3(iw))*(w3(iw)/w4(iw))**(1.0_rp+b_m(iw))
1505  end do
1506  !-------------------------------------------------------
1507  !
1508  ah_vent(i_mp_qc,1:2) = (/1.0000_rp,1.0000_rp/) ! no effect
1509  ah_vent(i_mp_qr,1:2) = (/1.0000_rp,0.780_rp/)
1510  ah_vent(i_mp_qi,1:2) = (/1.0000_rp,0.860_rp/)
1511  ah_vent(i_mp_qs,1:2) = (/1.0000_rp,0.780_rp/)
1512  ah_vent(i_mp_qg,1:2) = (/1.0000_rp,0.780_rp/)
1513  bh_vent(i_mp_qc,1:2) = (/0.0000_rp,0.0000_rp/)
1514  bh_vent(i_mp_qr,1:2) = (/0.108_rp,0.308_rp/)
1515  bh_vent(i_mp_qi,1:2) = (/0.140_rp,0.280_rp/)
1516  bh_vent(i_mp_qs,1:2) = (/0.108_rp,0.308_rp/)
1517  bh_vent(i_mp_qg,1:2) = (/0.108_rp,0.308_rp/)
1518  !
1519  do iw=1, hydro_max
1520  n = 0
1521  if( (nu(iw) + b_m(iw) + n) > eps_gamma )then
1522  w1(iw) = gammafunc( (nu(iw) + b_m(iw) + n)/mu(iw) )
1523  w2(iw) = gammafunc( (nu(iw) + 1.0_rp)/mu(iw) )
1524  w3(iw) = gammafunc( (nu(iw) + 2.0_rp)/mu(iw) )
1525  ah_vent0(iw,1)= ah_vent(iw,1)*(w1(iw)/w2(iw))*(w2(iw)/w3(iw))**(b_m(iw)+n-1.0_rp)
1526  ah_vent0(iw,2)= ah_vent(iw,2)*(w1(iw)/w2(iw))*(w2(iw)/w3(iw))**(b_m(iw)+n-1.0_rp)
1527  flag_vent0(iw)=.true.
1528  else
1529  ah_vent0(iw,1)= 1.0_rp
1530  ah_vent0(iw,2)= 1.0_rp
1531  flag_vent0(iw)=.false.
1532  end if
1533  n = 1
1534  if( (nu(iw) + b_m(iw) + n) > eps_gamma )then
1535  w1(iw) = gammafunc( (nu(iw) + b_m(iw) + n)/mu(iw) )
1536  w2(iw) = gammafunc( (nu(iw) + 1.0_rp)/mu(iw) )
1537  w3(iw) = gammafunc( (nu(iw) + 2.0_rp)/mu(iw) )
1538  ah_vent1(iw,1)= ah_vent(iw,1)*(w1(iw)/w2(iw))*(w2(iw)/w3(iw))**(b_m(iw)+n-1.0_rp)
1539  ah_vent1(iw,2)= ah_vent(iw,2)*(w1(iw)/w2(iw))*(w2(iw)/w3(iw))**(b_m(iw)+n-1.0_rp)
1540  flag_vent1(iw)=.true.
1541  else
1542  ah_vent1(iw,1)= 1.0_rp
1543  ah_vent1(iw,2)= 1.0_rp
1544  flag_vent1(iw)=.true.
1545  end if
1546  end do
1547  do iw=1, hydro_max
1548  n = 0
1549  if( (nu(iw) + 1.5_rp*b_m(iw) + 0.5_rp*beta_v(iw,1) + n) < eps_gamma )then
1550  flag_vent0(iw)=.false.
1551  end if
1552  if(flag_vent0(iw))then
1553  w1(iw) = gammafunc( (nu(iw) + 1.5_rp*b_m(iw) + 0.5_rp*beta_v(iw,1) + n)/mu(iw) )
1554  w2(iw) = gammafunc( (nu(iw) + 1.0_rp)/mu(iw) )
1555  w3(iw) = gammafunc( (nu(iw) + 2.0_rp)/mu(iw) )
1556  ! [Add] 11/08/30 T.Mitsui
1557  w4(iw) = gammafunc( (nu(iw) + 2.0_rp*b_m(iw) + beta_v(iw,1) + n)/mu(iw) )
1558  bh_vent0(iw,1)=bh_vent(iw,1)*(w4(iw)/w2(iw))*(w2(iw)/w3(iw))**(2.00_rp*b_m(iw)+beta_v(iw,1)+n-1.0_rp)
1559  w5(iw) = gammafunc( (nu(iw) + 1.5_rp*b_m(iw) + 0.5_rp*beta_v(iw,2) + n)/mu(iw) )
1560  bh_vent0(iw,2)=bh_vent(iw,2)*(w5(iw)/w2(iw))*(w2(iw)/w3(iw))**(1.5_rp*b_m(iw)+0.5_rp*beta_v(iw,2)+n-1.0_rp)
1561  else
1562  bh_vent0(iw,1) = 0.0_rp
1563  bh_vent0(iw,2) = 0.0_rp
1564  end if
1565  !
1566  n = 1
1567  if( (nu(iw) + 1.5_rp*b_m(iw) + 0.5_rp*beta_v(iw,1) + n) < eps_gamma )then
1568  flag_vent1(iw)=.false.
1569  end if
1570  if(flag_vent1(iw))then
1571  w1(iw) = gammafunc( (nu(iw) + 1.5_rp*b_m(iw) + 0.5_rp*beta_v(iw,1) + n)/mu(iw) )
1572  w2(iw) = gammafunc( (nu(iw) + 1.0_rp)/mu(iw) )
1573  w3(iw) = gammafunc( (nu(iw) + 2.0_rp)/mu(iw) )
1574  ! [Add] 11/08/30 T.Mitsui
1575  w4(iw) = gammafunc( (nu(iw) + 2.0_rp*b_m(iw) + beta_v(iw,1) + n)/mu(iw) )
1576  bh_vent1(iw,1)=bh_vent(iw,1)*(w4(iw)/w2(iw))*(w2(iw)/w3(iw))**(2.00_rp*b_m(iw)+beta_v(iw,1)+n-1.0_rp)
1577  !
1578  w5(iw) = gammafunc( (nu(iw) + 1.5_rp*b_m(iw) + 0.5_rp*beta_v(iw,2) + n)/mu(iw) )
1579  bh_vent1(iw,2)=bh_vent(iw,2)*(w5(iw)/w2(iw))*(w2(iw)/w3(iw))**(1.5_rp*b_m(iw)+0.5_rp*beta_v(iw,2)+n-1.0_rp)
1580  else
1581  bh_vent1(iw,1) = 0.0_rp
1582  bh_vent1(iw,2) = 0.0_rp
1583  end if
1584  end do
1585  !-------------------------------------------------------
1586  ! coefficient for collision process
1587  ! stochastic coefficient for collision cross section
1588  ! sb06 (90) -- self collection
1589  do iw=1, hydro_max
1590  n = 0
1591  w1(iw) = gammafunc( (2.0_rp*b_rea(iw) + nu(iw) + 1.0_rp + n)/mu(iw) )
1592  w2(iw) = gammafunc( (nu(iw) + 1.0_rp)/mu(iw) )
1593  w3(iw) = gammafunc( (nu(iw) + 2.0_rp)/mu(iw) )
1594  delta_b0(iw) = w1(iw)/w2(iw) &
1595  *( w2(iw)/w3(iw) )**(2.0_rp*b_rea(iw) + n)
1596  n = 1
1597  w1(iw) = gammafunc( (2.0_rp*b_rea(iw) + nu(iw) + 1.0_rp + n)/mu(iw) )
1598  delta_b1(iw) = w1(iw)/w2(iw) &
1599  *( w2(iw)/w3(iw) )**(2.0_rp*b_rea(iw) + n)
1600  end do
1601  ! stochastic coefficient for collision cross section
1602  ! sb06(91) -- riming( collide with others )
1603  do iw=1, hydro_max
1604  n = 0
1605  w1(iw) = gammafunc( (b_rea(iw) + nu(iw) + 1.0_rp + n)/mu(iw) )
1606  w2(iw) = gammafunc( (nu(iw) + 1.0_rp)/mu(iw) )
1607  w3(iw) = gammafunc( (nu(iw) + 2.0_rp)/mu(iw) )
1608  w4(iw) = gammafunc( (b_rea(iw) + nu(iw) + 1.0_rp )/mu(iw) )
1609  n = 1
1610  w5(iw) = gammafunc( (b_rea(iw) + nu(iw) + 1.0_rp + n)/mu(iw) )
1611  end do
1612  ! ia > ib ( larger particles "a" catch smaller particles "b" )
1613  do ia=1, hydro_max
1614  do ib=1, hydro_max
1615  n=0 !
1616  ! NOTE, collected particle has a moment of n.
1617  ! collecting particle has only number(n=0).
1618  delta_ab0(ia,ib) = 2.0_rp*(w1(ib)/w2(ib))*(w4(ia)/w2(ia)) &
1619  * ( w2(ib)/w3(ib) )**(b_rea(ib)+n) &
1620  * ( w2(ia)/w3(ia) )**(b_rea(ia) )
1621  n=1 !
1622  delta_ab1(ia,ib) = 2.0_rp*(w5(ib)/w2(ib))*(w4(ia)/w2(ia)) &
1623  * ( w2(ib)/w3(ib) )**(b_rea(ib)+n) &
1624  * ( w2(ia)/w3(ia) )**(b_rea(ia) )
1625  end do
1626  end do
1627  ! stochastic coefficient for terminal velocity
1628  ! sb06(92) -- self collection
1629  ! assuming equivalent area circle.
1630  do iw=1, hydro_max
1631  n = 0
1632  w1(iw) = gammafunc( (2.0_rp*beta_v(iw,2) + 2.0_rp*b_rea(iw) + nu(iw) + 1.0_rp + n)/mu(iw) )
1633  w2(iw) = gammafunc( ( 2.0_rp*b_rea(iw) + nu(iw) + 1.0_rp + n)/mu(iw) )
1634  w3(iw) = gammafunc( (nu(iw) + 1.0_rp)/mu(iw) )
1635  w4(iw) = gammafunc( (nu(iw) + 2.0_rp)/mu(iw) )
1636  theta_b0(iw) = w1(iw)/w2(iw) * ( w3(iw)/w4(iw) )**(2.0_rp*beta_v(iw,2))
1637  n = 1
1638  w1(iw) = gammafunc( (2.0_rp*beta_v(iw,2) + 2.0_rp*b_rea(iw) + nu(iw) + 1.0_rp + n)/mu(iw) )
1639  w2(iw) = gammafunc( ( 2.0_rp*b_rea(iw) + nu(iw) + 1.0_rp + n)/mu(iw) )
1640  theta_b1(iw) = w1(iw)/w2(iw) * ( w3(iw)/w4(iw) )**(2.0_rp*beta_v(iw,2))
1641  end do
1642  !
1643  ! stochastic coefficient for terminal velocity
1644  ! sb06(93) -- riming( collide with others )
1645  do iw=1, hydro_max
1646  n = 0
1647  w1(iw) = gammafunc( (beta_v(iw,2) + 2.0_rp*b_rea(iw) + nu(iw) + 1.0_rp + n)/mu(iw) )
1648  w2(iw) = gammafunc( ( 2.0_rp*b_rea(iw) + nu(iw) + 1.0_rp + n)/mu(iw) )
1649  w3(iw) = gammafunc( (beta_v(iw,2) + 2.0_rp*b_rea(iw) + nu(iw) + 1.0_rp )/mu(iw) )
1650  w4(iw) = gammafunc( ( 2.0_rp*b_rea(iw) + nu(iw) + 1.0_rp )/mu(iw) )
1651  !
1652  w5(iw) = gammafunc( (nu(iw) + 1.0_rp)/mu(iw) )
1653  w6(iw) = gammafunc( (nu(iw) + 2.0_rp)/mu(iw) )
1654  n = 1
1655  w7(iw) = gammafunc( (beta_v(iw,2) + b_rea(iw) + nu(iw) + 1.0_rp + n)/mu(iw) )
1656  w8(iw) = gammafunc( ( b_rea(iw) + nu(iw) + 1.0_rp + n)/mu(iw) )
1657  end do
1658  ! ia > ib ( larger particles "a" catch smaller particles "b" )
1659  do ia=1, hydro_max
1660  do ib=1, hydro_max
1661  theta_ab0(ia,ib) = 2.0_rp * (w1(ib)/w2(ib))*(w3(ia)/w4(ia)) &
1662  * (w5(ia)/w6(ia))**beta_v(ia,2) &
1663  * (w5(ib)/w6(ib))**beta_v(ib,2)
1664  theta_ab1(ia,ib) = 2.0_rp * (w7(ib)/w8(ib))*(w3(ia)/w4(ia)) &
1665  * (w5(ia)/w6(ia))**beta_v(ia,2) &
1666  * (w5(ib)/w6(ib))**beta_v(ib,2)
1667  end do
1668  end do
1669 
1670  log_info("ATMOS_PHY_MP_sn14_init",'(100a16)') "LABEL ",wlabel(:)
1671  log_info("ATMOS_PHY_MP_sn14_init",'(a,100ES16.6)') "capacity ",cap(:) ! [Add] 11/08/30 T.Mitsui
1672  log_info("ATMOS_PHY_MP_sn14_init",'(a,100ES16.6)') "coef_m2 ",coef_m2(:)
1673  log_info("ATMOS_PHY_MP_sn14_init",'(a,100ES16.6)') "coef_d ",coef_d(:)
1674  !
1675  log_info("ATMOS_PHY_MP_sn14_init",'(a,100ES16.6)') "coef_d3 ",coef_d3(:)
1676  log_info("ATMOS_PHY_MP_sn14_init",'(a,100ES16.6)') "coef_d6 ",coef_d6(:)
1677  log_info("ATMOS_PHY_MP_sn14_init",'(a,100ES16.6)') "coef_d2v ",coef_d2v(:)
1678  log_info("ATMOS_PHY_MP_sn14_init",'(a,100ES16.6)') "coef_md2v ",coef_md2v(:)
1679  log_info("ATMOS_PHY_MP_sn14_init",'(a,100ES16.6)') "a_d2vt ",a_d2vt(:)
1680  log_info("ATMOS_PHY_MP_sn14_init",'(a,100ES16.6)') "b_d2vt ",b_d2vt(:)
1681  !
1682  log_info("ATMOS_PHY_MP_sn14_init",'(a,100ES16.6)') "coef_r2 ",coef_r2(:)
1683  log_info("ATMOS_PHY_MP_sn14_init",'(a,100ES16.6)') "coef_r3 ",coef_r3(:)
1684  log_info("ATMOS_PHY_MP_sn14_init",'(a,100ES16.6)') "coef_re ",coef_re(:)
1685  !
1686  log_info("ATMOS_PHY_MP_sn14_init",'(a,100ES16.6)') "a_area ",a_area(:)
1687  log_info("ATMOS_PHY_MP_sn14_init",'(a,100ES16.6)') "b_area ",b_area(:)
1688  log_info("ATMOS_PHY_MP_sn14_init",'(a,100ES16.6)') "ax_area ",ax_area(:)
1689  log_info("ATMOS_PHY_MP_sn14_init",'(a,100ES16.6)') "bx_area ",bx_area(:)
1690  log_info("ATMOS_PHY_MP_sn14_init",'(a,100ES16.6)') "a_rea ",a_rea(:)
1691  log_info("ATMOS_PHY_MP_sn14_init",'(a,100ES16.6)') "b_rea ",b_rea(:)
1692  log_info("ATMOS_PHY_MP_sn14_init",'(a,100ES16.6)') "a_rea3 ",a_rea3(:)
1693  log_info("ATMOS_PHY_MP_sn14_init",'(a,100ES16.6)') "b_rea3 ",b_rea3(:)
1694  !
1695  log_info("ATMOS_PHY_MP_sn14_init",'(a,100ES16.6)') "coef_rea2 ",coef_rea2(:)
1696  log_info("ATMOS_PHY_MP_sn14_init",'(a,100ES16.6)') "coef_rea3 ",coef_rea3(:)
1697  log_info("ATMOS_PHY_MP_sn14_init",'(a,100ES16.6)') "coef_vt0 ",coef_vt0(:,1)
1698  log_info("ATMOS_PHY_MP_sn14_init",'(a,100ES16.6)') "coef_vt1 ",coef_vt1(:,1)
1699  log_info("ATMOS_PHY_MP_sn14_init",'(a,100ES16.6)') "coef_A ",coef_a(:)
1700  log_info("ATMOS_PHY_MP_sn14_init",'(a,100ES16.6)') "coef_lambda ",coef_lambda(:)
1701 
1702  log_info("ATMOS_PHY_MP_sn14_init",'(a,100ES16.6)') "ah_vent0 sml",ah_vent0(:,1)
1703  log_info("ATMOS_PHY_MP_sn14_init",'(a,100ES16.6)') "ah_vent0 lrg",ah_vent0(:,2)
1704  log_info("ATMOS_PHY_MP_sn14_init",'(a,100ES16.6)') "ah_vent1 sml",ah_vent1(:,1)
1705  log_info("ATMOS_PHY_MP_sn14_init",'(a,100ES16.6)') "ah_vent1 lrg",ah_vent1(:,2)
1706  log_info("ATMOS_PHY_MP_sn14_init",'(a,100ES16.6)') "bh_vent0 sml",bh_vent0(:,1)
1707  log_info("ATMOS_PHY_MP_sn14_init",'(a,100ES16.6)') "bh_vent0 lrg",bh_vent0(:,2)
1708  log_info("ATMOS_PHY_MP_sn14_init",'(a,100ES16.6)') "bh_vent1 sml",bh_vent1(:,1)
1709  log_info("ATMOS_PHY_MP_sn14_init",'(a,100ES16.6)') "bh_vent1 lrg",bh_vent1(:,2)
1710 
1711  log_info("ATMOS_PHY_MP_sn14_init",'(a,100ES16.6)') "delta_b0 ",delta_b0(:)
1712  log_info("ATMOS_PHY_MP_sn14_init",'(a,100ES16.6)') "delta_b1 ",delta_b1(:)
1713  log_info("ATMOS_PHY_MP_sn14_init",'(a,100ES16.6)') "theta_b0 ",theta_b0(:)
1714  log_info("ATMOS_PHY_MP_sn14_init",'(a,100ES16.6)') "theta_b1 ",theta_b1(:)
1715 
1716  do ia=1, hydro_max
1717  log_info("ATMOS_PHY_MP_sn14_init",'(a,a10,a,100ES16.6)') "delta0(a,b)=(",trim(wlabel(ia)),",b)=",(delta_ab0(ia,ib),ib=1,hydro_max)
1718  enddo
1719  do ia=1, hydro_max
1720  log_info("ATMOS_PHY_MP_sn14_init",'(a,a10,a,100ES16.6)') "delta1(a,b)=(",trim(wlabel(ia)),",b)=",(delta_ab1(ia,ib),ib=1,hydro_max)
1721  enddo
1722  do ia=1, hydro_max
1723  log_info("ATMOS_PHY_MP_sn14_init",'(a,a10,a,100ES16.6)') "theta0(a,b)=(",trim(wlabel(ia)),",b)=",(theta_ab0(ia,ib),ib=1,hydro_max)
1724  enddo
1725  do ia=1, hydro_max
1726  log_info("ATMOS_PHY_MP_sn14_init",'(a,a10,a,100ES16.6)') "theta1(a,b)=(",trim(wlabel(ia)),",b)=",(theta_ab1(ia,ib),ib=1,hydro_max)
1727  enddo
1728 
1729  return
1730  end subroutine mp_sn14_init
1731  !-----------------------------------------------------------------------------
1732  subroutine mp_sn14 ( &
1733  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
1734  DENS, &
1735  W, &
1736  QTRC, &
1737  PRES0, &
1738  TEMP0, &
1739  Qdry, &
1740  CPtot0, &
1741  CVtot0, &
1742  CCN, &
1743  dt, &
1744  cz, &
1745  fz, &
1746  RHOQ_t, &
1747  RHOE_t, &
1748  CPtot_t, &
1749  CVtot_t, &
1750  EVAPORATE )
1751 
1752  use scale_atmos_hydrometeor, only: &
1753  cp_vapor, &
1754  cp_water, &
1755  cp_ice, &
1756  cv_vapor, &
1757  cv_water, &
1758  cv_ice
1759  use scale_atmos_saturation, only: &
1760  moist_psat_liq => atmos_saturation_psat_liq, &
1761  moist_psat_ice => atmos_saturation_psat_ice
1762  implicit none
1763 
1764  integer, intent(in) :: KA, KS, KE
1765  integer, intent(in) :: IA, IS, IE
1766  integer, intent(in) :: JA, JS, JE
1767 
1768  real(RP), intent(in) :: DENS (ka,ia,ja)
1769  real(RP), intent(in) :: W (ka,ia,ja)
1770  real(RP), intent(in) :: QTRC (ka,ia,ja,qa_mp)
1771  real(RP), intent(in) :: PRES0 (ka,ia,ja)
1772  real(RP), intent(in) :: TEMP0 (ka,ia,ja)
1773  real(RP), intent(in) :: Qdry (ka,ia,ja)
1774  real(RP), intent(in) :: CPtot0(ka, ia, ja)
1775  real(RP), intent(in) :: CVtot0(ka, ia, ja)
1776  real(RP), intent(in) :: CCN (ka,ia,ja)
1777  real(RP), intent(in) :: dt
1778  real(RP), intent(in) :: cz( ka,ia,ja)
1779  real(RP), intent(in) :: fz(0:ka,ia,ja)
1780 
1781  real(RP),intent(out) :: RHOQ_t(ka, ia, ja, qa_mp)
1782  real(RP),intent(out) :: RHOE_t(ka, ia, ja)
1783  real(RP),intent(out) :: CPtot_t(ka, ia, ja)
1784  real(RP),intent(out) :: CVtot_t(ka, ia, ja)
1785 
1786  real(RP), intent(out) :: EVAPORATE(ka,ia,ja) ! number of evaporated cloud [/m3/s]
1787 
1788  real(RP) :: pres(ka,ia,ja)
1789  real(RP) :: temp(ka,ia,ja)
1790  real(RP) :: cva(ka,ia,ja)
1791  real(RP) :: cpa(ka,ia,ja)
1792  real(RP) :: RHOQ0_t(ka, ia, ja, qa_mp)
1793  real(RP) :: RHOE0_t(ka, ia, ja)
1794  real(RP) :: CPtot0_t(ka, ia, ja)
1795  real(RP) :: CVtot0_t(ka, ia, ja)
1796  !
1797  ! diagnostic variables
1798  !
1799  real(RP) :: rhoe(ka,ia,ja)
1800  real(RP) :: rhoq(i_qv:i_ng,ka,ia,ja)
1801  real(RP) :: rhoq2(i_qv:i_ng,ka,ia,ja)
1802  !
1803  real(RP) :: xq(hydro_max,ka,ia,ja)
1804  !
1805  real(RP) :: dq_xa(hydro_max,ka,ia,ja)
1806  real(RP) :: vt_xa(hydro_max,2,ka,ia,ja) ! terminal velocity
1807 
1808  real(RP) :: wtemp(ka,ia,ja) ! filtered temperature
1809  real(RP) :: esw(ka,ia,ja) ! saturated vapor pressure(water)
1810  real(RP) :: esi(ka,ia,ja) ! saturated vapor pressure(ice)
1811  !
1812  real(RP) :: rho_fac
1813  real(RP) :: rho_fac_q(hydro_max,ka,ia,ja) ! factor for tracers, 1:cloud, 2:rain, 3:ice, 4: snow, 5:graupel
1814  !
1815  real(RP) :: drhogqv ! d (rho*qv*gsgam2)
1816  real(RP) :: drhogqc, drhognc ! qc, nc
1817  real(RP) :: drhogqr, drhognr ! qr, nr
1818  real(RP) :: drhogqi, drhogni ! qi, ni
1819  real(RP) :: drhogqs, drhogns ! qs, ns
1820  real(RP) :: drhogqg, drhogng ! qg, ng
1821 
1822  ! production rate
1823  real(RP) :: PQ(pq_max,ka,ia,ja)
1824 
1825  real(RP) :: wrm_dqc, wrm_dnc
1826  real(RP) :: wrm_dqr, wrm_dnr
1827 
1828  ! production rate of mixed-phase collection process
1829  real(RP) :: Pac(pac_max,ka,ia,ja)
1830 
1831  real(RP) :: gc_dqc, gc_dnc
1832  real(RP) :: sc_dqc, sc_dnc
1833  real(RP) :: ic_dqc, ic_dnc
1834  real(RP) :: rg_dqg, rg_dng
1835  real(RP) :: rg_dqr, rg_dnr
1836  real(RP) :: rs_dqr, rs_dnr, rs_dqs, rs_dns
1837  real(RP) :: ri_dqr, ri_dnr
1838  real(RP) :: ri_dqi, ri_dni
1839  real(RP) :: ii_dqi, ii_dni
1840  real(RP) :: is_dqi, is_dni, ss_dns
1841  real(RP) :: gs_dqs, gs_dns, gg_dng
1842  ! mixed-phase collection process total plus(clp_), total minus(clm_)
1843  real(RP) :: clp_dqc, clp_dnc, clm_dqc, clm_dnc
1844  real(RP) :: clp_dqr, clp_dnr, clm_dqr, clm_dnr
1845  real(RP) :: clp_dqi, clp_dni, clm_dqi, clm_dni
1846  real(RP) :: clp_dqs, clp_dns, clm_dqs, clm_dns
1847  real(RP) :: clp_dqg, clp_dng, clm_dqg, clm_dng
1848  real(RP) :: fac1, fac3, fac4, fac6, fac7, fac9, fac10
1849  ! production rate of partial conversion(ice, snow => graupel)
1850  real(RP) :: pco_dqi, pco_dni
1851  real(RP) :: pco_dqs, pco_dns
1852  real(RP) :: pco_dqg, pco_dng
1853  ! production rate of enhanced melting due to
1854  real(RP) :: eml_dqc, eml_dnc
1855  real(RP) :: eml_dqr, eml_dnr
1856  real(RP) :: eml_dqi, eml_dni
1857  real(RP) :: eml_dqs, eml_dns
1858  real(RP) :: eml_dqg, eml_dng
1859  ! production rate of ice multiplication by splintering
1860  real(RP) :: spl_dqi, spl_dni
1861  real(RP) :: spl_dqg, spl_dqs
1862 
1863  real(RP) :: rrho(ka,ia,ja)
1864 
1865  !-----------------------------------------------
1866  ! work for explicit supersaturation modeling
1867  !-----------------------------------------------
1868  real(RP) :: dTdt_equiv_d(ka,ia,ja) !
1869  !--------------------------------------------------
1870  !
1871  ! variables for output
1872  !
1873  !--------------------------------------------------
1874  ! work for column production term
1875  real(RP) :: sl_PLCdep(ia,ja)
1876  real(RP) :: sl_PLRdep(ia,ja), sl_PNRdep(ia,ja) !
1877  !--------------------------------------------------
1878  real(RP) :: qke_d(ka,ia,ja)
1879 
1880  real(RP), parameter :: eps = 1.e-19_rp
1881  real(RP), parameter :: eps_qv = 1.e-19_rp
1882  real(RP), parameter :: eps_rhoge = 1.e-19_rp
1883  real(RP), parameter :: eps_rhog = 1.e-19_rp
1884  integer :: ntdiv
1885 
1886  real(RP) :: sw
1887 
1888  integer :: k, i, j, iq
1889 
1890  real(RP) :: dqv, dqc, dqr, dqi, dqs, dqg, dcv, dcp
1891 
1892  !---------------------------------------------------------------------------
1893 
1894  ! total tendency
1895  rhoq_t(:,:,:,:) = 0.0_rp
1896  rhoe_t(:,:,:) = 0.0_rp
1897  cptot_t(:,:,:) = 0.0_rp
1898  cvtot_t(:,:,:) = 0.0_rp
1899 
1900  ! intermidiate variable
1901  do j = js, je
1902  do i = is, ie
1903  do k = ks, ke
1904  cpa(k,i,j) = cptot0(k,i,j)
1905  cva(k,i,j) = cvtot0(k,i,j)
1906  pres(k,i,j) = pres0(k,i,j)
1907  temp(k,i,j) = temp0(k,i,j)
1908  enddo
1909  enddo
1910  enddo
1911 
1912  !============================================================================
1913  !
1914  !-- Each process is integrated sequentially.
1915  ! 1. Nucleation and filter
1916  ! 2. Phase change
1917  ! 3. Collection
1918  !
1919  !============================================================================
1920 
1921  !----------------------------------------------------------------------------
1922  !
1923  ! 1.Nucleation of cloud water and cloud ice
1924  !
1925  !----------------------------------------------------------------------------
1926  call prof_rapstart('MP_Preprocess', 3)
1927 
1928  do j = js, je
1929  do i = is, ie
1930  do k = ks, ke
1931  do iq = i_qv, i_ng
1932  rhoq(iq,k,i,j) = dens(k,i,j) * qtrc(k,i,j,iq)
1933  enddo
1934  rhoq2(i_qv,k,i,j) = dens(k,i,j)*qtrc(k,i,j,i_qv)
1935  rhoq2(i_ni,k,i,j) = max( 0.0_rp, dens(k,i,j)*qtrc(k,i,j,i_ni) )
1936  rhoq2(i_nc,k,i,j) = max( 0.0_rp, dens(k,i,j)*qtrc(k,i,j,i_nc) )
1937  enddo
1938  enddo
1939  enddo
1940 
1941  do j = js, je
1942  do i = is, ie
1943  do k = ks, ke
1944  rrho(k,i,j) = 1.0_rp / dens(k,i,j)
1945  rhoe(k,i,j) = dens(k,i,j) * temp(k,i,j) * cva(k,i,j)
1946  wtemp(k,i,j) = max(temp(k,i,j), tem_min)
1947  enddo
1948  enddo
1949  enddo
1950 
1951  if( opt_debug_tem ) call debug_tem_kij( ka, ks, ke, ia, is, ie, ja, js, je, &
1952  1, temp(:,:,:), dens(:,:,:), pres(:,:,:), qtrc(:,:,:,i_qv) )
1953 
1954  do j = js, je
1955  do i = is, ie
1956  do k = ks, ke
1957  rho_fac = rho_0 / max(dens(k,i,j),rho_min)
1958  rho_fac_q(i_mp_qc,k,i,j) = rho_fac**gamma_v(i_mp_qc)
1959  rho_fac_q(i_mp_qr,k,i,j) = rho_fac**gamma_v(i_mp_qr)
1960  rho_fac_q(i_mp_qi,k,i,j) = (pres(k,i,j)/pre0_vt)**a_pre0_vt * (temp(k,i,j)/tem0_vt)**a_tem0_vt
1961  rho_fac_q(i_mp_qs,k,i,j) = rho_fac_q(i_mp_qi,k,i,j)
1962  rho_fac_q(i_mp_qg,k,i,j) = rho_fac_q(i_mp_qi,k,i,j)
1963  enddo
1964  enddo
1965  enddo
1966 
1967 !OCL XFILL
1968  do j = js, je
1969  do i = is, ie
1970  sl_plcdep(i,j) = 0.0_rp
1971  sl_plrdep(i,j) = 0.0_rp
1972  sl_pnrdep(i,j) = 0.0_rp
1973  end do
1974  end do
1975 
1976 !OCL XFILL
1977  do j = js, je
1978  do i = is, ie
1979  do k = ks, ke
1980  qke_d(k,i,j) = 0.0_rp ! 2*TKE
1981  enddo
1982  enddo
1983  enddo
1984 
1985 !OCL XFILL
1986  do j = js, je
1987  do i = is, ie
1988  do k = ks, ke
1989  dtdt_equiv_d(k,i,j) = 0.0_rp
1990  enddo
1991  enddo
1992  enddo
1993 
1994  call prof_rapend ('MP_Preprocess', 3)
1995 
1996  call prof_rapstart('MP_Nucleation', 3)
1997 
1998  call nucleation_kij( &
1999  ka, ks, ke, ia, is, ie, ja, js, je, &
2000  cz(:,:,:), & ! (in)
2001  fz(:,:,:), & ! (in)
2002  w(:,:,:), & ! (in)
2003  dens(:,:,:), & ! (in)
2004  wtemp(:,:,:), & ! (in)
2005  pres(:,:,:), & ! (in)
2006  qdry(:,:,:), & ! (in)
2007  rhoq2(:,:,:,:), & ! (in)
2008  cpa(:,:,:), & ! (in)
2009  dtdt_equiv_d(:,:,:), & ! (in)
2010  qke_d(:,:,:), & ! (in)
2011  ccn(:,:,:), & ! (in)
2012  dt, & ! (in)
2013  pq(:,:,:,:) ) ! (out)
2014 
2015  ! tendency
2016  do j = js, je
2017  do i = is, ie
2018  do k = ks, ke
2019  ! nucleation
2020  drhogqc = dt * pq(i_lcccn,k,i,j)
2021  drhognc = dt * pq(i_ncccn,k,i,j)
2022  drhogqi = dt * pq(i_liccn,k,i,j)
2023  drhogni = dt * pq(i_niccn,k,i,j)
2024  drhogqv = max( -rhoq(i_qv,k,i,j), -drhogqc-drhogqi )
2025  fac1 = drhogqv / min( -drhogqc-drhogqi, -eps ) ! limiting coefficient
2026 
2027  dqv = drhogqv * rrho(k,i,j)
2028  dqc = drhogqc*fac1 * rrho(k,i,j)
2029  dqi = drhogqi*fac1 * rrho(k,i,j)
2030 
2031  dcv = cv_vapor * dqv + cv_water * dqc + cv_ice * dqi
2032  dcp = cp_vapor * dqv + cp_water * dqc + cp_ice * dqi
2033 
2034  rhoq0_t(k,i,j,i_qv) = drhogqv / dt
2035  rhoq0_t(k,i,j,i_qc) = drhogqc*fac1 / dt
2036  rhoq0_t(k,i,j,i_qi) = drhogqi*fac1 / dt
2037  rhoq0_t(k,i,j,i_nc) = drhognc / dt
2038  rhoq0_t(k,i,j,i_ni) = drhogni / dt
2039 
2040  rhoe0_t(k,i,j) = - lhv * drhogqv / dt + lhf * drhogqi*fac1 / dt
2041 
2042  cvtot0_t(k,i,j) = dcv/dt
2043  cptot0_t(k,i,j) = dcp/dt
2044 
2045  enddo
2046  enddo
2047  enddo
2048 
2049  ! total tendency
2050  do j = js, je
2051  do i = is, ie
2052  do k = ks, ke
2053 
2054  rhoe_t(k,i,j) = rhoe_t(k,i,j) + rhoe0_t(k,i,j)
2055  cvtot_t(k,i,j) = cvtot_t(k,i,j) + cvtot0_t(k,i,j)
2056  cptot_t(k,i,j) = cptot_t(k,i,j) + cptot0_t(k,i,j)
2057 
2058  enddo
2059  enddo
2060  enddo
2061 
2062  ! update intermidiate variable
2063  do j = js, je
2064  do i = is, ie
2065  do k = ks, ke
2066 
2067  rhoq(i_qv,k,i,j) = rhoq(i_qv,k,i,j) + rhoq0_t(k,i,j,i_qv)*dt
2068  rhoq(i_qc,k,i,j) = max(0.0_rp, rhoq(i_qc,k,i,j) + rhoq0_t(k,i,j,i_qc)*dt )
2069  rhoq(i_qi,k,i,j) = max(0.0_rp, rhoq(i_qi,k,i,j) + rhoq0_t(k,i,j,i_qi)*dt )
2070  rhoq(i_nc,k,i,j) = max(0.0_rp, rhoq(i_nc,k,i,j) + rhoq0_t(k,i,j,i_nc)*dt )
2071  rhoq(i_ni,k,i,j) = max(0.0_rp, rhoq(i_ni,k,i,j) + rhoq0_t(k,i,j,i_ni)*dt )
2072 
2073  ! cloud number concentration filter
2074  rhoq(i_nc,k,i,j) = min( rhoq(i_nc,k,i,j), nc_uplim_d(1,i,j) )
2075 
2076  rhoe(k,i,j) = rhoe(k,i,j) + rhoe0_t(k,i,j)*dt
2077  cva(k,i,j) = cva(k,i,j) + cvtot0_t(k,i,j)*dt
2078  cpa(k,i,j) = cpa(k,i,j) + cptot0_t(k,i,j)*dt
2079 
2080  temp(k,i,j) = rhoe(k,i,j) / ( dens(k,i,j) * cva(k,i,j) )
2081  pres(k,i,j) = dens(k,i,j) * (cpa(k,i,j)-cva(k,i,j)) * temp(k,i,j)
2082  wtemp(k,i,j) = max( temp(k,i,j), tem_min )
2083  enddo
2084  enddo
2085  enddo
2086 
2087 ! if( opt_debug ) call debugreport_nucleation
2088  if( opt_debug_tem ) call debug_tem_kij( ka, ks, ke, ia, is, ie, ja, js, je, &
2089  2, temp(:,:,:), dens(:,:,:), pres(:,:,:), qtrc(:,:,:,i_qv) )
2090 
2091 
2092  call prof_rapend ('MP_Nucleation', 3)
2093  !----------------------------------------------------------------------------
2094  !
2095  ! 2.Phase change: Freezing, Melting, Vapor deposition
2096  !
2097  !----------------------------------------------------------------------------
2098  call prof_rapstart('MP_Phase_change', 3)
2099 
2100  do j = js, je
2101  do i = is, ie
2102  do k = ks, ke
2103  rhoq2(i_qr,k,i,j) = rhoq(i_qr,k,i,j)
2104  rhoq2(i_nr,k,i,j) = rhoq(i_nr,k,i,j)
2105  xq(i_mp_qr,k,i,j) = max(xr_min, min(xr_max, rhoq2(i_qr,k,i,j)/(rhoq2(i_nr,k,i,j)+nr_min) ))
2106 
2107  dq_xa(i_mp_qr,k,i,j) = a_m(i_mp_qr)*xq(i_mp_qr,k,i,j)**b_m(i_mp_qr)
2108  vt_xa(i_mp_qr,1,k,i,j) = alpha_v(i_mp_qr,1)*(xq(i_mp_qr,k,i,j)**beta_v(i_mp_qr,1))*rho_fac_q(i_mp_qr,k,i,j)
2109  vt_xa(i_mp_qr,2,k,i,j) = vt_xa(i_mp_qr,1,k,i,j)
2110 
2111  !! Following values shoud be already filtered to be non-zero before sbroutine was called.
2112  ! Mass concentration [kg/m3]
2113  rhoq2(i_qv,k,i,j) = rhoq(i_qv,k,i,j)
2114  rhoq2(i_qc,k,i,j) = rhoq(i_qc,k,i,j)
2115  rhoq2(i_qi,k,i,j) = rhoq(i_qi,k,i,j)
2116  rhoq2(i_qs,k,i,j) = rhoq(i_qs,k,i,j)
2117  rhoq2(i_qg,k,i,j) = rhoq(i_qg,k,i,j)
2118  ! Number concentration[/m3] (should be filtered to avoid zero division.)
2119  rhoq2(i_nc,k,i,j) = rhoq(i_nc,k,i,j)
2120  rhoq2(i_ni,k,i,j) = rhoq(i_ni,k,i,j)
2121  rhoq2(i_ns,k,i,j) = rhoq(i_ns,k,i,j)
2122  rhoq2(i_ng,k,i,j) = rhoq(i_ng,k,i,j)
2123 
2124  ! Mass of mean particle [kg]
2125  ! SB06(94)
2126  !
2127  xq(i_mp_qc,k,i,j) = min(xc_max, max(xc_min, rhoq2(i_qc,k,i,j)/(rhoq2(i_nc,k,i,j)+nc_min) ))
2128  xq(i_mp_qi,k,i,j) = min(xi_max, max(xi_min, rhoq2(i_qi,k,i,j)/(rhoq2(i_ni,k,i,j)+ni_min) ))
2129  xq(i_mp_qs,k,i,j) = min(xs_max, max(xs_min, rhoq2(i_qs,k,i,j)/(rhoq2(i_ns,k,i,j)+ns_min) ))
2130  xq(i_mp_qg,k,i,j) = min(xg_max, max(xg_min, rhoq2(i_qg,k,i,j)/(rhoq2(i_ng,k,i,j)+ng_min) ))
2131  ! diamter of average mass
2132  ! SB06(32)
2133  dq_xa(i_mp_qc,k,i,j) = a_m(i_mp_qc)*xq(i_mp_qc,k,i,j)**b_m(i_mp_qc)
2134  dq_xa(i_mp_qi,k,i,j) = a_m(i_mp_qi)*xq(i_mp_qi,k,i,j)**b_m(i_mp_qi)
2135  dq_xa(i_mp_qs,k,i,j) = a_m(i_mp_qs)*xq(i_mp_qs,k,i,j)**b_m(i_mp_qs)
2136  dq_xa(i_mp_qg,k,i,j) = a_m(i_mp_qg)*xq(i_mp_qg,k,i,j)**b_m(i_mp_qg)
2137 
2138  ! terminal velocity of average mass
2139  vt_xa(i_mp_qc,1,k,i,j) = alpha_v(i_mp_qc,1)*(xq(i_mp_qc,k,i,j)**beta_v(i_mp_qc,1))*rho_fac_q(i_mp_qc,k,i,j)
2140  vt_xa(i_mp_qi,1,k,i,j) = alpha_v(i_mp_qi,1)*(xq(i_mp_qi,k,i,j)**beta_v(i_mp_qi,1))*rho_fac_q(i_mp_qi,k,i,j)
2141  vt_xa(i_mp_qs,1,k,i,j) = alpha_v(i_mp_qs,1)*(xq(i_mp_qs,k,i,j)**beta_v(i_mp_qs,1))*rho_fac_q(i_mp_qs,k,i,j)
2142  vt_xa(i_mp_qg,1,k,i,j) = alpha_v(i_mp_qg,1)*(xq(i_mp_qg,k,i,j)**beta_v(i_mp_qg,1))*rho_fac_q(i_mp_qg,k,i,j)
2143  vt_xa(i_mp_qc,2,k,i,j) = alpha_v(i_mp_qc,2)*(xq(i_mp_qc,k,i,j)**beta_v(i_mp_qc,2))*rho_fac_q(i_mp_qc,k,i,j)
2144  vt_xa(i_mp_qi,2,k,i,j) = alpha_v(i_mp_qi,2)*(xq(i_mp_qi,k,i,j)**beta_v(i_mp_qi,2))*rho_fac_q(i_mp_qi,k,i,j)
2145  vt_xa(i_mp_qs,2,k,i,j) = alpha_v(i_mp_qs,2)*(xq(i_mp_qs,k,i,j)**beta_v(i_mp_qs,2))*rho_fac_q(i_mp_qs,k,i,j)
2146  vt_xa(i_mp_qg,2,k,i,j) = alpha_v(i_mp_qg,2)*(xq(i_mp_qg,k,i,j)**beta_v(i_mp_qg,2))*rho_fac_q(i_mp_qg,k,i,j)
2147 
2148  end do
2149  end do
2150  end do
2151 
2152  call moist_psat_liq( ka, ks, ke, ia, is, ie, ja, js, je, &
2153  wtemp(:,:,:), esw(:,:,:) ) ! [IN], [OUT]
2154  call moist_psat_ice( ka, ks, ke, ia, is, ie, ja, js, je, &
2155  wtemp(:,:,:), esi(:,:,:) ) ! [IN], [OUT]
2156 
2157  call freezing_water_kij( &
2158  ka, ks, ke, ia, is, ie, ja, js, je, &
2159  dt, & ! (in)
2160  rhoq2(:,:,:,:), & ! (in)
2161  xq(:,:,:,:), & ! (in)
2162  temp(:,:,:), & ! (in)
2163  pq(:,:,:,:) ) ! (inout)
2164 
2165  call dep_vapor_melt_ice_kij( &
2166  ka, ks, ke, ia, is, ie, ja, js, je, &
2167  dens(:,:,:), & ! (in)
2168  wtemp(:,:,:), & ! (in)
2169  pres(:,:,:), & ! (in)
2170  qdry(:,:,:), & ! (in)
2171  rhoq2(:,:,:,:), & ! (in)
2172  esw(:,:,:), & ! (in)
2173  esi(:,:,:), & ! (in)
2174  xq(:,:,:,:), & ! (in)
2175  vt_xa(:,:,:,:,:), & ! (in)
2176  dq_xa(:,:,:,:), & ! (in)
2177  pq(:,:,:,:) ) ! (inout)
2178  !
2179  ! update subroutine
2180  !
2182  ka, ks, ke, ia, is, ie, ja, js, je, &
2183  ntdiv, &
2184  ntmax_phase_change, & ! (in)
2185  dt, & ! (in)
2186  !gsgam2_d(:,:,:), & ! (in)
2187  cz(:,:,:), fz(:,:,:), & ! (in)
2188  w(:,:,:), & ! (in)
2189  dtdt_equiv_d(:,:,:), & ! (in)
2190  dens(:,:,:), & ! (in)
2191  qdry(:,:,:), & ! (in)
2192  esw(:,:,:), & ! (in)
2193  esi(:,:,:), & ! (in)
2194  rhoq2(:,:,:,:), & ! (in)
2195  pres(:,:,:), & ! (in)
2196  temp(:,:,:), & ! (in)
2197  cpa(:,:,:), & ! (in)
2198  cva(:,:,:), & ! (in)
2199  pq(:,:,:,:), & ! (inout)
2200  sl_plcdep(:,:), & ! (inout)
2201  sl_plrdep(:,:), & ! (inout)
2202  sl_pnrdep(:,:), & ! (inout)
2203  rhoq0_t(:,:,:,:), & ! (out)
2204  rhoe0_t(:,:,:), & ! (out)
2205  cptot0_t(:,:,:), & ! (out)
2206  cvtot0_t(:,:,:), & ! (out)
2207  evaporate(:,:,:) ) ! (out)
2208 
2209  ! total tendency
2210  do j = js, je
2211  do i = is, ie
2212  do k = ks, ke
2213 
2214  rhoe_t(k,i,j) = rhoe_t(k,i,j) + rhoe0_t(k,i,j)
2215  cvtot_t(k,i,j) = cvtot_t(k,i,j) + cvtot0_t(k,i,j)
2216  cptot_t(k,i,j) = cptot_t(k,i,j) + cptot0_t(k,i,j)
2217 
2218  enddo
2219  enddo
2220  enddo
2221 
2222  ! update intermidiate variable
2223  do j = js, je
2224  do i = is, ie
2225  do k = ks, ke
2226  do iq = 1, qa_mp
2227 
2228  rhoq(iq,k,i,j) = max(0.0_rp, rhoq(iq,k,i,j) + rhoq0_t(k,i,j,iq)*dt )
2229 
2230  enddo
2231  enddo
2232  enddo
2233  enddo
2234 
2235  do j = js, je
2236  do i = is, ie
2237  do k = ks, ke
2238 
2239  rhoe(k,i,j) = rhoe(k,i,j) + rhoe0_t(k,i,j)*dt
2240  cva(k,i,j) = cva(k,i,j) + cvtot0_t(k,i,j)*dt
2241  cpa(k,i,j) = cpa(k,i,j) + cptot0_t(k,i,j)*dt
2242  temp(k,i,j) = rhoe(k,i,j) / ( dens(k,i,j) * cva(k,i,j) )
2243  pres(k,i,j) = dens(k,i,j) * (cpa(k,i,j) - cva(k,i,j)) * temp(k,i,j)
2244 
2245  enddo
2246  enddo
2247  enddo
2248 
2249 ! if( opt_debug ) call debugreport_phasechange
2250  if( opt_debug_tem ) call debug_tem_kij( ka, ks, ke, ia, is, ie, ja, js, je, &
2251  3, temp(:,:,:), dens(:,:,:), pres(:,:,:), qtrc(:,:,:,i_qv) )
2252 
2253  call prof_rapend ('MP_Phase_change', 3)
2254 
2255  !---------------------------------------------------------------------------
2256  !
2257  ! 3.Collection process
2258  !
2259  !---------------------------------------------------------------------------
2260  call prof_rapstart('MP_Collection', 3)
2261 
2262  ! parameter setting
2263  do j = js, je
2264  do i = is, ie
2265  do k = ks, ke
2266 
2267  ! Mass concentration [kg/m3]
2268  rhoq2(i_qv,k,i,j) = rhoq(i_qv,k,i,j)
2269  rhoq2(i_qc,k,i,j) = rhoq(i_qc,k,i,j)
2270  rhoq2(i_qr,k,i,j) = rhoq(i_qr,k,i,j)
2271  rhoq2(i_qi,k,i,j) = rhoq(i_qi,k,i,j)
2272  rhoq2(i_qs,k,i,j) = rhoq(i_qs,k,i,j)
2273  rhoq2(i_qg,k,i,j) = rhoq(i_qg,k,i,j)
2274  ! Number concentration[/m3]
2275  rhoq2(i_nc,k,i,j) = rhoq(i_nc,k,i,j)
2276  rhoq2(i_nr,k,i,j) = rhoq(i_nr,k,i,j)
2277  rhoq2(i_ni,k,i,j) = rhoq(i_ni,k,i,j)
2278  rhoq2(i_ns,k,i,j) = rhoq(i_ns,k,i,j)
2279  rhoq2(i_ng,k,i,j) = rhoq(i_ng,k,i,j)
2280 
2281  ! Mass of mean particle [kg]
2282  xq(i_mp_qc,k,i,j) = min(xc_max, max(xc_min, rhoq2(i_qc,k,i,j)/(rhoq2(i_nc,k,i,j)+nc_min) ) )
2283  xq(i_mp_qr,k,i,j) = min(xr_max, max(xr_min, rhoq2(i_qr,k,i,j)/(rhoq2(i_nr,k,i,j)+nr_min) ) )
2284  xq(i_mp_qi,k,i,j) = min(xi_max, max(xi_min, rhoq2(i_qi,k,i,j)/(rhoq2(i_ni,k,i,j)+ni_min) ) )
2285  xq(i_mp_qs,k,i,j) = min(xs_max, max(xs_min, rhoq2(i_qs,k,i,j)/(rhoq2(i_ns,k,i,j)+ns_min) ) )
2286  xq(i_mp_qg,k,i,j) = min(xg_max, max(xg_min, rhoq2(i_qg,k,i,j)/(rhoq2(i_ng,k,i,j)+ng_min) ) )
2287 
2288  ! effective cross section is assume as area equivalent circle
2289  dq_xa(i_mp_qc,k,i,j) = 2.0_rp*a_rea(i_mp_qc)*xq(i_mp_qc,k,i,j)**b_rea(i_mp_qc)
2290  dq_xa(i_mp_qr,k,i,j) = 2.0_rp*a_rea(i_mp_qr)*xq(i_mp_qr,k,i,j)**b_rea(i_mp_qr)
2291  dq_xa(i_mp_qi,k,i,j) = 2.0_rp*a_rea(i_mp_qi)*xq(i_mp_qi,k,i,j)**b_rea(i_mp_qi)
2292  dq_xa(i_mp_qs,k,i,j) = 2.0_rp*a_rea(i_mp_qs)*xq(i_mp_qs,k,i,j)**b_rea(i_mp_qs)
2293  dq_xa(i_mp_qg,k,i,j) = 2.0_rp*a_rea(i_mp_qg)*xq(i_mp_qg,k,i,j)**b_rea(i_mp_qg)
2294 
2295  ! terminal velocity of average mass
2296  ! SB06(33)
2297  vt_xa(i_mp_qc,2,k,i,j) = alpha_v(i_mp_qc,2)*(xq(i_mp_qc,k,i,j)**beta_v(i_mp_qc,2))*rho_fac_q(i_mp_qc,k,i,j)
2298  vt_xa(i_mp_qr,2,k,i,j) = alpha_v(i_mp_qr,2)*(xq(i_mp_qr,k,i,j)**beta_v(i_mp_qr,2))*rho_fac_q(i_mp_qr,k,i,j)
2299  vt_xa(i_mp_qi,2,k,i,j) = alpha_v(i_mp_qi,2)*(xq(i_mp_qi,k,i,j)**beta_v(i_mp_qi,2))*rho_fac_q(i_mp_qi,k,i,j)
2300  vt_xa(i_mp_qs,2,k,i,j) = alpha_v(i_mp_qs,2)*(xq(i_mp_qs,k,i,j)**beta_v(i_mp_qs,2))*rho_fac_q(i_mp_qs,k,i,j)
2301  vt_xa(i_mp_qg,2,k,i,j) = alpha_v(i_mp_qg,2)*(xq(i_mp_qg,k,i,j)**beta_v(i_mp_qg,2))*rho_fac_q(i_mp_qg,k,i,j)
2302  enddo
2303  enddo
2304  enddo
2305 
2306  ! Auto-conversion, Accretion, Self-collection, Break-up
2307  ! [Mod] T.Seiki
2308  if ( mp_doautoconversion ) then
2309  call aut_acc_slc_brk_kij( &
2310  ka, ks, ke, ia, is, ie, ja, js, je, &
2311  rhoq2(:,:,:,:), & ! (in)
2312  xq(:,:,:,:), & ! (in)
2313  dq_xa(:,:,:,:), & ! (in)
2314  dens(:,:,:), & ! (in)
2315  pq(:,:,:,:) ) ! (inout)
2316 
2317  else
2318 !OCL XFILL
2319  do j = js, je
2320  do i = is, ie
2321  do k = ks, ke
2322  pq(i_lcaut,k,i,j) = 0.0_rp
2323  pq(i_ncaut,k,i,j) = 0.0_rp
2324  pq(i_nraut,k,i,j) = 0.0_rp
2325  pq(i_lcacc,k,i,j) = 0.0_rp
2326  pq(i_ncacc,k,i,j) = 0.0_rp
2327  pq(i_nrslc,k,i,j) = 0.0_rp
2328  pq(i_nrbrk,k,i,j) = 0.0_rp
2329  end do
2330  end do
2331  end do
2332  endif
2333 
2335  ! collection process
2336  ka, ks, ke, ia, is, ie, ja, js, je, &
2337  temp(:,:,:), & ! (in)
2338  rhoq2(:,:,:,:), & ! (in)
2339  xq(:,:,:,:), & ! (in)
2340  dq_xa(:,:,:,:), & ! (in)
2341  vt_xa(:,:,:,:,:), & ! (in)
2342  pq(:,:,:,:), & ! (inout)
2343  pac(:,:,:,:) ) ! (out)
2344 
2345  call ice_multiplication_kij( &
2346  ka, ks, ke, ia, is, ie, ja, js, je, &
2347  pac(:,:,:,:), & ! (in)
2348  temp(:,:,:), & ! (in)
2349  rhoq2(:,:,:,:), & ! (in)
2350  xq(:,:,:,:), & ! (in)
2351  pq(:,:,:,:) ) ! (inout)
2352 
2353  !
2354  ! update
2355  ! rhogq = l*gsgam
2356  !
2357  profile_start("sn14_update_rhoq")
2358  do j = js, je
2359  do i = is, ie
2360  do k = ks, ke
2361 
2362  ! warm collection process
2363  wrm_dqc = max( dt*( pq(i_lcaut,k,i,j)+pq(i_lcacc,k,i,j) ), -rhoq2(i_qc,k,i,j) )
2364  wrm_dnc = max( dt*( pq(i_ncaut,k,i,j)+pq(i_ncacc,k,i,j) ), -rhoq2(i_nc,k,i,j) )
2365  wrm_dnr = max( dt*( pq(i_nraut,k,i,j)+pq(i_nrslc,k,i,j)+pq(i_nrbrk,k,i,j) ), -rhoq2(i_nr,k,i,j) )
2366  wrm_dqr = -wrm_dqc
2367  ! mixed phase collection
2368  ! Pxxacyy2zz xx and yy decrease and zz increase .
2369  !
2370  ! At first fixer is applied to decreasing particles.
2371  ! order of fixer: graupel-cloud, snow-cloud, ice-cloud, graupel-rain, snow-rain, ice-rain,
2372  ! snow-ice, ice-ice, graupel-snow, snow-snow
2373  ! cloud mass decrease
2374  gc_dqc = max( dt*pac(i_lgaclc2lg,k,i,j) , min(0.0_rp, -rhoq2(i_qc,k,i,j)-wrm_dqc )) ! => dqg
2375  sc_dqc = max( dt*pac(i_lsaclc2ls,k,i,j) , min(0.0_rp, -rhoq2(i_qc,k,i,j)-wrm_dqc-gc_dqc )) ! => dqs
2376  ic_dqc = max( dt*pac(i_liaclc2li,k,i,j) , min(0.0_rp, -rhoq2(i_qc,k,i,j)-wrm_dqc-gc_dqc-sc_dqc )) ! => dqi
2377  ! cloud num. decrease
2378  gc_dnc = max( dt*pac(i_ngacnc2ng,k,i,j) , min(0.0_rp, -rhoq2(i_nc,k,i,j)-wrm_dnc )) ! => dnc:minus
2379  sc_dnc = max( dt*pac(i_nsacnc2ns,k,i,j) , min(0.0_rp, -rhoq2(i_nc,k,i,j)-wrm_dnc-gc_dnc )) ! => dnc:minus
2380  ic_dnc = max( dt*pac(i_niacnc2ni,k,i,j) , min(0.0_rp, -rhoq2(i_nc,k,i,j)-wrm_dnc-gc_dnc-sc_dnc )) ! => dnc:minus
2381 
2382  ! rain mass decrease ( tem < 273.15K)
2383  sw = sign(0.5_rp, t00-temp(k,i,j)) + 0.5_rp ! if( temp(k,i,j) <= T00 )then sw=1, else sw=0
2384  rg_dqr = max( dt*pac(i_lraclg2lg, k,i,j), min(0.0_rp, -rhoq2(i_qr,k,i,j)-wrm_dqr )) * sw
2385  rg_dqg = max( dt*pac(i_lraclg2lg, k,i,j), min(0.0_rp, -rhoq2(i_qg,k,i,j) )) * ( 1.0_rp - sw )
2386  rs_dqr = max( dt*pac(i_lracls2lg_r,k,i,j), min(0.0_rp, -rhoq2(i_qr,k,i,j)-wrm_dqr-rg_dqr )) * sw
2387  ri_dqr = max( dt*pac(i_lracli2lg_r,k,i,j), min(0.0_rp, -rhoq2(i_qr,k,i,j)-wrm_dqr-rg_dqr-rs_dqr )) * sw
2388  ! rain num. decrease
2389  rg_dnr = max( dt*pac(i_nracng2ng, k,i,j), min(0.0_rp, -rhoq2(i_nr,k,i,j)-wrm_dnr )) * sw
2390  rg_dng = max( dt*pac(i_nracng2ng, k,i,j), min(0.0_rp, -rhoq2(i_ng,k,i,j) )) * ( 1.0_rp - sw )
2391  rs_dnr = max( dt*pac(i_nracns2ng_r,k,i,j), min(0.0_rp, -rhoq2(i_nr,k,i,j)-wrm_dnr-rg_dnr )) * sw
2392  ri_dnr = max( dt*pac(i_nracni2ng_r,k,i,j), min(0.0_rp, -rhoq2(i_nr,k,i,j)-wrm_dnr-rg_dnr-rs_dnr )) * sw
2393 
2394  ! ice mass decrease
2395  fac1 = (ri_dqr-eps)/ (dt*pac(i_lracli2lg_r,k,i,j)-eps) ! suppress factor by filter of rain
2396  ri_dqi = max( dt*pac(i_lracli2lg_i,k,i,j)*fac1, min(0.0_rp, -rhoq2(i_qi,k,i,j)+ic_dqc )) ! => dqg
2397  ii_dqi = max( dt*pac(i_liacli2ls,k,i,j) , min(0.0_rp, -rhoq2(i_qi,k,i,j)+ic_dqc-ri_dqi )) ! => dqs
2398  is_dqi = max( dt*pac(i_liacls2ls,k,i,j) , min(0.0_rp, -rhoq2(i_qi,k,i,j)+ic_dqc-ri_dqi-ii_dqi )) ! => dqs
2399  ! ice num. decrease
2400  fac4 = (ri_dnr-eps)/ (dt*pac(i_nracni2ng_r,k,i,j)-eps) ! suppress factor by filter of rain
2401  ri_dni = max( dt*pac(i_nracni2ng_i,k,i,j)*fac4, min(0.0_rp, -rhoq2(i_ni,k,i,j) )) ! => dni:minus
2402  ii_dni = max( dt*pac(i_niacni2ns,k,i,j) , min(0.0_rp, -rhoq2(i_ni,k,i,j)-ri_dni )) ! => dni:minus,dns:plus(*0.5)
2403  is_dni = max( dt*pac(i_niacns2ns,k,i,j) , min(0.0_rp, -rhoq2(i_ni,k,i,j)-ri_dni-ii_dni )) ! => dni:minus,dns:plus
2404  ! snow mass decrease
2405  fac3 = (rs_dqr-eps)/(dt*pac(i_lracls2lg_r,k,i,j)-eps) ! suppress factor by filter of rain
2406  rs_dqs = max( dt*pac(i_lracls2lg_s,k,i,j)*fac3, min(0.0_rp, -rhoq2(i_qs,k,i,j)+sc_dqc+ii_dqi+is_dqi )) ! => dqg
2407  gs_dqs = max( dt*pac(i_lgacls2lg,k,i,j) , min(0.0_rp, -rhoq2(i_qs,k,i,j)+sc_dqc+ii_dqi+is_dqi-rs_dqs )) ! => dqg
2408  ! snow num. decrease
2409  fac6 = (rs_dnr-eps)/(dt*pac(i_nracns2ng_r,k,i,j)-eps) ! suppress factor by filter of rain
2410 ! fac7 = (is_dni-eps)/(dt*Pac(I_NIacNS2NS, k,i,j)-eps) ! suppress factor by filter of ice
2411  rs_dns = max( dt*pac(i_nracns2ng_s,k,i,j)*fac6, min(0.0_rp, -rhoq2(i_ns,k,i,j)+0.50_rp*ii_dni+is_dni )) ! => dns:minus
2412  gs_dns = max( dt*pac(i_ngacns2ng,k,i,j) , min(0.0_rp, -rhoq2(i_ns,k,i,j)+0.50_rp*ii_dni+is_dni-rs_dns )) ! => dns:minus
2413  ss_dns = max( dt*pac(i_nsacns2ns,k,i,j) , min(0.0_rp, -rhoq2(i_ns,k,i,j)+0.50_rp*ii_dni+is_dni-rs_dns-gs_dns ))
2414  !
2415  gg_dng = max( dt*pac(i_ngacng2ng,k,i,j) , min(0.0_rp, -rhoq2(i_ng,k,i,j) ))
2416  !
2417  ! total plus in mixed phase collection(clp_)
2418  ! mass
2419  ! if( temp(k,i,j) <= T00 )then sw=1, else sw=0
2420  clp_dqc = 0.0_rp
2421  clp_dqr = (-rg_dqg-rs_dqs-ri_dqi) * (1.0_rp-sw)
2422  clp_dqi = -ic_dqc
2423  clp_dqs = -sc_dqc-ii_dqi-is_dqi
2424  clp_dqg = -gc_dqc -gs_dqs + (-rg_dqr-rs_dqr-rs_dqs-ri_dqr-ri_dqi) * sw
2425  ! num.( number only increase when a+b=>c, dnc=-dna)
2426  clp_dnc = 0.0_rp
2427  clp_dnr = 0.0_rp
2428  clp_dni = 0.0_rp
2429  clp_dns = -ii_dni*0.5_rp
2430  clp_dng = (-rs_dnr-ri_dnr) * sw
2431  ! total minus in mixed phase collection(clm_)
2432  ! mass
2433  clm_dqc = gc_dqc+sc_dqc+ic_dqc
2434  clm_dqr = (rg_dqr+rs_dqr+ri_dqr) * sw
2435  clm_dqi = ri_dqi+ii_dqi+is_dqi
2436  clm_dqs = rs_dqs+gs_dqs
2437  clm_dqg = rg_dqg * (1.0_rp-sw)
2438  ! num.
2439  clm_dnc = gc_dnc+sc_dnc+ic_dnc
2440  clm_dnr = (rg_dnr+rs_dnr+ri_dnr) * sw
2441  clm_dni = ri_dni+ii_dni+is_dni
2442  clm_dns = rs_dns+ss_dns+gs_dns
2443  clm_dng = gg_dng + rg_dng * (1.0_rp-sw)
2444 
2445  ! partial conversion
2446  ! 08/05/08 [Mod] T.Mitsui
2447  pco_dqi = max( dt*pq(i_licon,k,i,j), -clp_dqi )
2448  pco_dqs = max( dt*pq(i_lscon,k,i,j), -clp_dqs )
2449  pco_dqg = -pco_dqi-pco_dqs
2450  ! 08/05/08 [Mod] T.Mitsui
2451  pco_dni = max( dt*pq(i_nicon,k,i,j), -clp_dni )
2452  pco_dns = max( dt*pq(i_nscon,k,i,j), -clp_dns )
2453  pco_dng = -pco_dni-pco_dns
2454  ! enhanced melting ( always negative value )
2455  ! ice-cloud melting produces cloud, others produce rain
2456  eml_dqi = max( dt*pq(i_liacm,k,i,j), min(0.0_rp, -rhoq2(i_qi,k,i,j)-(clp_dqi+clm_dqi)-pco_dqi ))
2457  eml_dqs = max( dt*pq(i_lsacm,k,i,j), min(0.0_rp, -rhoq2(i_qs,k,i,j)-(clp_dqs+clm_dqs)-pco_dqs ))
2458  eml_dqg = max( dt*(pq(i_lgacm,k,i,j)+pq(i_lgarm,k,i,j)+pq(i_lsarm,k,i,j)+pq(i_liarm,k,i,j)), &
2459  min(0.0_rp, -rhoq2(i_qg,k,i,j)-(clp_dqg+clm_dqg)-pco_dqg ))
2460  eml_dqc = -eml_dqi
2461  eml_dqr = -eml_dqs-eml_dqg
2462  !
2463  eml_dni = max( dt*pq(i_niacm,k,i,j), min(0.0_rp, -rhoq2(i_ni,k,i,j)-(clp_dni+clm_dni)-pco_dni ))
2464  eml_dns = max( dt*pq(i_nsacm,k,i,j), min(0.0_rp, -rhoq2(i_ns,k,i,j)-(clp_dns+clm_dns)-pco_dns ))
2465  eml_dng = max( dt*(pq(i_ngacm,k,i,j)+pq(i_ngarm,k,i,j)+pq(i_nsarm,k,i,j)+pq(i_niarm,k,i,j)), &
2466  min(0.0_rp, -rhoq2(i_ng,k,i,j)-(clp_dng+clm_dng)-pco_dng ))
2467  eml_dnc = -eml_dni
2468  eml_dnr = -eml_dns-eml_dng
2469  !
2470  ! ice multiplication
2471  spl_dqg = max( dt*pq(i_lgspl,k,i,j), min(0.0_rp, -rhoq2(i_qg,k,i,j)-(clp_dqg+clm_dqg)-pco_dqg-eml_dqg ))
2472  spl_dqs = max( dt*pq(i_lsspl,k,i,j), min(0.0_rp, -rhoq2(i_qs,k,i,j)-(clp_dqs+clm_dqs)-pco_dqs-eml_dqs ))
2473  spl_dqi = -spl_dqg-spl_dqs
2474  fac9 = (spl_dqg-eps)/(dt*pq(i_lgspl,k,i,j)-eps)
2475  fac10 = (spl_dqs-eps)/(dt*pq(i_lsspl,k,i,j)-eps)
2476  spl_dni = dt*pq(i_nispl,k,i,j)*fac9*fac10
2477  !
2478  ! total cloud change
2479  drhogqc = (wrm_dqc + clp_dqc + clm_dqc + eml_dqc )
2480  drhognc = (wrm_dnc + clp_dnc + clm_dnc + eml_dnc )
2481  ! total rain change
2482  drhogqr = (wrm_dqr + clp_dqr + clm_dqr + eml_dqr )
2483  drhognr = (wrm_dnr + clp_dnr + clm_dnr + eml_dnr )
2484  ! total ice change
2485  drhogqi = ( clp_dqi + clm_dqi + pco_dqi + eml_dqi + spl_dqi)
2486  drhogni = ( clp_dni + clm_dni + pco_dni + eml_dni + spl_dni)
2487  ! total snow change
2488  drhogqs = ( clp_dqs + clm_dqs + pco_dqs + eml_dqs + spl_dqs)
2489  drhogns = ( clp_dns + clm_dns + pco_dns + eml_dns )
2490  ! total graupel change
2491  drhogqg = ( clp_dqg + clm_dqg + pco_dqg + eml_dqg + spl_dqg)
2492  drhogng = ( clp_dng + clm_dng + pco_dng + eml_dng )
2493  !
2494 
2495  ! tendency
2496  rhoq0_t(k,i,j,i_qc) = drhogqc / dt
2497  rhoq0_t(k,i,j,i_nc) = drhognc / dt
2498  rhoq0_t(k,i,j,i_qr) = drhogqr / dt
2499  rhoq0_t(k,i,j,i_nr) = drhognr / dt
2500  rhoq0_t(k,i,j,i_qi) = drhogqi / dt
2501  rhoq0_t(k,i,j,i_ni) = drhogni / dt
2502  rhoq0_t(k,i,j,i_qs) = drhogqs / dt
2503  rhoq0_t(k,i,j,i_ns) = drhogns / dt
2504  rhoq0_t(k,i,j,i_qg) = drhogqg / dt
2505  rhoq0_t(k,i,j,i_ng) = drhogng / dt
2506 
2507  rhoe0_t(k,i,j) = lhf * ( drhogqi + drhogqs + drhogqg ) / dt
2508 
2509  dqc = drhogqc * rrho(k,i,j)
2510  dqr = drhogqr * rrho(k,i,j)
2511  dqi = drhogqi * rrho(k,i,j)
2512  dqs = drhogqs * rrho(k,i,j)
2513  dqg = drhogqg * rrho(k,i,j)
2514 
2515  dcv = cv_water * ( dqc + dqr ) + cv_ice * ( dqi + dqs + dqg )
2516  dcp = cp_water * ( dqc + dqr ) + cp_ice * ( dqi + dqs + dqg )
2517 
2518  cvtot0_t(k,i,j) = dcv/dt
2519  cptot0_t(k,i,j) = dcp/dt
2520 
2521  enddo
2522  enddo
2523  enddo
2524  profile_stop("sn14_update_rhoq")
2525 
2526  call prof_rapend ('MP_Collection', 3)
2527 
2528  call prof_rapstart('MP_Postprocess', 3)
2529 
2530  ! total tendency
2531  do j = js, je
2532  do i = is, ie
2533  do k = ks, ke
2534 
2535  rhoe_t(k,i,j) = rhoe_t(k,i,j) + rhoe0_t(k,i,j)
2536  cvtot_t(k,i,j) = cvtot_t(k,i,j) + cvtot0_t(k,i,j)
2537  cptot_t(k,i,j) = cptot_t(k,i,j) + cptot0_t(k,i,j)
2538 
2539  enddo
2540  enddo
2541  enddo
2542 
2543  !--- update
2544 
2545  do j = js, je
2546  do i = is, ie
2547  do k = ks, ke
2548  do iq = i_qc, i_ng
2549 
2550  rhoq(iq,k,i,j) = max(0.0_rp, rhoq(iq,k,i,j) + rhoq0_t(k,i,j,iq)*dt )
2551 
2552  enddo
2553  enddo
2554  enddo
2555  enddo
2556 
2557  ! total tendency
2558  do iq = i_qv, i_ng
2559  do j = js, je
2560  do i = is, ie
2561  do k = ks, ke
2562 
2563  rhoq_t(k,i,j,iq) = ( rhoq(iq,k,i,j) - dens(k,i,j)*qtrc(k,i,j,iq) )/dt
2564 
2565  enddo
2566  enddo
2567  enddo
2568  enddo
2569 
2570 ! if( opt_debug ) call debugreport_collection
2571  if( opt_debug_tem ) call debug_tem_kij( ka, ks, ke, ia, is, ie, ja, js, je, &
2572  4, temp(:,:,:), dens(:,:,:), pres(:,:,:), qtrc(:,:,:,i_qv) )
2573 
2574  call prof_rapend ('MP_Postprocess', 3)
2575 
2576  return
2577  end subroutine mp_sn14
2578 
2579  !-----------------------------------------------------------------------------
2580  subroutine debug_tem_kij( &
2581  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
2582  point, &
2583  tem, &
2584  rho, &
2585  pre, &
2586  qv )
2587  use scale_prc, only: &
2588  prc_myrank
2589  implicit none
2590 
2591  integer, intent(in) :: KA, KS, KE
2592  integer, intent(in) :: IA, IS, IE
2593  integer, intent(in) :: JA, JS, JE
2594 
2595  integer, intent(in) :: point
2596  real(RP), intent(in) :: tem(ka,ia,ja)
2597  real(RP), intent(in) :: rho(ka,ia,ja)
2598  real(RP), intent(in) :: pre(ka,ia,ja)
2599  real(RP), intent(in) :: qv (ka,ia,ja)
2600 
2601  integer :: k ,i, j
2602  !---------------------------------------------------------------------------
2603 
2604  do j = js, je
2605  do i = is, ie
2606  do k = ks, ke
2607  if ( tem(k,i,j) < tem_min &
2608  .OR. rho(k,i,j) < rho_min &
2609  .OR. pre(k,i,j) < 1.0_rp ) then
2610 
2611  log_info("ATMOS_PHY_MP_SN14_debug_tem_kij",'(A,I3,A,4(F16.5),3(I6))') &
2612  "point: ", point, " low tem,rho,pre:", tem(k,i,j), rho(k,i,j), pre(k,i,j), qv(k,i,j), k, i, j, prc_myrank
2613  endif
2614  enddo
2615  enddo
2616  enddo
2617 
2618  return
2619  end subroutine debug_tem_kij
2620 
2621  subroutine nucleation_kij( &
2622  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
2623  cz, fz, w, &
2624  rho, tem, pre, qdry, &
2625  rhoq, &
2626  cpa, & ! in
2627  dTdt_rad, & ! in
2628  qke, & ! in
2629  CCN, & ! in
2630  dt, & ! in
2631  PQ )
2632  use scale_prc, only: &
2633  prc_abort
2634  use scale_atmos_saturation, only: &
2635  moist_psat_liq => atmos_saturation_psat_liq, &
2636  moist_psat_ice => atmos_saturation_psat_ice, &
2637  moist_pres2qsat_liq => atmos_saturation_pres2qsat_liq, &
2638  moist_pres2qsat_ice => atmos_saturation_pres2qsat_ice, &
2639  moist_dqsi_dtem_dens => atmos_saturation_dqs_dtem_dens_liq
2640  implicit none
2641 
2642  integer, intent(in) :: KA, KS, KE
2643  integer, intent(in) :: IA, IS, IE
2644  integer, intent(in) :: JA, JS, JE
2645 
2646  real(RP), intent(in) :: cz( ka,ia,ja) !
2647  real(RP), intent(in) :: fz(0:ka,ia,ja) !
2648  real(RP), intent(in) :: w (ka,ia,ja) ! w of full level
2649  real(RP), intent(in) :: rho(ka,ia,ja) ! [Add] 09/08/18 T.Mitsui
2650  real(RP), intent(in) :: tem(ka,ia,ja) ! [Add] 09/08/18 T.Mitsui
2651  real(RP), intent(in) :: pre(ka,ia,ja) ! [Add] 09/08/18 T.Mitsui
2652  real(RP), intent(in) :: qdry(ka,ia,ja)
2653  !
2654  real(RP), intent(in) :: rhoq(i_qv:i_ng,ka,ia,ja) !
2655  !
2656  real(RP), intent(in) :: cpa(ka,ia,ja) ! in 09/08/18 [Add] T.Mitsui
2657  real(RP), intent(in) :: dTdt_rad(ka,ia,ja) ! 09/08/18 T.Mitsui
2658  real(RP), intent(in) :: qke(ka,ia,ja) ! 09/08/18 T.Mitsui
2659  real(RP), intent(in) :: dt
2660  real(RP), intent(in) :: CCN(ka,ia,ja)
2661  !
2662  real(RP), intent(out) :: PQ(pq_max,ka,ia,ja)
2663  !
2664  ! namelist variables
2665  !
2666  ! total aerosol number concentration [/m3]
2667  real(RP), parameter :: c_ccn_ocean= 1.00e+8_rp
2668  real(RP), parameter :: c_ccn_land = 1.26e+9_rp
2669  real(RP), save :: c_ccn = 1.00e+8_rp
2670  ! aerosol activation factor
2671  real(RP), parameter :: kappa_ocean= 0.462_rp
2672  real(RP), parameter :: kappa_land = 0.308_rp
2673  real(RP), save :: kappa = 0.462_rp
2674  real(RP), save :: c_in = 1.0_rp
2675  ! SB06 (36)
2676  real(RP), save :: nm_M92 = 1.e+3_rp
2677  real(RP), save :: am_M92 = -0.639_rp
2678  real(RP), save :: bm_M92 = 12.96_rp
2679  !
2680  real(RP), save :: in_max = 1000.e+3_rp ! max num. of Ice-Nuclei [num/m3]
2681  real(RP), save :: ssi_max= 0.60_rp
2682  real(RP), save :: ssw_max= 1.1_rp ! [%]
2683  !
2684  logical, save :: flag_first = .true.
2685  real(RP), save :: qke_min = 0.03_rp ! sigma=0.1[m/s], 09/08/18 T.Mitsui
2686  real(RP), save :: tem_ccn_low=233.150_rp ! = -40 degC ! [Add] 10/08/03 T.Mitsui
2687  real(RP), save :: tem_in_low =173.150_rp ! = -100 degC ! [Add] 10/08/03 T.Mitsui
2688  logical, save :: nucl_twomey = .false.
2689  logical, save :: inucl_w = .false.
2690  !
2691  namelist / param_atmos_phy_mp_sn14_nucleation / &
2692  in_max, & !
2693  c_ccn, kappa, & ! cloud nucleation
2694  nm_m92, am_m92, bm_m92, & ! ice nucleation
2695  xc_ccn, xi_ccn, &
2696  tem_ccn_low, & ! [Add] 10/08/03 T.Mitsui
2697  tem_in_low, & ! [Add] 10/08/03 T.Mitsui
2698  ssw_max, ssi_max, &
2699  nucl_twomey, inucl_w ! [Add] 13/01/30 Y.Sato
2700  !
2701 ! real(RP) :: c_ccn_map(1,IA,JA) ! c_ccn horizontal distribution
2702 ! real(RP) :: kappa_map(1,IA,JA) ! kappa horizontal distribution
2703 ! real(RP) :: c_in_map(1,IA,JA) ! c_in horizontal distribution ! [Add] 11/08/30 T.Mitsui
2704  real(RP) :: esw(ka,ia,ja) ! saturation vapor pressure, water
2705  real(RP) :: esi(ka,ia,ja) ! ice
2706  real(RP) :: ssw(ka,ia,ja) ! super saturation (water)
2707  real(RP) :: ssi(ka,ia,ja) ! super saturation (ice)
2708 ! real(RP) :: w_dsswdz(KA,IA,JA) ! w*(d_ssw/ d_z) super saturation(water) flux
2709  real(RP) :: w_dssidz(ka,ia,ja) ! w*(d_ssi/ d_z), 09/04/14 T.Mitsui
2710 ! real(RP) :: ssw_below(KA,IA,JA)! ssw(k-1)
2711  real(RP) :: ssi_below(ka,ia,ja)! ssi(k-1), 09/04/14 T.Mitsui
2712  real(RP) :: z_below(ka,ia,ja) ! z(k-1)
2713  real(RP) :: dzh ! z(k)-z(k-1)
2714  real(RP) :: pv ! vapor pressure
2715  ! work variables for Twomey Equation.
2716  real(RP) :: qsw(ka,ia,ja)
2717  real(RP) :: qsi(ka,ia,ja)
2718  real(RP) :: dqsidtem_rho(ka,ia,ja)
2719  real(RP) :: dssidt_rad(ka,ia,ja)
2720 ! real(RP) :: dni_ratio(KA,IA,JA)
2721  real(RP) :: wssi, wdssi
2722  !
2723 ! real(RP) :: xi_nuc(1,IA,JA) ! xi use the value @ cloud base
2724 ! real(RP) :: alpha_nuc(1,IA,JA) ! alpha_nuc
2725 ! real(RP) :: eta_nuc(1,IA,JA) ! xi use the value @ cloud base
2726  !
2727  real(RP) :: sigma_w(ka,ia,ja)
2728  real(RP) :: weff(ka,ia,ja)
2729  real(RP) :: weff_max(ka,ia,ja)
2730  real(RP) :: velz(ka)
2731  !
2732  real(RP) :: coef_ccn(ia,ja)
2733  real(RP) :: slope_ccn(ia,ja)
2734  real(RP) :: nc_new(ka,ia,ja)
2735  real(RP) :: nc_new_below(ka,ia,ja)
2736  real(RP) :: dnc_new
2737  real(RP) :: nc_new_max ! Lohmann (2002),
2738  real(RP) :: a_max
2739  real(RP) :: b_max
2740  logical :: flag_nucleation(ka,ia,ja)
2741  !
2742  real(RP) :: r_gravity
2743  real(RP), parameter :: r_sqrt3=0.577350269_rp ! = sqrt(1.d0/3.d0)
2744  real(RP), parameter :: eps=1.e-30_rp
2745  !====> ! 09/08/18
2746  !
2747  real(RP) :: dlcdt_max, dli_max ! defined by supersaturation
2748  real(RP) :: dncdt_max, dni_max ! defined by supersaturation
2749  real(RP) :: rdt
2750  !
2751  integer :: i, j, k
2752  !
2753  if( flag_first )then
2754  rewind(io_fid_conf)
2755  read(io_fid_conf, nml=param_atmos_phy_mp_sn14_nucleation, end=100)
2756 100 log_nml(param_atmos_phy_mp_sn14_nucleation)
2757  flag_first=.false.
2758 
2759  if ( mp_couple_aerosol .AND. nucl_twomey ) then
2760  log_error("ATMOS_PHY_MP_SN14_nucleation_kij",*) "nucl_twomey should be false when MP_couple_aerosol is true, stop"
2761  call prc_abort
2762  endif
2763  endif
2764  !
2765 ! c_ccn_map(1,:,:) = c_ccn
2766 ! kappa_map(1,:,:) = kappa
2767 ! c_in_map(1,:,:) = c_in
2768  !
2769 ! nc_uplim_d(1,:,:) = c_ccn_map(1,:,:)*1.5_RP
2770  do j = js, je
2771  do i = is, ie
2772  nc_uplim_d(1,i,j) = c_ccn*1.5_rp
2773  end do
2774  end do
2775  !
2776  rdt = 1.0_rp/dt
2777  r_gravity = 1.0_rp/grav
2778  !
2779  call moist_psat_liq ( ka, ks, ke, ia, is, ie, ja, js, je, &
2780  tem(:,:,:), esw(:,:,:) ) ! [IN], [OUT]
2781  call moist_psat_ice ( ka, ks, ke, ia, is, ie, ja, js, je, &
2782  tem(:,:,:), esi(:,:,:) ) ! [IN], [OUT]
2783  call moist_pres2qsat_liq ( ka, ks, ke, ia, is, ie, ja, js, je, &
2784  tem(:,:,:), pre(:,:,:), qdry(:,:,:), & ! [IN]
2785  qsw(:,:,:) ) ! [OUT]
2786  call moist_pres2qsat_ice ( ka, ks, ke, ia, is, ie, ja, js, je, &
2787  tem(:,:,:), pre(:,:,:), qdry(:,:,:), & ! [IN]
2788  qsi(:,:,:) ) ! [OUT]
2789  call moist_dqsi_dtem_dens( ka, ks, ke, ia, is, ie, ja, js, je, &
2790  tem(:,:,:), rho(:,:,:), & ! [IN]
2791  dqsidtem_rho(:,:,:) ) ! [OUT]
2792  !
2793  ! Lohmann (2002),JAS, eq.(1) but changing unit [cm-3] => [m-3]
2794  a_max = 1.e+6_rp*0.1_rp*(1.e-6_rp)**1.27_rp
2795  b_max = 1.27_rp
2796  !
2797  ssi_max = 1.0_rp
2798 
2799  do j = js, je
2800  do i = is, ie
2801  do k = ks, ke
2802  pv = rhoq(i_qv,k,i,j)*rvap*tem(k,i,j)
2803  ssw(k,i,j) = min( mp_ssw_lim, ( pv/esw(k,i,j)-1.0_rp ) )*100.0_rp
2804  ssi(k,i,j) = (pv/esi(k,i,j) - 1.00_rp)
2805 ! ssw_below(k+1,i,j) = ssw(k,i,j)
2806  ssi_below(k+1,i,j) = ssi(k,i,j)
2807  z_below(k+1,i,j) = cz(k,i,j)
2808  end do
2809 ! ssw_below(KS,i,j) = ssw(KS,i,j)
2810  ssi_below(ks,i,j) = ssi(ks,i,j)
2811  z_below(ks,i,j) = cz(ks-1,i,j)
2812 
2813  ! dS/dz is evaluated by first order upstream difference
2814  !*** Solution for Twomey Equation ***
2815 ! coef_ccn(i,j) = 1.E+6_RP*0.88_RP*(c_ccn_map(1,i,j)*1.E-6_RP)**(2.0_RP/(kappa_map(1,i,j) + 2.0_RP)) * &
2816  coef_ccn(i,j) = 1.e+6_rp*0.88_rp*(c_ccn*1.e-6_rp)**(2.0_rp/(kappa + 2.0_rp)) &
2817 ! * (70.0_RP)**(kappa_map(1,i,j)/(kappa_map(1,i,j) + 2.0_RP))
2818  * (70.0_rp)**(kappa/(kappa + 2.0_rp))
2819 ! slope_ccn(i,j) = 1.5_RP*kappa_map(1,i,j)/(kappa_map(1,i,j) + 2.0_RP)
2820  slope_ccn(i,j) = 1.5_rp*kappa/(kappa + 2.0_rp)
2821  !
2822  do k=ks, ke
2823  sigma_w(k,i,j) = r_sqrt3*sqrt(max(qke(k,i,j),qke_min))
2824  end do
2825  sigma_w(ks-1,i,j) = sigma_w(ks,i,j)
2826  sigma_w(ke+1,i,j) = sigma_w(ke,i,j)
2827  ! effective vertical velocity
2828  do k=ks, ke
2829  weff(k,i,j) = w(k,i,j) - cpa(k,i,j)*r_gravity*dtdt_rad(k,i,j)
2830  end do
2831 
2832  end do
2833  end do
2834  !
2835  if( mp_couple_aerosol ) then
2836 
2837  do j = js, je
2838  do i = is, ie
2839  do k = ks, ke
2840  if( ssw(k,i,j) > 1.e-10_rp .AND. pre(k,i,j) > 300.e+2_rp ) then
2841  nc_new(k,i,j) = max( ccn(k,i,j), c_ccn )
2842  else
2843  nc_new(k,i,j) = 0.0_rp
2844  endif
2845  enddo
2846  enddo
2847  enddo
2848 
2849  else
2850 
2851  if( nucl_twomey ) then
2852  ! diagnose cloud condensation nuclei
2853 
2854  do j = js, je
2855  do i = is, ie
2856  do k = ks, ke
2857  ! effective vertical velocity (maximum vertical velocity in turbulent flow)
2858  weff_max(k,i,j) = weff(k,i,j) + sigma_w(k,i,j)
2859  ! large scale upward motion region and saturated
2860  if( (weff(k,i,j) > 1.e-8_rp) .AND. (ssw(k,i,j) > 1.e-10_rp) .AND. pre(k,i,j) > 300.e+2_rp )then
2861  ! Lohmann (2002), eq.(1)
2862  nc_new_max = coef_ccn(i,j)*weff_max(k,i,j)**slope_ccn(i,j)
2863  nc_new(k,i,j) = a_max*nc_new_max**b_max
2864  else
2865  nc_new(k,i,j) = 0.0_rp
2866  end if
2867  end do
2868  end do
2869  end do
2870  else
2871  ! calculate cloud condensation nuclei
2872  do j = js, je
2873  do i = is, ie
2874  do k = ks, ke
2875  if( ssw(k,i,j) > 1.e-10_rp .AND. pre(k,i,j) > 300.e+2_rp ) then
2876  nc_new(k,i,j) = c_ccn*ssw(k,i,j)**kappa
2877  else
2878  nc_new(k,i,j) = 0.0_rp
2879  endif
2880  enddo
2881  enddo
2882  enddo
2883  endif
2884 
2885  endif
2886 
2887  do j = js, je
2888  do i = is, ie
2889  do k = ks, ke
2890  ! nc_new is bound by upper limit
2891  if( nc_new(k,i,j) > nc_uplim_d(1,i,j) )then ! no more CCN
2892  flag_nucleation(k,i,j) = .false.
2893  nc_new_below(k+1,i,j) = 1.e+30_rp
2894  else if( nc_new(k,i,j) > eps )then ! nucleation can occur
2895  flag_nucleation(k,i,j) = .true.
2896  nc_new_below(k+1,i,j) = nc_new(k,i,j)
2897  else ! nucleation cannot occur(unsaturated or negative w)
2898  flag_nucleation(k,i,j) = .false.
2899  nc_new_below(k+1,i,j) = 0.0_rp
2900  end if
2901  end do
2902  nc_new_below(ks,i,j) = 0.0_rp
2903 ! do k=KS, KE
2904  ! search maximum value of nc_new
2905 ! if( ( nc_new(k,i,j) < nc_new_below(k,i,j) ) .OR. &
2906 ! ( nc_new_below(k,i,j) > c_ccn_map(1,i,j)*0.05_RP ) )then ! 5% of c_ccn
2907 ! ( nc_new_below(k,i,j) > c_ccn*0.05_RP ) )then ! 5% of c_ccn
2908 ! flag_nucleation(k,i,j) = .false.
2909 ! end if
2910 ! end do
2911 
2912  end do
2913  end do
2914 
2915  if( mp_couple_aerosol ) then
2916 
2917  do j = js, je
2918  do i = is, ie
2919  do k = ks, ke
2920  ! nucleation occurs at only cloud base.
2921  ! if CCN is more than below parcel, nucleation newly occurs
2922  ! effective vertical velocity
2923  if ( flag_nucleation(k,i,j) .AND. & ! large scale upward motion region and saturated
2924  tem(k,i,j) > tem_ccn_low ) then
2925  dlcdt_max = (rhoq(i_qv,k,i,j) - esw(k,i,j)/(rvap*tem(k,i,j)))*rdt
2926  dncdt_max = dlcdt_max/xc_min
2927 ! dnc_new = nc_new(k,i,j)-rhoq(I_NC,k,i,j)
2928  dnc_new = nc_new(k,i,j)
2929  pq(i_ncccn,k,i,j) = min( dncdt_max, dnc_new*rdt )
2930  pq(i_lcccn,k,i,j) = min( dlcdt_max, xc_min*pq(i_ncccn,k,i,j) )
2931  else
2932  pq(i_ncccn,k,i,j) = 0.0_rp
2933  pq(i_lcccn,k,i,j) = 0.0_rp
2934  end if
2935  end do
2936  end do
2937  end do
2938 
2939  else
2940 
2941  if( nucl_twomey ) then
2942  do j = js, je
2943  do i = is, ie
2944  do k = ks, ke
2945  ! nucleation occurs at only cloud base.
2946  ! if CCN is more than below parcel, nucleation newly occurs
2947  ! effective vertical velocity
2948  if ( flag_nucleation(k,i,j) .AND. & ! large scale upward motion region and saturated
2949  tem(k,i,j) > tem_ccn_low .AND. &
2950  nc_new(k,i,j) > rhoq(i_nc,k,i,j) ) then
2951  dlcdt_max = (rhoq(i_qv,k,i,j) - esw(k,i,j)/(rvap*tem(k,i,j)))*rdt
2952  dncdt_max = dlcdt_max/xc_min
2953  dnc_new = nc_new(k,i,j)-rhoq(i_nc,k,i,j)
2954  pq(i_ncccn,k,i,j) = min( dncdt_max, dnc_new*rdt )
2955  pq(i_lcccn,k,i,j) = min( dlcdt_max, xc_min*pq(i_ncccn,k,i,j) )
2956  else
2957  pq(i_ncccn,k,i,j) = 0.0_rp
2958  pq(i_lcccn,k,i,j) = 0.0_rp
2959  end if
2960  end do
2961  end do
2962  end do
2963  else
2964  do j = js, je
2965  do i = is, ie
2966  do k = ks, ke
2967  ! effective vertical velocity
2968  if( tem(k,i,j) > tem_ccn_low .AND. &
2969  nc_new(k,i,j) > rhoq(i_nc,k,i,j) ) then
2970  dlcdt_max = (rhoq(i_qv,k,i,j) - esw(k,i,j)/(rvap*tem(k,i,j)))*rdt
2971  dncdt_max = dlcdt_max/xc_min
2972  dnc_new = nc_new(k,i,j)-rhoq(i_nc,k,i,j)
2973  pq(i_ncccn,k,i,j) = min( dncdt_max, dnc_new*rdt )
2974  pq(i_lcccn,k,i,j) = min( dlcdt_max, xc_min*pq(i_ncccn,k,i,j) )
2975  else
2976  pq(i_ncccn,k,i,j) = 0.0_rp
2977  pq(i_lcccn,k,i,j) = 0.0_rp
2978  end if
2979  end do
2980  end do
2981  end do
2982  endif
2983  endif
2984 
2985  !
2986  ! ice nucleation
2987  !
2988  ! +++ NOTE ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
2989  ! Based on Phillips etal.(2006).
2990  ! However this approach doesn't diagnose Ni itself but diagnose tendency.
2991  ! Original approach adjust Ni instantaneously .
2992  ! +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
2993  do j = js, je
2994  do i = is, ie
2995  do k = ks, ke-1
2996  velz(k) = ( w(k,i,j) * ( cz(k+1,i,j) - fz(k,i,j) ) + w(k+1,i,j) * ( fz(k,i,j) - cz(k,i,j) ) ) / ( cz(k+1,i,j) - cz(k,i,j) ) ! @ half level
2997  end do
2998  velz(ke) = 0.0_rp
2999  do k = ks, ke
3000  dzh = cz(k,i,j) - z_below(k,i,j)
3001  w_dssidz(k,i,j) = velz(k) * (ssi(k,i,j) - ssi_below(k,i,j))/dzh ! 09/04/14 [Add] T.Mitsui
3002  dssidt_rad(k,i,j) = -rhoq(i_qv,k,i,j)/(rho(k,i,j)*qsi(k,i,j)*qsi(k,i,j))*dqsidtem_rho(k,i,j)*dtdt_rad(k,i,j)
3003  dli_max = (rhoq(i_qv,k,i,j) - esi(k,i,j)/(rvap*tem(k,i,j)))*rdt
3004  dni_max = min( dli_max/xi_ccn, (in_max-rhoq(i_ni,k,i,j))*rdt )
3005  wdssi = min( w_dssidz(k,i,j)+dssidt_rad(k,i,j), 0.01_rp)
3006  wssi = min( ssi(k,i,j), ssi_max)
3007  ! SB06(34),(35)
3008  if( ( wdssi > eps ) .AND. & !
3009  (tem(k,i,j) < 273.15_rp ) .AND. & !
3010  (rhoq(i_ni,k,i,j) < in_max ) .AND. &
3011  (wssi >= eps ) )then !
3012 ! PNIccn(k,i,j) = min(dni_max, c_in_map(1,i,j)*bm_M92*nm_M92*0.3_RP*exp(0.3_RP*bm_M92*(wssi-0.1_RP))*wdssi)
3013  if( inucl_w ) then
3014  pq(i_niccn,k,i,j) = min(dni_max, c_in*bm_m92*nm_m92*0.3_rp*exp(0.3_rp*bm_m92*(wssi-0.1_rp))*wdssi)
3015  else
3016  pq(i_niccn,k,i,j) = min(dni_max, max(c_in*nm_m92*exp(0.3_rp*bm_m92*(wssi-0.1_rp) )-rhoq(i_ni,k,i,j),0.0_rp )*rdt )
3017  endif
3018  pq(i_liccn,k,i,j) = min(dli_max, pq(i_niccn,k,i,j)*xi_ccn )
3019  ! only for output
3020 ! dni_ratio(k,i,j) = dssidt_rad(k,i,j)/( w_dssidz(k,i,j)+dssidt_rad(k,i,j) )
3021  else
3022  pq(i_niccn,k,i,j) = 0.0_rp
3023  pq(i_liccn,k,i,j) = 0.0_rp
3024  end if
3025  end do
3026  end do
3027  end do
3028 
3029  return
3030  end subroutine nucleation_kij
3031  !----------------------------
3032  subroutine ice_multiplication_kij( &
3033  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
3034  Pac, & ! in
3035  tem, rhoq, xq, & ! in
3036  PQ ) ! inout
3038  ! ice multiplication by splintering
3039  ! we consider Hallet-Mossop process
3040  use scale_specfunc, only: &
3041  gammafunc => sf_gamma
3042  implicit none
3043 
3044  integer, intent(in) :: KA, KS, KE
3045  integer, intent(in) :: IA, IS, IE
3046  integer, intent(in) :: JA, JS, JE
3047  !
3048  real(RP), intent(in) :: Pac(pac_max,ka,ia,ja)
3049  real(RP), intent(in) :: tem(ka,ia,ja)
3050  real(RP), intent(in) :: rhoq(i_qv:i_ng,ka,ia,ja)
3051  real(RP), intent(in) :: xq(hydro_max,ka,ia,ja)
3052  !
3053  real(RP), intent(inout):: PQ(pq_max,ka,ia,ja)
3054  !
3055  logical, save :: flag_first = .true.
3056  ! production of (350.d3 ice particles)/(cloud riming [g]) => 350*d6 [/kg]
3057  real(RP), parameter :: pice = 350.0e+6_rp
3058  ! production of (1 ice particle)/(250 cloud particles riming)
3059  real(RP), parameter :: pnc = 250.0_rp
3060  real(RP), parameter :: rc_cr= 12.e-6_rp ! critical size[micron]
3061  ! temperature factor
3062  real(RP) :: fp
3063  ! work for incomplete gamma function
3064  real(RP), save :: xc_cr ! mass[kg] of cloud with r=critical size[micron]
3065  real(RP), save :: alpha ! slope parameter of gamma function
3066  real(RP), save :: gm, lgm ! gamma(alpha), log(gamma(alpha))
3067  real(RP) :: igm ! in complete gamma(x,alpha)
3068  real(RP) :: x
3069  ! coefficient of expansion using in calculation of igm
3070  real(RP) :: a0,a1,a2,a3,a4,a5
3071  real(RP) :: a6,a7,a8,a9,a10
3072  real(RP) :: an1,an2,b0,b1,b2,c0,c1,c2
3073  real(RP) :: d0,d1,d2,e1,e2,h0,h1,h2
3074  real(RP), parameter :: eps=1.0e-30_rp
3075  ! number of cloud droplets larger than 12 micron(radius).
3076  real(RP) :: n12
3077  !
3078  real(RP) :: wn, wni, wns, wng
3079  integer :: i, j, k
3080  !
3081  if( flag_first )then
3082  flag_first = .false.
3083  ! work for Incomplete Gamma function
3084  xc_cr = (2.0_rp*rc_cr/a_m(i_mp_qc))**(1.0_rp/b_m(i_mp_qc))
3085  alpha = (nu(i_mp_qc)+1.0_rp)/mu(i_mp_qc)
3086  gm = gammafunc(alpha)
3087  lgm = log(gm)
3088  end if
3089  !
3090  do j = js, je
3091  do i = is, ie
3092  do k = ks, ke
3093  ! Here we assume particle temperature is same as environment temperature.
3094  ! If you want to treat in a better manner,
3095  ! you can diagnose with eq.(64) in CT(86)
3096  if (tem(k,i,j) > 270.16_rp)then
3097  fp = 0.0_rp
3098  else if(tem(k,i,j) >= 268.16_rp)then
3099  fp = (270.16_rp-tem(k,i,j))*0.5_rp
3100  else if(tem(k,i,j) >= 265.16_rp)then
3101  fp = (tem(k,i,j)-265.16_rp)*0.333333333_rp
3102  else
3103  fp = 0.0_rp
3104  end if
3105  ! Approximation of Incomplete Gamma function
3106  ! Here we calculate with algorithm by Numerical Recipes.
3107  ! This approach is based on eq.(78) in Cotton etal.(1986),
3108  ! but more accurate and expanded for Generalized Gamma Distribution.
3109  x = coef_lambda(i_mp_qc)*(xc_cr/xq(i_mp_qc,k,i,j))**mu(i_mp_qc)
3110  !
3111  if(x<1.e-2_rp*alpha)then ! negligible
3112  igm = 0.0_rp
3113  else if(x<alpha+1.0_rp)then ! series expansion
3114  ! 10th-truncation is enough for cloud droplet.
3115  a0 = 1.0_rp/alpha ! n=0
3116  a1 = a0*x/(alpha+1.0_rp) ! n=1
3117  a2 = a1*x/(alpha+2.0_rp) ! n=2
3118  a3 = a2*x/(alpha+3.0_rp) ! n=3
3119  a4 = a3*x/(alpha+4.0_rp) ! n=4
3120  a5 = a4*x/(alpha+5.0_rp) ! n=5
3121  a6 = a5*x/(alpha+6.0_rp) ! n=6
3122  a7 = a6*x/(alpha+7.0_rp) ! n=7
3123  a8 = a7*x/(alpha+8.0_rp) ! n=8
3124  a9 = a8*x/(alpha+9.0_rp) ! n=9
3125  a10 = a9*x/(alpha+10.0_rp) ! n=10
3126  igm = (a0+a1+a2+a3+a4+a5+a6+a7+a8+a9+a10)*exp( -x + alpha*log(x) - lgm )
3127  else if(x<alpha*1.d2) then ! continued fraction expansion
3128  ! 2nd-truncation is enough for cloud droplet.
3129  ! setup
3130  b0 = x+1.0_rp-alpha
3131  c0 = 1.0_rp/eps
3132  d0 = 1.0_rp/b0
3133  h0 = d0
3134  ! n=1
3135  an1 = -(1.0_rp-alpha)
3136  b1 = b0 + 2.0_rp
3137  d1 = 1.0_rp/(an1*d0+b1)
3138  c1 = b1+an1/c0
3139  e1 = d1*c1
3140  h1 = h0*e1
3141  ! n=2
3142  an2 = -2.0_rp*(2.0_rp-alpha)
3143  b2 = b1 + 2.0_rp
3144  d2 = 1.0_rp/(an2*d1+b2)
3145  c2 = b2+an2/c1
3146  e2 = d2*c2
3147  h2 = h1*e2
3148  !
3149  igm = 1.0_rp - exp( -x + alpha*log(x) - lgm )*h2
3150  else ! negligible
3151  igm = 1.0_rp
3152  end if
3153  ! n12 is number of cloud droplets larger than 12 micron.
3154  n12 = rhoq(i_nc,k,i,j)*(1.0_rp-igm)
3155  ! eq.(82) CT(86)
3156  wn = (pice + n12/((rhoq(i_qc,k,i,j)+xc_min)*pnc) )*fp ! filtered by xc_min
3157  wni = wn*(-pac(i_liaclc2li,k,i,j) ) ! riming production rate is all negative
3158  wns = wn*(-pac(i_lsaclc2ls,k,i,j) )
3159  wng = wn*(-pac(i_lgaclc2lg,k,i,j) )
3160  pq(i_nispl,k,i,j) = wni+wns+wng
3161  !
3162  pq(i_lsspl,k,i,j) = - wns*xq(i_mp_qi,k,i,j) ! snow => ice
3163  pq(i_lgspl,k,i,j) = - wng*xq(i_mp_qi,k,i,j) ! graupel => ice
3164  !
3165  end do
3166  end do
3167  end do
3168  !
3169  return
3170  end subroutine ice_multiplication_kij
3171  !----------------------------
3172  subroutine mixed_phase_collection_kij( &
3173  ! collection process
3174  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
3175  wtem, rhoq, & ! in
3176  xq, dq_xave, vt_xave, & ! in
3177  ! rho ! [Add] 11/08/30
3178  PQ, & ! inout
3179  Pac ) ! out
3181  moist_psat_ice => atmos_saturation_psat_ice
3182  implicit none
3183 
3184  integer, intent(in) :: KA, KS, KE
3185  integer, intent(in) :: IA, IS, IE
3186  integer, intent(in) :: JA, JS, JE
3187 
3188  !--- mixed-phase collection process
3189  ! And all we set all production term as a negative sign to avoid confusion.
3190  !
3191  real(RP), intent(in) :: wtem(ka,ia,ja)
3192  !--- mass/number concentration[kg/m3]
3193  real(RP), intent(in) :: rhoq(i_qv:i_ng,ka,ia,ja)
3194  ! necessary ?
3195  real(RP), intent(in) :: xq(hydro_max,ka,ia,ja)
3196  !--- diameter of averaged mass( D(ave x) )
3197  real(RP), intent(in) :: dq_xave(hydro_max,ka,ia,ja)
3198  !--- terminal velocity of averaged mass( vt(ave x) )
3199  real(RP), intent(in) :: vt_xave(hydro_max,2,ka,ia,ja)
3200  ! [Add] 11/08/30 T.Mitsui, for autoconversion of ice
3201  ! real(RP), intent(in) :: rho(KA,IA,JA)
3202  !--- partial conversion
3203  real(RP), intent(inout):: PQ(pq_max,ka,ia,ja)
3204  !
3205  real(RP), intent(out):: Pac(pac_max,ka,ia,ja)
3206  !
3207  ! namelist variables
3208  !=== for collection
3209  !--- threshold of diameters to collide with others
3210  real(RP), save :: dc0 = 15.0e-6_rp ! lower threshold of cloud
3211  real(RP), save :: dc1 = 40.0e-6_rp ! upper threshold of cloud
3212  real(RP), save :: di0 = 150.0e-6_rp ! lower threshold of cloud
3213  real(RP), save :: ds0 = 150.0e-6_rp ! lower threshold of cloud
3214  real(RP), save :: dg0 = 150.0e-6_rp ! lower threshold of cloud
3215  !--- standard deviation of terminal velocity[m/s]
3216  real(RP), save :: sigma_c=0.00_rp ! cloud
3217  real(RP), save :: sigma_r=0.00_rp ! rain
3218  real(RP), save :: sigma_i=0.2_rp ! ice
3219  real(RP), save :: sigma_s=0.2_rp ! snow
3220  real(RP), save :: sigma_g=0.00_rp ! graupel
3221  !--- max collection efficiency for cloud
3222  real(RP), save :: E_im = 0.80_rp ! ice max
3223  real(RP), save :: E_sm = 0.80_rp ! snow max
3224  real(RP), save :: E_gm = 1.00_rp ! graupel max
3225  !--- collection efficiency between 2 species
3226  real(RP), save :: E_ir=1.0_rp ! ice x rain
3227  real(RP), save :: E_sr=1.0_rp ! snow x rain
3228  real(RP), save :: E_gr=1.0_rp ! graupel x rain
3229  real(RP), save :: E_ii=1.0_rp ! ice x ice
3230  real(RP), save :: E_si=1.0_rp ! snow x ice
3231  real(RP), save :: E_gi=1.0_rp ! graupel x ice
3232  real(RP), save :: E_ss=1.0_rp ! snow x snow
3233  real(RP), save :: E_gs=1.0_rp ! graupel x snow
3234  real(RP), save :: E_gg=1.0_rp ! graupel x graupel
3235  !=== for partial conversion
3236  !--- flag: 1=> partial conversion to graupel, 0=> no conversion
3237  integer, save :: i_iconv2g=1 ! ice => graupel
3238  integer, save :: i_sconv2g=1 ! snow => graupel
3239  !--- bulk density of graupel
3240  real(RP), save :: rho_g = 900.0_rp ! [kg/m3]
3241  !--- space filling coefficient [%]
3242  real(RP), save :: cfill_i = 0.68_rp ! ice
3243  real(RP), save :: cfill_s = 0.01_rp ! snow
3244  !--- critical diameter for ice conversion
3245  real(RP), save :: di_cri = 500.e-6_rp ! [m]
3246  ! [Add] 10/08/03 T.Mitsui
3247  logical, save :: opt_stick_KS96=.false.
3248  logical, save :: opt_stick_CO86=.false.
3249  real(RP), parameter :: a_dec = 0.883_rp
3250  real(RP), parameter :: b_dec = 0.093_rp
3251  real(RP), parameter :: c_dec = 0.00348_rp
3252  real(RP), parameter :: d_dec = 4.5185e-5_rp
3253  !
3254  logical, save :: flag_first = .true.
3255  namelist / param_atmos_phy_mp_sn14_collection / &
3256  dc0, dc1, di0, ds0, dg0, &
3257  sigma_c, sigma_r, sigma_i, sigma_s, sigma_g, &
3258  opt_stick_ks96, &
3259  opt_stick_co86, &
3260  e_im, e_sm, e_gm, &
3261  e_ir, e_sr, e_gr, e_ii, e_si, e_gi, e_ss, e_gs, e_gg, &
3262  i_iconv2g, i_sconv2g, rho_g, cfill_i, cfill_s, di_cri
3263  !
3264  real(RP) :: tem(ka,ia,ja)
3265  !
3266  !--- collection efficency of each specie
3267  real(RP) :: E_c, E_r, E_i, E_s, E_g !
3268  real(RP) :: E_ic, E_sc, E_gc !
3269  !--- sticking efficiency
3270  real(RP) :: E_stick(ka,ia,ja)
3271  ! [Add] 10/08/03 T.Mitsui
3272  real(RP) :: temc, temc2, temc3
3273  real(RP) :: E_dec
3274  real(RP) :: esi_rat
3275  real(RP) :: esi(ka,ia,ja)
3276  !
3277  real(RP) :: temc_p, temc_m ! celcius tem.
3278  ! [Add] 11/08/30 T.Mitsui, estimation of autoconversion time
3279 ! real(RP) :: ci_aut(KA,IA,JA)
3280 ! real(RP) :: taui_aut(KA,IA,JA)
3281 ! real(RP) :: tau_sce(KA,IA,JA)
3282  !--- DSD averaged diameter for each species
3283  real(RP) :: ave_dc ! cloud
3284 ! real(RP) :: ave_dr ! rain
3285  real(RP) :: ave_di ! ice
3286  real(RP) :: ave_ds ! snow
3287  real(RP) :: ave_dg ! graupel
3288  !--- coefficient of collection equations(L:mass, N:number)
3289  real(RP) :: coef_acc_LCI, coef_acc_NCI ! cloud - cloud ice
3290  real(RP) :: coef_acc_LCS, coef_acc_NCS ! cloud - snow
3291  !
3292  real(RP) :: coef_acc_LCG, coef_acc_NCG ! cloud - graupel
3293  real(RP) :: coef_acc_LRI_I, coef_acc_NRI_I ! rain - cloud ice
3294  real(RP) :: coef_acc_LRI_R, coef_acc_NRI_R ! rain - cloud ice
3295  real(RP) :: coef_acc_LRS_S, coef_acc_NRS_S ! rain - snow
3296  real(RP) :: coef_acc_LRS_R, coef_acc_NRS_R ! rain - snow
3297  real(RP) :: coef_acc_LRG, coef_acc_NRG ! rain - graupel
3298  real(RP) :: coef_acc_LII, coef_acc_NII ! cloud ice - cloud ice
3299  real(RP) :: coef_acc_LIS, coef_acc_NIS ! cloud ice - snow
3300  real(RP) :: coef_acc_NSS ! snow - snow
3301  real(RP) :: coef_acc_NGG ! grauepl - graupel
3302  real(RP) :: coef_acc_LSG, coef_acc_NSG ! snow - graupel
3303  !--- (diameter) x (diameter)
3304  real(RP) :: dcdc, dcdi, dcds, dcdg
3305  real(RP) :: drdr, drdi, drds, drdg
3306  real(RP) :: didi, dids, didg
3307  real(RP) :: dsds, dsdg
3308  real(RP) :: dgdg
3309  !--- (terminal velocity) x (terminal velocity)
3310  real(RP) :: vcvc, vcvi, vcvs, vcvg
3311  real(RP) :: vrvr, vrvi, vrvs, vrvg
3312  real(RP) :: vivi, vivs, vivg
3313  real(RP) :: vsvs, vsvg
3314  real(RP) :: vgvg
3315  !
3316  real(RP) :: wx_cri, wx_crs
3317  real(RP) :: coef_emelt
3318  real(RP) :: w1
3319 
3320  real(RP) :: sw
3321  !
3322  integer :: i, j, k
3323  !
3324  if( flag_first )then
3325  rewind( io_fid_conf )
3326  read( io_fid_conf, nml=param_atmos_phy_mp_sn14_collection, end=100 )
3327 100 log_nml(param_atmos_phy_mp_sn14_collection)
3328  flag_first = .false.
3329  end if
3330  !
3331  ! [Add] 10/08/03 T.Mitsui
3332  do j = js, je
3333  do i = is, ie
3334  do k = ks, ke
3335  tem(k,i,j) = max( wtem(k,i,j), tem_min ) ! 11/08/30 T.Mitsui
3336  end do
3337  end do
3338  end do
3339 
3340  call moist_psat_ice( ka, ks, ke, ia, is, ie, ja, js, je, &
3341  tem(:,:,:), esi(:,:,:) ) ! [IN], [OUT]
3342 
3343  if( opt_stick_ks96 )then
3344  do j = js, je
3345  do i = is, ie
3346  do k = ks, ke
3347  ! Khain and Sednev (1996), eq.(3.15)
3348  temc = tem(k,i,j) - t00
3349  temc2 = temc*temc
3350  temc3 = temc2*temc
3351  e_dec = max(0.0_rp, a_dec + b_dec*temc + c_dec*temc2 + d_dec*temc3 )
3352  esi_rat = rhoq(i_qv,k,i,j)*rvap*tem(k,i,j)/esi(k,i,j)
3353  e_stick(k,i,j) = min(1.0_rp, e_dec*esi_rat)
3354  end do
3355  end do
3356  end do
3357  else if( opt_stick_co86 )then
3358  do j = js, je
3359  do i = is, ie
3360  do k = ks, ke
3361  ! [Add] 11/08/30 T.Mitsui, Cotton et al. (1986)
3362  temc = min(tem(k,i,j) - t00,0.0_rp)
3363  w1 = 0.035_rp*temc-0.7_rp
3364  e_stick(k,i,j) = 10._rp**w1
3365  end do
3366  end do
3367  end do
3368  else
3369  do j = js, je
3370  do i = is, ie
3371  do k = ks, ke
3372  ! Lin et al. (1983)
3373  temc_m = min(tem(k,i,j) - t00,0.0_rp) ! T < 273.15
3374  e_stick(k,i,j) = exp(0.09_rp*temc_m)
3375  end do
3376  end do
3377  end do
3378  end if
3379  !
3380  profile_start("sn14_collection")
3381 !OCL NOSIMD
3382  do j = js, je
3383  do i = is, ie
3384  do k = ks, ke
3385  !
3386 ! temc_m = min(tem(k,i,j) - T00,0.0_RP) ! T < 273.15
3387  temc_p = max(tem(k,i,j) - t00,0.0_rp) ! T > 273.15
3388  ! averaged diameter using SB06(82)
3389  ave_dc = coef_d(i_mp_qc)*xq(i_mp_qc,k,i,j)**b_m(i_mp_qc)
3390  ave_di = coef_d(i_mp_qi)*xq(i_mp_qi,k,i,j)**b_m(i_mp_qi)
3391  ave_ds = coef_d(i_mp_qs)*xq(i_mp_qs,k,i,j)**b_m(i_mp_qs)
3392  ave_dg = coef_d(i_mp_qg)*xq(i_mp_qg,k,i,j)**b_m(i_mp_qg)
3393  !------------------------------------------------------------------------
3394  ! coellection efficiency are given as follows
3395  e_c = max(0.0_rp, min(1.0_rp, (ave_dc-dc0)/(dc1-dc0) ))
3396  sw = 0.5_rp - sign(0.5_rp, di0-ave_di) ! if(ave_di>di0)then sw=1
3397  e_i = e_im * sw
3398  sw = 0.5_rp - sign(0.5_rp, ds0-ave_ds) ! if(ave_ds>ds0)then sw=1
3399  e_s = e_sm * sw
3400  sw = 0.5_rp - sign(0.5_rp, dg0-ave_dg) ! if(ave_dg>dg0)then sw=1
3401  e_g = e_gm * sw
3402  e_ic = e_i*e_c
3403  e_sc = e_s*e_c
3404  e_gc = e_g*e_c
3405  !------------------------------------------------------------------------
3406  ! Collection: a collects b ( assuming particle size a>b )
3407  dcdc = dq_xave(i_mp_qc,k,i,j) * dq_xave(i_mp_qc,k,i,j)
3408  drdr = dq_xave(i_mp_qr,k,i,j) * dq_xave(i_mp_qr,k,i,j)
3409  didi = dq_xave(i_mp_qi,k,i,j) * dq_xave(i_mp_qi,k,i,j)
3410  dsds = dq_xave(i_mp_qs,k,i,j) * dq_xave(i_mp_qs,k,i,j)
3411  dgdg = dq_xave(i_mp_qg,k,i,j) * dq_xave(i_mp_qg,k,i,j)
3412  dcdi = dq_xave(i_mp_qc,k,i,j) * dq_xave(i_mp_qi,k,i,j)
3413  dcds = dq_xave(i_mp_qc,k,i,j) * dq_xave(i_mp_qs,k,i,j)
3414  dcdg = dq_xave(i_mp_qc,k,i,j) * dq_xave(i_mp_qg,k,i,j)
3415  drdi = dq_xave(i_mp_qr,k,i,j) * dq_xave(i_mp_qi,k,i,j)
3416  drds = dq_xave(i_mp_qr,k,i,j) * dq_xave(i_mp_qs,k,i,j)
3417  drdg = dq_xave(i_mp_qr,k,i,j) * dq_xave(i_mp_qg,k,i,j)
3418  dids = dq_xave(i_mp_qi,k,i,j) * dq_xave(i_mp_qs,k,i,j)
3419 ! didg = dq_xave(I_mp_QI,k,i,j) * dq_xave(I_mp_QG,k,i,j)
3420  dsdg = dq_xave(i_mp_qs,k,i,j) * dq_xave(i_mp_qg,k,i,j)
3421  !
3422  vcvc = vt_xave(i_mp_qc,2,k,i,j)* vt_xave(i_mp_qc,2,k,i,j)
3423  vrvr = vt_xave(i_mp_qr,2,k,i,j)* vt_xave(i_mp_qr,2,k,i,j)
3424  vivi = vt_xave(i_mp_qi,2,k,i,j)* vt_xave(i_mp_qi,2,k,i,j)
3425  vsvs = vt_xave(i_mp_qs,2,k,i,j)* vt_xave(i_mp_qs,2,k,i,j)
3426  vgvg = vt_xave(i_mp_qg,2,k,i,j)* vt_xave(i_mp_qg,2,k,i,j)
3427  vcvi = vt_xave(i_mp_qc,2,k,i,j)* vt_xave(i_mp_qi,2,k,i,j)
3428  vcvs = vt_xave(i_mp_qc,2,k,i,j)* vt_xave(i_mp_qs,2,k,i,j)
3429  vcvg = vt_xave(i_mp_qc,2,k,i,j)* vt_xave(i_mp_qg,2,k,i,j)
3430  vrvi = vt_xave(i_mp_qr,2,k,i,j)* vt_xave(i_mp_qi,2,k,i,j)
3431  vrvs = vt_xave(i_mp_qr,2,k,i,j)* vt_xave(i_mp_qs,2,k,i,j)
3432  vrvg = vt_xave(i_mp_qr,2,k,i,j)* vt_xave(i_mp_qg,2,k,i,j)
3433  vivs = vt_xave(i_mp_qi,2,k,i,j)* vt_xave(i_mp_qs,2,k,i,j)
3434 ! vivg = vt_xave(I_mp_QI,2,k,i,j)* vt_xave(I_mp_QG,2,k,i,j)
3435  vsvg = vt_xave(i_mp_qs,2,k,i,j)* vt_xave(i_mp_qg,2,k,i,j)
3436  !------------------------------------------------------------------------
3437  !
3438  !+++ pattern 1: a + b => a (a>b)
3439  ! (i-c, s-c, g-c, s-i, g-r, s-g)
3440  !------------------------------------------------------------------------
3441  ! cloud-ice => ice
3442  ! reduction term of cloud
3443  coef_acc_lci = &
3444  ( delta_b1(i_mp_qc)*dcdc + delta_ab1(i_mp_qi,i_mp_qc)*dcdi + delta_b0(i_mp_qi)*didi ) &
3445  * ( theta_b1(i_mp_qc)*vcvc - theta_ab1(i_mp_qi,i_mp_qc)*vcvi + theta_b0(i_mp_qi)*vivi &
3446  + sigma_i + sigma_c )**0.5_rp
3447  coef_acc_nci = &
3448  ( delta_b0(i_mp_qc)*dcdc + delta_ab0(i_mp_qi,i_mp_qc)*dcdi + delta_b0(i_mp_qi)*didi ) &
3449  * ( theta_b0(i_mp_qc)*vcvc - theta_ab0(i_mp_qi,i_mp_qc)*vcvi + theta_b0(i_mp_qi)*vivi &
3450  + sigma_i + sigma_c )**0.5_rp
3451  pac(i_liaclc2li,k,i,j)= -0.25_rp*pi*e_ic*rhoq(i_ni,k,i,j)*rhoq(i_qc,k,i,j)*coef_acc_lci
3452  pac(i_niacnc2ni,k,i,j)= -0.25_rp*pi*e_ic*rhoq(i_ni,k,i,j)*rhoq(i_nc,k,i,j)*coef_acc_nci
3453  ! cloud-snow => snow
3454  ! reduction term of cloud
3455  coef_acc_lcs = &
3456  ( delta_b1(i_mp_qc)*dcdc + delta_ab1(i_mp_qs,i_mp_qc)*dcds + delta_b0(i_mp_qs)*dsds ) &
3457  * ( theta_b1(i_mp_qc)*vcvc - theta_ab1(i_mp_qs,i_mp_qc)*vcvs + theta_b0(i_mp_qs)*vsvs &
3458  + sigma_s + sigma_c )**0.5_rp
3459  coef_acc_ncs = &
3460  ( delta_b0(i_mp_qc)*dcdc + delta_ab0(i_mp_qs,i_mp_qc)*dcds + delta_b0(i_mp_qs)*dsds ) &
3461  * ( theta_b0(i_mp_qc)*vcvc - theta_ab0(i_mp_qs,i_mp_qc)*vcvs + theta_b0(i_mp_qs)*vsvs &
3462  + sigma_s + sigma_c )**0.5_rp
3463  pac(i_lsaclc2ls,k,i,j)= -0.25_rp*pi*e_sc*rhoq(i_ns,k,i,j)*rhoq(i_qc,k,i,j)*coef_acc_lcs
3464  pac(i_nsacnc2ns,k,i,j)= -0.25_rp*pi*e_sc*rhoq(i_ns,k,i,j)*rhoq(i_nc,k,i,j)*coef_acc_ncs
3465  ! cloud-graupel => graupel
3466  ! reduction term of cloud
3467  coef_acc_lcg = &
3468  ( delta_b1(i_mp_qc)*dcdc + delta_ab1(i_mp_qg,i_mp_qc)*dcdg + delta_b0(i_mp_qg)*dgdg ) &
3469  * ( theta_b1(i_mp_qc)*vcvc - theta_ab1(i_mp_qg,i_mp_qc)*vcvg + theta_b0(i_mp_qg)*vgvg &
3470  + sigma_g + sigma_c )**0.5_rp
3471  coef_acc_ncg = &
3472  ( delta_b0(i_mp_qc)*dcdc + delta_ab0(i_mp_qg,i_mp_qc)*dcdg + delta_b0(i_mp_qg)*dgdg ) &
3473  * ( theta_b0(i_mp_qc)*vcvc - theta_ab0(i_mp_qg,i_mp_qc)*vcvg + theta_b0(i_mp_qg)*vgvg &
3474  + sigma_g + sigma_c )**0.5_rp
3475  pac(i_lgaclc2lg,k,i,j)= -0.25_rp*pi*e_gc*rhoq(i_ng,k,i,j)*rhoq(i_qc,k,i,j)*coef_acc_lcg
3476  pac(i_ngacnc2ng,k,i,j)= -0.25_rp*pi*e_gc*rhoq(i_ng,k,i,j)*rhoq(i_nc,k,i,j)*coef_acc_ncg
3477  ! snow-graupel => graupel
3478  coef_acc_lsg = &
3479  ( delta_b1(i_mp_qs)*dsds + delta_ab1(i_mp_qg,i_mp_qs)*dsdg + delta_b0(i_mp_qg)*dgdg ) &
3480  * ( theta_b1(i_mp_qs)*vsvs - theta_ab1(i_mp_qg,i_mp_qs)*vsvg + theta_b0(i_mp_qg)*vgvg &
3481  + sigma_g + sigma_s )**0.5_rp
3482  coef_acc_nsg = &
3483  ( delta_b0(i_mp_qs)*dsds + delta_ab0(i_mp_qg,i_mp_qs)*dsdg + delta_b0(i_mp_qg)*dgdg ) &
3484  ! [fix] T.Mitsui 08/05/08
3485  * ( theta_b0(i_mp_qs)*vsvs - theta_ab0(i_mp_qg,i_mp_qs)*vsvg + theta_b0(i_mp_qg)*vgvg &
3486  + sigma_g + sigma_s )**0.5_rp
3487  pac(i_lgacls2lg,k,i,j)= -0.25_rp*pi*e_stick(k,i,j)*e_gs*rhoq(i_ng,k,i,j)*rhoq(i_qs,k,i,j)*coef_acc_lsg
3488  pac(i_ngacns2ng,k,i,j)= -0.25_rp*pi*e_stick(k,i,j)*e_gs*rhoq(i_ng,k,i,j)*rhoq(i_ns,k,i,j)*coef_acc_nsg
3489  !------------------------------------------------------------------------
3490  ! ice-snow => snow
3491  ! reduction term of ice
3492  coef_acc_lis = &
3493  ( delta_b1(i_mp_qi)*didi + delta_ab1(i_mp_qs,i_mp_qi)*dids + delta_b0(i_mp_qs)*dsds ) &
3494  * ( theta_b1(i_mp_qi)*vivi - theta_ab1(i_mp_qs,i_mp_qi)*vivs + theta_b0(i_mp_qs)*vsvs &
3495  + sigma_i + sigma_s )**0.5_rp
3496  coef_acc_nis = &
3497  ( delta_b0(i_mp_qi)*didi + delta_ab0(i_mp_qs,i_mp_qi)*dids + delta_b0(i_mp_qs)*dsds ) &
3498  * ( theta_b0(i_mp_qi)*vivi - theta_ab0(i_mp_qs,i_mp_qi)*vivs + theta_b0(i_mp_qs)*vsvs &
3499  + sigma_i + sigma_s )**0.5_rp
3500  pac(i_liacls2ls,k,i,j)= -0.25_rp*pi*e_stick(k,i,j)*e_si*rhoq(i_ns,k,i,j)*rhoq(i_qi,k,i,j)*coef_acc_lis
3501  pac(i_niacns2ns,k,i,j)= -0.25_rp*pi*e_stick(k,i,j)*e_si*rhoq(i_ns,k,i,j)*rhoq(i_ni,k,i,j)*coef_acc_nis
3502  !
3503  sw = sign(0.5_rp, t00-tem(k,i,j)) + 0.5_rp
3504  ! if ( tem(k,i,j) <= T00 )then
3505  ! rain-graupel => graupel
3506  ! reduction term of rain
3507  ! sw = 1
3508  ! else
3509  ! rain-graupel => rain
3510  ! reduction term of graupel
3511  ! sw = 0
3512  coef_acc_lrg = &
3513  ( ( delta_b1(i_mp_qr)*drdr + delta_ab1(i_mp_qg,i_mp_qr)*drdg + delta_b0(i_mp_qg)*dgdg ) * sw &
3514  + ( delta_b1(i_mp_qg)*dgdg + delta_ab1(i_mp_qr,i_mp_qg)*drdg + delta_b0(i_mp_qr)*drdr ) * (1.0_rp-sw) ) &
3515  * sqrt( ( theta_b1(i_mp_qr)*vrvr - theta_ab1(i_mp_qg,i_mp_qr)*vrvg + theta_b0(i_mp_qg)*vgvg ) * sw &
3516  + ( theta_b1(i_mp_qg)*vgvg - theta_ab1(i_mp_qr,i_mp_qg)*vrvg + theta_b0(i_mp_qr)*vrvr ) * (1.0_rp-sw) &
3517  + sigma_r + sigma_g )
3518  pac(i_lraclg2lg,k,i,j) = -0.25_rp*pi*e_gr*coef_acc_lrg &
3519  * ( rhoq(i_ng,k,i,j)*rhoq(i_qr,k,i,j) * sw &
3520  + rhoq(i_nr,k,i,j)*rhoq(i_qg,k,i,j) * (1.0_rp-sw) )
3521  coef_acc_nrg = &
3522  ( delta_b0(i_mp_qr)*drdr + delta_ab0(i_mp_qg,i_mp_qr)*drdg + delta_b0(i_mp_qg)*dgdg ) &
3523  * ( theta_b0(i_mp_qr)*vrvr - theta_ab0(i_mp_qg,i_mp_qr)*vrvg + theta_b0(i_mp_qg)*vgvg &
3524  + sigma_r + sigma_g )**0.5_rp
3525  pac(i_nracng2ng,k,i,j) = -0.25_rp*pi*e_gr*rhoq(i_ng,k,i,j)*rhoq(i_nr,k,i,j)*coef_acc_nrg
3526  !
3527  !------------------------------------------------------------------------
3528  !
3529  !+++ pattern 2: a + b => c (a>b)
3530  ! (r-i,r-s)
3531  !------------------------------------------------------------------------
3532  ! rain-ice => graupel
3533  ! reduction term of ice
3534  coef_acc_lri_i = &
3535  ( delta_b1(i_mp_qi)*didi + delta_ab1(i_mp_qr,i_mp_qi)*drdi + delta_b0(i_mp_qr)*drdr ) &
3536  * ( theta_b1(i_mp_qi)*vivi - theta_ab1(i_mp_qr,i_mp_qi)*vrvi + theta_b0(i_mp_qr)*vrvr &
3537  + sigma_r + sigma_i )**0.5_rp
3538  coef_acc_nri_i = &
3539  ( delta_b0(i_mp_qi)*didi + delta_ab0(i_mp_qr,i_mp_qi)*drdi + delta_b0(i_mp_qr)*drdr ) &
3540  * ( theta_b0(i_mp_qi)*vivi - theta_ab0(i_mp_qr,i_mp_qi)*vrvi + theta_b0(i_mp_qr)*vrvr &
3541  + sigma_r + sigma_i )**0.5_rp
3542  pac(i_lracli2lg_i,k,i,j)= -0.25_rp*pi*e_ir*rhoq(i_nr,k,i,j)*rhoq(i_qi,k,i,j)*coef_acc_lri_i
3543  pac(i_nracni2ng_i,k,i,j)= -0.25_rp*pi*e_ir*rhoq(i_nr,k,i,j)*rhoq(i_ni,k,i,j)*coef_acc_nri_i
3544  ! reduction term of rain
3545  coef_acc_lri_r = &
3546  ( delta_b1(i_mp_qr)*drdr + delta_ab1(i_mp_qi,i_mp_qr)*drdi + delta_b0(i_mp_qi)*didi ) &
3547  * ( theta_b1(i_mp_qr)*vrvr - theta_ab1(i_mp_qi,i_mp_qr)*vrvi + theta_b0(i_mp_qi)*vivi &
3548  + sigma_r + sigma_i )**0.5_rp
3549  coef_acc_nri_r = &
3550  ( delta_b0(i_mp_qr)*drdr + delta_ab0(i_mp_qi,i_mp_qr)*drdi + delta_b0(i_mp_qi)*didi ) &
3551  * ( theta_b0(i_mp_qr)*vrvr - theta_ab0(i_mp_qi,i_mp_qr)*vrvi + theta_b0(i_mp_qi)*vivi &
3552  + sigma_r + sigma_i )**0.5_rp
3553  pac(i_lracli2lg_r,k,i,j)= -0.25_rp*pi*e_ir*rhoq(i_ni,k,i,j)*rhoq(i_qr,k,i,j)*coef_acc_lri_r
3554  pac(i_nracni2ng_r,k,i,j)= -0.25_rp*pi*e_ir*rhoq(i_ni,k,i,j)*rhoq(i_nr,k,i,j)*coef_acc_nri_r
3555  ! rain-snow => graupel
3556  ! reduction term of snow
3557  coef_acc_lrs_s = &
3558  ( delta_b1(i_mp_qs)*dsds + delta_ab1(i_mp_qr,i_mp_qs)*drds + delta_b0(i_mp_qr)*drdr ) &
3559  * ( theta_b1(i_mp_qs)*vsvs - theta_ab1(i_mp_qr,i_mp_qs)*vrvs + theta_b0(i_mp_qr)*vrvr &
3560  + sigma_r + sigma_s )**0.5_rp
3561  coef_acc_nrs_s = &
3562  ( delta_b0(i_mp_qs)*dsds + delta_ab0(i_mp_qr,i_mp_qs)*drds + delta_b0(i_mp_qr)*drdr ) &
3563  * ( theta_b0(i_mp_qs)*vsvs - theta_ab0(i_mp_qr,i_mp_qs)*vrvs + theta_b0(i_mp_qr)*vrvr &
3564  + sigma_r + sigma_s )**0.5_rp
3565  pac(i_lracls2lg_s,k,i,j)= -0.25_rp*pi*e_sr*rhoq(i_nr,k,i,j)*rhoq(i_qs,k,i,j)*coef_acc_lrs_s
3566  pac(i_nracns2ng_s,k,i,j)= -0.25_rp*pi*e_sr*rhoq(i_nr,k,i,j)*rhoq(i_ns,k,i,j)*coef_acc_nrs_s
3567  ! reduction term of rain
3568  coef_acc_lrs_r = &
3569  ( delta_b1(i_mp_qr)*drdr + delta_ab1(i_mp_qs,i_mp_qr)*drds + delta_b0(i_mp_qs)*dsds ) &
3570  * ( theta_b1(i_mp_qr)*vrvr - theta_ab1(i_mp_qs,i_mp_qr)*vrvs + theta_b0(i_mp_qs)*vsvs &
3571  + sigma_r + sigma_s )**0.5_rp
3572  coef_acc_nrs_r = &
3573  ( delta_b0(i_mp_qr)*drdr + delta_ab0(i_mp_qs,i_mp_qr)*drds + delta_b0(i_mp_qs)*dsds ) &
3574  * ( theta_b0(i_mp_qr)*vrvr - theta_ab0(i_mp_qs,i_mp_qr)*vrvs + theta_b0(i_mp_qs)*vsvs &
3575  + sigma_r + sigma_s )**0.5_rp
3576  pac(i_lracls2lg_r,k,i,j)= -0.25_rp*pi*e_sr*rhoq(i_ns,k,i,j)*rhoq(i_qr,k,i,j)*coef_acc_lrs_r
3577  pac(i_nracns2ng_r,k,i,j)= -0.25_rp*pi*e_sr*rhoq(i_ns,k,i,j)*rhoq(i_nr,k,i,j)*coef_acc_nrs_r
3578  !------------------------------------------------------------------------
3579  !
3580  !+++ pattern 3: a + a => b (i-i)
3581  !
3582  !------------------------------------------------------------------------
3583  ! ice-ice ( reduction is double, but includes double-count)
3584  coef_acc_lii = &
3585  ( delta_b0(i_mp_qi)*didi + delta_ab1(i_mp_qi,i_mp_qi)*didi + delta_b1(i_mp_qi)*didi ) &
3586  * ( theta_b0(i_mp_qi)*vivi - theta_ab1(i_mp_qi,i_mp_qi)*vivi + theta_b1(i_mp_qi)*vivi &
3587  + sigma_i + sigma_i )**0.5_rp
3588  coef_acc_nii = &
3589  ( delta_b0(i_mp_qi)*didi + delta_ab0(i_mp_qi,i_mp_qi)*didi + delta_b0(i_mp_qi)*didi ) &
3590  * ( theta_b0(i_mp_qi)*vivi - theta_ab0(i_mp_qi,i_mp_qi)*vivi + theta_b0(i_mp_qi)*vivi &
3591  + sigma_i + sigma_i )**0.5_rp
3592  pac(i_liacli2ls,k,i,j)= -0.25_rp*pi*e_stick(k,i,j)*e_ii*rhoq(i_ni,k,i,j)*rhoq(i_qi,k,i,j)*coef_acc_lii
3593  pac(i_niacni2ns,k,i,j)= -0.25_rp*pi*e_stick(k,i,j)*e_ii*rhoq(i_ni,k,i,j)*rhoq(i_ni,k,i,j)*coef_acc_nii
3594  !
3595 ! ci_aut(k,i,j) = 0.25_RP*pi*E_ii*rhoq(I_NI,k,i,j)*coef_acc_LII
3596 ! taui_aut(k,i,j) = 1._RP/max(E_stick(k,i,j)*ci_aut(k,i,j),1.E-10_RP)
3597 ! tau_sce(k,i,j) = rhoq(I_QI,k,i,j)/max(rhoq(I_QI,k,i,j)+rhoq(I_QS,k,i,j),1.E-10_RP)
3598  !------------------------------------------------------------------------
3599  !
3600  !+++ pattern 4: a + a => a (s-s)
3601  !
3602  !------------------------------------------------------------------------
3603  ! snow-snow => snow
3604  coef_acc_nss = &
3605  ( delta_b0(i_mp_qs)*dsds + delta_ab0(i_mp_qs,i_mp_qs)*dsds + delta_b0(i_mp_qs)*dsds ) &
3606  * ( theta_b0(i_mp_qs)*vsvs - theta_ab0(i_mp_qs,i_mp_qs)*vsvs + theta_b0(i_mp_qs)*vsvs &
3607  + sigma_s + sigma_s )**0.5_rp
3608  pac(i_nsacns2ns,k,i,j)= -0.125_rp*pi*e_stick(k,i,j)*e_ss*rhoq(i_ns,k,i,j)*rhoq(i_ns,k,i,j)*coef_acc_nss
3609  !
3610  ! graupel-grauple => graupel
3611  coef_acc_ngg = &
3612  ( delta_b0(i_mp_qg)*dgdg + delta_ab0(i_mp_qg,i_mp_qg)*dgdg + delta_b0(i_mp_qg)*dgdg ) &
3613  * ( theta_b0(i_mp_qg)*vgvg - theta_ab0(i_mp_qg,i_mp_qg)*vgvg + theta_b0(i_mp_qg)*vgvg &
3614  + sigma_g + sigma_g )**0.5_rp
3615  pac(i_ngacng2ng,k,i,j)= -0.125_rp*pi*e_stick(k,i,j)*e_gg*rhoq(i_ng,k,i,j)*rhoq(i_ng,k,i,j)*coef_acc_ngg
3616  !
3617  !------------------------------------------------------------------------
3618  !--- Partial conversion
3619  ! SB06(70),(71)
3620  ! i_iconv2g: option whether partial conversions work or not
3621  ! ice-cloud => graupel
3622  sw = 0.5_rp - sign(0.5_rp,di_cri-ave_di) ! if( ave_di > di_cri )then sw=1
3623  wx_cri = cfill_i*dwatr/rho_g*( pi/6.0_rp*rho_g*ave_di*ave_di*ave_di/xq(i_mp_qi,k,i,j) - 1.0_rp ) * sw
3624  pq(i_licon,k,i,j) = i_iconv2g * pac(i_liaclc2li,k,i,j)/max(1.0_rp, wx_cri) * sw
3625  pq(i_nicon,k,i,j) = i_iconv2g * pq(i_licon,k,i,j)/xq(i_mp_qi,k,i,j) * sw
3626 
3627  ! snow-cloud => graupel
3628  wx_crs = cfill_s*dwatr/rho_g*( pi/6.0_rp*rho_g*ave_ds*ave_ds*ave_ds/xq(i_mp_qs,k,i,j) - 1.0_rp )
3629  pq(i_lscon,k,i,j) = i_sconv2g * (pac(i_lsaclc2ls,k,i,j))/max(1.0_rp, wx_crs)
3630  pq(i_nscon,k,i,j) = i_sconv2g * pq(i_lscon,k,i,j)/xq(i_mp_qs,k,i,j)
3631  !------------------------------------------------------------------------
3632  !--- enhanced melting( due to collection-freezing of water droplets )
3633  ! originally from Rutledge and Hobbs(1984). eq.(A.21)
3634  ! if T > 273.15 then temc_p=T-273.15, else temc_p=0
3635  ! 08/05/08 [fix] T.Mitsui LHF00 => LHF0
3636  ! melting occurs around T=273K, so LHF0 is suitable both SIMPLE and EXACT,
3637  ! otherwise LHF can have sign both negative(EXACT) and positive(SIMPLE).
3638 !!$ coef_emelt = -CL/LHF00*temc_p
3639  coef_emelt = cl/lhf0*temc_p
3640  ! cloud-graupel
3641  pq(i_lgacm,k,i,j) = coef_emelt*pac(i_lgaclc2lg,k,i,j)
3642  pq(i_ngacm,k,i,j) = pq(i_lgacm,k,i,j)/xq(i_mp_qg,k,i,j)
3643  ! rain-graupel
3644  pq(i_lgarm,k,i,j) = coef_emelt*pac(i_lraclg2lg,k,i,j)
3645  pq(i_ngarm,k,i,j) = pq(i_lgarm,k,i,j)/xq(i_mp_qg,k,i,j)
3646  ! cloud-snow
3647  pq(i_lsacm,k,i,j) = coef_emelt*(pac(i_lsaclc2ls,k,i,j))
3648  pq(i_nsacm,k,i,j) = pq(i_lsacm,k,i,j)/xq(i_mp_qs,k,i,j)
3649  ! rain-snow
3650  pq(i_lsarm,k,i,j) = coef_emelt*(pac(i_lracls2lg_r,k,i,j)+pac(i_lracls2lg_s,k,i,j))
3651  pq(i_nsarm,k,i,j) = pq(i_lsarm,k,i,j)/xq(i_mp_qg,k,i,j) ! collect? might be I_mp_QS
3652  ! cloud-ice
3653  pq(i_liacm,k,i,j) = coef_emelt*pac(i_liaclc2li,k,i,j)
3654  pq(i_niacm,k,i,j) = pq(i_liacm,k,i,j)/xq(i_mp_qi,k,i,j)
3655  ! rain-ice
3656  pq(i_liarm,k,i,j) = coef_emelt*(pac(i_lracli2lg_r,k,i,j)+pac(i_lracli2lg_i,k,i,j))
3657  pq(i_niarm,k,i,j) = pq(i_liarm,k,i,j)/xq(i_mp_qg,k,i,j) ! collect? might be I_mp_QI
3658  end do
3659  end do
3660  end do
3661  profile_stop("sn14_collection")
3662 
3663  !
3664  return
3665  end subroutine mixed_phase_collection_kij
3666  !----------------------------
3667  ! Auto-conversion, Accretion, Self-collection, Break-up
3668  subroutine aut_acc_slc_brk_kij( &
3669  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
3670  rhoq, xq, dq_xave, &
3671  rho, &
3672  PQ )
3673  implicit none
3674 
3675  integer, intent(in) :: KA, KS, KE
3676  integer, intent(in) :: IA, IS, IE
3677  integer, intent(in) :: JA, JS, JE
3678  !
3679  real(RP), intent(in) :: rhoq(i_qv:i_ng,ka,ia,ja)
3680  real(RP), intent(in) :: xq(hydro_max,ka,ia,ja)
3681  real(RP), intent(in) :: dq_xave(hydro_max,ka,ia,ja)
3682  real(RP), intent(in) :: rho(ka,ia,ja)
3683  !
3684  real(RP), intent(inout) :: PQ(pq_max,ka,ia,ja)
3685  !
3686  ! parameter for autoconversion
3687  real(RP), parameter :: kcc = 4.44e+9_rp ! collision efficiency [m3/kg2/sec]
3688  real(RP), parameter :: tau_min = 1.e-20_rp ! empirical filter by T.Mitsui
3689  real(RP), parameter :: rx_sep = 1.0_rp/x_sep ! 1/x_sep, 10/08/03 [Add] T.Mitsui
3690  !
3691  ! parameter for accretion
3692  real(RP), parameter :: kcr = 5.8_rp ! collision efficiency [m3/kg2/sec]
3693  real(RP), parameter :: thr_acc = 5.e-5_rp ! threshold for universal function original
3694  !
3695  ! parameter for self collection and collison break-up
3696  real(RP), parameter :: krr = 4.33_rp ! k_rr, S08 (35)
3697  real(RP), parameter :: kaprr = 60.7_rp ! kappa_rr, SB06(11)
3698  real(RP), parameter :: kbr = 1000._rp ! k_br, SB06(14)
3699  real(RP), parameter :: kapbr = 2.3e+3_rp ! kappa_br, SB06(15)
3700  real(RP), parameter :: dr_min = 0.35e-3_rp ! minimum diameter, SB06(13)-(15)
3701  !
3702  ! work variables
3703  real(RP) :: coef_nuc0 ! coefficient of number for Auto-conversion
3704  real(RP) :: coef_nuc1 ! mass
3705  real(RP) :: coef_aut0 ! number
3706  real(RP) :: coef_aut1 ! mass
3707  real(RP) :: lwc ! lc+lr
3708  real(RP) :: tau ! conversion ratio: qr/(qc+qr) ranges [0:1]
3709  real(RP) :: rho_fac ! factor of air density
3710  real(RP) :: psi_aut ! Universal function of Auto-conversion
3711  real(RP) :: psi_acc ! Universal function of Accretion
3712  real(RP) :: psi_brk ! Universal function of Breakup
3713  real(RP) :: ddr ! diameter difference from equilibrium
3714  !
3715  integer :: i, j, k
3716  !
3717  coef_nuc0 = (nu(i_mp_qc)+2.0_rp)/(nu(i_mp_qc)+1.0_rp)
3718  coef_nuc1 = (nu(i_mp_qc)+2.0_rp)*(nu(i_mp_qc)+4.0_rp)/(nu(i_mp_qc)+1.0_rp)/(nu(i_mp_qc)+1.0_rp)
3719  coef_aut0 = -kcc*coef_nuc0
3720  coef_aut1 = -kcc/x_sep/20._rp*coef_nuc1
3721  !
3722  do j = js, je
3723  do i = is, ie
3724  do k = ks, ke
3725  lwc = rhoq(i_qr,k,i,j)+rhoq(i_qc,k,i,j)
3726  if( lwc > xc_min )then
3727  tau = max(tau_min, rhoq(i_qr,k,i,j)/lwc)
3728  else
3729  tau = tau_min
3730  end if
3731  rho_fac = sqrt(rho_0/max(rho(k,i,j),rho_min))
3732  !
3733  ! Auto-conversion ( cloud-cloud => rain )
3734  psi_aut = 400._rp*(tau**0.7_rp)*(1.0_rp - (tau**0.7_rp))**3 ! (6) SB06
3735  pq(i_ncaut,k,i,j) = coef_aut0*rhoq(i_qc,k,i,j)*rhoq(i_qc,k,i,j)*rho_fac*rho_fac ! (9) SB06 sc+aut
3736  ! lc = lwc*(1-tau), lr = lwc*tau
3737  pq(i_lcaut,k,i,j) = coef_aut1*lwc*lwc*xq(i_mp_qc,k,i,j)*xq(i_mp_qc,k,i,j) & ! (4) SB06
3738  *((1.0_rp-tau)*(1.0_rp-tau) + psi_aut)*rho_fac*rho_fac !
3739  pq(i_nraut,k,i,j) = -rx_sep*pq(i_lcaut,k,i,j) ! (A7) SB01
3740  !
3741  ! Accretion ( cloud-rain => rain )
3742  psi_acc =(tau/(tau+thr_acc))**4 ! (8) SB06
3743  pq(i_lcacc,k,i,j) = -kcr*rhoq(i_qc,k,i,j)*rhoq(i_qr,k,i,j)*rho_fac*psi_acc ! (7) SB06
3744  pq(i_ncacc,k,i,j) = -kcr*rhoq(i_nc,k,i,j)*rhoq(i_qr,k,i,j)*rho_fac*psi_acc ! (A6) SB01
3745  !
3746  ! Self-collection ( rain-rain => rain )
3747  pq(i_nrslc,k,i,j) = -krr*rhoq(i_nr,k,i,j)*rhoq(i_qr,k,i,j)*rho_fac ! (A.8) SB01
3748  !
3749  ! Collisional breakup of rain
3750  ddr = min(1.e-3_rp, dq_xave(i_mp_qr,k,i,j) - dr_eq )
3751  if ( dq_xave(i_mp_qr,k,i,j) < dr_min ) then ! negligible
3752  psi_brk = -1.0_rp
3753  else if ( dq_xave(i_mp_qr,k,i,j) <= dr_eq ) then
3754  psi_brk = kbr*ddr ! (14) SB06
3755  else
3756  psi_brk = exp(kapbr*ddr) - 1.0_rp ! (15) SB06 (SB06 has a typo)
3757  end if
3758  pq(i_nrbrk,k,i,j) = - (psi_brk + 1.0_rp)*pq(i_nrslc,k,i,j) ! (13) SB06
3759  !
3760  end do
3761  end do
3762  end do
3763  !
3764  return
3765  end subroutine aut_acc_slc_brk_kij
3766  ! Vapor Deposition, Ice Melting
3767  subroutine dep_vapor_melt_ice_kij( &
3768  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
3769  rho, tem, pre, qd, & ! in
3770  rhoq, & ! in
3771  esw, esi, & ! in
3772  xq, vt_xave, dq_xave, & ! in
3773  PQ ) ! inout
3774  use scale_const, only: &
3775  eps => const_eps
3776  implicit none
3777 
3778  integer, intent(in) :: KA, KS, KE
3779  integer, intent(in) :: IA, IS, IE
3780  integer, intent(in) :: JA, JS, JE
3781 
3782  ! Diffusion growth or Evaporation, Sublimation
3783  real(RP), intent(inout) :: PQ(pq_max,ka,ia,ja) ! mass change for cloud, [Add] 09/08/18 T.Mitsui
3784 
3785  real(RP), intent(in) :: rho(ka,ia,ja) ! air density
3786  real(RP), intent(in) :: tem(ka,ia,ja) ! air temperature
3787  real(RP), intent(in) :: pre(ka,ia,ja) ! air pressure
3788  real(RP), intent(in) :: qd(ka,ia,ja) ! mixing ratio of dry air
3789  real(RP), intent(in) :: esw(ka,ia,ja) ! saturation vapor pressure(liquid water)
3790  real(RP), intent(in) :: esi(ka,ia,ja) ! saturation vapor pressure(solid water)
3791  real(RP), intent(in) :: rhoq(i_qv:i_ng,ka,ia,ja)
3792  real(RP), intent(in) :: xq(hydro_max,ka,ia,ja) ! mean mass
3793  ! Notice following values differ from mean terminal velocity or diameter.
3794  ! mean(vt(x)) /= vt(mean(x)) and mean(D(x)) /= D(mean(x))
3795  ! Following ones are vt(mean(x)) and D(mean(x)).
3796  real(RP), intent(in) :: vt_xave(hydro_max,2,ka,ia,ja) ! terminal velocity of mean cloud 09/08/18 [Add], T.Mitsui
3797  !
3798  real(RP), intent(in) :: dq_xave(hydro_max,ka,ia,ja) ! diameter
3799  !
3800  real(RP) :: rho_lim ! limited density 09/08/18 T.Mitsui
3801  real(RP) :: temc_lim ! limited temperature[celsius] 09/08/18 T.Mitsui
3802  real(RP) :: pre_lim ! limited density 09/08/18 T.Mitsui
3803  real(RP) :: temc ! temperature[celsius]
3804 ! real(RP) :: pv ! vapor pressure
3805  real(RP) :: qv ! mixing ratio of water vapor [Add] 09/08/18
3806 ! real(RP) :: ssw ! super saturation ratio(liquid water)
3807 ! real(RP) :: ssi ! super saturation ratio(ice water)
3808  real(RP) :: nua, r_nua ! kinematic viscosity of air
3809  real(RP) :: mua ! viscosity of air
3810  real(RP) :: Kalfa ! thermal conductance
3811  real(RP) :: Dw ! diffusivity of water vapor
3812  real(RP) :: Dt ! diffusivity of heat
3813  real(RP) :: Gw, Gi ! diffusion factor by balance between heat and vapor
3814  real(RP) :: Gwr, Gii, Gis, Gig ! for rain, ice, snow and graupel.
3815  real(RP) :: Gm ! melting factor by balance between heat and vapor
3816  real(RP) :: Nsc_r3 !
3817  ! [Mod] 11/08/30 T.Mitsui, considering large and small branches
3818 ! real(RP) :: Nrecs_r2 ! 09/08/18 [Add] T.Mitsui
3819  real(RP) :: Nrers_r2, Nreis_r2 !
3820  real(RP) :: Nress_r2, Nregs_r2 !
3821 ! real(RP) :: Nrecl_r2 ! 09/08/18 [Add] T.Mitsui
3822  real(RP) :: Nrerl_r2, Nreil_r2 !
3823  real(RP) :: Nresl_r2, Nregl_r2 !
3824  real(RP) :: NscNrer_s, NscNrer_l
3825  real(RP) :: NscNrei_s, NscNrei_l
3826  real(RP) :: NscNres_s, NscNres_l
3827  real(RP) :: NscNreg_s, NscNreg_l
3828  real(RP) :: ventLR_s, ventLR_l
3829  real(RP) :: ventNI_s, ventNI_l, ventLI_s, ventLI_l
3830  real(RP) :: ventNS_s, ventNS_l, ventLS_s, ventLS_l
3831  real(RP) :: ventNG_s, ventNG_l, ventLG_s, ventLG_l
3832  !
3833  real(RP) :: wtr, wti, wts, wtg
3834  real(RP), parameter :: r_14=1.0_rp/1.4_rp
3835  real(RP), parameter :: r_15=1.0_rp/1.5_rp
3836  !
3837  real(RP) :: ventLR
3838  real(RP) :: ventNI, ventLI
3839  real(RP) :: ventNS, ventLS
3840  real(RP) :: ventNG, ventLG
3841  !
3842  real(RP), parameter :: Re_max=1.e+3_rp
3843  real(RP), parameter :: Re_min=1.e-4_rp
3844 
3845  real(RP) :: sw
3846  !
3847  integer :: i, j, k
3848  !
3849  ! Notice,T.Mitsui
3850  ! Vapor deposition and melting would not be solved iteratively to reach equilibrium.
3851  ! Because following phenomena are not adjustment but transition.
3852  ! Just time-scales differ among them.
3853  ! If we would treat more appropreately, there would be time-splitting method to solve each ones.
3854 
3855  profile_start("sn14_dep_vapor")
3856  do j = js, je
3857  do i = is, ie
3858  do k = ks, ke
3859  temc = tem(k,i,j) - t00 ! degC
3860  temc_lim= max(temc, -40._rp ) ! [Add] 09/08/18 T.Mitsui, Pruppacher and Klett(1997),(13-3)
3861  rho_lim = max(rho(k,i,j),rho_min) ! [Add] 09/08/18 T.Mitsui
3862  qv = rhoq(i_qv,k,i,j)/rho_lim
3863  pre_lim = rho_lim*(qd(k,i,j)*rdry + qv*rvap)*(temc_lim+t00) ![Add] 09/08/18 T.Mitsui
3864  !--------------------------------------------------------------------
3865  ! Diffusion growth part is described in detail
3866  ! by Pruppacher and Klett (1997) Sec. 13.2(liquid) and 13.3(solid)
3867  !
3868  ! G:factor of thermal diffusion(1st.term) and vapor diffusion(2nd. term)
3869  ! SB06(23),(38), Lin et al(31),(52) or others
3870  ! Dw is introduced by Pruppacher and Klett(1997),(13-3)
3871  dw = 0.211e-4_rp* (((temc_lim+t00)/t00)**1.94_rp) *(p00/pre_lim)
3872  kalfa = ka0 + temc_lim*dka_dt
3873  mua = mua0 + temc_lim*dmua_dt
3874  nua = mua/rho_lim
3875  r_nua = 1.0_rp/nua
3876  gw = (lhv0/kalfa/tem(k,i,j))*(lhv0/rvap/tem(k,i,j)-1.0_rp)+(rvap*tem(k,i,j)/dw/esw(k,i,j))
3877  gi = (lhs0/kalfa/tem(k,i,j))*(lhs0/rvap/tem(k,i,j)-1.0_rp)+(rvap*tem(k,i,j)/dw/esi(k,i,j))
3878  ! capacities account for their surface geometries
3879  gwr = 4.0_rp*pi/cap(i_mp_qr)/gw
3880  gii = 4.0_rp*pi/cap(i_mp_qi)/gi
3881  gis = 4.0_rp*pi/cap(i_mp_qs)/gi
3882  gig = 4.0_rp*pi/cap(i_mp_qg)/gi
3883  ! vent: ventilation effect( asymmetry vapor field around particles due to aerodynamic )
3884  ! SB06 (30),(31) and each coefficient is by (88),(89)
3885  nsc_r3 = (nua/dw)**(0.33333333_rp) ! (Schmidt number )^(1/3)
3886  !
3887 ! Nrecs_r2 = sqrt(max(Re_min,min(Re_max,vt_xave(I_mp_QC,1,k,i,j)*dq_xave(I_mp_QC,k,i,j)*r_nua))) ! (Reynolds number)^(1/2) cloud
3888  nrers_r2 = sqrt(max(re_min,min(re_max,vt_xave(i_mp_qr,1,k,i,j)*dq_xave(i_mp_qr,k,i,j)*r_nua))) ! (Reynolds number)^(1/2) rain
3889  nreis_r2 = sqrt(max(re_min,min(re_max,vt_xave(i_mp_qi,1,k,i,j)*dq_xave(i_mp_qi,k,i,j)*r_nua))) ! (Reynolds number)^(1/2) cloud ice
3890  nress_r2 = sqrt(max(re_min,min(re_max,vt_xave(i_mp_qs,1,k,i,j)*dq_xave(i_mp_qs,k,i,j)*r_nua))) ! (Reynolds number)^(1/2) snow
3891  nregs_r2 = sqrt(max(re_min,min(re_max,vt_xave(i_mp_qg,1,k,i,j)*dq_xave(i_mp_qg,k,i,j)*r_nua))) ! (Reynolds number)^(1/2) graupel
3892  !
3893 ! Nrecl_r2 = sqrt(max(Re_min,min(Re_max,vt_xave(I_mp_QC,2,k,i,j)*dq_xave(I_mp_QC,k,i,j)*r_nua))) ! (Reynolds number)^(1/2) cloud
3894  nrerl_r2 = sqrt(max(re_min,min(re_max,vt_xave(i_mp_qr,2,k,i,j)*dq_xave(i_mp_qr,k,i,j)*r_nua))) ! (Reynolds number)^(1/2) rain
3895  nreil_r2 = sqrt(max(re_min,min(re_max,vt_xave(i_mp_qi,2,k,i,j)*dq_xave(i_mp_qi,k,i,j)*r_nua))) ! (Reynolds number)^(1/2) cloud ice
3896  nresl_r2 = sqrt(max(re_min,min(re_max,vt_xave(i_mp_qs,2,k,i,j)*dq_xave(i_mp_qs,k,i,j)*r_nua))) ! (Reynolds number)^(1/2) snow
3897  nregl_r2 = sqrt(max(re_min,min(re_max,vt_xave(i_mp_qg,2,k,i,j)*dq_xave(i_mp_qg,k,i,j)*r_nua))) ! (Reynolds number)^(1/2) graupel
3898  nscnrer_s=nsc_r3*nrers_r2 ! small rain
3899  nscnrer_l=nsc_r3*nrerl_r2 ! large rain
3900  !
3901  nscnrei_s=nsc_r3*nreis_r2 ! small ice
3902  nscnrei_l=nsc_r3*nreil_r2 ! large ice
3903  !
3904  nscnres_s=nsc_r3*nress_r2 ! small snow
3905  nscnres_l=nsc_r3*nresl_r2 ! large snow
3906  !
3907  nscnreg_s=nsc_r3*nregs_r2 ! small snow
3908  nscnreg_l=nsc_r3*nregl_r2 ! large snow
3909  !
3910  ventlr_s = ah_vent1(i_mp_qr,1) + bh_vent1(i_mp_qr,1)*nscnrer_s
3911  ventlr_l = ah_vent1(i_mp_qr,2) + bh_vent1(i_mp_qr,2)*nscnrer_l
3912  !
3913  ventni_s = ah_vent0(i_mp_qi,1) + bh_vent0(i_mp_qi,1)*nscnrei_s
3914  ventni_l = ah_vent0(i_mp_qi,2) + bh_vent0(i_mp_qi,2)*nscnrei_l
3915  ventli_s = ah_vent1(i_mp_qi,1) + bh_vent1(i_mp_qi,1)*nscnrei_s
3916  ventli_l = ah_vent1(i_mp_qi,2) + bh_vent1(i_mp_qi,2)*nscnrei_l
3917  !
3918  ventns_s = ah_vent0(i_mp_qs,1) + bh_vent0(i_mp_qs,1)*nscnres_s
3919  ventns_l = ah_vent0(i_mp_qs,2) + bh_vent0(i_mp_qs,2)*nscnres_l
3920  ventls_s = ah_vent1(i_mp_qs,1) + bh_vent1(i_mp_qs,1)*nscnres_s
3921  ventls_l = ah_vent1(i_mp_qs,2) + bh_vent1(i_mp_qs,2)*nscnres_l
3922  !
3923  ventng_s = ah_vent0(i_mp_qg,1) + bh_vent0(i_mp_qg,1)*nscnreg_s
3924  ventng_l = ah_vent0(i_mp_qg,2) + bh_vent0(i_mp_qg,2)*nscnreg_l
3925  ventlg_s = ah_vent1(i_mp_qg,1) + bh_vent1(i_mp_qg,1)*nscnreg_s
3926  ventlg_l = ah_vent1(i_mp_qg,2) + bh_vent1(i_mp_qg,2)*nscnreg_l
3927  !
3928  ! branch is 1.4 for rain, snow, graupel; is 1.0 for ice (PK97, 13-60,-61,-88,-89).
3929  !
3930  wtr = ( min(max( nscnrer_s*r_14, 0.5_rp), 2.0_rp) -0.5_rp )*r_15 ! weighting between 1.4*0.5 and 1.4*2
3931  wti = ( min(max( nscnrei_s , 0.5_rp), 2.0_rp) -0.5_rp )*r_15 ! weighting between 1.0*0.5 and 1.0*2
3932  wts = ( min(max( nscnres_s*r_14, 0.5_rp), 2.0_rp) -0.5_rp )*r_15 ! weighting between 1.4*0.5 and 1.4*2
3933  wtg = ( min(max( nscnreg_s*r_14, 0.5_rp), 2.0_rp) -0.5_rp )*r_15 ! weighting between 1.4*0.5 and 1.4*2
3934  ! interpolation between two branches
3935  ventni = (1.0_rp-wti)*ventni_s + wti*ventni_l
3936  ventns = (1.0_rp-wts)*ventns_s + wts*ventns_l
3937  ventng = (1.0_rp-wtg)*ventng_s + wtg*ventng_l
3938  !
3939  ventlr = (1.0_rp-wtr)*ventlr_s + wtr*ventlr_l
3940  ventli = (1.0_rp-wti)*ventli_s + wti*ventli_l
3941  ventls = (1.0_rp-wts)*ventls_s + wts*ventls_l
3942  ventlg = (1.0_rp-wtg)*ventlg_s + wtg*ventlg_l
3943  !
3944  ! SB06(29)
3945  ! [Mod] 08/05/08 T.Mitsui, recover PNXdep, and rain is only evaporation.
3946  ! Ni, Ns, Ng should decrease in nature so we add this term.
3947  ! And vapor deposition never occur unless number exist.
3948  ! [Add comment] 09/08/18
3949  ! recover condensation/evaporation of rain,
3950  ! and ventilation effects are not taken into account for cloud.
3951  !
3952 !!$***************************************************************************
3953 !!$ NOTICE:
3954 !!$ 09/08/18 [Mod] Hereafter PLxdep means inverse of timescale.
3955 !!$***************************************************************************
3956 !!$ PQ(I_LCdep,k,i,j) = Gwr*ssw*rhoq(I_NC,k,i,j)*dq_xave(I_mp_QC,k,i,j)*coef_deplc
3957 !!$ PQ(I_LRdep,k,i,j) = Gwr*ssw*rhoq(I_NR,k,i,j)*dq_xave(I_mp_QR,k,i,j)*ventLR
3958 !!$ PQ(I_LIdep,k,i,j) = Gii*ssi*rhoq(I_NI,k,i,j)*dq_xave(I_mp_QI,k,i,j)*ventLI
3959 !!$ PQ(I_LSdep,k,i,j) = Gis*ssi*rhoq(I_NS,k,i,j)*dq_xave(I_mp_QS,k,i,j)*ventLS
3960 !!$ PQ(I_LGdep,k,i,j) = Gig*ssi*rhoq(I_NG,k,i,j)*dq_xave(I_mp_QG,k,i,j)*ventLG
3961  pq(i_lcdep,k,i,j) = gwr*rhoq(i_nc,k,i,j)*dq_xave(i_mp_qc,k,i,j)*coef_deplc
3962  pq(i_lrdep,k,i,j) = gwr*rhoq(i_nr,k,i,j)*dq_xave(i_mp_qr,k,i,j)*ventlr
3963  pq(i_lidep,k,i,j) = gii*rhoq(i_ni,k,i,j)*dq_xave(i_mp_qi,k,i,j)*ventli
3964  pq(i_lsdep,k,i,j) = gis*rhoq(i_ns,k,i,j)*dq_xave(i_mp_qs,k,i,j)*ventls
3965  pq(i_lgdep,k,i,j) = gig*rhoq(i_ng,k,i,j)*dq_xave(i_mp_qg,k,i,j)*ventlg
3966  pq(i_nrdep,k,i,j) = pq(i_lrdep,k,i,j)/xq(i_mp_qr,k,i,j)
3967  pq(i_nidep,k,i,j) = 0.0_rp
3968  pq(i_nsdep,k,i,j) = pq(i_lsdep,k,i,j)/xq(i_mp_qs,k,i,j)
3969  pq(i_ngdep,k,i,j) = pq(i_lgdep,k,i,j)/xq(i_mp_qg,k,i,j)
3970  !
3971  !------------------------------------------------------------------------
3972  ! Melting part is described by Pruppacher and Klett (1997) Sec.16.3.1
3973  ! Here we omit "Shedding" of snow-flakes and ice-particles.
3974  ! "Shedding" may be applicative if you refer
3975  ! eq.(38) in Cotton etal.(1986) Jour. Clim. Appl. Meteor. p.1658-1680.
3976  ! SB06(73)
3977  dt = kalfa/(cpvap*rho_0)
3978  ! Gm: factor caused by balance between
3979  ! "water evaporation cooling(1st.)" and "fusion heating(2nd.)"
3980  ! SB06(76)
3981  ! [fix] 08/05/08 T.Mitsui LHF00 => EMELT and esw => PSAT0
3982  ! LHS0 is more suitable than LHS because melting occurs around 273.15 K.
3983  gm = 2.0_rp*pi/emelt&
3984  * ( (kalfa*dt/dw)*(temc) + (dw*lhs0/rvap)*(esi(k,i,j)/tem(k,i,j)-psat0/t00) )
3985  ! SB06(76)
3986  ! Notice! melting only occurs where T > 273.15 K else doesn't.
3987  ! [fix] 08/05/08 T.Mitsui, Gm could be both positive and negative value.
3988  ! See Pruppacher and Klett(1997) eq.(16-79) or Rasmussen and Pruppacher(1982)
3989  sw = ( sign(0.5_rp,temc) + 0.5_rp ) * ( sign(0.5_rp,gm-eps) + 0.5_rp ) ! sw = 1 if( (temc>=0.0_RP) .AND. (Gm>0.0_RP) ), otherwise sw = 0
3990  ! if Gm==0 then rh and tem is critical value for melting process.
3991  ! 08/05/16 [Mod] T.Mitsui, change term of PLimlt. N_i => L_i/ (limited x_i)
3992  ! because melting never occur when N_i=0.
3993  pq(i_limlt,k,i,j) = - gm * rhoq(i_qi,k,i,j)*dq_xave(i_mp_qi,k,i,j)*ventli/xq(i_mp_qi,k,i,j) * sw
3994  ! [Mod] 08/08/23 T.Mitsui for Seifert(2008)
3995  pq(i_nimlt,k,i,j) = - gm * rhoq(i_ni,k,i,j)*dq_xave(i_mp_qi,k,i,j)*ventni/xq(i_mp_qi,k,i,j) * sw ! 09/08/18 [Mod] recover, T.Mitsui
3996  pq(i_lsmlt,k,i,j) = - gm * rhoq(i_qs,k,i,j)*dq_xave(i_mp_qs,k,i,j)*ventls/xq(i_mp_qs,k,i,j) * sw
3997  ! [Mod] 08/08/23 T.Mitsui for Seifert(2008)
3998  pq(i_nsmlt,k,i,j) = - gm * rhoq(i_ns,k,i,j)*dq_xave(i_mp_qs,k,i,j)*ventns/xq(i_mp_qs,k,i,j) * sw ! 09/08/18 [Mod] recover, T.Mitsui
3999  pq(i_lgmlt,k,i,j) = - gm * rhoq(i_qg,k,i,j)*dq_xave(i_mp_qg,k,i,j)*ventlg/xq(i_mp_qg,k,i,j) * sw
4000  ! [Mod] 08/08/23 T.Mitsui for Seifert(2008)
4001  pq(i_ngmlt,k,i,j) = - gm * rhoq(i_ng,k,i,j)*dq_xave(i_mp_qg,k,i,j)*ventng/xq(i_mp_qg,k,i,j) * sw ! 09/08/18 [Mod] recover, T.Mitsui
4002 
4003  end do
4004  end do
4005  end do
4006  profile_stop("sn14_dep_vapor")
4007  !
4008  return
4009  end subroutine dep_vapor_melt_ice_kij
4010  !-----------------------------------------------------------------------------
4011  subroutine freezing_water_kij( &
4012  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
4013  dt, &
4014  rhoq, xq, tem, &
4015  PQ )
4016  implicit none
4017  !
4018  ! In this subroutine,
4019  ! We assumed surface temperature of droplets are same as environment.
4020 
4021  integer, intent(in) :: KA, KS, KE
4022  integer, intent(in) :: IA, IS, IE
4023  integer, intent(in) :: JA, JS, JE
4024 
4025  real(RP), intent(in) :: dt
4026  !
4027  real(RP), intent(in) :: tem(ka,ia,ja)
4028  !
4029  real(RP), intent(in) :: rhoq(i_qv:i_ng,ka,ia,ja)
4030  real(RP), intent(in) :: xq(hydro_max,ka,ia,ja)
4031  !
4032  real(RP), intent(inout):: PQ(pq_max,ka,ia,ja)
4033  !
4034  real(RP), parameter :: temc_min = -65.0_rp
4035  real(RP), parameter :: a_het = 0.2_rp ! SB06 (44)
4036  real(RP), parameter :: b_het = 0.65_rp ! SB06 (44)
4037  !
4038  real(RP) :: coef_m2_c
4039  real(RP) :: coef_m2_r
4040  ! temperature [celsius]
4041  real(RP) :: temc, temc2, temc3, temc4
4042  ! temperature function of homegenous/heterogenous freezing
4043  real(RP) :: Jhom, Jhet
4044  real(RP) :: rdt
4045  real(RP) :: tmp
4046  !
4047  integer :: i,j,k
4048  !
4049  rdt = 1.0_rp/dt
4050  !
4051  coef_m2_c = coef_m2(i_mp_qc)
4052  coef_m2_r = coef_m2(i_mp_qr)
4053  !
4054  profile_start("sn14_freezing")
4055  do j = js, je
4056  do i = is, ie
4057  do k = ks, ke
4058  temc = max( tem(k,i,j) - t00, temc_min )
4059  ! These cause from aerosol-droplet interaction.
4060  ! Bigg(1953) formula, Khain etal.(2000) eq.(4.5), Pruppacher and Klett(1997) eq.(9-48)
4061  jhet = a_het*exp( -b_het*temc - 1.0_rp )
4062  ! These cause in nature.
4063  ! Cotton and Field 2002, QJRMS. (12)
4064  if( temc < -65.0_rp )then
4065  jhom = 10.0_rp**(24.37236_rp)*1.e+3_rp
4066  jhet = a_het*exp( 65.0_rp*b_het - 1.0_rp ) ! 09/04/14 [Add], fixer T.Mitsui
4067  else if( temc < -30.0_rp ) then
4068  temc2 = temc*temc
4069  temc3 = temc*temc2
4070  temc4 = temc2*temc2
4071  jhom = 10.0_rp**(&
4072  - 243.40_rp - 14.75_rp*temc - 0.307_rp*temc2 &
4073  - 0.00287_rp*temc3 - 0.0000102_rp*temc4 ) *1.e+3_rp
4074  else if( temc < 0.0_rp) then
4075  jhom = 10._rp**(-7.63_rp-2.996_rp*(temc+30.0_rp))*1.e+3_rp
4076  else
4077  jhom = 0.0_rp
4078  jhet = 0.0_rp
4079  end if
4080  ! Note, xc should be limited in range[xc_min:xc_max].
4081  ! and PNChom need to be calculated by NC
4082  ! because reduction rate of Nc need to be bound by NC.
4083  ! For the same reason PLChom also bound by LC and xc.
4084  ! Basically L and N should be independent
4085  ! but average particle mass x should be in suitable range.
4086  ! Homogenous Freezing
4087  pq(i_lchom,k,i,j) = 0.0_rp
4088  pq(i_nchom,k,i,j) = 0.0_rp
4089  ! Heterogenous Freezing
4090 #if defined(PGI) || defined(SX)
4091  tmp = min( xq(i_mp_qc,k,i,j)*(jhet+jhom)*dt, 1.e+3_rp) ! apply exp limiter
4092  pq(i_lchet,k,i,j) = -rdt*rhoq(i_qc,k,i,j)*( 1.0_rp - exp( -coef_m2_c*tmp ) )
4093  pq(i_nchet,k,i,j) = -rdt*rhoq(i_nc,k,i,j)*( 1.0_rp - exp( - tmp ) )
4094 
4095  tmp = min( xq(i_mp_qr,k,i,j)*(jhet+jhom)*dt, 1.e+3_rp) ! apply exp limiter
4096  pq(i_lrhet,k,i,j) = -rdt*rhoq(i_qr,k,i,j)*( 1.0_rp - exp( -coef_m2_r*tmp ) )
4097  pq(i_nrhet,k,i,j) = -rdt*rhoq(i_nr,k,i,j)*( 1.0_rp - exp( - tmp ) )
4098 #else
4099  pq(i_lchet,k,i,j) = -rdt*rhoq(i_qc,k,i,j)*( 1.0_rp - exp( -coef_m2_c*xq(i_mp_qc,k,i,j)*(jhet+jhom)*dt ) )
4100  pq(i_nchet,k,i,j) = -rdt*rhoq(i_nc,k,i,j)*( 1.0_rp - exp( - xq(i_mp_qc,k,i,j)*(jhet+jhom)*dt ) )
4101  pq(i_lrhet,k,i,j) = -rdt*rhoq(i_qr,k,i,j)*( 1.0_rp - exp( -coef_m2_r*xq(i_mp_qr,k,i,j)*(jhet+jhom)*dt ) )
4102  pq(i_nrhet,k,i,j) = -rdt*rhoq(i_nr,k,i,j)*( 1.0_rp - exp( - xq(i_mp_qr,k,i,j)*(jhet+jhom)*dt ) )
4103 #endif
4104  end do
4105  end do
4106  end do
4107  profile_stop("sn14_freezing")
4108  !
4109  return
4110  end subroutine freezing_water_kij
4111 
4112  !-----------------------------------------------------------------------------
4116 !OCL SERIAL
4118  KA, KS, KE, &
4119  DENS, &
4120  TEMP, &
4121  RHOQ, &
4122  PRES, &
4123  vterm )
4124  use scale_const, only: &
4125  const_undef
4126  implicit none
4127 
4128  integer, intent(in) :: KA, KS, KE
4129 
4130  real(RP), intent(in) :: RHOQ(ka,i_qc:i_ng) ! rho * q
4131  real(RP), intent(in) :: DENS(ka) ! rho
4132  real(RP), intent(in) :: TEMP(ka) ! temperature
4133  real(RP), intent(in) :: PRES(ka) ! pressure
4134 
4135  real(RP), intent(out) :: vterm(ka,qa_mp-1) ! terminal velocity of cloud mass
4136 
4137  real(RP) :: xq ! average mass of 1 particle( mass/number )
4138 
4139  real(RP) :: rhofac ! density factor for terminal velocity( air friction )
4140  real(RP) :: rhofac_q(hydro_max)
4141 
4142  real(RP) :: rlambdar ! work for diagnosis of Rain DSD ( Seifert, 2008 )
4143  real(RP) :: mud_r
4144  real(RP) :: dq, dql ! weigthed diameter. Improved Rogers etal. (1993) formula by T.Mitsui
4145 
4146 
4147  real(RP) :: weight ! weighting coefficient for 2-branches is determined by ratio between 0.745mm and weighted diameter. SB06 Table.1
4148  real(RP) :: velq_s ! terminal velocity for small branch of Rogers formula
4149  real(RP) :: velq_l ! terminal velocity for large branch of Rogers formula
4150 
4151  real(RP) :: tmp
4152  integer :: k, i, j, iq
4153  !---------------------------------------------------------------------------
4154 
4155  profile_start("sn14_terminal_vel")
4156 
4157  mud_r = 3.0_rp * nu(i_mp_qr) + 2.0_rp
4158 
4159 
4160  do k = ks, ke
4161  rhofac = rho_0 / max( dens(k), rho_min )
4162 
4163  ! QC
4164  rhofac_q(i_mp_qc) = rhofac ** gamma_v(i_mp_qc)
4165  xq = max( xqmin(i_mp_qc), min( xqmax(i_mp_qc), rhoq(k,i_qc) / ( rhoq(k,i_nc) + nqmin(i_mp_qc) ) ) )
4166 
4167  vterm(k,i_mp_qc) = -rhofac_q(i_mp_qc) * coef_vt1(i_mp_qc,1) * xq**beta_v(i_mp_qc,1)
4168  ! NC
4169  vterm(k,i_mp_nc) = -rhofac_q(i_mp_qc) * coef_vt0(i_mp_qc,1) * xq**beta_vn(i_mp_qc,1)
4170 
4171  ! QR
4172  rhofac_q(i_mp_qr) = rhofac ** gamma_v(i_mp_qr)
4173  xq = max( xqmin(i_mp_qr), min( xqmax(i_mp_qr), rhoq(k,i_qr) / ( rhoq(k,i_nr) + nqmin(i_mp_qr) ) ) )
4174 
4175  rlambdar = a_m(i_mp_qr) * xq**b_m(i_mp_qr) &
4176  * ( (mud_r+3.0_rp) * (mud_r+2.0_rp) * (mud_r+1.0_rp) )**(-0.333333333_rp)
4177  dq = ( 4.0_rp + mud_r ) * rlambdar ! D^(3)+mu weighted mean diameter
4178  dql = dq
4179  weight = min( 1.0_rp, max( 0.0_rp, 0.5_rp * ( 1.0_rp + tanh( pi * log( dq/d_vtr_branch ) ) ) ) )
4180 
4181  velq_s = coef_vtr_ar2 * dq &
4182  * ( 1.0_rp - ( 1.0_rp + coef_vtr_br2*rlambdar )**(-5.0_rp-mud_r) )
4183  velq_l = coef_vtr_ar1 - coef_vtr_br1 &
4184  * ( 1.0_rp + coef_vtr_cr1*rlambdar )**(-4.0_rp-mud_r)
4185  vterm(k,i_mp_qr) = -rhofac_q(i_mp_qr) &
4186  * ( velq_l * ( weight ) &
4187  + velq_s * ( 1.0_rp - weight ) )
4188  ! NR
4189  dq = ( 1.0_rp + mud_r ) * rlambdar
4190  weight = min( 1.0_rp, max( 0.0_rp, 0.5_rp * ( 1.0_rp + tanh( pi * log( dq/d_vtr_branch ) ) ) ) )
4191 
4192  velq_s = coef_vtr_ar2 * dql &
4193  * ( 1.0_rp - ( 1.0_rp + coef_vtr_br2*rlambdar )**(-2.0_rp-mud_r) )
4194  velq_l = coef_vtr_ar1 - coef_vtr_br1 &
4195  * ( 1.0_rp + coef_vtr_cr1*rlambdar )**(-1.0_rp-mud_r)
4196  vterm(k,i_mp_nr) = -rhofac_q(i_mp_qr) &
4197  * ( velq_l * ( weight ) &
4198  + velq_s * ( 1.0_rp - weight ) )
4199 
4200  ! QI
4201  rhofac_q(i_mp_qi) = ( pres(k)/pre0_vt )**a_pre0_vt * ( temp(k)/tem0_vt )**a_tem0_vt
4202  xq = max( xqmin(i_mp_qi), min( xqmax(i_mp_qi), rhoq(k,i_qi) / ( rhoq(k,i_ni) + nqmin(i_mp_qi) ) ) )
4203 
4204  tmp = a_m(i_mp_qi) * xq**b_m(i_mp_qi)
4205  dq = coef_dave_l(i_mp_qi) * tmp
4206  weight = min( 1.0_rp, max( 0.0_rp, 0.5_rp * ( 1.0_rp + log( dq/d0_li ) ) ) )
4207 
4208  velq_s = coef_vt1(i_mp_qi,1) * xq**beta_v(i_mp_qi,1)
4209  velq_l = coef_vt1(i_mp_qi,2) * xq**beta_v(i_mp_qi,2)
4210  vterm(k,i_mp_qi) = -rhofac_q(i_mp_qi) &
4211  * ( velq_l * ( weight ) &
4212  + velq_s * ( 1.0_rp - weight ) )
4213  ! NI
4214  dq = coef_dave_n(i_mp_qi) * tmp
4215  weight = min( 1.0_rp, max( 0.0_rp, 0.5_rp * ( 1.0_rp + log( dq/d0_ni ) ) ) )
4216 
4217  velq_s = coef_vt0(i_mp_qi,1) * xq**beta_vn(i_mp_qi,1)
4218  velq_l = coef_vt0(i_mp_qi,2) * xq**beta_vn(i_mp_qi,2)
4219  vterm(k,i_mp_ni) = -rhofac_q(i_mp_qi) &
4220  * ( velq_l * ( weight ) &
4221  + velq_s * ( 1.0_rp - weight ) )
4222 
4223  ! QS
4224  rhofac_q(i_mp_qs) = rhofac_q(i_mp_qi)
4225  xq = max( xqmin(i_mp_qs), min( xqmax(i_mp_qs), rhoq(k,i_qs) / ( rhoq(k,i_ns) + nqmin(i_mp_qs) ) ) )
4226 
4227  tmp = a_m(i_mp_qs) * xq**b_m(i_mp_qs)
4228  dq = coef_dave_l(i_mp_qs) * tmp
4229  weight = min( 1.0_rp, max( 0.0_rp, 0.5_rp * ( 1.0_rp + log( dq/d0_ls ) ) ) )
4230 
4231  velq_s = coef_vt1(i_mp_qs,1) * xq**beta_v(i_mp_qs,1)
4232  velq_l = coef_vt1(i_mp_qs,2) * xq**beta_v(i_mp_qs,2)
4233  vterm(k,i_mp_qs) = -rhofac_q(i_mp_qs) &
4234  * ( velq_l * ( weight ) &
4235  + velq_s * ( 1.0_rp - weight ) )
4236  ! NS
4237  dq = coef_dave_n(i_mp_qs) * tmp
4238  weight = min( 1.0_rp, max( 0.0_rp, 0.5_rp * ( 1.0_rp + log( dq/d0_ns ) ) ) )
4239 
4240  velq_s = coef_vt0(i_mp_qs,1) * xq**beta_vn(i_mp_qs,1)
4241  velq_l = coef_vt0(i_mp_qs,2) * xq**beta_vn(i_mp_qs,2)
4242  vterm(k,i_mp_ns) = -rhofac_q(i_mp_qs) &
4243  * ( velq_l * ( weight ) &
4244  + velq_s * ( 1.0_rp - weight ) )
4245 
4246  ! QG
4247  rhofac_q(i_mp_qg) = rhofac_q(i_mp_qi)
4248  xq = max( xqmin(i_mp_qg), min( xqmax(i_mp_qg), rhoq(k,i_qg) / ( rhoq(k,i_ng) + nqmin(i_mp_qg) ) ) )
4249 
4250  tmp = a_m(i_mp_qg) * xq**b_m(i_mp_qg)
4251  dq = coef_dave_l(i_mp_qg) * tmp
4252  weight = min( 1.0_rp, max( 0.0_rp, 0.5_rp * ( 1.0_rp + log( dq/d0_lg ) ) ) )
4253 
4254  velq_s = coef_vt1(i_mp_qg,1) * xq**beta_v(i_mp_qg,1)
4255  velq_l = coef_vt1(i_mp_qg,2) * xq**beta_v(i_mp_qg,2)
4256  vterm(k,i_mp_qg) = -rhofac_q(i_mp_qg) &
4257  * ( velq_l * ( weight ) &
4258  + velq_s * ( 1.0_rp - weight ) )
4259  ! NG
4260  dq = coef_dave_n(i_mp_qg) * tmp
4261  weight = min( 1.0_rp, max( 0.0_rp, 0.5_rp * ( 1.0_rp + log( dq/d0_ng ) ) ) )
4262 
4263  velq_s = coef_vt0(i_mp_qg,1) * xq**beta_vn(i_mp_qg,1)
4264  velq_l = coef_vt0(i_mp_qg,2) * xq**beta_vn(i_mp_qg,2)
4265  vterm(k,i_mp_ng) = -rhofac_q(i_mp_qg) &
4266  * ( velq_l * ( weight ) &
4267  + velq_s * ( 1.0_rp - weight ) )
4268  enddo
4269 
4270  do iq = 1, qa_mp-1
4271  vterm(1:ks-2,iq) = const_undef
4272  vterm(ks-1,iq) = vterm(ks,iq)
4273  vterm(ke+1:ka,iq) = const_undef
4274  enddo
4275 
4276  profile_stop("sn14_terminal_vel")
4277 
4278  return
4280  !----------------------------------------------------------------
4281  subroutine update_by_phase_change_kij( &
4282  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
4283  ntdiv, ntmax, & ! in [Add] 10/08/03
4284  dt, & ! in
4285  !gsgam2, & ! in
4286  cz, & ! in
4287  fz, & ! in
4288  w, & ! in
4289  dTdt_rad, & ! in
4290  rho, & ! in
4291  qdry, & ! in
4292  esw, esi, rhoq2, & ! in
4293  pre, tem, & ! in
4294  cpa,cva, & ! in
4295  PQ, & ! inout
4296  sl_PLCdep, & ! inout
4297  sl_PLRdep, sl_PNRdep, & ! inout
4298  RHOQ_t, & ! out
4299  RHOE_t, & ! out
4300  CPtot_t, & ! out
4301  CVtot_t, & ! out
4302  qc_evaporate ) ! out
4304  use scale_atmos_hydrometeor, only: &
4305  cp_vapor, &
4306  cp_water, &
4307  cp_ice, &
4308  cv_vapor, &
4309  cv_water, &
4310  cv_ice
4311  use scale_atmos_saturation, only: &
4312  moist_pres2qsat_liq => atmos_saturation_pres2qsat_liq, &
4313  moist_pres2qsat_ice => atmos_saturation_pres2qsat_ice, &
4314  moist_dqs_dtem_dens_liq => atmos_saturation_dqs_dtem_dens_liq, &
4315  moist_dqs_dtem_dens_ice => atmos_saturation_dqs_dtem_dens_ice, &
4316  moist_dqs_dtem_dpre_liq => atmos_saturation_dqs_dtem_dpre_liq, &
4317  moist_dqs_dtem_dpre_ice => atmos_saturation_dqs_dtem_dpre_ice
4318  implicit none
4319 
4320  integer, intent(in) :: KA, KS, KE
4321  integer, intent(in) :: IA, IS, IE
4322  integer, intent(in) :: JA, JS, JE
4323 
4324  integer, intent(in) :: ntdiv ! [Add] 10/08/03
4325  integer, intent(in) :: ntmax ! [Add] 10/08/03
4326  !
4327  real(RP), intent(in) :: dt ! time step[s]
4328  !real(RP), intent(in) :: gsgam2(KA,IA,JA) ! metric
4329  real(RP), intent(in) :: cz(ka,ia,ja) ! altitude [m]
4330  real(RP), intent(in) :: fz(ka,ia,ja) ! altitude difference [m]
4331  real(RP), intent(in) :: w(ka,ia,ja) ! vertical velocity @ full level [m/s]
4332  real(RP), intent(in) :: dTdt_rad(ka,ia,ja) ! temperture tendency by radiation[K/s]
4333  real(RP), intent(in) :: rho(ka,ia,ja) ! density[kg/m3]
4334  real(RP), intent(in) :: qdry(ka,ia,ja) ! dry air mass ratio [kg/kg]
4335  real(RP), intent(in) :: esw(ka,ia,ja) ! saturated vapor pressure for liquid
4336  real(RP), intent(in) :: esi(ka,ia,ja) ! for ice
4337  real(RP), intent(in) :: rhoq2(i_qv:i_ng,ka,ia,ja)
4338 
4339  real(RP), intent(in) :: tem(ka,ia,ja) ! temperature[K]
4340  real(RP), intent(in) :: pre(ka,ia,ja) ! pressure[Pa]
4341  real(RP), intent(in) :: cva(ka,ia,ja) ! specific heat at constant volume
4342  real(RP), intent(in) :: cpa(ka,ia,ja) !
4343 
4344  !+++ tendency[kg/m3/s]
4345  real(RP), intent(inout) :: PQ(pq_max,ka,ia,ja)
4346  !+++ Column integrated tendency[kg/m2/s]
4347  real(RP), intent(inout) :: sl_PLCdep(ia,ja)
4348  real(RP), intent(inout) :: sl_PLRdep(ia,ja), sl_PNRdep(ia,ja)
4349 
4350  real(RP),intent(out) :: RHOQ_t(ka, ia, ja, qa_mp)
4351  real(RP),intent(out) :: RHOE_t(ka, ia, ja)
4352  real(RP),intent(out) :: CPtot_t(ka,ia,ja)
4353  real(RP),intent(out) :: CVtot_t(ka,ia,ja)
4354 
4355  !+++ tendency[kg/m3/s]
4356  real(RP), intent(out) :: qc_evaporate(ka,ia,ja)
4357 
4358  real(RP) :: xi ! mean mass of ice particles
4359  real(RP) :: rrho ! 1/rho
4360  real(RP) :: wtem(ka,ia,ja) ! temperature[K]
4361  !
4362  real(RP) :: r_cva ! specific heat at constant volume
4363  !real(RP) :: cpa ! specific heat at constant pressure
4364  real(RP) :: r_cpa ! specific heat at constant pressure
4365  real(RP) :: qsw(ka,ia,ja) ! saturated mixing ratio for liquid
4366  real(RP) :: qsi(ka,ia,ja) ! saturated mixing ratio for solid
4367  real(RP) :: dqswdtem_rho(ka,ia,ja) ! (dqsw/dtem)_rho
4368  real(RP) :: dqsidtem_rho(ka,ia,ja) ! (dqsi/dtem)_rho
4369  real(RP) :: dqswdtem_pre(ka,ia,ja) ! (dqsw/dtem)_pre
4370  real(RP) :: dqsidtem_pre(ka,ia,ja) ! (dqsi/dtem)_pre
4371  real(RP) :: dqswdpre_tem(ka,ia,ja) ! (dqsw/dpre)_tem
4372  real(RP) :: dqsidpre_tem(ka,ia,ja) ! (dqsi/dpre)_tem
4373  !
4374  real(RP) :: w2 ! vetical velocity[m/s]
4375  real(RP) :: Acnd ! Pdynliq + Bergeron-Findeisen
4376  real(RP) :: Adep ! Pdyndep + Bergeron-Findeisen
4377  real(RP) :: aliqliq, asolliq
4378  real(RP) :: aliqsol, asolsol
4379  real(RP) :: Pdynliq ! production term of ssw by vertical motion
4380  real(RP) :: Pdynsol ! production term of ssi by vertical motion
4381  real(RP) :: Pradliq ! production term of ssw by radiation
4382  real(RP) :: Pradsol ! production term of ssi by radiation
4383  real(RP) :: taucnd, r_taucnd ! time scale of ssw change by MP
4384  real(RP) :: taudep, r_taudep ! time scale of ssi change by MP
4385  real(RP) :: taucnd_c(ka,ia,ja), r_taucnd_c ! by cloud
4386  real(RP) :: taucnd_r(ka,ia,ja), r_taucnd_r ! by rain
4387  real(RP) :: taudep_i(ka,ia,ja), r_taudep_i ! by ice
4388  real(RP) :: taudep_s(ka,ia,ja), r_taudep_s ! by snow
4389  real(RP) :: taudep_g(ka,ia,ja), r_taudep_g ! by graupel
4390  ! alternative tendency through changing ssw and ssi
4391  real(RP) :: PNCdep ! [Add] 11/08/30 T.Mitsui
4392  real(RP) :: PLR2NR, PLI2NI, PLS2NS, PLG2NG
4393  real(RP) :: coef_a_cnd, coef_b_cnd
4394  real(RP) :: coef_a_dep, coef_b_dep
4395  !
4396  real(RP) :: frz_dqc
4397  real(RP) :: frz_dnc
4398  real(RP) :: frz_dqr
4399  real(RP) :: frz_dnr
4400  real(RP) :: mlt_dqi
4401  real(RP) :: mlt_dni
4402  real(RP) :: mlt_dqs
4403  real(RP) :: mlt_dns
4404  real(RP) :: mlt_dqg
4405  real(RP) :: mlt_dng
4406  real(RP) :: dep_dqi
4407  real(RP) :: dep_dni
4408  real(RP) :: dep_dqs
4409  real(RP) :: dep_dns
4410  real(RP) :: dep_dqg
4411  real(RP) :: dep_dng
4412  real(RP) :: dep_dqr
4413  real(RP) :: dep_dnr
4414  real(RP) :: dep_dqc
4415  real(RP) :: dep_dnc ! 11/08/30 [Add] T.Mitsui, dep_dnc
4416  real(RP) :: r_xc_ccn, r_xi_ccn ! 11/08/30 [Add] T.Mitsui
4417  !
4418  real(RP) :: drhoqv
4419  real(RP) :: drhoqc, drhoqr, drhoqi, drhoqs, drhoqg
4420  real(RP) :: drhonc, drhonr, drhoni, drhons, drhong
4421  !
4422  real(RP) :: fac1, fac2, fac3, fac4, fac5, fac6
4423  real(RP) :: r_rvaptem ! 1/(Rvap*tem)
4424  real(RP) :: pv ! vapor pressure
4425  real(RP) :: lvsw, lvsi ! saturated vapor density
4426  real(RP) :: dlvsw, dlvsi !
4427  ! [Add] 11/08/30 T.Mitsui
4428  real(RP) :: dcnd, ddep ! total cndensation/deposition
4429  real(RP) :: uplim_cnd ! upper limit of condensation
4430  real(RP) :: lowlim_cnd ! lower limit of evaporation
4431  ! [Add] 11/08/30 T.Mitsui
4432  real(RP) :: uplim_dep ! upper limit of condensation
4433  real(RP) :: lowlim_dep ! lower limit of evaporation
4434  real(RP) :: ssw, ssi ! supersaturation ratio
4435  real(RP) :: r_esw, r_esi ! 1/esw, 1/esi
4436  real(RP) :: r_lvsw, r_lvsi ! 1/(lvsw*ssw), 1/(lvsi*ssi)
4437  real(RP) :: r_dt ! 1/dt
4438  real(RP) :: ssw_o, ssi_o
4439 ! real(RP) :: dt_dyn
4440 ! real(RP) :: dt_mp
4441  !
4442 ! real(RP) :: tem_lh(KA,IA,JA)
4443 ! real(RP) :: dtemdt_lh(KA,IA,JA)
4444  real(RP), save :: fac_cndc = 1.0_rp
4445  logical, save :: opt_fix_taucnd_c=.false.
4446  logical, save :: flag_first =.true.
4447  !
4448  namelist / param_atmos_phy_mp_sn14_condensation / &
4449  opt_fix_taucnd_c, fac_cndc
4450 
4451  real(RP) :: fac_cndc_wrk
4452  !
4453  real(RP), parameter :: tau100day = 1.e+7_rp
4454  real(RP), parameter :: r_tau100day = 1.e-7_rp
4455  real(RP), parameter :: eps=1.e-30_rp
4456  !
4457  real(RP) :: dz
4458  !
4459  integer :: i,j,k,iqw
4460  real(RP) :: sw
4461  real(RP) :: dqv, dqc, dqr, dqi, dqs, dqg, dcv, dcp
4462  !
4463 
4464  ! [Add] 11/08/30 T.Mitsui
4465  if( flag_first )then
4466  flag_first = .false.
4467  rewind(io_fid_conf)
4468  read (io_fid_conf,nml=param_atmos_phy_mp_sn14_condensation, end=100)
4469 100 log_nml(param_atmos_phy_mp_sn14_condensation)
4470  end if
4471  !
4472 ! dt_dyn = dt*ntmax
4473 ! dt_mp = dt*(ntdiv-1)
4474  !
4475  r_dt = 1.0_rp/dt
4476  !
4477  r_xc_ccn=1.0_rp/xc_ccn
4478 ! r_xi_ccn=1.0_RP/xi_ccn
4479  !
4480  if( opt_fix_taucnd_c )then
4481  fac_cndc_wrk = fac_cndc**(1.0_rp-b_m(i_mp_qc))
4482  do j = js, je
4483  do i = is, ie
4484  do k = ks, ke
4485  pq(i_lcdep,k,i,j) = pq(i_lcdep,k,i,j)*fac_cndc_wrk
4486  end do
4487  end do
4488  end do
4489  log_info("ATMOS_PHY_MP_SN14_update_by_phase_change_kij",*) "taucnd:fac_cndc_wrk=",fac_cndc_wrk
4490  end if
4491 
4492 !OCL XFILL
4493  do j = js, je
4494  do i = is, ie
4495  do k = ks, ke
4496  ! Temperature lower limit is only used for saturation condition.
4497  ! On the other hand original "tem" is used for calculation of latent heat or energy equation.
4498  wtem(k,i,j) = max( tem(k,i,j), tem_min )
4499  end do
4500  end do
4501  end do
4502 
4503  call moist_pres2qsat_liq( ka, ks, ke, ia, is, ie, ja, js, je, &
4504  wtem(:,:,:), pre(:,:,:), qdry(:,:,:), & ! [IN]
4505  qsw(:,:,:) ) ! [OUT]
4506  call moist_pres2qsat_ice( ka, ks, ke, ia, is, ie, ja, js, je, &
4507  wtem(:,:,:), pre(:,:,:), qdry(:,:,:), & ! [IN]
4508  qsi(:,:,:) ) ! [OUT]
4509  call moist_dqs_dtem_dens_liq( ka, ks, ke, ia, is, ie, ja, js, je, &
4510  wtem(:,:,:), rho(:,:,:), & ! [IN]
4511  dqswdtem_rho(:,:,:) ) ! [OUT]
4512  call moist_dqs_dtem_dens_ice( ka, ks, ke, ia, is, ie, ja, js, je, &
4513  wtem(:,:,:), rho(:,:,:), & ! [IN]
4514  dqsidtem_rho(:,:,:) ) ! [OUT]
4515  call moist_dqs_dtem_dpre_liq( ka, ks, ke, ia, is, ie, ja, js, je, &
4516  wtem(:,:,:), pre(:,:,:), qdry(:,:,:), & ! [IN]
4517  dqswdtem_pre(:,:,:), dqswdpre_tem(:,:,:) ) ! [OUT]
4518  call moist_dqs_dtem_dpre_ice( ka, ks, ke, ia, is, ie, ja, js, je, &
4519  wtem(:,:,:), pre(:,:,:), qdry(:,:,:), & ! [IN]
4520  dqsidtem_pre(:,:,:), dqsidpre_tem(:,:,:) ) ! [OUT]
4521 
4522  profile_start("sn14_update")
4523  do j = js, je
4524  do i = is, ie
4525  do k = ks, ke
4526  if( cz(k,i,j) <= 25000.0_rp )then
4527  w2 = w(k,i,j)
4528  else
4529  w2 = 0.0_rp
4530  end if
4531  if( pre(k,i,j) < esw(k,i,j)+1.e-10_rp )then
4532  qsw(k,i,j) = 1.0_rp
4533  dqswdtem_rho(k,i,j) = 0.0_rp
4534  dqswdtem_pre(k,i,j) = 0.0_rp
4535  dqswdpre_tem(k,i,j) = 0.0_rp
4536  end if
4537  if( pre(k,i,j) < esi(k,i,j)+1.e-10_rp )then
4538  qsi(k,i,j) = 1.0_rp
4539  dqsidtem_rho(k,i,j) = 0.0_rp
4540  dqsidtem_pre(k,i,j) = 0.0_rp
4541  dqsidpre_tem(k,i,j) = 0.0_rp
4542  end if
4543 
4544  r_rvaptem = 1.0_rp/(rvap*wtem(k,i,j))
4545  lvsw = esw(k,i,j)*r_rvaptem ! rho=p/(Rv*T)
4546  lvsi = esi(k,i,j)*r_rvaptem !
4547  pv = rhoq2(i_qv,k,i,j)*rvap*tem(k,i,j)
4548  r_esw = 1.0_rp/esw(k,i,j)
4549  r_esi = 1.0_rp/esi(k,i,j)
4550  ssw = min( mp_ssw_lim, ( pv*r_esw-1.0_rp ) )
4551  ssi = pv*r_esi - 1.0_rp
4552  r_lvsw = 1.0_rp/lvsw
4553  r_lvsi = 1.0_rp/lvsi
4554  r_taucnd_c = pq(i_lcdep,k,i,j)*r_lvsw
4555  r_taucnd_r = pq(i_lrdep,k,i,j)*r_lvsw
4556  r_taudep_i = pq(i_lidep,k,i,j)*r_lvsi
4557  r_taudep_s = pq(i_lsdep,k,i,j)*r_lvsi
4558  r_taudep_g = pq(i_lgdep,k,i,j)*r_lvsi
4559 ! taucnd_c(k,i,j) = 1.0_RP/(r_taucnd_c+r_tau100day)
4560 ! taucnd_r(k,i,j) = 1.0_RP/(r_taucnd_r+r_tau100day)
4561 ! taudep_i(k,i,j) = 1.0_RP/(r_taudep_i+r_tau100day)
4562 ! taudep_s(k,i,j) = 1.0_RP/(r_taudep_s+r_tau100day)
4563 ! taudep_g(k,i,j) = 1.0_RP/(r_taudep_g+r_tau100day)
4564 
4565  r_cva = 1.0_rp / cva(k,i,j)
4566  r_cpa = 1.0_rp / cpa(k,i,j)
4567 
4568  ! Coefficient of latent heat release for ssw change by PLCdep and PLRdep
4569  aliqliq = 1.0_rp &
4570  + r_cva*( lhv00 + (cvvap-cl)*tem(k,i,j) )*dqswdtem_rho(k,i,j)
4571  ! Coefficient of latent heat release for ssw change by PLIdep, PLSdep and PLGdep
4572  asolliq = 1.0_rp &
4573  + r_cva*( lhv00 + lhf00 + (cvvap-ci)*tem(k,i,j) )*dqswdtem_rho(k,i,j)
4574  ! Coefficient of latent heat release for ssi change by PLCdep and PLRdep
4575  aliqsol = 1.0_rp &
4576  + r_cva*( lhv00 + (cvvap-cl)*tem(k,i,j) )*dqsidtem_rho(k,i,j)
4577  ! Coefficient of latent heat release for ssi change by PLIdep, PLSdep and PLGdep
4578  asolsol = 1.0_rp &
4579  + r_cva*( lhv00 + lhf00 + (cvvap-ci)*tem(k,i,j) )*dqsidtem_rho(k,i,j)
4580  pdynliq = w2 * grav * ( r_cpa*dqswdtem_pre(k,i,j) + rho(k,i,j)*dqswdpre_tem(k,i,j) )
4581  pdynsol = w2 * grav * ( r_cpa*dqsidtem_pre(k,i,j) + rho(k,i,j)*dqsidpre_tem(k,i,j) )
4582  pradliq = -dtdt_rad(k,i,j) * dqswdtem_rho(k,i,j)
4583  pradsol = -dtdt_rad(k,i,j) * dqsidtem_rho(k,i,j)
4584 
4585  ssw_o = ssw
4586  ssi_o = ssi
4587 ! ssw_o = ssw - Pdynliq*(dt_dyn-dt_mp)/qsw(k,i,j) + Pradliq*r_qsw*dt_mp
4588 ! ssi_o = ssi - Pdynsol*(dt_dyn-dt_mp)/qsi(k,i,j) + Pradsol*r_qsi*dt_mp
4589 
4590  acnd = pdynliq + pradliq &
4591  - ( r_taudep_i+r_taudep_s+r_taudep_g ) * ( qsw(k,i,j) - qsi(k,i,j) )
4592  adep = pdynsol + pradsol &
4593  + ( r_taucnd_c+r_taucnd_r ) * ( qsw(k,i,j) - qsi(k,i,j) )
4594  r_taucnd = &
4595  + aliqliq*( r_taucnd_c+r_taucnd_r ) &
4596  + asolliq*( r_taudep_i+r_taudep_s+r_taudep_g )
4597  r_taudep = &
4598  + aliqsol*( r_taucnd_c+r_taucnd_r )&
4599  + asolsol*( r_taudep_i+r_taudep_s+r_taudep_g )
4600 
4601  uplim_cnd = max( rho(k,i,j)*ssw_o*qsw(k,i,j)*r_dt, 0.0_rp )
4602  lowlim_cnd = min( rho(k,i,j)*ssw_o*qsw(k,i,j)*r_dt, 0.0_rp )
4603  if( r_taucnd < r_tau100day )then
4604 ! taucnd = tau100day
4605  pq(i_lcdep,k,i,j) = max(lowlim_cnd, min(uplim_cnd, pq(i_lcdep,k,i,j)*ssw_o ))
4606  pq(i_lrdep,k,i,j) = max(lowlim_cnd, min(uplim_cnd, pq(i_lrdep,k,i,j)*ssw_o ))
4607  pq(i_nrdep,k,i,j) = min(0.0_rp, pq(i_nrdep,k,i,j)*ssw_o )
4608 ! PLR2NR = 0.0_RP
4609  else
4610  taucnd = 1.0_rp/r_taucnd
4611  ! Production term for liquid water content
4612  coef_a_cnd = rho(k,i,j)*acnd*taucnd
4613  coef_b_cnd = rho(k,i,j)*taucnd*r_dt*(ssw_o*qsw(k,i,j)-acnd*taucnd) * ( exp(-dt*r_taucnd) - 1.0_rp )
4614  pq(i_lcdep,k,i,j) = coef_a_cnd*r_taucnd_c - coef_b_cnd*r_taucnd_c
4615  plr2nr = pq(i_nrdep,k,i,j)/(pq(i_lrdep,k,i,j)+1.e-30_rp)
4616  pq(i_lrdep,k,i,j) = coef_a_cnd*r_taucnd_r - coef_b_cnd*r_taucnd_r
4617  pq(i_nrdep,k,i,j) = min(0.0_rp, pq(i_lrdep,k,i,j)*plr2nr )
4618  end if
4619 
4620  uplim_dep = max( rho(k,i,j)*ssi_o*qsi(k,i,j)*r_dt, 0.0_rp )
4621  lowlim_dep = min( rho(k,i,j)*ssi_o*qsi(k,i,j)*r_dt, 0.0_rp )
4622  if( r_taudep < r_tau100day )then
4623 ! taudep = tau100day
4624  pq(i_lidep,k,i,j) = max(lowlim_dep, min(uplim_dep, pq(i_lidep,k,i,j)*ssi_o ))
4625  pq(i_lsdep,k,i,j) = max(lowlim_dep, min(uplim_dep, pq(i_lsdep,k,i,j)*ssi_o ))
4626  pq(i_lgdep,k,i,j) = max(lowlim_dep, min(uplim_dep, pq(i_lgdep,k,i,j)*ssi_o ))
4627  pq(i_nidep,k,i,j) = min(0.0_rp, pq(i_nidep,k,i,j)*ssi_o )
4628  pq(i_nsdep,k,i,j) = min(0.0_rp, pq(i_nsdep,k,i,j)*ssi_o )
4629  pq(i_ngdep,k,i,j) = min(0.0_rp, pq(i_ngdep,k,i,j)*ssi_o )
4630  else
4631  taudep = 1.0_rp/r_taudep
4632  ! Production term for ice water content
4633  coef_a_dep = rho(k,i,j)*adep*taudep
4634  coef_b_dep = rho(k,i,j)*taudep*r_dt*(ssi_o*qsi(k,i,j)-adep*taudep) * ( exp(-dt*r_taudep) - 1.0_rp )
4635  pli2ni = pq(i_nidep,k,i,j)/max(pq(i_lidep,k,i,j),1.e-30_rp)
4636  pls2ns = pq(i_nsdep,k,i,j)/max(pq(i_lsdep,k,i,j),1.e-30_rp)
4637  plg2ng = pq(i_ngdep,k,i,j)/max(pq(i_lgdep,k,i,j),1.e-30_rp)
4638  pq(i_lidep,k,i,j) = coef_a_dep*r_taudep_i - coef_b_dep*r_taudep_i
4639  pq(i_lsdep,k,i,j) = coef_a_dep*r_taudep_s - coef_b_dep*r_taudep_s
4640  pq(i_lgdep,k,i,j) = coef_a_dep*r_taudep_g - coef_b_dep*r_taudep_g
4641  pq(i_nidep,k,i,j) = min(0.0_rp, pq(i_lidep,k,i,j)*pli2ni )
4642  pq(i_nsdep,k,i,j) = min(0.0_rp, pq(i_lsdep,k,i,j)*pls2ns )
4643  pq(i_ngdep,k,i,j) = min(0.0_rp, pq(i_lgdep,k,i,j)*plg2ng )
4644  end if
4645 
4646  sw = 0.5_rp - sign(0.5_rp, pq(i_lcdep,k,i,j)+eps) != 1 for PLCdep<=-eps
4647  pncdep = min(0.0_rp, ((rhoq2(i_qc,k,i,j)+pq(i_lcdep,k,i,j)*dt)*r_xc_ccn - rhoq2(i_nc,k,i,j))*r_dt ) * sw
4648 ! if( PQ(I_LCdep,k,i,j) < -eps )then
4649 ! PNCdep = min(0.0_RP, ((rhoq2(I_QC,k,i,j)+PQ(I_LCdep,k,i,j)*dt)*r_xc_ccn - rhoq2(I_NC,k,i,j))*r_dt )
4650 ! else
4651 ! PNCdep = 0.0_RP
4652 ! end if
4653 ! if( PQ(I_LIdep,k,i,j) < -eps )then
4654 ! PQ(I_NIdep,k,i,j) = min(0.0_RP, ((li(k,i,j)+PQ(I_LIdep,k,i,j)*dt)*r_xi_ccn - rhoq2(I_NI,k,i,j))*r_dt )
4655 ! else
4656 ! PQ(I_NIdep,k,i,j) = 0.0_RP
4657 ! end if
4658 
4659  !--- evaporation/condensation
4660  r_rvaptem = 1.0_rp/(rvap*wtem(k,i,j))
4661  lvsw = esw(k,i,j)*r_rvaptem
4662  dlvsw = rhoq2(i_qv,k,i,j)-lvsw
4663  dcnd = dt*(pq(i_lcdep,k,i,j)+pq(i_lrdep,k,i,j))
4664 
4665  sw = ( sign(0.5_rp,dcnd) + sign(0.5_rp,dlvsw) ) &
4666  * ( 0.5_rp + sign(0.5_rp,abs(dcnd)-eps) ) ! to avoid zero division
4667  ! sw= 1: always supersaturated
4668  ! sw=-1: always unsaturated
4669  ! sw= 0: partially unsaturated during timestep
4670  fac1 = min(dlvsw*sw,dcnd*sw)*sw / (abs(sw)-1.0_rp+dcnd) & ! sw=1,-1
4671  + 1.0_rp - abs(sw) ! sw=0
4672  dep_dqc = max( dt*pq(i_lcdep,k,i,j)*fac1, &
4673  -rhoq2(i_qc,k,i,j) - 1e30_rp*(sw+1.0_rp) )*abs(sw) != -lc for sw=-1, -inf for sw=1
4674  dep_dqr = max( dt*pq(i_lrdep,k,i,j)*fac1, &
4675  -rhoq2(i_qr,k,i,j) - 1e30_rp*(sw+1.0_rp) )*abs(sw) != -lr for sw=-1, -inf for sw=1
4676 ! if ( (dcnd > eps) .AND. (dlvsw > eps) )then
4677 ! ! always supersaturated
4678 ! fac1 = min(dlvsw,dcnd)/dcnd
4679 ! dep_dqc = dt*PQ(I_LCdep,k,i,j)*fac1
4680 ! dep_dqr = dt*PQ(I_LRdep,k,i,j)*fac1
4681 ! else if( (dcnd < -eps) .AND. (dlvsw < -eps) )then
4682 ! ! always unsaturated
4683 ! fac1 = max( dlvsw,dcnd )/dcnd
4684 ! dep_dqc = max( dt*PQ(I_LCdep,k,i,j)*fac1, -rhoq2(I_QC,k,i,j) )
4685 ! dep_dqr = max( dt*PQ(I_LRdep,k,i,j)*fac1, -rhoq2(I_QR,k,i,j) )
4686 ! else
4687 ! ! partially unsaturated during timestep
4688 ! fac1 = 1.0_RP
4689 ! dep_dqc = 0.0_RP
4690 ! dep_dqr = 0.0_RP
4691 ! end if
4692 
4693  ! evaporation always lose number(always negative).
4694  dep_dnc = max( dt*pncdep*fac1, -rhoq2(i_nc,k,i,j) ) ! ss>0 dep=0, ss<0 dep<0 ! [Add] 11/08/30 T.Mitsui
4695  dep_dnr = max( dt*pq(i_nrdep,k,i,j)*fac1, -rhoq2(i_nr,k,i,j) ) ! ss>0 dep=0, ss<0 dep<0
4696 
4697  qc_evaporate(k,i,j) = - dep_dnc ! [Add] Y.Sato 15/09/08
4698 
4699  !--- deposition/sublimation
4700  lvsi = esi(k,i,j)*r_rvaptem
4701  ddep = dt*(pq(i_lidep,k,i,j)+pq(i_lsdep,k,i,j)+pq(i_lgdep,k,i,j))
4702  dlvsi = rhoq2(i_qv,k,i,j)-lvsi ! limiter for esi>1.d0
4703 
4704  sw = ( sign(0.5_rp,ddep) + sign(0.5_rp,dlvsi) ) &
4705  * ( 0.5_rp + sign(0.5_rp,abs(ddep)-eps) ) ! to avoid zero division
4706  ! sw= 1: always supersaturated
4707  ! sw=-1: always unsaturated
4708  ! sw= 0: partially unsaturated during timestep
4709  fac2 = min(dlvsi*sw,ddep*sw)*sw / (abs(sw)-1.0_rp+ddep) & ! sw=1,-1
4710  + 1.0_rp - abs(sw) ! sw=0
4711  dep_dqi = max( dt*pq(i_lidep,k,i,j) &
4712  * ( 1.0_rp-abs(sw) + fac2*abs(sw) ), & != fac2 for sw=-1,1, 1 for sw=0
4713  -rhoq2(i_qi,k,i,j) - 1e30_rp*(sw+1.0_rp) ) != -li for sw=-1, -inf for sw=0,1
4714  dep_dqs = max( dt*pq(i_lsdep,k,i,j) &
4715  * ( 1.0_rp-abs(sw) + fac2*abs(sw) ), & != fac2 for sw=-1,1, 1 for sw=0
4716  -rhoq2(i_qs,k,i,j) - 1e30_rp*(sw+1.0_rp) ) != -ls for sw=-1, -inf for sw=0,1
4717  dep_dqg = max( dt*pq(i_lgdep,k,i,j) &
4718  * ( 1.0_rp-abs(sw) + fac2*abs(sw) ), & != fac2 for sw=-1,1, 1 for sw=0
4719  -rhoq2(i_qg,k,i,j) - 1e30_rp*(sw+1.0_rp) ) != -lg for sw=-1, -inf for sw=0,1
4720 ! if ( (ddep > eps) .AND. (dlvsi > eps) )then
4721 ! ! always supersaturated
4722 ! fac2 = min(dlvsi,ddep)/ddep
4723 ! dep_dqi = dt*PQ(I_LIdep,k,i,j)*fac2
4724 ! dep_dqs = dt*PQ(I_LSdep,k,i,j)*fac2
4725 ! dep_dqg = dt*PQ(I_LGdep,k,i,j)*fac2
4726 ! else if ( (ddep < -eps) .AND. (dlvsi < -eps) )then
4727 ! ! always unsaturated
4728 ! fac2 = max(dlvsi,ddep)/ddep
4729 ! dep_dqi = max(dt*PQ(I_LIdep,k,i,j)*fac2, -rhoq2(I_QI,k,i,j) )
4730 ! dep_dqs = max(dt*PQ(I_LSdep,k,i,j)*fac2, -rhoq2(I_QS,k,i,j) )
4731 ! dep_dqg = max(dt*PQ(I_LGdep,k,i,j)*fac2, -rhoq2(I_QG,k,i,j) )
4732 ! else
4733 ! ! partially unsaturated during timestep
4734 ! fac2 = 1.0_RP
4735 ! dep_dqi = dt*PQ(I_LIdep,k,i,j)
4736 ! dep_dqs = dt*PQ(I_LSdep,k,i,j)
4737 ! dep_dqg = dt*PQ(I_LGdep,k,i,j)
4738 ! end if
4739 
4740  ! evaporation always lose number(always negative).
4741  dep_dni = max( dt*pq(i_nidep,k,i,j)*fac2, -rhoq2(i_ni,k,i,j) ) ! ss>0 dep=0, ss<0 dep<0
4742  dep_dns = max( dt*pq(i_nsdep,k,i,j)*fac2, -rhoq2(i_ns,k,i,j) ) ! ss>0 dep=0, ss<0 dep<0
4743  dep_dng = max( dt*pq(i_ngdep,k,i,j)*fac2, -rhoq2(i_ng,k,i,j) ) ! ss>0 dep=0, ss<0 dep<0
4744 
4745  !--- freezing of cloud drop
4746  frz_dqc = max( dt*(pq(i_lchom,k,i,j)+pq(i_lchet,k,i,j)), -rhoq2(i_qc,k,i,j)-dep_dqc ) ! negative value
4747  frz_dnc = max( dt*(pq(i_nchom,k,i,j)+pq(i_nchet,k,i,j)), -rhoq2(i_nc,k,i,j)-dep_dnc ) ! negative value
4748  fac3 = ( frz_dqc-eps )/( dt*(pq(i_lchom,k,i,j)+pq(i_lchet,k,i,j))-eps )
4749  fac4 = ( frz_dnc-eps )/( dt*(pq(i_nchom,k,i,j)+pq(i_nchet,k,i,j))-eps )
4750  pq(i_lchom,k,i,j) = fac3*pq(i_lchom,k,i,j)
4751  pq(i_lchet,k,i,j) = fac3*pq(i_lchet,k,i,j)
4752  pq(i_nchom,k,i,j) = fac4*pq(i_nchom,k,i,j)
4753  pq(i_nchet,k,i,j) = fac4*pq(i_nchet,k,i,j)
4754 
4755  !--- melting
4756  ! ice change
4757  mlt_dqi = max( dt*pq(i_limlt,k,i,j), -rhoq2(i_qi,k,i,j)-dep_dqi ) ! negative value
4758  mlt_dni = max( dt*pq(i_nimlt,k,i,j), -rhoq2(i_ni,k,i,j)-dep_dni ) ! negative value
4759  ! snow change
4760  mlt_dqs = max( dt*pq(i_lsmlt,k,i,j), -rhoq2(i_qs,k,i,j)-dep_dqs ) ! negative value
4761  mlt_dns = max( dt*pq(i_nsmlt,k,i,j), -rhoq2(i_ns,k,i,j)-dep_dns ) ! negative value
4762  ! graupel change
4763  mlt_dqg = max( dt*pq(i_lgmlt,k,i,j), -rhoq2(i_qg,k,i,j)-dep_dqg ) ! negative value
4764  mlt_dng = max( dt*pq(i_ngmlt,k,i,j), -rhoq2(i_ng,k,i,j)-dep_dng ) ! negative value
4765 
4766  !--- freezing of larger droplets
4767  frz_dqr = max( dt*(pq(i_lrhet,k,i,j)), min(0.0_rp, -rhoq2(i_qr,k,i,j)-dep_dqr) ) ! negative value
4768  frz_dnr = max( dt*(pq(i_nrhet,k,i,j)), min(0.0_rp, -rhoq2(i_nr,k,i,j)-dep_dnr) ) ! negative value
4769 
4770  fac5 = ( frz_dqr-eps )/( dt*pq(i_lrhet,k,i,j)-eps )
4771  pq(i_lrhet,k,i,j) = fac5*pq(i_lrhet,k,i,j)
4772  fac6 = ( frz_dnr-eps )/( dt*pq(i_nrhet,k,i,j)-eps )
4773  pq(i_nrhet,k,i,j) = fac6*pq(i_nrhet,k,i,j)
4774 
4775  ! water vapor change
4776  drhoqv = -(dep_dqc+dep_dqi+dep_dqs+dep_dqg+dep_dqr)
4777 
4778  xi = min(xi_max, max(xi_min, rhoq2(i_qi,k,i,j)/(rhoq2(i_ni,k,i,j)+ni_min) ))
4779  sw = 0.5_rp + sign(0.5_rp,xi-x_sep) ! if (xi>=x_sep) then sw=1 else sw=0
4780  ! sw=1: large ice crystals turn into rain by melting
4781 
4782  ! total cloud change
4783  drhoqc = ( frz_dqc - mlt_dqi*(1.0_rp-sw) + dep_dqc )
4784  drhonc = ( frz_dnc - mlt_dni*(1.0_rp-sw) + dep_dnc )
4785  ! total rain change
4786  drhoqr = ( frz_dqr - mlt_dqg - mlt_dqs - mlt_dqi*sw + dep_dqr )
4787  drhonr = ( frz_dnr - mlt_dng - mlt_dns - mlt_dni*sw + dep_dnr )
4788 
4789  ! total ice change
4790  drhoqi = (-frz_dqc + mlt_dqi + dep_dqi )
4791  drhoni = (-frz_dnc + mlt_dni + dep_dni )
4792 
4793  ! total snow change
4794  drhoqs = ( mlt_dqs + dep_dqs )
4795  drhons = ( mlt_dns + dep_dns )
4796 
4797  ! total graupel change
4798  drhoqg = (-frz_dqr + mlt_dqg + dep_dqg )
4799  drhong = (-frz_dnr + mlt_dng + dep_dng )
4800 
4801  ! tendency
4802  rhoq_t(k,i,j,i_qv) = drhoqv / dt
4803  rhoq_t(k,i,j,i_qc) = drhoqc / dt
4804  rhoq_t(k,i,j,i_nc) = drhonc / dt
4805  rhoq_t(k,i,j,i_qr) = drhoqr / dt
4806  rhoq_t(k,i,j,i_nr) = drhonr / dt
4807  rhoq_t(k,i,j,i_qi) = drhoqi / dt
4808  rhoq_t(k,i,j,i_ni) = drhoni / dt
4809  rhoq_t(k,i,j,i_qs) = drhoqs / dt
4810  rhoq_t(k,i,j,i_ns) = drhons / dt
4811  rhoq_t(k,i,j,i_qg) = drhoqg / dt
4812  rhoq_t(k,i,j,i_ng) = drhong / dt
4813 
4814  rhoe_t(k,i,j) = ( - lhv * drhoqv + lhf * ( drhoqi+ drhoqs + drhoqg ) ) / dt
4815 
4816  rrho = 1.0_rp/rho(k,i,j)
4817  dqv = drhoqv * rrho
4818  dqc = drhoqc * rrho
4819  dqr = drhoqr * rrho
4820  dqi = drhoqi * rrho
4821  dqs = drhoqs * rrho
4822  dqg = drhoqg * rrho
4823 
4824  dcv = cv_vapor * dqv + cv_water * ( dqc + dqr ) + cv_ice * ( dqi + dqs + dqg )
4825  dcp = cp_vapor * dqv + cp_water * ( dqc + dqr ) + cp_ice * ( dqi + dqs + dqg )
4826 
4827  cvtot_t(k,i,j) = dcv/dt
4828  cptot_t(k,i,j) = dcp/dt
4829 
4830  dz = fz(k,i,j) - fz(k,i,j)
4831  sl_plcdep(i,j) = sl_plcdep(i,j) + dep_dqc*dz
4832  sl_plrdep(i,j) = sl_plrdep(i,j) + dep_dqr*dz
4833  sl_pnrdep(i,j) = sl_pnrdep(i,j) + dep_dnr*dz
4834  end do
4835  end do
4836  end do
4837  profile_stop("sn14_update")
4838 
4839  return
4840  end subroutine update_by_phase_change_kij
4841 
4842 end module scale_atmos_phy_mp_sn14
character(len=h_short), dimension(qa_mp), parameter, public atmos_phy_mp_sn14_tracer_names
subroutine, public atmos_phy_mp_sn14_terminal_velocity(KA, KS, KE, DENS, TEMP, RHOQ, PRES, vterm)
ATMOS_PHY_MP_sn14_terminal_velocity Calculate terminal velocity.
real(rp), public const_cvdry
specific heat (dry air,constant volume) [J/kg/K]
Definition: scale_const.F90:57
real(rp), parameter, public const_psat0
saturate pressure of water vapor at 0C [Pa]
Definition: scale_const.F90:81
real(rp), public const_cpdry
specific heat (dry air,constant pressure) [J/kg/K]
Definition: scale_const.F90:56
module atmosphere / saturation
module ATMOSPHERE / Physics Cloud Microphysics
real(rp), public cv_ice
CV for ice [J/kg/K].
real(rp), parameter, public const_ci
specific heat (ice) [J/kg/K]
Definition: scale_const.F90:67
subroutine ice_multiplication_kij(KA, KS, KE, IA, IS, IE, JA, JS, JE, Pac, tem, rhoq, xq, PQ)
real(rp), parameter, public const_dwatr
density of water [kg/m3]
Definition: scale_const.F90:82
real(rp), public cp_ice
CP for ice [J/kg/K].
integer, parameter, public i_hs
snow
real(rp) function, public sf_gamma(x)
Gamma function.
subroutine, public atmos_phy_mp_sn14_cloud_fraction(KA, KS, KE, IA, IS, IE, JA, JS, JE, QTRC, mask_criterion, cldfrac)
ATMOS_PHY_MP_sn14_cloud_fraction Calculate Cloud Fraction.
real(rp), parameter, public const_cl
specific heat (liquid water) [J/kg/K]
Definition: scale_const.F90:66
subroutine mixed_phase_collection_kij(KA, KS, KE, IA, IS, IE, JA, JS, JE, wtem, rhoq, xq, dq_xave, vt_xave, PQ, Pac)
integer, parameter, public i_hr
liquid water rain
integer, parameter, public i_hi
ice water cloud
subroutine, public atmos_phy_mp_sn14_setup(KA, IA, JA)
ATMOS_PHY_MP_sn14_setup setup.
real(rp), parameter, public const_tem00
temperature reference (0C) [K]
Definition: scale_const.F90:90
integer, public io_fid_conf
Config file ID.
Definition: scale_io.F90:55
real(rp), public const_cvvap
specific heat (water vapor, constant volume) [J/kg/K]
Definition: scale_const.F90:65
integer, parameter, public i_hh
hail
real(rp), public const_lhf0
latent heat of fusion at 0C [J/kg]
Definition: scale_const.F90:79
subroutine, public atmos_phy_mp_sn14_tendency(KA, KS, KE, IA, IS, IE, JA, JS, JE, DENS, W, QTRC, PRES, TEMP, Qdry, CPtot, CVtot, CCN, dt, cz, fz, RHOQ_t, RHOE_t, CPtot_t, CVtot_t, EVAPORATE)
ATMOS_PHY_MP_sn14_tendency calculate tendency.
subroutine, public atmos_phy_mp_sn14_qtrc2nhyd(KA, KS, KE, IA, IS, IE, JA, JS, JE, QTRC0, Ne)
Calculate number concentration of each category.
real(rp), public const_rdry
specific gas constant (dry air) [J/kg/K]
Definition: scale_const.F90:55
real(rp), public const_undef
Definition: scale_const.F90:41
character(len=h_short), dimension(qa_mp), parameter, public atmos_phy_mp_sn14_tracer_units
subroutine, public atmos_phy_mp_sn14_qtrc2qhyd(KA, KS, KE, IA, IS, IE, JA, JS, JE, QTRC0, Qe)
ATMOS_PHY_MP_sn14_qtrc2qhyd Calculate mass ratio of each category.
real(rp), parameter, public const_lhs0
latent heat of sublimation at 0C [J/kg]
Definition: scale_const.F90:77
real(rp), public cv_vapor
CV for vapor [J/kg/K].
module atmosphere / hydrometeor
integer, parameter, public atmos_phy_mp_sn14_ntracers
character(len=h_mid), dimension(qa_mp), parameter, public atmos_phy_mp_sn14_tracer_descriptions
real(rp), parameter, public const_lhv0
latent heat of vaporizaion at 0C [J/kg]
Definition: scale_const.F90:75
real(rp), public const_pre00
pressure reference [Pa]
Definition: scale_const.F90:88
integer, parameter, public atmos_phy_mp_sn14_nwaters
module PROCESS
Definition: scale_prc.F90:11
integer, parameter, public qa_mp
subroutine nucleation_kij(KA, KS, KE, IA, IS, IE, JA, JS, JE, cz, fz, w, rho, tem, pre, qdry, rhoq, cpa, dTdt_rad, qke, CCN, dt, PQ)
real(rp), public const_lhf00
latent heat of fusion at 0K [J/kg]
Definition: scale_const.F90:80
real(rp), public const_grav
standard acceleration of gravity [m/s2]
Definition: scale_const.F90:46
subroutine freezing_water_kij(KA, KS, KE, IA, IS, IE, JA, JS, JE, dt, rhoq, xq, tem, PQ)
real(rp), public const_lhv00
latent heat of vaporizaion at 0K [J/kg]
Definition: scale_const.F90:76
integer, public prc_myrank
process num in local communicator
Definition: scale_prc.F90:89
module SPECFUNC
real(rp), parameter, public const_rvap
specific gas constant (water vapor) [J/kg/K]
Definition: scale_const.F90:63
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:338
module CONSTANT
Definition: scale_const.F90:11
integer, parameter, public i_hc
liquid water cloud
integer, parameter, public atmos_phy_mp_sn14_nices
subroutine debug_tem_kij(KA, KS, KE, IA, IS, IE, JA, JS, JE, point, tem, rho, pre, qv)
real(rp), parameter, public const_emelt
Definition: scale_const.F90:72
subroutine, public atmos_phy_mp_sn14_qhyd2qtrc(KA, KS, KE, IA, IS, IE, JA, JS, JE, Qe, QTRC, QNUM)
subroutine aut_acc_slc_brk_kij(KA, KS, KE, IA, IS, IE, JA, JS, JE, rhoq, xq, dq_xave, rho, PQ)
subroutine, public prof_rapstart(rapname_base, level)
Start raptime.
Definition: scale_prof.F90:157
real(dp), parameter, public const_undef8
undefined value (REAL8)
Definition: scale_const.F90:40
module profiler
Definition: scale_prof.F90:11
real(rp), public const_eps
small number
Definition: scale_const.F90:33
subroutine, public atmos_phy_mp_sn14_effective_radius(KA, KS, KE, IA, IS, IE, JA, JS, JE, DENS0, TEMP0, QTRC0, Re)
ATMOS_PHY_MP_sn14_effective_radius Calculate Effective Radius.
module PRECISION
real(rp), public const_pi
pi
Definition: scale_const.F90:31
module STDIO
Definition: scale_io.F90:10
integer, parameter, public n_hyd
real(rp), parameter, public const_cpvap
specific heat (water vapor, constant pressure) [J/kg/K]
Definition: scale_const.F90:64
subroutine, public prof_rapend(rapname_base, level)
Save raptime.
Definition: scale_prof.F90:210
real(rp), public cp_vapor
CP for vapor [J/kg/K].
real(rp), public cp_water
CP for water [J/kg/K].
integer, parameter, public i_hg
graupel
subroutine update_by_phase_change_kij(KA, KS, KE, IA, IS, IE, JA, JS, JE, ntdiv, ntmax, dt, cz, fz, w, dTdt_rad, rho, qdry, esw, esi, rhoq2, pre, tem, cpa, cva, PQ, sl_PLCdep, sl_PLRdep, sl_PNRdep, RHOQ_t, RHOE_t, CPtot_t, CVtot_t, qc_evaporate)
real(rp), public cv_water
CV for water [J/kg/K].