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