SCALE-RM
scale_atmos_phy_mp_suzuki10.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
21 !-------------------------------------------------------------------------------
22 #include "scalelib.h"
24  !-----------------------------------------------------------------------------
25  !
26  !++ Used modules
27  !
28  use scale_precision
29  use scale_io
30  use scale_prof
31 
32  !-----------------------------------------------------------------------------
33  implicit none
34  private
35  !-----------------------------------------------------------------------------
36  !
37  !++ Public procedure
38  !
50 
51  !-----------------------------------------------------------------------------
52  !
53  !++ Public parameters & variables
54  !
55  !-----------------------------------------------------------------------------
56  integer, public :: atmos_phy_mp_suzuki10_ntracers
57  integer, public :: atmos_phy_mp_suzuki10_nwaters
58  integer, public :: atmos_phy_mp_suzuki10_nices
59  integer, public :: atmos_phy_mp_suzuki10_nccn
60  integer, public :: atmos_phy_mp_suzuki10_nbnd
61 
62 
63  character(len=H_SHORT), public, allocatable :: atmos_phy_mp_suzuki10_tracer_names(:)
64  character(len=H_MID) , public, allocatable :: atmos_phy_mp_suzuki10_tracer_descriptions(:)
65  character(len=H_SHORT), public, allocatable :: atmos_phy_mp_suzuki10_tracer_units(:)
66 
67 # include "kernels.h"
68  !-----------------------------------------------------------------------------
69  !
70  !++ Private procedure
71  !
72  private :: mp_suzuki10
73  private :: nucleat
74  private :: nucleata
75  private :: cndevpsbl
76  private :: cndevpsbla
77  private :: liqphase
78  private :: icephase
79  private :: mixphase
80  private :: ice_nucleat
81  private :: freezing
82  private :: melting
83  private :: collmain
84  private :: collmainf
85  private :: collcoag
86  private :: getrule
87  private :: faero
88  private :: random_setup
89  private :: r_collcoag
90  private :: mkpara
91  private :: rdkdat
92  private :: sdfgrid
93  private :: getcp
94  private :: fcpc
95  private :: getck
96  private :: fckrn
97  private :: getvt
98  private :: fvterm
99  private :: getbr
100  private :: fbulkrad
101  private :: paraout
102  private :: tinvss
103  private :: getknot
104  private :: fbspl
105 ! private :: fpb
106  private :: getmatrx
107  private :: getcoef
108  private :: fspline
109  private :: getcoef2
110  private :: fspline2
111 
112  !-----------------------------------------------------------------------------
113  !
114  !++ Private parameters
115  !
116  integer :: qa
117  !$acc declare create(QA)
118  integer, parameter :: i_qv = 1
119 
120  character(len=3) :: namspc (8) = (/ 'Qcl', &
121  'Qic', &
122  'Qip', &
123  'Qid', &
124  'Qis', &
125  'Qig', &
126  'Qih', &
127  'Qae' /)
128 
129  character(len=27) :: lnamspc(8) = (/ 'Mixing ratio of cloud bin', &
130  'Mixing ratio of colum bin', &
131  'Mixing ratio of plate bin', &
132  'Mixing ratio of dendrit bin', &
133  'Mixing ratio of snow bin', &
134  'Mixing ratio of graupel bin', &
135  'Mixing ratio of hail bin', &
136  'Mixing ratio of aerosol bin' /)
137 
138  integer, public :: nbin = 33 ! tentatively public
139  integer :: nspc = 7
140  integer, public :: nccn = 0 ! tentatively public
141  integer :: kphase = 0
142  integer :: iceflg = 1
143  integer :: num_hyd = 0
144 
145  integer, parameter :: i_mp_qc = 1
146  integer, parameter :: i_mp_qcl = 2
147  integer, parameter :: i_mp_qp = 3
148  integer, parameter :: i_mp_qd = 4
149  integer, parameter :: i_mp_qs = 5
150  integer, parameter :: i_mp_qg = 6
151  integer, parameter :: i_mp_qh = 7
152 
153  logical :: mp_doautoconversion = .true. ! apply collision process ?
154  logical :: mp_couple_aerosol = .false. ! apply CCN effect?
155 
156  !--- array definition
157  integer :: num_start_waters
158  integer :: num_end_waters
159  integer :: num_start_ices
160  integer :: num_end_ices
161 
162  !--- Indeces for determining species of cloud particle
163  integer, parameter :: il = 1 !--- index for liquid water
164  integer, parameter :: ic = 2 !--- index for columnar ice
165  integer, parameter :: ip = 3 !--- index for plate ice
166  integer, parameter :: id = 4 !--- index for dendrite ice
167  integer, parameter :: iss = 5 !--- index for snow
168  integer, parameter :: ig = 6 !--- index for graupel
169  integer, parameter :: ih = 7 !--- index for hail
170 
171  !--- bin information of hydrometeors
172  real(rp) :: dxmic !--- d( log(m) ) of hydrometeor bin
173  real(rp), allocatable :: xctr( : ) !--- log( m ) value of bin center [xctr = exp( 4/3*pi*DWATER*radc^3 )]
174  real(rp), allocatable :: xbnd( : ) !--- log( m ) value of bin boundary
175  real(rp), allocatable :: radc( : ) !--- radius of hydrometeor at bin center [m]
176  real(rp), allocatable :: cctr( :,: ) !--- capacitance of hydrometeor at bin center (C of equation A.17 in Suzuki 2004)
177  real(rp), allocatable :: cbnd( :,: ) !--- capacitance of hydrometeor at bin boundary (C of equation A.17 in Suzuki 2004)
178  real(rp), allocatable :: ck( :,:,:,: ) !-- collection kernel (K of equation A.20 in Suzuki 2004)
179  real(rp), allocatable :: vt( :,: ) !--- terminal velocity of hydrometeor [m/s]
180  real(rp), allocatable :: br( :,: ) !--- bulk density of hydrometeor [kg/m^3]
181  integer, allocatable :: ifrsl( :,:,: ) !--- type of species after collision
182  !--- bin information of aerosol (not supported)
183  real(rp), allocatable :: xactr( : ) !--- log( ma ) value of bin center ( ma is mass of aerosol )
184  real(rp), allocatable :: xabnd( : ) !--- log( ma ) value of bin boundary ( ma is mass of aerosol )
185  real(rp), allocatable :: rada( : ) !--- radius of aerosol at bin center [m]
186 
187  real(rp), allocatable :: expxctr( : ) !--- exp( xctr )
188  real(rp), allocatable :: expxbnd( : ) !--- exp( xbnd )
189  real(rp), allocatable :: expxactr( : ) !--- exp( xactr )
190  real(rp), allocatable :: expxabnd( : ) !--- exp( xabnd )
191  real(rp), allocatable :: rexpxctr( : ) !--- 1/exp( xctr )
192  real(rp), allocatable :: rexpxbnd( : ) !--- 1/exp( xbnd )
193  real(rp), allocatable :: rexpxactr( : ) !--- 1/exp( xactr )
194  real(rp), allocatable :: rexpxabnd( : ) !--- 1/exp( xabnd1 )
195 
196  real(rp) :: dxaer !--- d( log(ma) ) of aerosol bin
197  real(rp) :: xasta !--- exponential of mass of aerosol for smallest aerosol bin
198  real(rp) :: xaend !--- exponential of mass of aerosol for largest aerosol bin
199 
200  real(rp), allocatable, save :: vterm(:) !--- terminal velocity
201  !$acc declare create(vterm)
202 
203  !--- constant for bin
204  real(rp), parameter :: cldmin = 1.0e-10_rp !--- threshould for cloud is regarded existing
205  real(rp), parameter :: oneovthird = 1.0_rp / 3.0_rp
206  real(rp), parameter :: thirdovforth = 3.0_rp / 4.0_rp
207  real(rp), parameter :: twoovthird = 2.0_rp / 3.0_rp
208 
209  real(rp) :: rbnd = 40.e-06_rp ! boundary radius of cloud and rian
210  integer :: nbnd ! boundary bin number corresponding to rbnd
211 
212  !--- constant for aerosol
213  real(rp) :: rhoa = 2.25e+03_rp ! density of aerosol ( NaCl )
214  real(rp) :: emaer = 58.0_rp ! molecular weight of aerosol ( salt )
215  real(rp) :: emwtr = 18.0_rp ! molecular weight of water
216  real(rp) :: rasta = 1.e-08_rp ! minimum radius of aerosol (m)
217  real(rp) :: raend = 1.e-06_rp ! maximum radius of aerosol (m)
218  real(rp) :: r0a = 1.e-07_rp ! average radius of aerosol (m)
219 
220  logical :: flg_regeneration ! flag regeneration of aerosol
221  logical :: flg_nucl ! flag nucleated cloud move into smallest bin
222  logical :: flg_icenucl ! flag ice nucleation
223  logical :: flg_sf_aero ! flag surface flux of aerosol
224  logical :: flg_rndm ! flag for sthastic integration for coll.-coag.
225 
226  real(rp), allocatable :: marate( : ) ! mass rate of each aerosol bin to total aerosol mass
227 ! integer, save :: K10_1, K10_2 ! scaling factor for 10m value (momentum)
228 ! real(RP) :: R10M1, R10M2 ! scaling factor for 10m value (momentum)
229 ! real(RP) :: R10H1, R10H2 ! scaling factor for 10m value (heat)
230 ! real(RP) :: R10E1, R10E2 ! scaling factor for 10m value (tracer)
231 
232  character(len=11), parameter :: fname_micpara="micpara.dat" !--- file name
233 
234  integer(4) :: fid_micpara
235 
236  !--- Use for stochastic method
237  integer, allocatable :: blrg( :,: )
238  integer, allocatable :: bsml( :,: )
239  real(rp) :: wgtbin
240  integer :: mspc, mbin
241 
242  !--- use for ice_nucleation
243  real(rp) :: n0_icenucl = 1.e+3_rp ! N0 of Meyer et al. (1992)
244 
245  !--- use for model without aerosol
246  real(rp) :: c_ccn = 100.e+6_rp ! N0 of Nc = N0*s^kappa
247  real(rp) :: kappa = 0.462_rp ! kappa of Nc = N0*s^kappa
248 
249  !--- use for aerosol coupled model
250  real(rp) :: sigma = 7.5e-02_rp ! water surface tension [ N/m2 ] (sigma in eq. (A.11) of Suzuki (2004) )
251  real(rp) :: vhfct = 2.0_rp ! van't hoff factor (i in eq.(A.11) of Suzuki (2004))
252 
253  real(rp), parameter :: tcrit = 271.15_rp
254  integer, allocatable :: kindx( :,: )
255 
256  !--- for creating micpara.dat (mkpara)
257  integer, parameter :: ndat = 33, icemax = 3
258  integer, parameter :: kdeg = 4, ldeg = 4, nspc_mk = 7
259 
260  real(dp) :: dxmic_mk
261 
262  real(dp), allocatable :: radc_mk( : ), xctr_mk( : ), xbnd_mk( : )
263  real(dp), allocatable :: cctr_mk( :,: ), cbnd_mk( :,: )
264  real(dp), allocatable :: ck_mk( :,:,:,: )
265  real(dp), allocatable :: vt_mk( :,: )
266  real(dp), allocatable :: br_mk( :,: )
267 
268  real(dp) :: xmss( nspc_mk,ndat ), zcap( nspc_mk,ndat ), vtrm( nspc_mk,ndat )
269  real(dp) :: blkr( nspc_mk,ndat ), blkd( nspc_mk,ndat ), ykrn( nspc_mk,nspc_mk,ndat,ndat )
270 
271  real(dp) :: ywll( ndat,ndat ), ywli( ndat,ndat,icemax ), ywls( ndat,ndat )
272  real(dp) :: ywlg( ndat,ndat ), ywlh( ndat,ndat )
273 
274  real(dp) :: ywil( ndat,ndat,icemax ), ywii( ndat,ndat,icemax,icemax )
275  real(dp) :: ywis( ndat,ndat,icemax ), ywig( ndat,ndat,icemax )
276  real(dp) :: ywih( ndat,ndat,icemax )
277 
278  real(dp) :: ywsl( ndat,ndat ), ywsi( ndat,ndat,icemax ), ywss( ndat,ndat )
279  real(dp) :: ywsg( ndat,ndat ), ywsh( ndat,ndat )
280 
281  real(dp) :: ywgl( ndat,ndat ), ywgi( ndat,ndat,icemax ), ywgs( ndat,ndat )
282  real(dp) :: ywgg( ndat,ndat ), ywgh( ndat,ndat )
283 
284  real(dp) :: ywhl( ndat,ndat ), ywhi( ndat,ndat,icemax ), ywhs( ndat,ndat )
285  real(dp) :: ywhg( ndat,ndat ), ywhh( ndat,ndat )
286 
287  ! for qhyd2qtrc
288  real(rp) :: sigma_sdf(5), r0_sdf(5), n0_sdf(5), rho_sdf(5)
289 
290  ! for lightning
291  real(rp), allocatable :: flg_noninduct(:,:)
292  real(rp), allocatable :: ecoll(:,:,:,:), rcoll(:,:,:,:)
293  real(rp), allocatable, save :: num_end(:,:,:,:)
294  real(rp), save :: ecoal_gsi = 0.0_rp
295  !----------------------------------------------------------------------------
296 contains
297  !-----------------------------------------------------------------------------
300  use scale_prc, only: &
301  prc_abort
302  implicit none
303 
304  namelist / param_atmos_phy_mp_suzuki10_bin / &
305  nbin, &
306  nccn, &
307  iceflg, &
308  kphase
309 
310  integer :: m, n, ierr
311  !---------------------------------------------------------------------------
312 
313  log_newline
314  log_info("ATMOS_PHY_MP_suzuki10_tracer_setup",*) 'Setup'
315  log_info("ATMOS_PHY_MP_suzuki10_tracer_setup",*) 'Tracers setup for Suzuki (2010) Spectral BIN model'
316 
317  log_newline
318  log_info("ATMOS_PHY_MP_suzuki10_tracer_setup",*) 'READ BIN NUMBER'
319 
320  rewind(io_fid_conf)
321  read(io_fid_conf,nml=param_atmos_phy_mp_suzuki10_bin,iostat=ierr)
322 
323  if( ierr < 0 ) then !--- missing
324  log_info("ATMOS_PHY_MP_suzuki10_tracer_setup",*) 'Not found namelist. Default used.'
325  elseif( ierr > 0 ) then !--- fatal error
326  log_error("ATMOS_PHY_MP_suzuki10_tracer_setup",*) 'Not appropriate names in namelist PARAM_ATMOS_PHY_MP_SUZUKI10_bin, Check!'
327  call prc_abort
328  end if
329 
330  log_nml(param_atmos_phy_mp_suzuki10_bin)
331 
332  if( iceflg == 0 ) then
333  nspc = 1
334  elseif( iceflg == 1 ) then
335  nspc = 7
336  else
337  log_error("ATMOS_PHY_MP_suzuki10_tracer_setup",*) "ICEFLG should be 0 (warm rain) or 1 (mixed rain) check!!"
338  call prc_abort
339  endif
340 
341  atmos_phy_mp_suzuki10_ntracers = 1 + nbin*nspc + nccn ! number of total tracers
342  atmos_phy_mp_suzuki10_nwaters = nbin ! number of liquid water
343  atmos_phy_mp_suzuki10_nices = nbin * ( nspc - 1 ) ! number of ice water
344  atmos_phy_mp_suzuki10_nccn = nccn ! number of ccn
345 
346  num_hyd = nbin * nspc
347 
348  num_start_waters = i_qv + 1
349  num_end_waters = i_qv + atmos_phy_mp_suzuki10_nwaters
350  num_start_ices = num_end_waters + 1
351  num_end_ices = num_end_waters + atmos_phy_mp_suzuki10_nices
352 
354 
355  allocate( atmos_phy_mp_suzuki10_tracer_names(qa) )
357  allocate( atmos_phy_mp_suzuki10_tracer_units(qa) )
358 
359  !---------------------------------------------------------------------------
360  !
361  !++ calculate each category and aerosol
362  !
363  !---------------------------------------------------------------------------
364 
365  do n = 1, qa
366  write(atmos_phy_mp_suzuki10_tracer_units(n),'(a)') 'kg/kg'
367  enddo
368 
369  write(atmos_phy_mp_suzuki10_tracer_names(1),'(a)') 'QV'
370  write(atmos_phy_mp_suzuki10_tracer_descriptions(1),'(a)') 'Water Vapor mixing ratio'
371 
372  do m = 1, nspc
373  do n = 1, nbin
374  write(atmos_phy_mp_suzuki10_tracer_names(1+nbin*(m-1)+n),'(a,i0)') trim(namspc(m)), n
375  write(atmos_phy_mp_suzuki10_tracer_descriptions(1+nbin*(m-1)+n),'(a,i0)') trim(lnamspc(m)), n
376  enddo
377  enddo
378 
379  do n = 1, nccn
380  write(atmos_phy_mp_suzuki10_tracer_names(1+nbin*nspc+n),'(a,i0)') trim(namspc(8)), n
381  write(atmos_phy_mp_suzuki10_tracer_descriptions(1+nbin*nspc+n),'(a,i0)') trim(lnamspc(8)), n
382  enddo
383 
384  !$acc update device(QA)
385 
386  return
388 
389  !-----------------------------------------------------------------------------
391  subroutine atmos_phy_mp_suzuki10_setup( &
392  KA, IA, JA, &
393  flg_lt )
394  use scale_prc, only: &
395  prc_abort, &
396  prc_masterrank, &
398  use scale_const, only: &
399  pi => const_pi, &
400  dwatr => const_dwatr, &
401  dice => const_dice
402  use scale_comm_cartesc, only: &
403  comm_bcast
404  use scale_atmos_hydrometeor, only: &
405  i_hc, &
406  i_hr, &
407  i_hi, &
408  i_hs, &
409  i_hg, &
410  i_hh
411  implicit none
412 
413  integer, intent(in) :: ka
414  integer, intent(in) :: ia
415  integer, intent(in) :: ja
416  logical, intent(in), optional :: flg_lt
417 
418  real(rp) :: rho_aero !--- density of aerosol
419  real(rp) :: r0_aero !--- center radius of aerosol (um)
420  real(rp) :: r_min !--- minimum radius of aerosol (um)
421  real(rp) :: r_max !--- maximum radius of aerosol (um)
422  real(rp) :: s10_emaer !--- moleculer weight of aerosol
423 
424  logical :: s10_flag_regene = .false. !--- flag of regeneration
425  logical :: s10_flag_nucleat = .false. !--- flag of regeneration
426  logical :: s10_flag_icenucleat = .false. !--- flag of ice nucleation
427  logical :: s10_flag_sfaero = .false. !--- flag of surface flux of aeorol
428  logical :: s10_flag_rndm = .false. !--- flag for sthastic integration for coll.-coag.
429  integer :: s10_rndm_mspc
430  integer :: s10_rndm_mbin
431 
432  namelist / param_atmos_phy_mp_suzuki10 / &
433  mp_doautoconversion, &
434  mp_couple_aerosol, &
435  rho_aero, &
436  r_min, &
437  r_max, &
438  r0_aero, &
439  s10_emaer, &
440  s10_flag_regene, &
441  s10_flag_nucleat, &
442  s10_flag_icenucleat, &
443  s10_flag_sfaero, &
444  s10_flag_rndm, &
445  s10_rndm_mspc, &
446  s10_rndm_mbin, &
447  c_ccn, kappa, &
448  n0_icenucl, &
449  sigma, vhfct, &
450  ecoal_gsi
451 
452  real(rp), parameter :: max_term_vel = 10.0_rp !-- terminal velocity for calculate dt of sedimentation
453 
454  logical :: flg_lt_
455  integer :: nnspc, nnbin
456  integer :: nn, mm, mmyu, nnyu
457  integer :: myu, nyu, i, j, k, n, ierr
458  !---------------------------------------------------------------------------
459 
460  log_newline
461  log_info("ATMOS_PHY_MP_suzuki10_setup",*) 'Setup'
462  log_info("ATMOS_PHY_MP_suzuki10_setup",*) 'Suzuki (2010) Spectral BIN model'
463 
464  !--- allocation
465  allocate( xctr( nbin ) )
466  allocate( xbnd( nbin+1 ) )
467  allocate( radc( nbin ) )
468  allocate( cctr( nbin,nspc_mk ) )
469  allocate( cbnd( nbin+1,nspc_mk ) )
470  allocate( ck( nspc_mk,nspc_mk,nbin,nbin ) )
471  allocate( vt( nspc_mk,nbin ) )
472  allocate( br( nspc_mk,nbin ) )
473  allocate( ifrsl( 2,nspc_mk,nspc_mk ) )
474  allocate( expxctr( nbin ) )
475  allocate( expxbnd( nbin+1 ) )
476  allocate( rexpxctr( nbin ) )
477  allocate( rexpxbnd( nbin+1 ) )
478  if ( nccn /= 0 ) then
479  allocate( xactr( nccn ) )
480  allocate( xabnd( nccn+1 ) )
481  allocate( rada( nccn ) )
482  allocate( expxactr( nccn ) )
483  allocate( expxabnd( nccn+1 ) )
484  allocate( rexpxactr( nccn ) )
485  allocate( rexpxabnd( nccn+1 ) )
486  endif
487 
488  mbin = nbin/2
489  mspc = nspc_mk*nspc_mk
490 
491  rho_aero = rhoa
492  s10_emaer = emaer
493  r_min = rasta
494  r_max = raend
495  r0_aero = r0a
496  s10_rndm_mspc = mspc
497  s10_rndm_mbin = mbin
498 
499  rewind(io_fid_conf)
500  read(io_fid_conf,nml=param_atmos_phy_mp_suzuki10,iostat=ierr)
501 
502  if ( ierr < 0 ) then !--- missing
503  log_info("ATMOS_PHY_MP_suzuki10_setup",*) 'Not found namelist. Default used.'
504  elseif( ierr > 0 ) then !--- fatal error
505  log_error("ATMOS_PHY_MP_suzuki10_setup",*) 'Not appropriate names in namelist PARAM_ATMOS_PHY_MP_SUZUKI10, Check!'
506  call prc_abort
507  endif
508  log_nml(param_atmos_phy_mp_suzuki10)
509 
510  if ( nspc /= 1 .AND. nspc /= 7 ) then
511  log_error("ATMOS_PHY_MP_suzuki10_setup",*) 'nspc should be set as 1 (warm rain) or 7 (mixed phase) check!'
512  call prc_abort
513  endif
514 
515  rhoa = rho_aero
516  emaer = s10_emaer
517  rasta = r_min
518  raend = r_max
519  r0a = r0_aero
520  flg_regeneration = s10_flag_regene
521  flg_nucl = s10_flag_nucleat
522  flg_icenucl = s10_flag_icenucleat
523  flg_sf_aero = s10_flag_sfaero
524  flg_rndm = s10_flag_rndm
525  mspc = s10_rndm_mspc
526  mbin = s10_rndm_mbin
527 
528  !--- read micpara.dat (microphysical parameter) and broad cast
529  if ( prc_ismaster ) then
530 
531  fid_micpara = io_get_available_fid()
532  !--- open parameter of cloud microphysics
533  open ( fid_micpara, file = fname_micpara, form = 'formatted', status = 'old', iostat=ierr )
534 
535  !--- micpara.dat does not exist
536  if ( ierr == 0 ) then
537 
538  read( fid_micpara,* ) nnspc, nnbin
539 
540  if ( nnbin /= nbin ) then
541  log_error("ATMOS_PHY_MP_suzuki10_setup",*) 'nbin in namelist and nbin in micpara.dat is different check!'
542  call prc_abort
543  endif
544 
545  ! grid parameter
546  log_info("ATMOS_PHY_MP_suzuki10_setup",*) 'Radius of cloud *'
547  do n = 1, nbin
548  read( fid_micpara,* ) nn, xctr( n ), radc( n )
549  log_info("ATMOS_PHY_MP_suzuki10_setup",'(A,1x,I3,1x,A,1x,ES15.7,1x,A)') &
550  "Radius of ", n, "th cloud bin (bin center)= ", radc( n ) , "[m]"
551  enddo
552  do n = 1, nbin+1
553  read( fid_micpara,* ) nn, xbnd( n )
554  enddo
555  read( fid_micpara,* ) dxmic
556  log_info("ATMOS_PHY_MP_suzuki10_setup",*) 'Width of Cloud SDF= ', dxmic
557 
558  ! capacity
559  do myu = 1, nspc_mk
560  do n = 1, nbin
561 ! read( fid_micpara,* ) mmyu, nn, cctr( myu,n )
562  read( fid_micpara,* ) mmyu, nn, cctr( n,myu )
563  enddo
564  do n = 1, nbin+1
565 ! read( fid_micpara,* ) mmyu, nn, cbnd( myu,n )
566  read( fid_micpara,* ) mmyu, nn, cbnd( n,myu )
567  enddo
568  enddo
569 
570  ! collection kernel
571  do myu = 1, nspc_mk
572  do nyu = 1, nspc_mk
573  do i = 1, nbin
574  do j = 1, nbin
575  read( fid_micpara,* ) mmyu, nnyu, mm, nn, ck( myu,nyu,i,j )
576  enddo
577  enddo
578  enddo
579  enddo
580 
581  ! terminal velocity
582  do myu = 1, nspc_mk
583  do n = 1, nbin
584  read( fid_micpara,* ) mmyu, nn, vt( myu,n )
585  enddo
586  enddo
587 
588  ! bulk density
589  do myu = 1, nspc_mk
590  do n = 1, nbin
591  read( fid_micpara,* ) mmyu, nn, br( myu,n )
592  enddo
593  enddo
594 
595  close ( fid_micpara )
596 
597  !--- micpara.dat does not exist
598  else
599 
600  log_info("ATMOS_PHY_MP_suzuki10_setup",*) 'micpara.dat is created'
601  call mkpara
602 
603  fid_micpara = io_get_available_fid()
604  !--- open parameter of cloud microphysics
605  open ( fid_micpara, file = fname_micpara, form = 'formatted', status = 'old', iostat=ierr )
606 
607  read( fid_micpara,* ) nnspc, nnbin
608 
609  if ( nnbin /= nbin ) then
610  log_error("ATMOS_PHY_MP_suzuki10_setup",*) 'nbin in inc_tracer and nbin in micpara.dat is different check!'
611  call prc_abort
612  endif
613 
614  ! grid parameter
615  log_info("ATMOS_PHY_MP_suzuki10_setup",*) 'Radius of cloud *'
616  do n = 1, nbin
617  read( fid_micpara,* ) nn, xctr( n ), radc( n )
618  log_info("ATMOS_PHY_MP_suzuki10_setup",'(A,1x,I3,1x,A,1x,ES15.7,1x,A)') &
619  "Radius of ", n, "th cloud bin (bin center)= ", radc( n ) , "[m]"
620  enddo
621  do n = 1, nbin+1
622  read( fid_micpara,* ) nn, xbnd( n )
623  enddo
624  read( fid_micpara,* ) dxmic
625  log_info("ATMOS_PHY_MP_suzuki10_setup",*) 'Width of Cloud SDF= ', dxmic
626 
627  ! capacity
628  do myu = 1, nspc_mk
629  do n = 1, nbin
630 ! read( fid_micpara,* ) mmyu, nn, cctr( myu,n )
631  read( fid_micpara,* ) mmyu, nn, cctr( n,myu )
632  enddo
633  do n = 1, nbin+1
634 ! read( fid_micpara,* ) mmyu, nn, cbnd( myu,n )
635  read( fid_micpara,* ) mmyu, nn, cbnd( n,myu )
636  enddo
637  enddo
638 
639  ! collection kernel
640  do myu = 1, nspc_mk
641  do nyu = 1, nspc_mk
642  do i = 1, nbin
643  do j = 1, nbin
644  read( fid_micpara,* ) mmyu, nnyu, mm, nn, ck( myu,nyu,i,j )
645  enddo
646  enddo
647  enddo
648  enddo
649 
650  ! terminal velocity
651  do myu = 1, nspc_mk
652  do n = 1, nbin
653  read( fid_micpara,* ) mmyu, nn, vt( myu,n )
654  enddo
655  enddo
656 
657  ! bulk density
658  do myu = 1, nspc_mk
659  do n = 1, nbin
660  read( fid_micpara,* ) mmyu, nn, br( myu,n )
661  enddo
662  enddo
663 
664  close ( fid_micpara )
665 
666  endif
667 
668  endif
669 
670  call comm_bcast( nbin, radc(:) )
671  call comm_bcast( nbin, xctr(:) )
672  call comm_bcast( dxmic )
673  call comm_bcast( nbin+1, xbnd(:) )
674  call comm_bcast( nbin, nspc_mk, cctr(:,:) )
675  call comm_bcast( nbin+1, nspc_mk, cbnd(:,:) )
676  call comm_bcast( nspc_mk, nspc_mk, nbin, nbin, ck(:,:,:,:) )
677  call comm_bcast( nspc_mk, nbin, br(:,:) )
678  call comm_bcast( nspc_mk, nbin, vt(:,:) )
679 
680  allocate( flg_noninduct( nspc,nspc ) )
681  allocate( ecoll( nspc,nspc,nbin,nbin ) )
682  allocate( rcoll( nspc,nspc,nbin,nbin ) )
683  !$omp workshare
684  flg_noninduct(:,:) = 0.0_rp
685  ecoll( :,:,:,: ) = 0.0_rp
686  rcoll( :,:,:,: ) = 0.0_rp
687  !$omp end workshare
688 
689  if ( present(flg_lt) ) then
690  flg_lt_ = flg_lt
691  else
692  flg_lt_ = .false.
693  end if
694  if( flg_lt_ ) then
695 
696  !$omp parallel do
697  do nyu = 1, nspc
698  do myu = 1, nspc
699  if( ( myu >= ic .and. myu <= iss ) .and. ( nyu == ig .or. nyu == ih ) ) then
700  flg_noninduct( myu,nyu ) = 1.0_rp
701  endif
702  enddo
703  enddo
704 
705  !$omp parallel do collapse(2)
706  do j = 1, nbin
707  do i = 1, nbin
708  do nyu = 1, nspc
709  do myu = 1, nspc
710  if( vt( myu,i ) /= vt( nyu,j ) ) then
711  ecoll( myu,nyu,i,j ) = ck( myu,nyu,i,j ) &
712  / ( pi*( radc( i )+radc( j ) )**2 * abs( vt( myu,i )-vt( nyu,j ) ) )
713  ecoll( myu,nyu,i,j ) = max( min( 1.0_rp, ecoll( myu,nyu,i,j ) ),0.0_rp )
714  else
715  ecoll( myu,nyu,i,j ) = 0.0_rp
716  endif
717 
718  if( ecoal_gsi /= 0.0_rp ) then
719  ecoll( myu,nyu,i,j ) = ecoal_gsi
720  endif
721 
722  if( ecoll( myu,nyu,i,j ) /= 0.0_rp ) then
723  rcoll( myu,nyu,i,j ) = ( 1.0_rp-ecoll( myu,nyu,i,j ) ) / ecoll( myu,nyu,i,j )
724  elseif( ecoll( myu,nyu,i,j ) == 0.0_rp ) then
725  rcoll( myu,nyu,i,j ) = 1.0_rp
726  endif
727  enddo
728  enddo
729  enddo
730  enddo
731 
732  endif
733 
734  !--- aerosol ( CCN ) (not supported)
735  if ( nccn /= 0 ) then
736 
737  xasta = log( rhoa*4.0_rp/3.0_rp*pi * ( rasta )**3 )
738  xaend = log( rhoa*4.0_rp/3.0_rp*pi * ( raend )**3 )
739 
740  dxaer = ( xaend-xasta )/nccn
741 
742  !$omp parallel do
743  do n = 1, nccn+1
744  xabnd( n ) = xasta + dxaer*( n-1 )
745  enddo
746  !$omp parallel do
747  do n = 1, nccn
748  xactr( n ) = ( xabnd( n )+xabnd( n+1 ) )*0.50_rp
749  rada( n ) = ( exp( xactr( n ) )*thirdovforth/pi/rhoa )**( oneovthird )
750  end do
751  do n = 1, nccn
752  log_info("ATMOS_PHY_MP_suzuki10_setup",'(A,1x,I3,1x,A,1x,ES15.7,1x,A)') &
753  "Radius of ", n, "th aerosol bin (bin center)= ", rada( n ) , "[m]"
754  enddo
755 
756  if ( flg_sf_aero ) then
757  log_error("ATMOS_PHY_MP_suzuki10_setup",*) "flg_sf_aero=true is not supported stop!! "
758  call prc_abort
759 ! if ( CZ(KS) >= 10.0_RP ) then
760 ! R10M1 = 10.0_RP / CZ(KS) * 0.50_RP ! scale with height
761 ! R10M2 = 10.0_RP / CZ(KS) * 0.50_RP ! scale with height
762 ! R10H1 = 1.0_RP * 0.50_RP
763 ! R10H2 = 1.0_RP * 0.0_RP
764 ! R10E1 = 1.0_RP * 0.50_RP
765 ! R10E2 = 1.0_RP * 0.50_RP
766 ! K10_1 = KS
767 ! K10_2 = KS
768 ! else
769 ! k = 1
770 ! do while ( CZ(k) < 10.0_RP )
771 ! k = k + 1
772 ! K10_1 = k
773 ! K10_2 = k + 1
774 ! R10M1 = ( CZ(k+1) - 10.0_RP ) / CDZ(k)
775 ! R10M2 = ( 10.0_RP - CZ(k) ) / CDZ(k)
776 ! R10H1 = ( CZ(k+1) - 10.0_RP ) / CDZ(k)
777 ! R10H2 = ( 10.0_RP - CZ(k) ) / CDZ(k)
778 ! R10E1 = ( CZ(k+1) - 10.0_RP ) / CDZ(k)
779 ! R10E2 = ( 10.0_RP - CZ(k) ) / CDZ(k)
780 ! enddo
781 ! endif
782  endif
783 
784  endif
785 
786  !--- determine nbnd
787  do n = 1, nbin
788  if( radc( n ) > rbnd ) then
789  nbnd = n
790  exit
791  endif
792  enddo
793  log_info("ATMOS_PHY_MP_suzuki10_setup",'(A,ES15.7,A)') 'Radius between cloud and rain is ', radc(nbnd), '[m]'
795 
796  !--- random number setup for stochastic method
797  if ( flg_rndm ) then
798  call random_setup( ia*ja*ka )
799  endif
800 
801  if ( mp_couple_aerosol .AND. nccn /=0 ) then
802  log_error("ATMOS_PHY_MP_suzuki10_setup",*) 'nccn should be 0 when MP_couple_aerosol = .true. !! stop'
803  call prc_abort
804  endif
805 
806  if ( nccn /= 0 ) then
807  !$omp parallel do
808  do n = 1, nccn
809  expxactr( n ) = exp( xactr( n ) )
810  rexpxactr( n ) = 1.0_rp / exp( xactr( n ) )
811  enddo
812  !$omp parallel do
813  do n = 1, nccn+1
814  expxabnd( n ) = exp( xabnd( n ) )
815  rexpxabnd( n ) = 1.0_rp / exp( xabnd( n ) )
816  enddo
817  endif
818 
819  allocate( vterm(qa-1) )
820  !$omp workshare
821  vterm(:) = 0.0_rp
822  !$omp end workshare
823  !$omp parallel do
824  do myu = 1, nspc
825  do n = 1, nbin
826  vterm((myu-1)*nbin+n) = -vt( myu,n )
827  enddo
828  enddo
829  !$omp parallel do
830  do n = 1, nbin
831  expxctr( n ) = exp( xctr( n ) )
832  rexpxctr( n ) = 1.0_rp / exp( xctr( n ) )
833  enddo
834  !$omp parallel do
835  do n = 1, nbin+1
836  expxbnd( n ) = exp( xbnd( n ) )
837  rexpxbnd( n ) = 1.0_rp / exp( xbnd( n ) )
838  enddo
839 
840  allocate( kindx(nbin,nbin) )
841  call getrule( ifrsl,kindx )
842 
843  !--- determine the parameters for interpolating SDF from qxx, Nxx of parent domain
844  sigma_sdf(1) = 0.2_rp
845  sigma_sdf(2) = 0.35_rp
846  sigma_sdf(3) = 0.35_rp
847  sigma_sdf(4) = 0.35_rp
848  sigma_sdf(5) = 0.35_rp
849  r0_sdf(1) = 5.e-6_rp
850  r0_sdf(2) = 2.61e-6_rp
851  r0_sdf(3) = 5.e-6_rp
852  r0_sdf(4) = 2.61e-6_rp
853  r0_sdf(5) = 2.61e-6_rp ! to be corrected
854  n0_sdf(1) = 8.0e+6_rp
855  n0_sdf(2) = 0.0_rp
856  n0_sdf(3) = 3.0e+6_rp
857  n0_sdf(4) = 4.0e+6_rp
858  n0_sdf(5) = 4.0e+6_rp ! to be corrected
859  rho_sdf(1) = dwatr
860  rho_sdf(2) = dice
861  rho_sdf(3) = 100.0_rp
862  rho_sdf(4) = 400.0_rp
863  rho_sdf(5) = 400.0_rp ! to be corrected
864 
865  !--- determine the parameters for interpolating SDF from qxx, Nxx of parent domain
866  sigma_sdf(1) = 0.2_rp
867  sigma_sdf(2) = 0.35_rp
868  sigma_sdf(3) = 0.35_rp
869  sigma_sdf(4) = 0.35_rp
870  sigma_sdf(5) = 0.35_rp
871  r0_sdf(1) = 5.e-6_rp
872  r0_sdf(2) = 2.61e-6_rp
873  r0_sdf(3) = 5.e-6_rp
874  r0_sdf(4) = 2.61e-6_rp
875  r0_sdf(5) = 2.61e-6_rp ! to be corrected
876  n0_sdf(1) = 8.0e+6_rp
877  n0_sdf(2) = 0.0_rp
878  n0_sdf(3) = 3.0e+6_rp
879  n0_sdf(4) = 4.0e+6_rp
880  n0_sdf(5) = 4.0e+6_rp ! to be corrected
881  rho_sdf(1) = dwatr
882  rho_sdf(2) = dice
883  rho_sdf(3) = 100.0_rp
884  rho_sdf(4) = 400.0_rp
885  rho_sdf(5) = 400.0_rp ! to be corrected
886 
887  !$acc update device(vterm)
888 
889  return
890  end subroutine atmos_phy_mp_suzuki10_setup
891 
892  !-----------------------------------------------------------------------------
895 
899 
900  deallocate( xctr )
901  deallocate( xbnd )
902  deallocate( radc )
903  deallocate( cctr )
904  deallocate( cbnd )
905  deallocate( ck )
906  deallocate( vt )
907  deallocate( br )
908  deallocate( ifrsl )
909  deallocate( expxctr )
910  deallocate( expxbnd )
911  deallocate( rexpxctr )
912  deallocate( rexpxbnd )
913  if ( nccn /= 0 ) then
914  deallocate( xactr )
915  deallocate( xabnd )
916  deallocate( rada )
917  deallocate( expxactr )
918  deallocate( expxabnd )
919  deallocate( rexpxactr )
920  deallocate( rexpxabnd )
921  endif
922 
923  deallocate( flg_noninduct )
924  deallocate( ecoll )
925  deallocate( rcoll )
926 
927  deallocate( vterm )
928 
929  deallocate( kindx )
930 
931 
932  return
933  end subroutine atmos_phy_mp_suzuki10_finalize
934  !-----------------------------------------------------------------------------
936  subroutine atmos_phy_mp_suzuki10_tendency( &
937  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
938  KIJMAX, &
939  dt, &
940  DENS, PRES, TEMP, &
941  QTRC, QDRY, &
942  CPtot, CVtot, &
943  CCN, &
944  RHOQ_t, RHOE_t, &
945  CPtot_t, CVtot_t, &
946  EVAPORATE, &
947  flg_lt, &
948  d0_crg, v0_crg, &
949  dqcrg, &
950  beta_crg, &
951  QTRC_crg, &
952  QSPLT_in, Sarea, &
953  RHOC_t_mp )
954  use scale_const, only: &
955  tem00 => const_tem00, &
956  pi => const_pi
957  use scale_atmos_saturation, only: &
958  atmos_saturation_pres2qsat_liq, &
959  atmos_saturation_pres2qsat_ice
960  implicit none
961 
962  integer, intent(in) :: ka, ks, ke
963  integer, intent(in) :: ia, is, ie
964  integer, intent(in) :: ja, js, je
965  integer, intent(in) :: kijmax
966 
967  real(dp), intent(in) :: dt
968  real(rp), intent(in) :: dens (ka,ia,ja)
969  real(rp), intent(in) :: pres (ka,ia,ja)
970  real(rp), intent(in) :: temp (ka,ia,ja)
971  real(rp), intent(in) :: qtrc (ka,ia,ja,qa)
972  real(rp), intent(in) :: qdry (ka,ia,ja)
973  real(rp), intent(in) :: cptot(ka,ia,ja)
974  real(rp), intent(in) :: cvtot(ka,ia,ja)
975  real(rp), intent(in) :: ccn (ka,ia,ja)
976 
977  real(rp), intent(out) :: rhoq_t (ka,ia,ja,qa)
978  real(rp), intent(out) :: rhoe_t (ka,ia,ja)
979  real(rp), intent(out) :: cptot_t(ka,ia,ja)
980  real(rp), intent(out) :: cvtot_t(ka,ia,ja)
981  real(rp), intent(out) :: evaporate(ka,ia,ja) !--- number of evaporated cloud [/m3]
982 
983  ! Optional for Lightning
984  logical, intent(in), optional :: flg_lt
985  real(rp), intent(in), optional :: d0_crg, v0_crg
986  real(rp), intent(in), optional :: dqcrg(ka,ia,ja)
987  real(rp), intent(in), optional :: beta_crg(ka,ia,ja)
988  real(rp), intent(in), optional :: qtrc_crg(ka,ia,ja,num_hyd)
989  real(rp), intent(out), optional :: qsplt_in(ka,ia,ja,3)
990  real(rp), intent(out), optional :: sarea(ka,ia,ja,num_hyd)
991  real(rp), intent(out), optional :: rhoc_t_mp(ka,ia,ja,num_hyd)
992 
993  real(rp) :: qsat_l(ka,ia,ja)
994  real(rp) :: qsat_i(ka,ia,ja)
995  real(rp) :: ssliq(ka,ia,ja)
996  real(rp) :: ssice(ka,ia,ja)
997 
998  integer :: ijk_index (kijmax,3)
999  integer :: index_cld (kijmax)
1000  integer :: index_cold(kijmax)
1001  integer :: index_warm(kijmax)
1002  integer :: ijkcount, ijkcount_cold, ijkcount_warm
1003  integer :: ijk, indirect
1004 
1005  real(rp) :: dens_ijk(kijmax)
1006  real(rp) :: pres_ijk(kijmax)
1007  real(rp) :: temp_ijk(kijmax)
1008  real(rp) :: qdry_ijk(kijmax)
1009  real(rp) :: qvap_ijk(kijmax)
1010  real(rp) :: ccn_ijk(kijmax)
1011  real(rp) :: cp_ijk(kijmax)
1012  real(rp) :: cv_ijk(kijmax)
1013  real(rp) :: evaporate_ijk(kijmax)
1014  real(rp) :: ghyd_ijk(nbin,nspc,kijmax)
1015  real(rp) :: gaer_ijk(max(nccn,1),kijmax)
1016  real(rp) :: cldsum
1017  integer :: countbin
1018  real(rp) :: rhoq_new
1019 
1020  !--- for lithgning
1021  logical :: flg_lt_l
1022  real(rp) :: d0_crg_l, v0_crg_l
1023  real(rp) :: gcrg_ijk(nbin,nspc,kijmax)
1024  real(rp) :: crg_sep_ijk(nspc,kijmax)
1025  real(rp) :: dqcrg_ijk(kijmax)
1026  real(rp) :: beta_crg_ijk(kijmax)
1027 
1028  integer :: step
1029  integer :: k, i, j, m, n, iq
1030  !---------------------------------------------------------------------------
1031 
1032  if ( nspc == 1 ) then
1033  log_progress(*) 'atmosphere / physics / microphysics / SBM (Liquid water only)'
1034  elseif( nspc > 1 ) then
1035  log_progress(*) 'atmosphere / physics / microphysics / SBM (Mixed phase)'
1036  endif
1037 
1038  flg_lt_l = .false.
1039  if ( present(flg_lt) ) then
1040  flg_lt_l = flg_lt
1041  end if
1042 
1043  if ( flg_lt_l ) then
1044  !$omp workshare
1045  crg_sep_ijk(:,:) = 0.0_rp
1046  gcrg_ijk(:,:,:) = 0.0_rp
1047  qsplt_in(:,:,:,:) = 0.0_rp
1048  !$omp end workshare
1049  d0_crg_l = d0_crg
1050  v0_crg_l = v0_crg
1051  end if
1052 
1053  call atmos_saturation_pres2qsat_liq( ka, ks, ke, & ! [IN]
1054  ia, is, ie, & ! [IN]
1055  ja, js, je, & ! [IN]
1056  temp(:,:,:), pres(:,:,:), qdry(:,:,:), & ! [IN]
1057  qsat_l(:,:,:) ) ! [OUT]
1058 
1059  call atmos_saturation_pres2qsat_ice( ka, ks, ke, & ! [IN]
1060  ia, is, ie, & ! [IN]
1061  ja, js, je, & ! [IN]
1062  temp(:,:,:), pres(:,:,:), qdry(:,:,:), & ! [IN]
1063  qsat_i(:,:,:) ) ! [OUT]
1064 
1065  !$omp parallel do
1066  do j = js, je
1067  do i = is, ie
1068  do k = ks, ke
1069  ssliq(k,i,j) = qtrc(k,i,j,i_qv) / qsat_l(k,i,j) - 1.0_rp
1070  ssice(k,i,j) = qtrc(k,i,j,i_qv) / qsat_i(k,i,j) - 1.0_rp
1071  enddo
1072  enddo
1073  enddo
1074 
1075  if ( nspc == 1 ) then
1076  !$omp workshare
1077  ssice(:,:,:) = 0.0_rp
1078  !$omp end workshare
1079  endif
1080 
1081 !--- store initial SDF of aerosol
1082 !--- this option is not supported
1083 ! if ( ofirst_sdfa ) then
1084 ! allocate( marate( nccn ) )
1085 ! do j = JS, JE
1086 ! do i = IS, IE
1087 ! do k = KS, KE
1088 ! sum2 = 0.0_RP
1089 ! do n = 1, nccn
1090 ! marate( n ) = gdga(k,i,j,n)*rexpxactr( n )
1091 ! sum2 = sum2 + gdga(k,i,j,n)*rexpxactr( n )
1092 ! enddo
1093 ! enddo
1094 ! enddo
1095 ! enddo
1096 ! if ( sum2 /= 0.0_RP ) then
1097 ! marate( 1:nccn ) = marate( 1:nccn )/sum2
1098 ! ofirst_sdfa = .false.
1099 ! endif
1100 ! endif
1101 
1102  !--- Arrange array for microphysics
1103 
1104  call prof_rapstart('MP_ijkconvert', 3)
1105 
1106  !$omp parallel do private(ijk)
1107  do j = js, je
1108  do i = is, ie
1109  do k = ks, ke
1110  ijk = 1 + (k-ks) + (ke-ks+1) * ( (i-is) + (ie-is+1) * (j-js) )
1111  ijk_index(ijk,1) = i
1112  ijk_index(ijk,2) = j
1113  ijk_index(ijk,3) = k
1114  end do
1115  end do
1116  end do
1117 
1118  ijkcount = 0
1119  do j = js, je
1120  do i = is, ie
1121  do k = ks, ke
1122  ! calc total hydrometeors
1123  cldsum = 0.0_rp
1124  countbin = i_qv + 1
1125  do m = 1, nspc
1126  do n = 1, nbin
1127  cldsum = cldsum + qtrc(k,i,j,countbin) * dens(k,i,j) / dxmic
1128  countbin = countbin + 1
1129  enddo
1130  enddo
1131 
1132  if ( cldsum > cldmin &
1133  .OR. ssliq(k,i,j) > 0.0_rp &
1134  .OR. ssice(k,i,j) > 0.0_rp ) then
1135  ijkcount = ijkcount + 1
1136  ijk = 1 + (k-ks) + (ke-ks+1) * ( (i-is) + (ie-is+1) * (j-js) )
1137  index_cld(ijkcount) = ijk
1138 
1139  else
1140 
1141  ! no hudrometeors and undersaturation (no microphysical process occcurs)
1142  do iq = 1, qa
1143  rhoq_t(k,i,j,iq) = 0.0_rp
1144  end do
1145  rhoe_t(k,i,j) = 0.0_rp
1146  cptot_t(k,i,j) = 0.0_rp
1147  cvtot_t(k,i,j) = 0.0_rp
1148  evaporate(k,i,j) = 0.0_rp
1149 
1150  end if
1151 
1152  end do
1153  end do
1154  end do
1155 
1156  !$omp parallel do private(indirect,i,j,k,countbin)
1157  do ijk = 1, ijkcount
1158  indirect = index_cld(ijk)
1159  i = ijk_index(indirect,1)
1160  j = ijk_index(indirect,2)
1161  k = ijk_index(indirect,3)
1162 
1163  dens_ijk(ijk) = dens(k,i,j)
1164  pres_ijk(ijk) = pres(k,i,j)
1165  temp_ijk(ijk) = temp(k,i,j)
1166  qdry_ijk(ijk) = qdry(k,i,j)
1167  cp_ijk(ijk) = cptot(k,i,j)
1168  cv_ijk(ijk) = cvtot(k,i,j)
1169  ccn_ijk(ijk) = ccn(k,i,j)
1170  qvap_ijk(ijk) = qtrc(k,i,j,i_qv)
1171 
1172  countbin = i_qv + 1
1173  do m = 1, nspc
1174  do n = 1, nbin
1175  ghyd_ijk(n,m,ijk) = qtrc(k,i,j,countbin) * dens(k,i,j) / dxmic
1176  countbin = countbin + 1
1177  enddo
1178  enddo
1179 
1180  do n = 1, nccn
1181  gaer_ijk(n,ijk) = qtrc(k,i,j,countbin) * dens(k,i,j) / dxaer
1182  countbin = countbin + 1
1183  enddo
1184 
1185  end do
1186 
1187  ijkcount_cold = 0
1188  ijkcount_warm = 0
1189  if ( nspc > 1 ) then
1190  do ijk = 1, ijkcount
1191  if ( temp_ijk(ijk) < tem00 ) then ! cold
1192  ijkcount_cold = ijkcount_cold + 1
1193  index_cold(ijkcount_cold) = ijk
1194  else ! warm
1195  ijkcount_warm = ijkcount_warm + 1
1196  index_warm(ijkcount_warm) = ijk
1197  endif
1198  end do
1199  else ! warn
1200  !$omp parallel do
1201  do ijk = 1, ijkcount
1202  index_warm(ijk) = ijk
1203  end do
1204  ijkcount_warm = ijkcount
1205  end if
1206 
1207  if ( flg_lt_l ) then
1208  !$omp parallel do private(indirect,i,j,k,countbin)
1209  do ijk = 1, ijkcount
1210  indirect = index_cld(ijk)
1211  i = ijk_index(indirect,1)
1212  j = ijk_index(indirect,2)
1213  k = ijk_index(indirect,3)
1214 
1215  countbin = 1
1216  do m = 1, nspc
1217  do n = 1, nbin
1218  gcrg_ijk(n,m,ijk) = qtrc_crg(k,i,j,countbin) * dens(k,i,j)
1219  countbin = countbin + 1
1220  enddo
1221  enddo
1222  beta_crg_ijk(ijk) = beta_crg(k,i,j)
1223  dqcrg_ijk(ijk) = dqcrg(k,i,j)
1224  end do
1225  end if
1226 
1227  call prof_rapend ('MP_ijkconvert', 3)
1228 
1229  ! tentative timername registration
1230  call prof_rapstart('MP_suzuki10', 3)
1231  call prof_rapend ('MP_suzuki10', 3)
1232  call prof_rapstart('_SBM_Nucleat', 3)
1233  call prof_rapend ('_SBM_Nucleat', 3)
1234 ! call PROF_rapstart('_SBM_NucleatA', 3)
1235 ! call PROF_rapend ('_SBM_NucleatA', 3)
1236  call prof_rapstart('_SBM_Liqphase', 3)
1237  call prof_rapend ('_SBM_Liqphase', 3)
1238  call prof_rapstart('_SBM_Icephase', 3)
1239  call prof_rapend ('_SBM_Icephase', 3)
1240  call prof_rapstart('_SBM_Mixphase', 3)
1241  call prof_rapend ('_SBM_Mixphase', 3)
1242  call prof_rapstart('_SBM_AdvLiq', 3)
1243  call prof_rapend ('_SBM_AdvLiq', 3)
1244  call prof_rapstart('_SBM_AdvIce', 3)
1245  call prof_rapend ('_SBM_AdvIce', 3)
1246  call prof_rapstart('_SBM_AdvMix', 3)
1247  call prof_rapend ('_SBM_AdvMix', 3)
1248 ! call PROF_rapstart('_SBM_FAero', 3)
1249 ! call PROF_rapend ('_SBM_FAero', 3)
1250  call prof_rapstart('_SBM_Freezing', 3)
1251  call prof_rapend ('_SBM_Freezing', 3)
1252  call prof_rapstart('_SBM_IceNucleat', 3)
1253  call prof_rapend ('_SBM_IceNucleat', 3)
1254  call prof_rapstart('_SBM_Melting', 3)
1255  call prof_rapend ('_SBM_Melting', 3)
1256  call prof_rapstart('_SBM_CollCoag', 3)
1257  call prof_rapend ('_SBM_CollCoag', 3)
1258 ! call PROF_rapstart('_SBM_CollCoagR', 3)
1259 ! call PROF_rapend ('_SBM_CollCoagR', 3)
1260 
1261  if ( ijkcount > 0 ) then
1262 
1263  call prof_rapstart('MP_suzuki10', 3)
1264 
1265  call mp_suzuki10( ka, ia, ja, & ! [IN]
1266  ijkcount, & ! [IN]
1267  ijkcount_cold, & ! [IN]
1268  ijkcount_warm, & ! [IN]
1269  index_cold( 1:ijkcount), & ! [IN]
1270  index_warm( 1:ijkcount), & ! [IN]
1271  dens_ijk( 1:ijkcount), & ! [IN]
1272  pres_ijk( 1:ijkcount), & ! [IN]
1273  qdry_ijk( 1:ijkcount), & ! [IN]
1274  ccn_ijk( 1:ijkcount), & ! [IN]
1275  temp_ijk( 1:ijkcount), & ! [INOUT]
1276  qvap_ijk( 1:ijkcount), & ! [INOUT]
1277  ghyd_ijk(:,:,1:ijkcount), & ! [INOUT]
1278  gaer_ijk(:, 1:ijkcount), & ! [INOUT]
1279  cp_ijk( 1:ijkcount), & ! [INOUT]
1280  cv_ijk( 1:ijkcount), & ! [INOUT]
1281  evaporate_ijk(1:ijkcount), & ! [OUT]
1282  dt, & ! [IN]
1283  flg_lt_l, & ! [IN]
1284  d0_crg_l, v0_crg_l, & ! [IN]
1285  dqcrg_ijk( 1:ijkcount), & ! [IN]
1286  beta_crg_ijk( 1:ijkcount), & ! [IN]
1287  gcrg_ijk(:,:,1:ijkcount), & ! [INOUT]
1288  crg_sep_ijk(:,1:ijkcount) ) ! [OUT]
1289 
1290  call prof_rapend ('MP_suzuki10', 3)
1291 
1292 ! if ( flg_sf_aero ) then
1293 ! do j = JS-2, JE+2
1294 ! do i = IS-2, IE+1
1295 ! VELX(i,j) = MOMX(K10_1,i,j) / ( DENS(K10_1,i+1,j)+DENS(K10_1,i,j) ) * R10M1 &
1296 ! + MOMX(K10_2,i,j) / ( DENS(K10_2,i+1,j)+DENS(K10_2,i,j) ) * R10M2
1297 ! enddo
1298 ! enddo
1299 !
1300 ! do j = JS-2, JE+1
1301 ! do i = IS-2, IE+2
1302 ! VELY(i,j) = MOMY(K10_1,i,j) / ( DENS(K10_1,i,j+1)+DENS(K10_1,i,j) ) * R10M1 &
1303 ! + MOMY(K10_2,i,j) / ( DENS(K10_2,i,j+1)+DENS(K10_2,i,j) ) * R10M2
1304 ! enddo
1305 ! enddo
1306 ! endif
1307 !
1308 ! !--- SURFACE FLUX by Monahan et al. (1986)
1309 ! if ( flg_sf_aero .AND. nccn /= 0 ) then
1310 ! do j = JS, JE
1311 ! do i = IS, IE
1312 ! ijk = ( j - JS ) * KMAX * IMAX &
1313 ! + ( i - IS ) * KMAX &
1314 ! + ( KS - KS ) &
1315 ! + 1
1316 ! Uabs = sqrt( ( ( VELX(i,j) + VELX(i-1,j ) ) * 0.50_RP )**2 &
1317 ! + ( ( VELY(i,j) + VELY(i ,j-1) ) * 0.50_RP )**2 )
1318 ! do n = 1, nccn
1319 ! if ( rada( n ) <= 2.0E-5_RP .AND. rada( n ) >= 3.0E-7_RP ) then
1320 ! bparam = ( 0.38_RP - log( rada( n ) ) )/0.65_RP
1321 ! SFLX_AERO(i,j,n) = 1.373_RP * Uabs**( 3.41_RP ) * rada( n )**( -3.0_RP ) &
1322 ! * ( 1.0_RP + 0.057_RP * rada( n )**( 1.05_RP ) ) &
1323 ! * 10.0_RP**( 1.19_RP * exp( -bparam*bparam ) )
1324 ! ! convert from [#/m^2/um/s] -> [kg/m^3/unit log (m)]
1325 ! SFLX_AERO(i,j,n) = SFLX_AERO(i,j,n) / DENS(KS,i,j) &
1326 ! / CDZ(KS) * rada( n ) / 3.0_RP * dt * expxactr( n )
1327 ! Gaer_ijk(n,ijk) = Gaer_ijk(n,ijk) + SFLX_AERO(i,j,n)/dxaer
1328 ! endif
1329 ! enddo
1330 ! enddo
1331 ! enddo
1332 ! endif
1333 
1334  call prof_rapstart('MP_ijkconvert', 3)
1335 
1336  !---- return original array
1337  !$omp parallel do private(indirect,i,j,k,countbin,rhoq_new)
1338  do ijk = 1, ijkcount
1339  indirect = index_cld(ijk)
1340  i = ijk_index(indirect,1)
1341  j = ijk_index(indirect,2)
1342  k = ijk_index(indirect,3)
1343 
1344  rhoe_t(k,i,j) = ( temp_ijk(ijk) * cv_ijk(ijk) - temp(k,i,j) * cvtot(k,i,j) ) * dens(k,i,j) / dt
1345  cptot_t(k,i,j) = ( cp_ijk(ijk) - cptot(k,i,j) ) / dt
1346  cvtot_t(k,i,j) = ( cv_ijk(ijk) - cvtot(k,i,j) ) / dt
1347  evaporate(k,i,j) = evaporate_ijk(ijk) / dt ! [#/m3/s]
1348 
1349  rhoq_t(k,i,j,i_qv) = ( qvap_ijk(ijk) - qtrc(k,i,j,i_qv) ) * dens(k,i,j) / dt
1350 
1351  countbin = i_qv + 1
1352  do m = 1, nspc
1353  do n = 1, nbin
1354  rhoq_new = ghyd_ijk(n,m,ijk) * dxmic
1355  rhoq_t(k,i,j,countbin) = ( rhoq_new - qtrc(k,i,j,countbin)*dens(k,i,j) ) / dt
1356  countbin = countbin + 1
1357  enddo
1358  enddo
1359 
1360  do n = 1, nccn
1361  rhoq_new = gaer_ijk(n,ijk) * dxaer
1362  rhoq_t(k,i,j,countbin) = ( rhoq_new - qtrc(k,i,j,countbin)*dens(k,i,j) ) / dt
1363  countbin = countbin + 1
1364  enddo
1365 
1366  if( flg_lt_l ) then
1367  countbin = 1
1368  do m = 1, nspc
1369  do n = 1, nbin
1370  rhoq_new = gcrg_ijk(n,m,ijk)
1371  rhoc_t_mp(k,i,j,countbin) = ( rhoq_new - qtrc_crg(k,i,j,countbin)*dens(k,i,j) ) / dt
1372  countbin = countbin + 1
1373  enddo
1374  enddo
1375 
1376  sarea(k,i,j,:) = 0.0_rp
1377  countbin = 1
1378  do m = 1, nspc
1379  do n = 1, nbin
1380  rhoq_new = ghyd_ijk(n,m,ijk) * dxmic
1381  sarea(k,i,j,countbin) = 4.0_rp*pi*radc(n)**2*rhoq_new*rexpxctr( n )
1382  countbin = countbin + 1
1383  enddo
1384  enddo
1385 
1386  qsplt_in(k,i,j,1) = crg_sep_ijk(ig,ijk) * dens(k,i,j) / dt
1387  qsplt_in(k,i,j,2) = ( crg_sep_ijk(ic,ijk) &
1388  + crg_sep_ijk(ip,ijk) &
1389  + crg_sep_ijk(id,ijk) ) * dens(k,i,j) / dt
1390  qsplt_in(k,i,j,3) = crg_sep_ijk(iss,ijk) * dens(k,i,j) / dt
1391 
1392  endif
1393 
1394  enddo
1395 
1396 
1397 ! if ( nccn /= 0 ) then
1398 ! AMR(:,:,:) = 0.0_RP
1399 ! do j = JS, JE
1400 ! do i = IS, IE
1401 ! do k = KS, KE
1402 ! do n = 1, nccn
1403 ! AMR(k,i,j) = AMR(k,i,j) + QTRC(k,i,j,QQE-1+n)
1404 ! enddo
1405 ! enddo
1406 ! enddo
1407 ! enddo
1408 ! endif
1409 
1410  call prof_rapend ('MP_ijkconvert', 3)
1411 
1412  endif
1413 
1414  return
1415  end subroutine atmos_phy_mp_suzuki10_tendency
1416 
1417  !-----------------------------------------------------------------------------
1419 !OCL SERIAL
1421  KA, &
1422  vterm_o )
1423  !$acc routine vector
1424  implicit none
1425 
1426  integer, intent(in) :: ka
1427 
1428  real(rp), intent(out) :: vterm_o(ka,qa-1)
1429 
1430  integer :: iq
1431  !---------------------------------------------------------------------------
1432 
1433  do iq = 1, qa-1
1434  vterm_o(:,iq) = vterm(iq)
1435  end do
1436 
1437  return
1439 
1440  !-----------------------------------------------------------------------------
1443  KA, KS, KE, &
1444  IA, IS, IE, &
1445  JA, JS, JE, &
1446  QTRC0, &
1447  mask_criterion, &
1448  cldfrac )
1449  implicit none
1450 
1451  integer, intent(in) :: ka, ks, ke
1452  integer, intent(in) :: ia, is, ie
1453  integer, intent(in) :: ja, js, je
1454 
1455  real(rp), intent(in) :: qtrc0 (ka,ia,ja,num_hyd)
1456  real(rp), intent(in) :: mask_criterion
1457  real(rp), intent(out) :: cldfrac(ka,ia,ja)
1458 
1459  real(rp) :: qhydro
1460  integer :: k, i, j, iq, ihydro
1461  !---------------------------------------------------------------------------
1462 
1463  if( nspc > 1 ) then
1464  !$omp parallel do private(qhydro)
1465  do j = js, je
1466  do i = is, ie
1467  do k = ks, ke
1468  qhydro = 0.0_rp
1469  do ihydro = 1, nspc
1470  do iq = nbin*(ihydro-1)+1, nbin*ihydro
1471  qhydro = qhydro + qtrc0(k,i,j,iq)
1472  enddo
1473  enddo
1474  cldfrac(k,i,j) = 0.5_rp + sign(0.5_rp,qhydro-mask_criterion)
1475  enddo
1476  enddo
1477  enddo
1478  elseif( nspc == 1 ) then
1479  !$omp parallel do private(qhydro)
1480  do j = js, je
1481  do i = is, ie
1482  do k = ks, ke
1483  qhydro = 0.0_rp
1484  do ihydro = 1, i_mp_qc
1485  do iq = nbin*(ihydro-1)+1, nbin*ihydro
1486  qhydro = qhydro + qtrc0(k,i,j,iq)
1487  enddo
1488  enddo
1489  cldfrac(k,i,j) = 0.5_rp + sign(0.5_rp,qhydro-mask_criterion)
1490  enddo
1491  enddo
1492  enddo
1493  endif
1494 
1495  return
1497 
1498  !-----------------------------------------------------------------------------
1501  KA, KS, KE, &
1502  IA, IS, IE, &
1503  JA, JS, JE, &
1504  DENS0, &
1505  TEMP0, &
1506  QTRC0, &
1507  Re )
1508  use scale_const, only: &
1509  eps => const_eps
1510  use scale_atmos_hydrometeor, only: &
1511  n_hyd, &
1512  i_qv, &
1513  i_hc, &
1514  i_hr, &
1515  i_hi, &
1516  i_hs, &
1517  i_hg, &
1518  i_hh
1519  implicit none
1520 
1521  integer, intent(in) :: ka, ks, ke
1522  integer, intent(in) :: ia, is, ie
1523  integer, intent(in) :: ja, js, je
1524 
1525  real(rp), intent(in) :: dens0(ka,ia,ja) ! density [kg/m3]
1526  real(rp), intent(in) :: temp0(ka,ia,ja) ! temperature [K]
1527  real(rp), intent(in) :: qtrc0(ka,ia,ja,num_hyd) ! tracer mass concentration [kg/kg]
1528  real(rp), intent(out) :: re (ka,ia,ja,n_hyd) ! effective radius [cm]
1529 
1530  real(rp), parameter :: um2cm = 100.0_rp
1531 
1532  real(rp) :: sum0(nspc), sum2, sum3, re_tmp(nspc)
1533  integer :: i, j, k, iq, ihydro
1534  !---------------------------------------------------------------------------
1535 
1536  !$omp parallel do private(sum2,sum3,ihydro)
1537  do k = ks, ke
1538  do j = js, je
1539  do i = is, ie
1540  re(k,i,j,:) = 0.0_rp
1541 
1542  ! HC
1543  sum3 = 0.0_rp
1544  sum2 = 0.0_rp
1545  ihydro = i_mp_qc
1546  do iq = 1, nbnd
1547  sum3 = sum3 &
1548  + ( ( qtrc0(k,i,j,iq) * dens0(k,i,j) ) & !--- [kg/kg] -> [kg/m3]
1549  * rexpxctr( iq-nbin*(ihydro-1) ) & !--- mass -> number
1550  * radc( iq-nbin*(ihydro-1) )**3.0_rp )
1551  sum2 = sum2 &
1552  + ( ( qtrc0(k,i,j,iq) * dens0(k,i,j) ) & !--- [kg/kg] -> [kg/m3]
1553  * rexpxctr( iq-nbin*(ihydro-1) ) & !--- mass -> number
1554  * radc( iq-nbin*(ihydro-1) )**2.0_rp )
1555  enddo
1556  sum3 = max( sum3, 0.0_rp )
1557  sum2 = max( sum2, 0.0_rp )
1558  if ( sum2 /= 0.0_rp ) then
1559  re(k,i,j,i_hc) = sum3 / sum2 * um2cm
1560  else
1561  re(k,i,j,i_hc) = 0.0_rp
1562  endif
1563 
1564  ! HR
1565  sum3 = 0.0_rp
1566  sum2 = 0.0_rp
1567  ihydro = i_mp_qc
1568  do iq = nbnd+1, nbin
1569  sum3 = sum3 &
1570  + ( ( qtrc0(k,i,j,iq) * dens0(k,i,j) ) & !--- [kg/kg] -> [kg/m3]
1571  * rexpxctr( iq-nbin*(ihydro-1) ) & !--- mass -> number
1572  * radc( iq-nbin*(ihydro-1) )**3.0_rp )
1573  sum2 = sum2 &
1574  + ( ( qtrc0(k,i,j,iq) * dens0(k,i,j) ) & !--- [kg/kg] -> [kg/m3]
1575  * rexpxctr( iq-nbin*(ihydro-1) ) & !--- mass -> number
1576  * radc( iq-nbin*(ihydro-1) )**2.0_rp )
1577  enddo
1578  sum3 = max( sum3, 0.0_rp )
1579  sum2 = max( sum2, 0.0_rp )
1580  if ( sum2 /= 0.0_rp ) then
1581  re(k,i,j,i_hr) = sum3 / sum2 * um2cm
1582  else
1583  re(k,i,j,i_hr) = 0.0_rp
1584  endif
1585 
1586  enddo
1587  enddo
1588  enddo
1589 
1590  ! other hydrometeors
1591  if ( nspc > 1 ) then
1592  !$omp parallel do private(sum0,sum2,sum3,re_tmp)
1593  do k = ks, ke
1594  do j = js, je
1595  do i = is, ie
1596  do ihydro = 2, nspc
1597  sum0(ihydro) = 0.0_rp
1598  sum2 = 0.0_rp
1599  sum3 = 0.0_rp
1600  do iq = nbin*(ihydro-1)+1, nbin*ihydro
1601  sum0(ihydro) = sum0(ihydro) &
1602  + ( qtrc0(k,i,j,iq) * dens0(k,i,j) ) !--- [kg/kg] -> [kg/m3]
1603  sum3 = sum3 &
1604  + ( ( qtrc0(k,i,j,iq) * dens0(k,i,j) ) & !--- [kg/kg] -> [kg/m3]
1605  * rexpxctr( iq-nbin*(ihydro-1) ) & !--- mass -> number
1606  * radc( iq-nbin*(ihydro-1) )**3.0_rp )
1607  sum2 = sum2 &
1608  + ( ( qtrc0(k,i,j,iq) * dens0(k,i,j) ) & !--- [kg/kg] -> [kg/m3]
1609  * rexpxctr( iq-nbin*(ihydro-1) ) & !--- mass -> number
1610  * radc( iq-nbin*(ihydro-1) )**2.0_rp )
1611  enddo
1612  sum3 = max( sum3, 0.0_rp )
1613  sum2 = max( sum2, 0.0_rp )
1614  if ( sum2 == 0.0_rp ) then
1615  re_tmp(ihydro) = 0.0_rp
1616  else
1617  re_tmp(ihydro) = sum3 / sum2 * um2cm
1618  end if
1619  end do
1620 
1621  re(k,i,j,i_hi) = ( re_tmp(i_mp_qcl) * sum0(i_mp_qcl) &
1622  + re_tmp(i_mp_qp ) * sum0(i_mp_qp ) &
1623  + re_tmp(i_mp_qd ) * sum0(i_mp_qd ) ) &
1624  / ( sum0(i_mp_qcl) + sum0(i_mp_qp) + sum0(i_mp_qd) + eps )
1625  re(k,i,j,i_hs) = re_tmp(i_mp_qs)
1626  re(k,i,j,i_hg) = re_tmp(i_mp_qg)
1627  re(k,i,j,i_hh) = re_tmp(i_mp_qh)
1628 
1629  enddo
1630  enddo
1631  enddo
1632 
1633  end if
1634 
1635  return
1637 
1638  !-----------------------------------------------------------------------------
1640  subroutine atmos_phy_mp_suzuki10_qtrc2qhyd( &
1641  KA, KS, KE, &
1642  IA, IS, IE, &
1643  JA, JS, JE, &
1644  QTRC0, &
1645  Qe )
1647  n_hyd, &
1648  i_hc, &
1649  i_hr, &
1650  i_hi, &
1651  i_hs, &
1652  i_hg, &
1653  i_hh
1654  implicit none
1655 
1656  integer, intent(in) :: ka, ks, ke
1657  integer, intent(in) :: ia, is, ie
1658  integer, intent(in) :: ja, js, je
1659  real(rp), intent(in) :: qtrc0(ka,ia,ja,num_hyd) ! tracer mass concentration [kg/kg]
1660  real(rp), intent(out) :: qe (ka,ia,ja,n_hyd) ! mixing ratio of each cateory [kg/kg]
1661 
1662  integer :: ihydro, ibin, iq, icateg
1663  integer :: k, i, j
1664  !---------------------------------------------------------------------------
1665 
1666 !OCL XFILL
1667  !$omp workshare
1668  qe(:,:,:,:) = 0.0_rp
1669  !$omp end workshare
1670 
1671  do ihydro = 1, nspc
1672  do ibin = 1, nbin
1673  iq = nbin*(ihydro-1) + ibin
1674 
1675  if ( iq > 0 .AND. iq <= nbin*(i_mp_qc ) ) then ! liquid
1676  if ( iq > 0 .AND. iq <= nbnd ) then ! cloud
1677  icateg = i_hc
1678  elseif( iq > nbnd .AND. iq <= nbin ) then ! rain
1679  icateg = i_hr
1680  endif
1681  elseif( iq > nbin*(i_mp_qcl-1) .AND. iq <= nbin*(i_mp_qcl) ) then ! ice (column)
1682  icateg = i_hi
1683  elseif( iq > nbin*(i_mp_qp -1) .AND. iq <= nbin*(i_mp_qp ) ) then ! ice (plate)
1684  icateg = i_hi
1685  elseif( iq > nbin*(i_mp_qd -1) .AND. iq <= nbin*(i_mp_qd ) ) then ! ice (dendrite)
1686  icateg = i_hi
1687  elseif( iq > nbin*(i_mp_qs -1) .AND. iq <= nbin*(i_mp_qs ) ) then ! snow
1688  icateg = i_hs
1689  elseif( iq > nbin*(i_mp_qg -1) .AND. iq <= nbin*(i_mp_qg ) ) then ! graupel
1690  icateg = i_hg
1691  elseif( iq > nbin*(i_mp_qh -1) .AND. iq <= num_hyd ) then ! hail
1692  icateg = i_hh
1693  endif
1694 
1695  !$omp parallel do collapse(2)
1696  do j = js, je
1697  do i = is, ie
1698  do k = ks, ke
1699  qe(k,i,j,icateg) = qe(k,i,j,icateg) + qtrc0(k,i,j,iq)
1700  enddo
1701  enddo
1702  enddo
1703  enddo
1704  enddo
1705 
1706  return
1707  end subroutine atmos_phy_mp_suzuki10_qtrc2qhyd
1708 
1710  subroutine atmos_phy_mp_suzuki10_qtrc2nhyd( &
1711  KA, KS, KE, &
1712  IA, IS, IE, &
1713  JA, JS, JE, &
1714  DENS, &
1715  QTRC0, &
1716  Ne )
1718  n_hyd, &
1719  i_hc, &
1720  i_hr, &
1721  i_hi, &
1722  i_hs, &
1723  i_hg, &
1724  i_hh
1725  implicit none
1726 
1727  integer, intent(in) :: ka, ks, ke
1728  integer, intent(in) :: ia, is, ie
1729  integer, intent(in) :: ja, js, je
1730  real(rp), intent(in) :: dens (ka,ia,ja) ! density [kg/m3]
1731  real(rp), intent(in) :: qtrc0(ka,ia,ja,num_hyd) ! tracer mass concentration [kg/kg]
1732  real(rp), intent(out) :: ne (ka,ia,ja,n_hyd) ! number concentration of each cateory [1/m3]
1733 
1734  integer :: ihydro, ibin, iq, icateg
1735  integer :: k, i, j
1736  !---------------------------------------------------------------------------
1737 
1738 !OCL XFILL
1739  !$omp workshare
1740  ne(:,:,:,:) = 0.0_rp
1741  !$omp end workshare
1742 
1743  do ihydro = 1, nspc
1744  do ibin = 1, nbin
1745  iq = nbin*(ihydro-1) + ibin
1746 
1747  if ( iq > 0 .AND. iq <= nbin*(i_mp_qc ) ) then ! liquid
1748  if ( iq > 0 .AND. iq <= nbnd ) then ! cloud
1749  icateg = i_hc
1750  elseif( iq > nbnd .AND. iq <= nbin ) then ! rain
1751  icateg = i_hr
1752  endif
1753  elseif( iq > nbin*(i_mp_qcl-1) .AND. iq <= nbin*(i_mp_qcl) ) then ! ice (column)
1754  icateg = i_hi
1755  elseif( iq > nbin*(i_mp_qp -1) .AND. iq <= nbin*(i_mp_qp ) ) then ! ice (plate)
1756  icateg = i_hi
1757  elseif( iq > nbin*(i_mp_qd -1) .AND. iq <= nbin*(i_mp_qd ) ) then ! ice (dendrite)
1758  icateg = i_hi
1759  elseif( iq > nbin*(i_mp_qs -1) .AND. iq <= nbin*(i_mp_qs ) ) then ! snow
1760  icateg = i_hs
1761  elseif( iq > nbin*(i_mp_qg -1) .AND. iq <= nbin*(i_mp_qg ) ) then ! graupel
1762  icateg = i_hg
1763  elseif( iq > nbin*(i_mp_qh -1) .AND. iq <= num_hyd ) then ! hail
1764  icateg = i_hh
1765  endif
1766 
1767  !$omp parallel do collapse(2)
1768  do j = js, je
1769  do i = is, ie
1770  do k = ks, ke
1771  ne(k,i,j,icateg) = ne(k,i,j,icateg) + dens(k,i,j) * qtrc0(k,i,j,iq) * rexpxctr(ibin)
1772  enddo
1773  enddo
1774  enddo
1775  enddo
1776  enddo
1777 
1778  return
1779  end subroutine atmos_phy_mp_suzuki10_qtrc2nhyd
1780 
1781  !-----------------------------------------------------------------------------
1783  subroutine atmos_phy_mp_suzuki10_qhyd2qtrc( &
1784  KA, KS, KE, IA, IS, IE, JA, JS, JE, &
1785  Qe, &
1786  QTRC, &
1787  QNUM )
1788  use scale_const, only: &
1789  pi => const_pi, &
1790  eps => const_eps
1791  use scale_atmos_hydrometeor, only: &
1792  n_hyd, &
1793  i_hc, &
1794  i_hr, &
1795  i_hi, &
1796  i_hs, &
1797  i_hg, &
1798  i_hh
1799  use scale_specfunc, only: &
1800  sf_gamma
1801  implicit none
1802  integer, intent(in) :: ka, ks, ke
1803  integer, intent(in) :: ia, is, ie
1804  integer, intent(in) :: ja, js, je
1805 
1806  real(rp), intent(in) :: qe(ka,ia,ja,n_hyd) ! mass ratio of each cateory [kg/kg]
1807 
1808  real(rp), intent(out) :: qtrc(ka,ia,ja,qa-1)
1809 
1810  real(rp), intent(in), optional :: qnum(ka,ia,ja,n_hyd) ! number concentratio
1811 
1812  real(rp) :: coef0, coef1, coef2
1813  real(rp) :: dummy(nbin)
1814  real(rp) :: tmp_hyd, num_hyd_l, lambda_hyd
1815 
1816  integer :: k, i, j, iq
1817 
1818  if ( present(qnum) ) then
1819  log_warn("ATMOS_PHY_MP_suzuki10_qhyd2qtrc",*) 'At this moment, number concentratio is ignored'
1820  end if
1821 
1822  !--- define coefficients
1823  coef0 = 4.0_rp/3.0_rp*pi
1824  coef1 = 4.0_rp/3.0_rp*sqrt(pi/2.0_rp)
1825 
1826  if( nspc == 1 ) then !--- put all hydrometeors to liquid (warm bin)
1827 
1828  !$omp parallel do private(tmp_hyd,dummy,coef2)
1829  do j = js, je
1830  do i = is, ie
1831  do k = ks, ke
1832 
1833  tmp_hyd = 0.0_rp
1834  do iq = 1, nbin
1835  dummy(iq) = coef1 / sigma_sdf(1) * rho_sdf(1) * radc( iq )**3 &
1836  * exp( &
1837  - ( log( radc(iq) )-log( r0_sdf(1) ) )**2*0.5_rp &
1838  / sigma_sdf(1) / sigma_sdf(1) &
1839  )
1840  tmp_hyd = tmp_hyd + dummy(iq)
1841  enddo
1842 
1843  coef2 = ( qe(k,i,j,i_hc) + qe(k,i,j,i_hr) &
1844  + qe(k,i,j,i_hi) + qe(k,i,j,i_hs) + qe(k,i,j,i_hg) ) &
1845  / ( tmp_hyd + ( 0.50_rp - sign(0.50_rp,tmp_hyd-eps) ) )
1846 
1847  do iq = 1, nbin
1848  qtrc(k+2,i,j,(il-1)*nbin+iq) = coef2 * dummy(iq)
1849  enddo
1850 
1851  enddo
1852  enddo
1853  enddo
1854 
1855  elseif( nspc > 1 ) then !--- put each hydrometeor to each category (ice bin)
1856 
1857  !$omp parallel do private(tmp_hyd,dummy,coef2,num_hyd_l,lambda_hyd)
1858  do j = js, je
1859  do i = is, ie
1860  do k = ks, ke
1861 
1862  !--- Rain and Cloud put into liquid bin (log-normal)
1863  tmp_hyd = 0.0_rp
1864  do iq = 1, nbin
1865  dummy(iq) = coef1 / sigma_sdf(1) * rho_sdf(1) * radc( iq )**3 &
1866  * exp( &
1867  - ( log( radc(iq) )-log( r0_sdf(1) ) )**2*0.5_rp &
1868  / sigma_sdf(1) / sigma_sdf(1) &
1869  )
1870  tmp_hyd = tmp_hyd + dummy(iq)
1871  enddo
1872 
1873  coef2 = ( qe(k,i,j,i_hc) + qe(k,i,j,i_hr) ) &
1874  / ( tmp_hyd + ( 0.50_rp - sign(0.50_rp,tmp_hyd-eps) ) )
1875 
1876  do iq = 1, nbin
1877  qtrc(k,i,j,(il-1)*nbin+iq) = coef2 * dummy(iq)
1878  enddo
1879 
1880  !--- Ice put into plate bin (log-normal)
1881  tmp_hyd = 0.0_rp
1882  do iq = 1, nbin
1883  dummy(iq) = coef1 / sigma_sdf(2) * rho_sdf(2) * radc( iq )**3 &
1884  * exp( &
1885  - ( log( radc(iq) )-log( r0_sdf(2) ) )**2*0.5_rp &
1886  / sigma_sdf(2) / sigma_sdf(2) &
1887  )
1888  tmp_hyd = tmp_hyd + dummy(iq)
1889  enddo
1890 
1891  coef2 = ( qe(k,i,j,i_hi) ) &
1892  / ( tmp_hyd + ( 0.50_rp - sign(0.50_rp,tmp_hyd-eps) ) )
1893 
1894  do iq = 1, nbin
1895  qtrc(k,i,j,(ic-1)*nbin+iq) = 0.0_rp
1896  qtrc(k,i,j,(ip-1)*nbin+iq) = coef2 * dummy(iq)
1897  qtrc(k,i,j,(id-1)*nbin+iq) = 0.0_rp
1898  enddo
1899 
1900  !--- Snow put into snow bin (gamma)
1901  num_hyd_l = coef0 * n0_sdf(3) * rho_sdf(3)
1902  lambda_hyd = ( pi * rho_sdf(3) / 6.0_rp *n0_sdf(3) * sf_gamma(4.0_rp) &
1903  / ( qe(k,i,j,i_hs) &
1904  + (0.50_rp-sign(0.50_rp,qe(k,i,j,i_hs)-eps)) &
1905  ) )**(0.25_rp)
1906 
1907  tmp_hyd = 0.0_rp
1908  do iq = 1, nbin
1909  dummy(iq) = num_hyd_l * radc( iq )**3 &
1910  * exp( -lambda_hyd * 0.5_rp * radc( iq ) )
1911  tmp_hyd = tmp_hyd + dummy(iq)
1912  enddo
1913 
1914  coef2 = ( qe(k,i,j,i_hs) ) &
1915  / ( tmp_hyd + ( 0.50_rp - sign(0.50_rp,tmp_hyd-eps) ) )
1916 
1917  do iq = 1, nbin
1918  qtrc(k,i,j,(iss-1)*nbin+iq) = coef2 * dummy(iq)
1919  enddo
1920 
1921  !--- Graupel put into Graupel bin (gamma)
1922  num_hyd_l = coef0 * n0_sdf(4) * rho_sdf(4)
1923  lambda_hyd = ( pi * rho_sdf(4) / 6.0_rp *n0_sdf(4) * sf_gamma(4.0_rp) &
1924  / ( qe(k,i,j,i_hg) &
1925  + (0.50_rp-sign(0.50_rp,qe(k,i,j,i_hg)-eps)) &
1926  ) )**(0.25_rp)
1927 
1928  tmp_hyd = 0.0_rp
1929  do iq = 1, nbin
1930  dummy(iq) = num_hyd_l * radc( iq )**3 &
1931  * exp( -lambda_hyd * 0.5_rp * radc( iq ) )
1932  tmp_hyd = tmp_hyd + dummy(iq)
1933  enddo
1934 
1935  coef2 = ( qe(k,i,j,i_hg) ) &
1936  / ( tmp_hyd + ( 0.50_rp - sign(0.50_rp,tmp_hyd-eps) ) )
1937 
1938  do iq = 1, nbin
1939  qtrc(k,i,j,(ig-1)*nbin+iq) = coef2 * dummy(iq)
1940  enddo
1941 
1942  !--- Hail put into Hail bin (gamma)
1943  num_hyd_l = coef0 * n0_sdf(5) * rho_sdf(5)
1944  lambda_hyd = ( pi * rho_sdf(5) / 6.0_rp *n0_sdf(5) * sf_gamma(4.0_rp) &
1945  / ( qe(k,i,j,i_hh) &
1946  + (0.50_rp-sign(0.50_rp,qe(k,i,j,i_hh)-eps)) &
1947  ) )**(0.25_rp)
1948 
1949  tmp_hyd = 0.0_rp
1950  do iq = 1, nbin
1951  dummy(iq) = num_hyd_l * radc( iq )**3 &
1952  * exp( -lambda_hyd * 0.5_rp * radc( iq ) )
1953  tmp_hyd = tmp_hyd + dummy(iq)
1954  enddo
1955 
1956  coef2 = ( qe(k,i,j,i_hh) ) &
1957  / ( tmp_hyd + ( 0.50_rp - sign(0.50_rp,tmp_hyd-eps) ) )
1958 
1959  do iq = 1, nbin
1960  qtrc(k,i,j,(ih-1)*nbin+iq) = coef2 * dummy(iq)
1961  enddo
1962 
1963  enddo
1964  enddo
1965  enddo
1966 
1967  endif
1968 
1969  !$omp parallel do collapse(3)
1970  do iq = num_hyd+1, qa-1
1971  do j = js, je
1972  do i = is, ie
1973  do k = ks, ke
1974  qtrc(k,i,j,iq) = 0.0_rp
1975  end do
1976  end do
1977  end do
1978  end do
1979 
1980  return
1981  end subroutine atmos_phy_mp_suzuki10_qhyd2qtrc
1982 
1983  !-----------------------------------------------------------------------------
1986  KA, KS, KE, &
1987  IA, IS, IE, &
1988  JA, JS, JE, &
1989  QTRC0, &
1990  Qecrg )
1992  n_hyd, &
1993  i_hc, &
1994  i_hr, &
1995  i_hi, &
1996  i_hs, &
1997  i_hg, &
1998  i_hh
1999  implicit none
2000 
2001  integer, intent(in) :: ka, ks, ke
2002  integer, intent(in) :: ia, is, ie
2003  integer, intent(in) :: ja, js, je
2004  real(rp), intent(in) :: qtrc0(ka,ia,ja,num_hyd) ! tracer charge density [fC/kg]
2005  real(rp), intent(out) :: qecrg(ka,ia,ja,n_hyd) ! charge density ratio of each cateory [fC/kg]
2006 
2007  integer :: ihydro, ibin, iq, icateg
2008  integer :: k, i, j
2009  !---------------------------------------------------------------------------
2010 
2011 !OCL XFILL
2012  !$omp workshare
2013  qecrg(:,:,:,:) = 0.0_rp
2014  !$omp end workshare
2015 
2016  do ihydro = 1, nspc
2017  do ibin = 1, nbin
2018  iq = nbin*(ihydro-1) + ibin
2019 
2020  if ( iq > 0 .AND. iq <= nbin*(i_mp_qc ) ) then ! liquid
2021  if ( iq > 0 .AND. iq <= nbnd ) then ! cloud
2022  icateg = i_hc
2023  elseif( iq > nbnd .AND. iq <= nbin ) then ! rain
2024  icateg = i_hr
2025  endif
2026  elseif( iq > nbin*(i_mp_qcl-1) .AND. iq <= nbin*(i_mp_qcl) ) then ! ice (column)
2027  icateg = i_hi
2028  elseif( iq > nbin*(i_mp_qp -1) .AND. iq <= nbin*(i_mp_qp ) ) then ! ice (plate)
2029  icateg = i_hi
2030  elseif( iq > nbin*(i_mp_qd -1) .AND. iq <= nbin*(i_mp_qd ) ) then ! ice (dendrite)
2031  icateg = i_hi
2032  elseif( iq > nbin*(i_mp_qs -1) .AND. iq <= nbin*(i_mp_qs ) ) then ! snow
2033  icateg = i_hs
2034  elseif( iq > nbin*(i_mp_qg -1) .AND. iq <= nbin*(i_mp_qg ) ) then ! graupel
2035  icateg = i_hg
2036  elseif( iq > nbin*(i_mp_qh -1) .AND. iq <= num_hyd ) then ! hail
2037  icateg = i_hh
2038  endif
2039 
2040  !$omp parallel do collapse(2)
2041  do j = js, je
2042  do i = is, ie
2043  do k = ks, ke
2044  qecrg(k,i,j,icateg) = qecrg(k,i,j,icateg) + qtrc0(k,i,j,iq)
2045  enddo
2046  enddo
2047  enddo
2048  enddo
2049  enddo
2050 
2051  return
2053 
2054  !-----------------------------------------------------------------------------
2055  subroutine mp_suzuki10( &
2056  KA, IA, JA, &
2057  ijkmax, &
2058  ijkmax_cold, &
2059  ijkmax_warm, &
2060  index_cold, &
2061  index_warm, &
2062  dens, &
2063  pres, &
2064  qdry, &
2065  ccn, &
2066  temp, &
2067  qvap, &
2068  ghyd, &
2069  gaer, &
2070  cp, &
2071  cv, &
2072  evaporate, &
2073  dt, &
2074  flg_lt, &
2075  d0_crg, &
2076  v0_crg, &
2077  dqcrg, &
2078  beta_crg, &
2079  gcrg, &
2080  crg_sep )
2081  use scale_const, only: &
2082  temp00 => const_tem00
2083  implicit none
2084 
2085  integer, intent(in) :: ka
2086  integer, intent(in) :: ia
2087  integer, intent(in) :: ja
2088 
2089  integer, intent(in) :: ijkmax
2090  integer, intent(in) :: ijkmax_cold
2091  integer, intent(in) :: ijkmax_warm
2092  integer , intent(in) :: index_cold(ijkmax)
2093  integer , intent(in) :: index_warm(ijkmax)
2094  real(rp), intent(in) :: dens (ijkmax) ! Density [kg/m3]
2095  real(rp), intent(in) :: pres (ijkmax) ! Pressure [Pa]
2096  real(rp), intent(in) :: qdry (ijkmax) ! dry air mass ratio [kg/kg]
2097  real(rp), intent(in) :: ccn (ijkmax) ! Number concentration of CCN [#/m3]
2098  real(rp), intent(inout) :: temp (ijkmax) ! Temperature [K]
2099  real(rp), intent(inout) :: qvap (ijkmax) ! Specific humidity [kg/kg]
2100  real(rp), intent(inout) :: ghyd (nbin,nspc,ijkmax) ! Mass size distribution function of hydrometeor
2101  real(rp), intent(inout) :: gaer (max(nccn,1),ijkmax) ! Mass size distribution function of aerosol
2102  real(rp), intent(inout) :: cp (ijkmax) ! specific heat
2103  real(rp), intent(inout) :: cv (ijkmax) ! specific heat
2104  real(rp), intent(out) :: evaporate (ijkmax) ! Number concentration of evaporated cloud [/m3/s]
2105  real(dp), intent(in) :: dt ! Time step interval
2106 
2107  ! Optional for Lightning
2108  logical, intent(in) :: flg_lt
2109  real(rp), intent(in) :: d0_crg, v0_crg
2110  real(rp), intent(in) :: dqcrg(ijkmax), beta_crg(ijkmax)
2111  real(rp), intent(inout) :: gcrg(nbin,nspc,ijkmax)
2112  real(rp), intent(out) :: crg_sep(nspc,ijkmax)
2113  !--- local
2114  integer :: m, n
2115  !---------------------------------------------------------------------------
2116 
2117  if( flg_lt ) then
2118  !$omp workshare
2119  crg_sep(:,:) = 0.0_rp
2120  !$omp end workshare
2121  endif
2122 
2123  if ( nccn /= 0 ) then
2124  if ( nspc == 1 ) then
2125  !---< warm rain only with aerosol tracer >---
2126 
2127  ! nucleation from aerosol
2128  call nucleata( ijkmax, & ! [IN]
2129  dens(:), & ! [IN]
2130  pres(:), & ! [IN]
2131  qdry(:), & ! [IN]
2132  temp(:), & ! [INOUT]
2133  qvap(:), & ! [INOUT]
2134  ghyd(:,:,:), & ! [INOUT]
2135  gaer(:,:), & ! [INOUT]
2136  cp(:), & ! [INOUT]
2137  cv(:), & ! [INOUT]
2138  dt ) ! [IN]
2139 
2140  ! condensation / evaporation
2141  call cndevpsbla( ijkmax, & ! [IN]
2142  flg_lt, & ! [IN]
2143  dens(:), & ! [IN]
2144  pres(:), & ! [IN]
2145  qdry(:), & ! [IN]
2146  temp(:), & ! [INOUT]
2147  qvap(:), & ! [INOUT]
2148  ghyd(:,:,:), & ! [INOUT]
2149  gaer(:,:), & ! [INOUT]
2150  cp(:), & ! [INOUT]
2151  cv(:), & ! [INOUT]
2152  evaporate(:), & ! [OUT]
2153  dt, & ! [IN]
2154  gcrg(:,:,:) ) ! [INOUT]
2155 
2156  if ( mp_doautoconversion ) then
2157  ! collision-coagulation
2158 
2159  call collmain( ka, ia, ja, & ! [IN]
2160  ijkmax, & ! [IN]
2161  flg_lt, & ! [IN]
2162  d0_crg, & ! [IN]
2163  v0_crg, & ! [IN]
2164  dqcrg(:), & ! [IN]
2165  beta_crg(:), & ! [IN]
2166  temp(:), & ! [IN]
2167  ghyd(:,:,:), & ! [INOUT]
2168  gcrg(:,:,:), & ! [INOUT]
2169  crg_sep(:,:), & ! [OUT]
2170  dt ) ! [IN]
2171 
2172  endif
2173 
2174  elseif( nspc > 1 ) then
2175  !---< mixed phase with aerosol tracer >---
2176 
2177  ! nucleation from aerosol
2178  call nucleata( ijkmax, & ! [IN]
2179  dens(:), & ! [IN]
2180  pres(:), & ! [IN]
2181  qdry(:), & ! [IN]
2182  temp(:), & ! [INOUT]
2183  qvap(:), & ! [INOUT]
2184  ghyd(:,:,:), & ! [INOUT]
2185  gaer(:,:), & ! [INOUT]
2186  cp(:), & ! [INOUT]
2187  cv(:), & ! [INOUT]
2188  dt ) ! [IN]
2189 
2190  ! freezing / melting
2191  call freezing( ijkmax, & ! [IN]
2192  ijkmax_cold, & ! [IN]
2193  index_cold(:), & ! [IN]
2194  flg_lt, & ! [IN]
2195  dens(:), & ! [IN]
2196  temp(:), & ! [INOUT]
2197  ghyd(:,:,:), & ! [INOUT]
2198  gcrg(:,:,:), & ! [INOUT]
2199  cp(:), & ! [INOUT]
2200  cv(:), & ! [INOUT]
2201  dt ) ! [IN]
2202 
2203  if( flg_icenucl ) then
2204  call ice_nucleat( ijkmax, & ! [IN]
2205  ijkmax_cold, & ! [IN]
2206  index_cold(:), & ! [IN]
2207  dens(:), & ! [IN]
2208  pres(:), & ! [IN]
2209  qdry(:), & ! [IN]
2210  temp(:), & ! [INOUT]
2211  qvap(:), & ! [INOUT]
2212  ghyd(:,:,:), & ! [INOUT]
2213  cp(:), & ! [INOUT]
2214  cv(:), & ! [INOUT]
2215  dt ) ! [IN]
2216  endif
2217 
2218  call melting( ijkmax, & ! [IN]
2219  ijkmax_warm, & ! [IN]
2220  index_warm(:), & ! [IN]
2221  flg_lt, & ! [IN]
2222  dens(:), & ! [IN]
2223  temp(:), & ! [INOUT]
2224  ghyd(:,:,:), & ! [INOUT]
2225  gcrg(:,:,:), & ! [INOUT]
2226  cp(:), & ! [INOUT]
2227  cv(:), & ! [INOUT]
2228  dt ) ! [IN]
2229 
2230  ! condensation / evaporation
2231  call cndevpsbla( ijkmax, & ! [IN]
2232  flg_lt, & ! [IN]
2233  dens(:), & ! [IN]
2234  pres(:), & ! [IN]
2235  qdry(:), & ! [IN]
2236  temp(:), & ! [INOUT]
2237  qvap(:), & ! [INOUT]
2238  ghyd(:,:,:), & ! [INOUT]
2239  gaer(:,:), & ! [INOUT]
2240  cp(:), & ! [INOUT]
2241  cv(:), & ! [INOUT]
2242  evaporate(:), & ! [OUT]
2243  dt, & ! [IN]
2244  gcrg(:,:,:) ) ! [INOUT]
2245 
2246  if ( mp_doautoconversion ) then
2247  ! collision-coagulation
2248  call collmainf( ka, ia, ja, & ! [IN]
2249  ijkmax, & ! [IN]
2250  flg_lt, & ! [IN]
2251  d0_crg, & ! [IN]
2252  v0_crg, & ! [IN]
2253  dqcrg(:), & ! [IN]
2254  beta_crg(:), & ! [IN]
2255  temp(:), & ! [IN]
2256  ghyd(:,:,:), & ! [INOUT]
2257  gcrg(:,:,:), & ! [INOUT]
2258  crg_sep(:,:), & ! [OUT]
2259  dt ) ! [IN]
2260 
2261  endif
2262 
2263  endif
2264 
2265  elseif( nccn == 0 ) then
2266 
2267  if ( nspc == 1 ) then
2268  !---< warm rain only without aerosol tracer >---
2269 
2270  ! nucleation
2271  call nucleat( ijkmax, & ! [IN]
2272  dens(:), & ! [IN]
2273  pres(:), & ! [IN]
2274  qdry(:), & ! [IN]
2275  ccn(:), & ! [IN]
2276  temp(:), & ! [INOUT]
2277  qvap(:), & ! [INOUT]
2278  ghyd(:,:,:), & ! [INOUT]
2279  cp(:), & ! [INOUT]
2280  cv(:), & ! [INOUT]
2281  dt ) ! [IN]
2282 
2283  ! condensation / evaporation
2284  call cndevpsbl( ijkmax, & ! [IN]
2285  flg_lt, & ! [IN]
2286  dens(:), & ! [IN]
2287  pres(:), & ! [IN]
2288  qdry(:), & ! [IN]
2289  temp(:), & ! [INOUT]
2290  qvap(:), & ! [INOUT]
2291  ghyd(:,:,:), & ! [INOUT]
2292  cp(:), & ! [INOUT]
2293  cv(:), & ! [INOUT]
2294  evaporate(:), & ! [OUT]
2295  dt, & ! [IN]
2296  gcrg(:,:,:) ) ! [INOUT]
2297 
2298  if ( mp_doautoconversion ) then
2299  ! collision-coagulation
2300  call collmain( ka, ia, ja, & ! [IN]
2301  ijkmax, & ! [IN]
2302  flg_lt, & ! [IN]
2303  d0_crg, & ! [IN]
2304  v0_crg, & ! [IN]
2305  dqcrg(:), & ! [IN]
2306  beta_crg(:), & ! [IN]
2307  temp(:), & ! [IN]
2308  ghyd(:,:,:), & ! [INOUT]
2309  gcrg(:,:,:), & ! [INOUT]
2310  crg_sep(:,:), & ! [OUT]
2311  dt ) ! [IN]
2312 
2313  endif
2314 
2315  elseif( nspc > 1 ) then
2316  !---< mixed phase without aerosol tracer >---
2317 
2318  ! nucleation
2319  call nucleat( ijkmax, & ! [IN]
2320  dens(:), & ! [IN]
2321  pres(:), & ! [IN]
2322  qdry(:), & ! [IN]
2323  ccn(:), & ! [IN]
2324  temp(:), & ! [INOUT]
2325  qvap(:), & ! [INOUT]
2326  ghyd(:,:,:), & ! [INOUT]
2327  cp(:), & ! [INOUT]
2328  cv(:), & ! [INOUT]
2329  dt ) ! [IN]
2330 
2331  ! freezing / melting
2332  call freezing( ijkmax, & ! [IN]
2333  ijkmax_cold, & ! [IN]
2334  index_cold(:), & ! [IN]
2335  flg_lt, & ! [IN]
2336  dens(:), & ! [IN]
2337  temp(:), & ! [INOUT]
2338  ghyd(:,:,:), & ! [INOUT]
2339  gcrg(:,:,:), & ! [INOUT]
2340  cp(:), & ! [INOUT]
2341  cv(:), & ! [INOUT]
2342  dt ) ! [IN]
2343 
2344  if( flg_icenucl ) then
2345  call ice_nucleat( ijkmax, & ! [IN]
2346  ijkmax_cold, & ! [IN]
2347  index_cold(:), & ! [IN]
2348  dens(:), & ! [IN]
2349  pres(:), & ! [IN]
2350  qdry(:), & ! [IN]
2351  temp(:), & ! [INOUT]
2352  qvap(:), & ! [INOUT]
2353  ghyd(:,:,:), & ! [INOUT]
2354  cp(:), & ! [INOUT]
2355  cv(:), & ! [INOUT]
2356  dt ) ! [IN]
2357  endif
2358 
2359  call melting( ijkmax, & ! [IN]
2360  ijkmax_warm, & ! [IN]
2361  index_warm(:), & ! [IN]
2362  flg_lt, & ! [IN]
2363  dens(:), & ! [IN]
2364  temp(:), & ! [INOUT]
2365  ghyd(:,:,:), & ! [INOUT]
2366  gcrg(:,:,:), & ! [INOUT]
2367  cp(:), & ! [INOUT]
2368  cv(:), & ! [INOUT]
2369  dt ) ! [IN]
2370 
2371  ! condensation / evaporation
2372  call cndevpsbl( ijkmax, & ! [IN]
2373  flg_lt, & ! [IN]
2374  dens(:), & ! [IN]
2375  pres(:), & ! [IN]
2376  qdry(:), & ! [IN]
2377  temp(:), & ! [INOUT]
2378  qvap(:), & ! [INOUT]
2379  ghyd(:,:,:), & ! [INOUT]
2380  cp(:), & ! [INOUT]
2381  cv(:), & ! [INOUT]
2382  evaporate(:), & ! [OUT]
2383  dt, & ! [IN]
2384  gcrg(:,:,:) ) ! [INOUT]
2385 
2386  if ( mp_doautoconversion ) then
2387  ! collision-coagulation
2388  call collmainf( ka, ia, ja, & ! [IN]
2389  ijkmax, & ! [IN]
2390  flg_lt, & ! [IN]
2391  d0_crg, & ! [IN]
2392  v0_crg, & ! [IN]
2393  dqcrg(:), & ! [IN]
2394  beta_crg(:), & ! [IN]
2395  temp(:), & ! [IN]
2396  ghyd(:,:,:), & ! [INOUT]
2397  gcrg(:,:,:), & ! [INOUT]
2398  crg_sep(:,:), & ! [OUT]
2399  dt ) ! [IN]
2400 
2401  endif
2402 
2403  endif
2404 
2405  endif
2406 
2407  return
2408  end subroutine mp_suzuki10
2409 
2410  !-----------------------------------------------------------------------------
2411  subroutine nucleat( &
2412  ijkmax, &
2413  dens, &
2414  pres, &
2415  qdry, &
2416  ccn, &
2417  temp, &
2418  qvap, &
2419  gc, &
2420  cp, &
2421  cv, &
2422  dtime )
2423  use scale_atmos_hydrometeor, only: &
2424  atmos_hydrometeor_lhv_para, &
2425  cp_vapor, &
2426  cp_water, &
2427  cv_vapor, &
2428  cv_water
2429  use scale_atmos_saturation, only: &
2430  atmos_saturation_pres2qsat_liq_para
2431  implicit none
2432 
2433  integer, intent(in) :: ijkmax
2434  real(rp), intent(in) :: dens(ijkmax) ! Density [kg/m3]
2435  real(rp), intent(in) :: pres(ijkmax) ! Pressure [Pa]
2436  real(rp), intent(in) :: qdry(ijkmax) ! dry air mass ratio [kg/kg]
2437  real(rp), intent(in) :: ccn(ijkmax) ! CCN number concentration [#/m3]
2438  real(rp), intent(inout) :: temp(ijkmax) ! Temperature [K]
2439  real(rp), intent(inout) :: qvap(ijkmax) ! Specific humidity [kg/kg]
2440  real(rp), intent(inout) :: gc (nbin,nspc,ijkmax) ! Mass size distribution function of hydrometeor
2441  real(rp), intent(inout) :: cp(ijkmax) ! specific heat
2442  real(rp), intent(inout) :: cv(ijkmax) ! specific heat
2443  real(dp), intent(in) :: dtime ! Time step interval
2444 
2445  real(rp) :: ssliq(ijkmax)
2446  real(rp) :: qlevp(ijkmax) ! LH
2447  real(rp) :: dmp
2448  real(rp) :: dqv
2449  integer :: n
2450  !
2451  real(rp) :: n_c
2452  real(rp) :: sumnum(ijkmax)
2453  real(rp) :: gcn( nbin,ijkmax ) ! number of cloud particles in each bin (=gc/exp(xctr))
2454  real(rp) :: qsat(ijkmax)
2455  integer :: ijk
2456 
2457  call prof_rapstart('_SBM_Nucleat', 3)
2458 
2459  ! lhv
2460  call atmos_hydrometeor_lhv_para( ijkmax, 1, ijkmax, temp(:), qlevp(:) )
2461 
2462  if( mp_couple_aerosol ) then
2463 
2464 
2465  !$omp parallel do private(dmp,dqv)
2466  do ijk = 1, ijkmax
2467  dmp = ccn(ijk) * expxctr( 1 )
2468  dmp = min( dmp,qvap(ijk)*dens(ijk) )
2469  gc( 1,il,ijk ) = gc( 1,il,ijk ) + dmp/dxmic
2470  dqv = dmp/dens(ijk)
2471  qvap(ijk) = qvap(ijk) - dqv
2472  temp(ijk) = temp(ijk) + dqv*qlevp(ijk)/cv(ijk)
2473  cp(ijk) = cp(ijk) + ( cp_water - cp_vapor ) * dqv
2474  cv(ijk) = cv(ijk) + ( cv_water - cv_vapor ) * dqv
2475  enddo
2476 
2477  else
2478 
2479  call atmos_saturation_pres2qsat_liq_para( ijkmax, 1, ijkmax, &
2480  temp(:), pres(:), qdry(:), & ! [IN]
2481  qsat(:) ) ! [OUT]
2482  !--- supersaturation
2483  !$omp parallel do
2484  do ijk = 1, ijkmax
2485  ssliq(ijk) = qvap(ijk)/qsat(ijk) - 1.0_rp
2486  enddo
2487 
2488  !$omp workshare
2489  sumnum(:) = 0.0_rp
2490  !$omp end workshare
2491  !$omp parallel do private(n_c,dmp,dqv)
2492  do ijk = 1, ijkmax
2493  ! if ( ssliq <= 0.0_RP ) cycle
2494  if ( ssliq(ijk) > 0.0_rp ) then
2495 
2496  !--- use for aerosol coupled model
2497  !--- mass -> number
2498  do n = 1, nbin
2499  gcn( n,ijk ) = gc( n,il,ijk )*rexpxctr( n )
2500  enddo
2501 
2502  do n = 1, nbin
2503  sumnum(ijk) = sumnum(ijk) + gcn( n,ijk )*dxmic
2504  enddo
2505  n_c = c_ccn * ( ssliq(ijk) * 1.e+2_rp )**( kappa )
2506  if ( n_c > sumnum(ijk) ) then
2507  dmp = ( n_c - sumnum(ijk) ) * expxctr( 1 )
2508  dmp = min( dmp,qvap(ijk)*dens(ijk) )
2509  gc( 1,il,ijk ) = gc( 1,il,ijk ) + dmp/dxmic
2510  dqv = dmp/dens(ijk)
2511  qvap(ijk) = qvap(ijk) - dqv
2512  qvap(ijk) = max( qvap(ijk),0.0_rp )
2513  temp(ijk) = temp(ijk) + dqv*qlevp(ijk)/cv(ijk)
2514  cp(ijk) = cp(ijk) + ( cp_water - cp_vapor ) * dqv
2515  cv(ijk) = cv(ijk) + ( cv_water - cv_vapor ) * dqv
2516  endif
2517  endif
2518 
2519  enddo
2520 
2521  endif
2522 
2523  call prof_rapend ('_SBM_Nucleat', 3)
2524 
2525  return
2526  end subroutine nucleat
2527 
2528  !-----------------------------------------------------------------------------
2529  subroutine nucleata( &
2530  ijkmax, &
2531  dens, &
2532  pres, &
2533  qdry, &
2534  temp, &
2535  qvap, &
2536  gc, &
2537  ga, &
2538  cp, &
2539  cv, &
2540  dtime )
2541  use scale_const, only: &
2542  pi => const_pi, &
2543  rvap => const_rvap, &
2544  rhow => const_dwatr
2545  use scale_atmos_hydrometeor, only: &
2546  atmos_hydrometeor_lhv_para, &
2547  cp_vapor, &
2548  cp_water, &
2549  cv_vapor, &
2550  cv_water
2551  use scale_atmos_saturation, only: &
2552  atmos_saturation_pres2qsat_liq_para
2553  implicit none
2554 
2555  integer, intent(in) :: ijkmax
2556  real(rp), intent(in) :: dens(ijkmax) ! Density [kg/m3]
2557  real(rp), intent(in) :: pres(ijkmax) ! Pressure [Pa]
2558  real(rp), intent(in) :: qdry(ijkmax) ! dry air mass ratio [kg/kg]
2559  real(rp), intent(inout) :: temp(ijkmax) ! Temperature [K]
2560  real(rp), intent(inout) :: qvap(ijkmax) ! Specific humidity [kg/kg]
2561  real(rp), intent(inout) :: gc (nbin,nspc,ijkmax) ! Mass size distribution function of hydrometeor
2562  real(rp), intent(inout) :: ga (nccn ,ijkmax) ! Mass size distribution function of aerosol
2563  real(rp), intent(inout) :: cp(ijkmax) ! specific heat
2564  real(rp), intent(inout) :: cv(ijkmax) ! specific heat
2565  real(dp), intent(in) :: dtime ! Time step interval
2566 
2567  real(rp) :: gan( nccn ) ! size distribution function ( aerosol ) : number ( gan = ga/exp( xactr ) )
2568  real(rp) :: ssliq ! supersaturatioin of liq. and ice, and LH
2569  real(rp) :: acoef, bcoef ! A and B in eq. (A.11) of Suzuki (2004)
2570  real(rp) :: rcrit ! critical radius (rcrit, r_N,crit of (A.11) of Suzuki (2004))
2571  real(rp) :: xcrit ! exp of hydrometeror whose radi is corresponding to rcrit (xcrit)
2572  real(rp) :: ractr, rcld, xcld, part, dmp
2573  integer :: n, nc, ncrit
2574  integer :: ncld(nccn)
2575 ! logical, save :: ofirst(1:ijkmax) = .true.
2576  !
2577  real(rp) :: qlevp(ijkmax)
2578  real(rp) :: qsatl(ijkmax)
2579  real(rp) :: dqv
2580  integer :: ijk
2581 
2582  call prof_rapstart('_SBM_NucleatA', 3)
2583 
2584  ! lhv
2585  call atmos_hydrometeor_lhv_para( ijkmax, 1, ijkmax, temp(:), qlevp(:) )
2586 
2587  call atmos_saturation_pres2qsat_liq_para( ijkmax, 1, ijkmax, &
2588  temp(:), pres(:), qdry(:), & ! [IN]
2589  qsatl(:) ) ! [OUT]
2590  !$omp parallel do private(ssliq,gan,ncld,acoef,bcoef,ractr,rcld,xcld,rcrit,xcrit,ncrit,part,nc,dmp,dqv)
2591  do ijk = 1, ijkmax
2592  !--- supersaturation
2593  ssliq = qvap(ijk)/qsatl(ijk) - 1.0_rp
2594 
2595  if ( ssliq <= 0.0_rp ) cycle
2596 
2597  !--- use for aerosol coupled model
2598  !--- mass -> number
2599  do n = 1, nccn
2600  gan( n ) = ga( n,ijk )*rexpxactr( n )
2601  enddo
2602 
2603  acoef = 2.0_rp*sigma/rvap/rhow/temp(ijk) ! A in (A.11) of Suzuki (2004)
2604  bcoef = vhfct* rhoa/rhow * emwtr/emaer ! B in (A.11) of Suzuki (2004)
2605 
2606  !--- relationship of bin number
2607  if ( flg_nucl ) then
2608  do n = 1, nccn
2609  ncld( n ) = 1
2610  end do
2611  else
2612  do n = 1, nccn
2613  ractr = ( expxactr( n )*thirdovforth/pi/rhoa )**( oneovthird )
2614  rcld = sqrt( 3.0_rp*bcoef*ractr*ractr*ractr / acoef )
2615  xcld = log( rhow * 4.0_rp*pi*oneovthird*rcld*rcld*rcld )
2616  ncld( n ) = int( ( xcld-xctr( 1 ) )/dxmic ) + 1
2617  ncld( n ) = min( max( ncld( n ),1 ),nbin )
2618  enddo
2619  endif
2620 
2621  !--- nucleation
2622  do n = nccn, 1, -1
2623  if ( ssliq <= 0.0_rp ) exit
2624  !--- use for aerosol coupled model
2625  acoef = 2.0_rp*sigma/rvap/rhow/temp(ijk) ! A in (A.11) of Suzuki (2004)
2626  rcrit = acoef*oneovthird * ( 4.0_rp/bcoef )**( oneovthird ) / ssliq**( twoovthird ) ! r_{N,crit} in (A.11) of Suzuki (2004)
2627  xcrit = log( rhoa * 4.0_rp*pi*oneovthird * rcrit*rcrit*rcrit )
2628  ncrit = int( ( xcrit-xabnd( 1 ) )/dxaer ) + 1
2629 
2630  if ( n == ncrit ) then
2631  part = ( xabnd( ncrit+1 )-xcrit )/dxaer
2632  elseif ( n > ncrit ) then
2633  part = 1.0_rp
2634  else
2635  exit
2636  endif
2637 
2638  !--- calculate mass change
2639  nc = ncld( n )
2640  dmp = part*gan( n )*dxaer*expxctr( nc )
2641  dmp = min( dmp,qvap(ijk)*dens(ijk) )
2642  gc( nc,il,ijk ) = gc( nc,il,ijk ) + dmp/dxmic
2643  gan( n ) = gan( n ) - dmp/dxaer*rexpxctr( nc )
2644  gan( n ) = max( gan( n ), 0.0_rp )
2645  dqv = dmp/dens(ijk)
2646  qvap(ijk) = qvap(ijk) - dqv
2647  qvap(ijk) = max( qvap(ijk),0.0_rp )
2648  temp(ijk) = temp(ijk) + dqv*qlevp(ijk)/cv(ijk)
2649  cp(ijk) = cp(ijk) + ( cp_water - cp_vapor ) * dqv
2650  cv(ijk) = cv(ijk) + ( cv_water - cv_vapor ) * dqv
2651  enddo
2652 
2653  !--- number -> mass
2654  do n = 1, nccn
2655  ga( n,ijk ) = gan( n )*expxactr( n )
2656  enddo
2657 
2658  enddo
2659 
2660  call prof_rapend ('_SBM_NucleatA', 3)
2661 
2662  return
2663  end subroutine nucleata
2664 
2665  !-----------------------------------------------------------------------------
2666  subroutine cndevpsbl( &
2667  ijkmax, &
2668  flg_lt, &
2669  dens, &
2670  pres, &
2671  qdry, &
2672  temp, &
2673  qvap, &
2674  gc, &
2675  cp, &
2676  cv, &
2677  evaporate, &
2678  dtime, &
2679  gcrg )
2680  implicit none
2681 
2682  integer, intent(in) :: ijkmax
2683  logical, intent(in) :: flg_lt
2684  real(rp), intent(in) :: dens(ijkmax) ! Density [kg/m3]
2685  real(rp), intent(in) :: pres(ijkmax) ! Pressure [Pa]
2686  real(rp), intent(in) :: qdry(ijkmax) ! dry air mass ratio [kg/kg]
2687  real(rp), intent(inout) :: temp(ijkmax) ! Temperature [K]
2688  real(rp), intent(inout) :: qvap(ijkmax) ! Specific humidity [kg/kg]
2689  real(rp), intent(inout) :: gc (nbin,nspc,ijkmax) ! Mass size distribution function of hydrometeor
2690  real(rp), intent(inout) :: cp(ijkmax) ! specific heat
2691  real(rp), intent(inout) :: cv(ijkmax) ! specific heat
2692  real(rp), intent(out) :: evaporate(ijkmax) ! Number concentration of evaporated cloud [/m3]
2693  real(dp), intent(in) :: dtime ! Time step interval
2694  real(rp), intent(inout) :: gcrg(nbin,nspc,ijkmax) ! Mass size distribution function of hydrometeor
2695  !---------------------------------------------------------------------------
2696 
2697  call liqphase( ijkmax, & ! [IN]
2698  flg_lt, & ! [IN]
2699  dens(:), & ! [IN]
2700  pres(:), & ! [IN]
2701  qdry(:), & ! [IN]
2702  temp(:), & ! [INOUT]
2703  qvap(:), & ! [INOUT]
2704  gc(:,:,:), & ! [INOUT]
2705  gcrg(:,:,:), & ! [INOUT]
2706  cp(:), & ! [INOUT]
2707  cv(:), & ! [INOUT]
2708  evaporate(:), & ! [OUT]
2709  dtime ) ! [IN]
2710 
2711  if( nspc > 1 ) then
2712 
2713  call icephase( ijkmax, & ! [IN]
2714  flg_lt, & ! [IN]
2715  dens(:), & ! [IN]
2716  pres(:), & ! [IN]
2717  qdry(:), & ! [IN]
2718  temp(:), & ! [INOUT]
2719  qvap(:), & ! [INOUT]
2720  gc(:,:,:), & ! [INOUT]
2721  gcrg(:,:,:), & ! [INOUT]
2722  cp(:), & ! [INOUT]
2723  cv(:), & ! [INOUT]
2724  dtime ) ! [IN]
2725 
2726  call mixphase( ijkmax, & ! [IN]
2727  flg_lt, & ! [IN]
2728  dens(:), & ! [IN]
2729  pres(:), & ! [IN]
2730  qdry(:), & ! [IN]
2731  temp(:), & ! [INOUT]
2732  qvap(:), & ! [INOUT]
2733  gc(:,:,:), & ! [INOUT]
2734  gcrg(:,:,:), & ! [INOUT]
2735  cp(:), & ! [INOUT]
2736  cv(:), & ! [INOUT]
2737  dtime ) ! [IN]
2738 
2739  endif
2740 
2741  return
2742  end subroutine cndevpsbl
2743 
2744  !-----------------------------------------------------------------------------
2745  subroutine cndevpsbla( &
2746  ijkmax, &
2747  flg_lt, &
2748  dens, &
2749  pres, &
2750  qdry, &
2751  temp, &
2752  qvap, &
2753  gc, &
2754  ga, &
2755  cp, &
2756  cv, &
2757  evaporate, &
2758  dtime, &
2759  gcrg )
2760  implicit none
2761 
2762  integer, intent(in) :: ijkmax
2763  logical, intent(in) :: flg_lt
2764  real(rp), intent(in) :: dens(ijkmax) ! Density [kg/m3]
2765  real(rp), intent(in) :: pres(ijkmax) ! Pressure [Pa]
2766  real(rp), intent(in) :: qdry(ijkmax) ! dry air mass ratio [kg/kg]
2767  real(rp), intent(inout) :: temp(ijkmax) ! Temperature [K]
2768  real(rp), intent(inout) :: qvap(ijkmax) ! Specific humidity [kg/kg]
2769  real(rp), intent(inout) :: gc (nbin,nspc,ijkmax) ! Mass size distribution function of hydrometeor
2770  real(rp), intent(inout) :: ga (nccn ,ijkmax) ! Mass size distribution function of aerosol
2771  real(rp), intent(inout) :: cp(ijkmax) ! specific heat
2772  real(rp), intent(inout) :: cv(ijkmax) ! specific heat
2773  real(rp), intent(out) :: evaporate(ijkmax) ! Number concentration of evaporated cloud [/m3]
2774  real(dp), intent(in) :: dtime ! Time step interval
2775  real(rp), intent(inout) :: gcrg(nbin,nspc,ijkmax) ! Mass size distribution function of hydrometeor
2776  !---------------------------------------------------------------------------
2777 
2778  call liqphase( ijkmax, & ! [IN]
2779  flg_lt, & ! [IN]
2780  dens(:), & ! [IN]
2781  pres(:), & ! [IN]
2782  qdry(:), & ! [IN]
2783  temp(:), & ! [INOUT]
2784  qvap(:), & ! [INOUT]
2785  gc(:,:,:), & ! [INOUT]
2786  gcrg(:,:,:), & ! [INOUT]
2787  cp(:), & ! [INOUT]
2788  cv(:), & ! [INOUT]
2789  evaporate(:), & ! [OUT]
2790  dtime ) ! [IN]
2791 
2792  ! regeneration of aerosol
2793  if ( flg_regeneration ) then
2794  call faero( ijkmax, & ! [IN]
2795 ! regene_gcn(:), & ! [IN]
2796  evaporate(:), & ! [IN]
2797  ga(:,:) ) ! [INOUT]
2798  endif
2799 
2800  if( nspc > 1 ) then
2801  call icephase( ijkmax, & ! [IN]
2802  flg_lt, & ! [IN]
2803  dens(:), & ! [IN]
2804  pres(:), & ! [IN]
2805  qdry(:), & ! [IN]
2806  temp(:), & ! [INOUT]
2807  qvap(:), & ! [INOUT]
2808  gc(:,:,:), & ! [INOUT]
2809  gcrg(:,:,:), & ! [INOUT]
2810  cp(:), & ! [INOUT]
2811  cv(:), & ! [INOUT]
2812  dtime ) ! [IN]
2813 
2814  call mixphase( ijkmax, & ! [IN]
2815  flg_lt, & ! [IN]
2816  dens(:), & ! [IN]
2817  pres(:), & ! [IN]
2818  qdry(:), & ! [IN]
2819  temp(:), & ! [INOUT]
2820  qvap(:), & ! [INOUT]
2821  gc(:,:,:), & ! [INOUT]
2822  gcrg(:,:,:), & ! [INOUT]
2823  cp(:), & ! [INOUT]
2824  cv(:), & ! [INOUT]
2825  dtime ) ! [IN]
2826  endif
2827 
2828  return
2829  end subroutine cndevpsbla
2830 
2831  !-----------------------------------------------------------------------------
2832  subroutine liqphase( &
2833  ijkmax, &
2834  flg_lt, &
2835  dens, &
2836  pres, &
2837  qdry, &
2838  temp, &
2839  qvap, &
2840  gc, &
2841  gcrg, &
2842  cp, &
2843  cv, &
2844  regene_gcn, &
2845  dtime )
2846  use scale_const, only: &
2847  pi => const_pi, &
2848  eps => const_eps, &
2849  esat0 => const_psat0, &
2850  rvap => const_rvap, &
2851  tmlt => const_tmelt, &
2852  tem00 => const_tem00
2853  use scale_atmos_hydrometeor, only: &
2854  atmos_hydrometeor_lhv_para, &
2855  cp_vapor, &
2856  cp_water, &
2857  cv_vapor, &
2858  cv_water
2859  use scale_atmos_saturation, only: &
2860  atmos_saturation_pres2qsat_liq_para
2861  implicit none
2862 
2863  integer, intent(in) :: ijkmax
2864  logical, intent(in) :: flg_lt
2865  real(rp), intent(in) :: dens(ijkmax) ! Density [kg/m3]
2866  real(rp), intent(in) :: pres(ijkmax) ! Pressure [Pa]
2867  real(rp), intent(in) :: qdry(ijkmax) ! dry air mass ratio [kg/kg]
2868  real(rp), intent(inout) :: temp(ijkmax) ! Temperature [K]
2869  real(rp), intent(inout) :: qvap(ijkmax) ! Specific humidity [kg/kg]
2870  real(rp), intent(inout) :: gc (nbin,nspc,ijkmax) ! Mass size distribution function of hydrometeor
2871  real(rp), intent(inout) :: gcrg(nbin,nspc,ijkmax)
2872  real(rp), intent(inout) :: cp(ijkmax) ! specific heat
2873  real(rp), intent(inout) :: cv(ijkmax) ! specific heat
2874  real(rp), intent(out) :: regene_gcn(ijkmax) ! mass of regenerated aerosol
2875  real(dp), intent(in) :: dtime ! Time step interval
2876 
2877  integer :: n, myu, ncount
2878  integer :: nloop(ijkmax) ! number of fractional step for condensation
2879  real(rp) :: gclold(ijkmax)
2880  real(rp) :: gclnew(ijkmax)
2881  real(rp) :: cndmss
2882  real(rp) :: dtcnd(ijkmax) ! dt for condensation with considering CFL condition
2883  real(rp) :: gtliq(ijkmax) ! G of eq. (A.17) of Suzuki (2004)
2884  real(rp) :: qlevp(ijkmax) ! LH
2885  real(rp) :: cefliq, a, sliqtnd
2886  real(rp) :: sumliq(ijkmax), umax(ijkmax)
2887  real(rp) :: ssliq ! super saturation
2888  real(rp) :: gcn( -1:nbin+2,nspc,ijkmax ) ! size distribution function(hydrometeor): number
2889  ! real(RP) :: gcnold( nbin,ijkmax )
2890  real(rp), parameter :: cflfct = 0.50_rp ! CFL limiter
2891  ! real(RP) :: old_sum_gcn, new_sum_gcn
2892  integer :: iflg( nspc,ijkmax ) ! flag whether calculation is conduct or not
2893  real(rp) :: csum
2894  real(rp) :: f1, f2, emu, cefd, cefk, festl
2895  real(rp) :: qsatl(ijkmax)
2896  real(rp), parameter :: afmyu = 1.72e-05_rp, bfmyu = 3.93e+2_rp
2897  real(rp), parameter :: cfmyu = 1.2e+02_rp, fct = 1.4e+3_rp
2898  real(rp) :: zerosw, qvtmp
2899  real(rp) :: dqv
2900  integer :: ijk, nn, mm, loopflg(ijkmax)
2901 
2902  !--- local for advection
2903  real(rp) :: uadv ( 0:nbin+2,nspc,ijkmax )
2904  real(rp) :: crn ( 0:nbin+2,nspc,ijkmax )
2905  real(rp) :: flq ( 1:nbin+1 )
2906  real(rp) :: acoef( 0:nbin+1,0:2 )
2907  real(rp) :: aip ( 0:nbin+1 )
2908  real(rp) :: aim ( 0:nbin+1 )
2909  real(rp) :: ai ( 0:nbin+1 )
2910  real(rp) :: cmins, cplus
2911  integer :: nloopmax
2912 
2913  !--- for lithgning component (if not use lightning component, arrays shown below work dummy array)
2914  real(rp) :: gcrgn( -1:nbin+2,nspc,ijkmax )
2915  real(rp) :: flq_c( 1:nbin+1 )
2916  real(rp) :: acoef_c( 0:nbin+1,0:2 )
2917  real(rp) :: aip_c ( 0:nbin+1 )
2918  real(rp) :: aim_c ( 0:nbin+1 )
2919  real(rp) :: ai_c ( 0:nbin+1 )
2920 
2921  call prof_rapstart('_SBM_Liqphase', 3)
2922 
2923  !$omp workshare
2924  iflg(:,:) = 0
2925  !$omp end workshare
2926  !$omp parallel do private(csum)
2927  do ijk = 1, ijkmax
2928  csum = 0.0_rp
2929  do n = 1, nbin
2930  csum = csum + gc( n,il,ijk )*dxmic
2931  enddo
2932  if( csum > cldmin ) iflg( il,ijk ) = 1
2933  enddo
2934 
2935  ! lhv
2936  call atmos_hydrometeor_lhv_para( ijkmax, 1, ijkmax, temp(:), qlevp(:) )
2937 
2938  call atmos_saturation_pres2qsat_liq_para( ijkmax, 1, ijkmax, &
2939  temp(:), pres(:), qdry(:), & ! [IN]
2940  qsatl(:) ) ! [OUT]
2941 
2942  !$omp workshare
2943  nloop(:) = 0
2944  gclold(:) = 0.0_rp
2945  !$omp end workshare
2946  !$omp parallel do private(ssliq,emu,cefd,cefk,festl,f1,f2)
2947  do ijk = 1, ijkmax
2948 
2949  do n = 1, nbin
2950  gclold(ijk) = gclold(ijk) + gc( n,il,ijk ) * dxmic
2951  enddo
2952  !
2953  !------- mass -> number
2954  do n = 1, nbin
2955  gcn( n,il,ijk ) = gc( n,il,ijk ) * rexpxctr( n )
2956  gcrgn( n,il,ijk ) = gcrg( n,il,ijk )
2957  enddo
2958  gcn( -1,il,ijk ) = 0.0_rp
2959  gcn( 0,il,ijk ) = 0.0_rp
2960  gcn( nbin+1,il,ijk ) = 0.0_rp
2961  gcn( nbin+2,il,ijk ) = 0.0_rp
2962  gcrgn( -1,il,ijk ) = 0.0_rp
2963  gcrgn( 0,il,ijk ) = 0.0_rp
2964  gcrgn( nbin+1,il,ijk ) = 0.0_rp
2965  gcrgn( nbin+2,il,ijk ) = 0.0_rp
2966  !
2967  !--- super saturation
2968  ssliq = qvap(ijk)/qsatl(ijk) - 1.0_rp
2969 
2970  emu = afmyu*( bfmyu/( temp(ijk)+cfmyu ) )*( temp(ijk)/tmlt )**1.50_rp
2971  cefd = emu/dens(ijk)
2972  cefk = fct*emu
2973 
2974  festl = esat0*exp( qlevp(ijk)/rvap*( 1.0_rp/tem00 - 1.0_rp/temp(ijk) ) )
2975  f1 = rvap*temp(ijk)/festl/cefd
2976  f2 = qlevp(ijk)/cefk/temp(ijk)*( qlevp(ijk)/rvap/temp(ijk) - 1.0_rp )
2977  gtliq(ijk) = 4.0_rp*pi/( f1+f2 ) !--- G of eq. (A.17) of Suzuki (2004)
2978  !------- CFL condition
2979  umax(ijk) = cbnd( 1,il )*rexpxbnd( 1 )*gtliq(ijk)*abs( ssliq )
2980  dtcnd(ijk) = cflfct*dxmic/umax(ijk)
2981  nloop(ijk) = int( dtime/dtcnd(ijk) ) + 1
2982  dtcnd(ijk) = dtime / nloop(ijk)
2983 
2984  nloop(ijk) = nloop(ijk) * iflg( il,ijk ) !--- for determing trivial loop
2985  enddo
2986  nloopmax = maxval(nloop,1)
2987 
2988  !
2989  !$omp workshare
2990  regene_gcn(:) = 0.0_rp
2991  !$omp end workshare
2992 !OCL LOOP_NOFISSION
2993 !OCL LOOP_NOINTERCHANGE
2994  do ncount = 1, nloopmax
2995 
2996  !$omp parallel do
2997  do ijk = 1, ijkmax
2998  loopflg(ijk) = min( 1, int(nloop(ijk)/ncount) ) ! 0 or 1
2999  enddo
3000 
3001  ! lhv
3002  call atmos_hydrometeor_lhv_para( ijkmax, 1, ijkmax, temp(:), qlevp(:) )
3003 
3004 !OCL LOOP_NOFUSION
3005  !$omp parallel do private(emu,cefd,cefk,festl,f1,f2)
3006  do ijk = 1, ijkmax
3007  do nn = 1, loopflg(ijk)
3008 
3009  emu = afmyu*( bfmyu/( temp(ijk)+cfmyu ) )*( temp(ijk)/tmlt )**1.50_rp
3010  cefd = emu/dens(ijk)
3011  cefk = fct*emu
3012  festl = esat0*exp( qlevp(ijk)/rvap*( 1.0_rp/tem00 - 1.0_rp/temp(ijk) ) )
3013  f1 = rvap*temp(ijk)/festl/cefd
3014  f2 = qlevp(ijk)/cefk/temp(ijk)*( qlevp(ijk)/rvap/temp(ijk) - 1.0_rp )
3015  gtliq(ijk) = 4.0_rp*pi/( f1+f2 ) !--- G of eq. (A.17) of Suzuki (2004)
3016 
3017  sumliq(ijk) = 0.0_rp
3018  ! old_sum_gcn = 0.0_RP
3019  do n = 1, nbin
3020  sumliq(ijk) = sumliq(ijk) + gcn( n,il,ijk )*cctr( n,il )*dxmic
3021  enddo
3022 
3023  enddo
3024  enddo
3025 
3026  !--- super saturation
3027  call atmos_saturation_pres2qsat_liq_para( ijkmax, 1, ijkmax, &
3028  temp(:), pres(:), qdry(:), & ! [IN]
3029  qsatl(:) ) ! [OUT]
3030 !OCL LOOP_NOFUSION
3031  !$omp parallel do private(ssliq,zerosw,qvtmp,cefliq,a,sliqtnd)
3032  do ijk = 1, ijkmax
3033  do nn = 1, loopflg(ijk)
3034 
3035  ssliq = qvap(ijk)/qsatl(ijk) - 1.0_rp
3036 
3037  !----- supersaturation tendency
3038  zerosw = 0.5_rp + sign( 0.5_rp,qvap(ijk)-eps ) !--- zerosw = 1 (qv>0), zerosw=0 (qv=0)
3039  qvtmp = qvap(ijk) * zerosw + ( qvap(ijk)+eps ) * ( 1.0_rp-zerosw )
3040  cefliq = ( ssliq+1.0_rp )*( 1.0_rp/qvtmp + qlevp(ijk)*qlevp(ijk)/cv(ijk)/rvap/temp(ijk)/temp(ijk) )
3041  a = - cefliq*sumliq(ijk)*gtliq(ijk)/dens(ijk) ! a of eq. (A.19) of Suzuki (2004)
3042  a = a + eps * ( 1.0_rp - zerosw ) !--- avoiding division by zero when qv = 0
3043 
3044  sliqtnd = zerosw * &
3045  ( ssliq * ( exp( a*dtcnd(ijk) )-1.0_rp )/( a*dtcnd(ijk) ) &
3046  * ( 0.5_rp + sign( 0.5_rp,abs(a*dtcnd(ijk)-0.1_rp) ) ) &
3047  + ssliq &
3048  * ( 0.5_rp - sign( 0.5_rp,abs(a*dtcnd(ijk)-0.1_rp) ) ) &
3049  ) &
3050  + ssliq * ( 1.0_rp - zerosw )
3051  !
3052  !----- change of SDF
3053  do mm = 1, iflg( il,ijk )
3054  !--- advection speed
3055  do n = 1, nbin+1
3056  uadv( n,il,ijk ) = cbnd( n,il )*rexpxbnd( n )*gtliq(ijk)*sliqtnd ! U of eq. (A.18) of Suzuki (2004)
3057  enddo
3058  uadv( 0 ,il,ijk ) = 0.0_rp
3059  uadv( nbin+2,il,ijk ) = 0.0_rp
3060 
3061 ! do n = 1, nbin
3062 ! gcnold( n ) = gcn( n,il,ijk )
3063 ! enddo
3064  enddo
3065 
3066  enddo
3067  enddo
3068 
3069  call prof_rapstart('_SBM_AdvLiq', 3)
3070 
3071  !$omp parallel do private(acoef,aip,aim,ai,flq,cplus,cmins)
3072  do ijk = 1, ijkmax
3073  do nn = 1, loopflg(ijk)
3074  do myu = 1, il
3075  do mm = 1, iflg( myu,ijk )
3076  do n = 0, nbin+2
3077  crn( n,myu,ijk ) = uadv( n,myu,ijk )*dtcnd(ijk)/dxmic
3078  enddo
3079 
3080  do n = 0, nbin+1
3081  acoef(n,0) = - ( gcn( n+1,myu,ijk )-26.0_rp*gcn( n,myu,ijk )+gcn( n-1,myu,ijk ) ) / 48.0_rp
3082  acoef(n,1) = ( gcn( n+1,myu,ijk ) -gcn( n-1,myu,ijk ) ) / 16.0_rp
3083  acoef(n,2) = ( gcn( n+1,myu,ijk )- 2.0_rp*gcn( n,myu,ijk )+gcn( n-1,myu,ijk ) ) / 48.0_rp
3084  end do
3085 
3086  do n = 0, nbin+1
3087  cplus = 1.0_rp - ( crn(n+1,myu,ijk) + abs(crn(n+1,myu,ijk)) )
3088 
3089  aip(n) = acoef(n,0) * ( 1.0_rp-cplus**1 ) &
3090  + acoef(n,1) * ( 1.0_rp-cplus**2 ) &
3091  + acoef(n,2) * ( 1.0_rp-cplus**3 )
3092  aip(n) = max( aip(n), 0.0_rp )
3093  end do
3094 
3095  do n = 0, nbin+1
3096  cmins = 1.0_rp - ( abs(crn(n,myu,ijk)) - crn(n,myu,ijk) )
3097 
3098  aim(n) = acoef(n,0) * ( 1.0_rp-cmins**1 ) &
3099  - acoef(n,1) * ( 1.0_rp-cmins**2 ) &
3100  + acoef(n,2) * ( 1.0_rp-cmins**3 )
3101  aim(n) = max( aim(n), 0.0_rp )
3102  end do
3103 
3104  do n = 0, nbin+1
3105  ai(n) = acoef(n,0) * 2.0_rp &
3106  + acoef(n,2) * 2.0_rp
3107  ai(n) = max( ai(n), aip(n)+aim(n)+cldmin )
3108  end do
3109 
3110  do n = 1, nbin+1
3111  flq(n) = ( aip(n-1)/ai(n-1)*gcn( n-1,myu,ijk ) &
3112  - aim(n )/ai(n )*gcn( n ,myu,ijk ) )*dxmic/dtcnd(ijk)
3113  end do
3114 
3115  do n = 1, nbin
3116  gcn( n,myu,ijk ) = gcn( n,myu,ijk ) - ( flq(n+1)-flq(n) )*dtcnd(ijk)/dxmic
3117  end do
3118 
3119  regene_gcn(ijk) = regene_gcn(ijk)+( -flq(1)*dtcnd(ijk)/dxmic ) &
3120  * min( uadv(1,myu,ijk),0.0_rp )/( uadv(1,myu,ijk) + eps )
3121  enddo
3122  enddo
3123  enddo
3124  enddo
3125 
3126 
3127  if ( flg_lt ) then
3128 
3129  !$omp parallel do private(acoef_c,aip_c,aim_c,ai_c,flq_c,cplus,cmins)
3130  do ijk = 1, ijkmax
3131  do nn = 1, loopflg(ijk)
3132  do myu = 1, il
3133  do mm = 1, iflg( myu,ijk )
3134  do n = 0, nbin+1
3135  acoef_c(n,0) = - ( gcrgn( n+1,myu,ijk )-26.0_rp*gcrgn( n,myu,ijk )+gcrgn( n-1,myu,ijk ) ) / 48.0_rp
3136  acoef_c(n,1) = ( gcrgn( n+1,myu,ijk ) -gcrgn( n-1,myu,ijk ) ) / 16.0_rp
3137  acoef_c(n,2) = ( gcrgn( n+1,myu,ijk )- 2.0_rp*gcrgn( n,myu,ijk )+gcrgn( n-1,myu,ijk ) ) / 48.0_rp
3138  enddo
3139 
3140  do n = 0, nbin+1
3141  cplus = 1.0_rp - ( crn(n+1,myu,ijk) + abs(crn(n+1,myu,ijk)) )
3142 
3143  aip_c(n) = acoef_c(n,0) * ( 1.0_rp-cplus**1 ) &
3144  + acoef_c(n,1) * ( 1.0_rp-cplus**2 ) &
3145  + acoef_c(n,2) * ( 1.0_rp-cplus**3 )
3146  aip_c(n) = max( aip_c(n), 0.0_rp )
3147  enddo
3148 
3149  do n = 0, nbin+1
3150  cmins = 1.0_rp - ( abs(crn(n,myu,ijk)) - crn(n,myu,ijk) )
3151 
3152  aim_c(n) = acoef_c(n,0) * ( 1.0_rp-cmins**1 ) &
3153  - acoef_c(n,1) * ( 1.0_rp-cmins**2 ) &
3154  + acoef_c(n,2) * ( 1.0_rp-cmins**3 )
3155  aim_c(n) = max( aim_c(n), 0.0_rp )
3156  enddo
3157 
3158  do n = 0, nbin+1
3159  ai_c(n) = acoef_c(n,0) * 2.0_rp &
3160  + acoef_c(n,2) * 2.0_rp
3161  ai_c(n) = max( ai_c(n), aip_c(n)+aim_c(n)+cldmin )
3162  enddo
3163 
3164  do n = 1, nbin+1
3165  flq_c(n) = ( aip_c(n-1)/ai_c(n-1)*gcrgn( n-1,myu,ijk ) &
3166  - aim_c(n )/ai_c(n )*gcrgn( n ,myu,ijk ) )/dtcnd(ijk)
3167  enddo
3168 
3169  do n = 1, nbin
3170  gcrgn( n,myu,ijk ) = gcrgn( n,myu,ijk ) - ( flq_c(n+1)-flq_c(n) )*dtcnd(ijk)
3171  enddo
3172  enddo
3173  enddo
3174  enddo
3175  enddo
3176 
3177  end if
3178 
3179  call prof_rapend ('_SBM_AdvLiq', 3)
3180 
3181 !OCL LOOP_NOFUSION
3182  !$omp parallel do private(cndmss,dqv)
3183  do ijk = 1, ijkmax
3184  do nn = 1, loopflg(ijk)
3185 
3186  !----- new mass
3187  gclnew(ijk) = 0.0_rp
3188 ! new_sum_gcn = 0.0_RP
3189  do n = 1, nbin
3190  gclnew(ijk) = gclnew(ijk) + gcn( n,il,ijk )*expxctr( n )
3191 ! old_sum_gcn = old_sum_gcn + gcnold( n )*dxmic
3192 ! new_sum_gcn = new_sum_gcn + gcn( n,ijk )*dxmic
3193  enddo
3194 
3195  gclnew(ijk) = gclnew(ijk)*dxmic
3196  !
3197  !----- change of humidity and temperature
3198  cndmss = gclnew(ijk) - gclold(ijk)
3199  dqv = cndmss/dens(ijk)
3200  qvap(ijk) = qvap(ijk) - dqv
3201  temp(ijk) = temp(ijk) + dqv*qlevp(ijk)/cv(ijk)
3202  cp(ijk) = cp(ijk) + ( cp_water - cp_vapor ) * dqv
3203  cv(ijk) = cv(ijk) + ( cv_water - cv_vapor ) * dqv
3204  !
3205  gclold(ijk) = gclnew(ijk)
3206  !
3207  !----- continue/end
3208  enddo
3209  enddo
3210 
3211  enddo !nloop
3212  !
3213  !$omp parallel do
3214  do ijk = 1, ijkmax
3215  !------- number -> mass
3216  do n = 1 , nbin
3217  gc( n,il,ijk ) = gcn( n,il,ijk )*expxctr( n )
3218  enddo
3219  enddo
3220  if ( flg_lt ) then
3221  !$omp parallel do
3222  do ijk = 1, ijkmax
3223  !------- number -> mass
3224  do n = 1 , nbin
3225  gcrg( n,il,ijk ) = gcrgn( n,il,ijk )
3226  enddo
3227  enddo
3228  end if
3229 
3230  ! lhv
3231  call atmos_hydrometeor_lhv_para( ijkmax, 1, ijkmax, temp(:), qlevp(:) )
3232 
3233  !$omp parallel do private(cndmss,dqv)
3234  do ijk = 1, ijkmax
3235  do n = 1 , nbin
3236  if ( gc( n,il,ijk ) < 0.0_rp ) then
3237  cndmss = -gc( n,il,ijk )
3238  gc( n,il,ijk ) = 0.0_rp
3239  gcrg( n,il,ijk ) = 0.0_rp
3240  dqv = cndmss/dens(ijk)
3241  qvap(ijk) = qvap(ijk) + dqv
3242  temp(ijk) = temp(ijk) - dqv*qlevp(ijk)/cv(ijk)
3243  cp(ijk) = cp(ijk) + ( cp_vapor - cp_water ) * dqv
3244  cv(ijk) = cv(ijk) + ( cv_vapor - cv_water ) * dqv
3245  endif
3246  enddo
3247  enddo
3248 
3249  call prof_rapend ('_SBM_Liqphase', 3)
3250 
3251  return
3252  end subroutine liqphase
3253 
3254  !-----------------------------------------------------------------------------
3255  subroutine icephase( &
3256  ijkmax, &
3257  flg_lt, &
3258  dens, &
3259  pres, &
3260  qdry, &
3261  temp, &
3262  qvap, &
3263  gc, &
3264  gcrg, &
3265  cp, &
3266  cv, &
3267  dtime )
3268  use scale_const, only: &
3269  pi => const_pi, &
3270  eps => const_eps, &
3271  esat0 => const_psat0, &
3272  tem00 => const_tem00, &
3273  rvap => const_rvap, &
3274  tmlt => const_tmelt
3275  use scale_atmos_hydrometeor, only: &
3276  atmos_hydrometeor_lhs_para, &
3277  cp_vapor, &
3278  cp_ice, &
3279  cv_vapor, &
3280  cv_ice
3281  use scale_atmos_saturation, only: &
3282  atmos_saturation_pres2qsat_ice_para
3283  implicit none
3284 
3285  integer, intent(in) :: ijkmax
3286  logical, intent(in) :: flg_lt
3287  real(rp), intent(in) :: dens(ijkmax) ! Density [kg/m3]
3288  real(rp), intent(in) :: pres(ijkmax) ! Pressure [Pa]
3289  real(rp), intent(in) :: qdry(ijkmax) ! dry air mass ratio [kg/kg]
3290  real(rp), intent(inout) :: temp(ijkmax) ! Temperature [K]
3291  real(rp), intent(inout) :: qvap(ijkmax) ! Specific humidity [kg/kg]
3292  real(rp), intent(inout) :: gc (nbin,nspc,ijkmax) ! Mass size distribution function of hydrometeor
3293  real(rp), intent(inout) :: gcrg(nbin,nspc,ijkmax)
3294  real(rp), intent(inout) :: cp(ijkmax) ! specific heat
3295  real(rp), intent(inout) :: cv(ijkmax) ! specific heat
3296  real(dp), intent(in) :: dtime ! Time step interval
3297 
3298  integer :: myu, n, ncount
3299  integer :: nloop(ijkmax) !number of fractional step for condensation
3300  real(rp) :: gciold(ijkmax), gcinew(ijkmax)
3301  real(rp) :: dtcnd(ijkmax) ! dt for condensation with considering CFL condition
3302  real(rp) :: sblmss
3303  real(rp) :: gtice(ijkmax), umax(ijkmax)
3304  real(rp) :: qlsbl(ijkmax)
3305  real(rp) :: cefice, d, uval, sicetnd
3306  real(rp) :: sumice(ijkmax)
3307  real(rp) :: ssice
3308  real(rp) :: gcn( -1:nbin+2,nspc,ijkmax ) ! size distribution function (Hydrometeor): number
3309  real(rp), parameter :: cflfct = 0.50_rp
3310 ! real(RP) :: dumm_regene(ijkmax)
3311  integer :: iflg( nspc,ijkmax )
3312  real(rp) :: csum
3313  real(rp) :: f1, f2, emu, cefd, cefk, festi
3314  real(rp) :: qsati(ijkmax)
3315  real(rp), parameter :: afmyu = 1.72e-05_rp, bfmyu = 3.93e+2_rp
3316  real(rp), parameter :: cfmyu = 1.2e+02_rp, fct = 1.4e+3_rp
3317  integer :: ijk, mm, nn, loopflg(ijkmax)
3318  real(rp) :: zerosw, qvtmp
3319  real(rp) :: dqv
3320 
3321  !--- local for advection
3322 ! real(RP) :: qadv( -1:nbin+2 ), uadv( 0:nbin+2 )
3323  real(rp) :: uadv ( 0:nbin+2,nspc,ijkmax )
3324  real(rp) :: crn ( 0:nbin+2,nspc,ijkmax )
3325  real(rp) :: flq ( 1:nbin+1 )
3326  real(rp) :: acoef( 0:nbin+1,0:2 )
3327  real(rp) :: aip ( 0:nbin+1 )
3328  real(rp) :: aim ( 0:nbin+1 )
3329  real(rp) :: ai ( 0:nbin+1 )
3330  real(rp) :: cmins, cplus
3331  integer :: nloopmax
3332 
3333  !--- for lithgning component (if not use lightning component, arrays shown below work dummy array)
3334  real(rp) :: gcrgn( -1:nbin+2,nspc,ijkmax )
3335  real(rp) :: flq_c( 1:nbin+1 )
3336  real(rp) :: acoef_c( 0:nbin+1,0:2 )
3337  real(rp) :: aip_c ( 0:nbin+1 )
3338  real(rp) :: aim_c ( 0:nbin+1 )
3339  real(rp) :: ai_c ( 0:nbin+1 )
3340  !---------------------------------------------------------------------------
3341 
3342  call prof_rapstart('_SBM_Icephase', 3)
3343 
3344  !$omp workshare
3345  iflg(:,:) = 0
3346  !$omp end workshare
3347  !$omp parallel do private(csum)
3348  do ijk = 1, ijkmax
3349 
3350  do myu = 2, nspc
3351  csum = 0.0_rp
3352  do n = 1, nbin
3353  csum = csum + gc( n,myu,ijk )*dxmic
3354  enddo
3355  if ( csum > cldmin ) iflg( myu,ijk ) = 1
3356  enddo
3357 
3358  enddo
3359 
3360  ! lhs
3361  call atmos_hydrometeor_lhs_para( ijkmax, 1, ijkmax, temp(:), qlsbl(:) )
3362 
3363  call atmos_saturation_pres2qsat_ice_para( ijkmax, 1, ijkmax, &
3364  temp(:), pres(:), qdry(:), & ! [IN]
3365  qsati(:) ) ! [OUT]
3366 
3367  !$omp workshare
3368  gciold(:) = 0.0_rp
3369  nloop(:) = 0
3370  !$omp end workshare
3371  !$omp parallel do private(ssice,emu,cefd,cefk,festi,f1,f2,uval)
3372  do ijk = 1, ijkmax
3373 
3374  !----- old mass
3375  do myu = 2, nspc
3376  do n = 1, nbin
3377  gciold(ijk) = gciold(ijk) + gc( n,myu,ijk )*dxmic
3378  enddo
3379  enddo
3380 
3381  !----- mass -> number
3382  do myu = 2, nspc
3383  do n = 1, nbin
3384  gcn( n,myu,ijk ) = gc( n,myu,ijk ) * rexpxctr( n )
3385  enddo
3386  gcn( -1,myu,ijk ) = 0.0_rp
3387  gcn( 0,myu,ijk ) = 0.0_rp
3388  gcn( nbin+1,myu,ijk ) = 0.0_rp
3389  gcn( nbin+2,myu,ijk ) = 0.0_rp
3390  end do
3391  if ( flg_lt ) then
3392  do myu = 2, nspc
3393  do n = 1, nbin
3394  gcrgn( n,myu,ijk ) = gcrg( n,myu,ijk )
3395  enddo
3396  gcrgn( -1,myu,ijk ) = 0.0_rp
3397  gcrgn( 0,myu,ijk ) = 0.0_rp
3398  gcrgn( nbin+1,myu,ijk ) = 0.0_rp
3399  gcrgn( nbin+2,myu,ijk ) = 0.0_rp
3400  enddo
3401  end if
3402 
3403 
3404  !--- supersaturation
3405  ssice = qvap(ijk)/qsati(ijk) - 1.0_rp
3406 
3407  emu = afmyu*( bfmyu/( temp(ijk)+cfmyu ) )*( temp(ijk)/tmlt )**1.50_rp
3408  cefd = emu/dens(ijk)
3409  cefk = fct*emu
3410  festi = esat0*exp( qlsbl(ijk)/rvap*( 1.0_rp/tem00 - 1.0_rp/temp(ijk) ) )
3411  f1 = rvap*temp(ijk)/festi/cefd
3412  f2 = qlsbl(ijk)/cefk/temp(ijk)*( qlsbl(ijk)/rvap/temp(ijk) - 1.0_rp )
3413  gtice(ijk) = 4.0_rp*pi/( f1+f2 )
3414  !----- CFL condition
3415  umax(ijk) = 0.0_rp
3416  do myu = 2, nspc
3417  uval = cbnd( 1,myu )*rexpxbnd( 1 )*gtice(ijk)*abs( ssice )
3418  umax(ijk) = max( umax(ijk),uval )
3419  enddo
3420 
3421  dtcnd(ijk) = cflfct*dxmic/umax(ijk)
3422  nloop(ijk) = int( dtime/dtcnd(ijk) ) + 1
3423  dtcnd(ijk) = dtime/nloop(ijk)
3424 
3425  nloop(ijk) = nloop(ijk) * maxval( iflg( 2:nspc,ijk ) ) !--- for determing trivial loop
3426 
3427  enddo
3428  nloopmax = maxval(nloop,1)
3429 
3430 
3431 !OCL LOOP_NOFISSION
3432 !OCL LOOP_NOINTERCHANGE
3433  do ncount = 1, nloopmax
3434 
3435  !$omp parallel do
3436  do ijk = 1, ijkmax
3437  loopflg(ijk) = min( 1, int(nloop(ijk)/ncount) ) ! 0 or 1
3438  enddo
3439 
3440  ! lhs
3441  call atmos_hydrometeor_lhs_para( ijkmax, 1, ijkmax, temp(:), qlsbl(:) )
3442 
3443 !OCL LOOP_NOFUSION
3444  !$omp parallel do private(emu,cefd,cefk,festi,f1,f2)
3445  do ijk = 1, ijkmax
3446  do nn = 1, loopflg(ijk)
3447  emu = afmyu*( bfmyu/( temp(ijk)+cfmyu ) )*( temp(ijk)/tmlt )**1.50_rp
3448  cefd = emu/dens(ijk)
3449  cefk = fct*emu
3450  festi = esat0*exp( qlsbl(ijk)/rvap*( 1.0_rp/tem00 - 1.0_rp/temp(ijk) ) )
3451  f1 = rvap*temp(ijk)/festi/cefd
3452  f2 = qlsbl(ijk)/cefk/temp(ijk)*( qlsbl(ijk)/rvap/temp(ijk) - 1.0_rp )
3453  gtice(ijk) = 4.0_rp*pi/( f1+f2 ) ! G of (A.17) of Suzuki (2004)
3454 
3455  sumice(ijk) = 0.0_rp
3456  do myu = 2, nspc
3457  do n = 1, nbin
3458  sumice(ijk) = sumice(ijk) + gcn( n,myu,ijk )*cctr( n,myu )*dxmic
3459  enddo
3460  enddo
3461 
3462  enddo
3463  enddo
3464 
3465 
3466  call atmos_saturation_pres2qsat_ice_para( ijkmax, 1, ijkmax, &
3467  temp(:), pres(:), qdry(:), & ! [IN]
3468  qsati(:) ) ! [OUT]
3469 
3470 !OCL LOOP_NOFUSION
3471  !$omp parallel do private(ssice,zerosw,qvtmp,cefice,d,sicetnd)
3472  do ijk = 1, ijkmax
3473  do nn = 1, loopflg(ijk)
3474 
3475  !--- supersaturation
3476  ssice = qvap(ijk)/qsati(ijk) - 1.0_rp
3477 
3478  !----- supersaturation tendency
3479  zerosw = 0.5_rp + sign( 0.5_rp,qvap(ijk)-eps ) !--- zerosw = 1 (qv>0), zerosw=0 (qv=0)
3480  qvtmp = qvap(ijk) * zerosw + ( qvap(ijk)+eps ) * ( 1.0_rp-zerosw )
3481 
3482  cefice = ( ssice+1.0_rp )*( 1.0_rp/qvtmp + qlsbl(ijk)*qlsbl(ijk)/cv(ijk)/rvap/temp(ijk)/temp(ijk) )
3483  d = - cefice*sumice(ijk)*gtice(ijk)/dens(ijk) ! d of (A.19) of Suzuki (2004)
3484  d = d + eps * ( 1.0_rp - zerosw ) !--- avoiding division by zero when qv = 0
3485 
3486  sicetnd = zerosw * &
3487  ( ssice * ( exp( d*dtcnd(ijk) )-1.0_rp )/( d*dtcnd(ijk) ) &
3488  * ( 0.5_rp + sign( 0.5_rp,abs(d*dtcnd(ijk)-0.1_rp) ) ) &
3489  + ssice &
3490  * ( 0.5_rp - sign( 0.5_rp,abs(d*dtcnd(ijk)-0.1_rp) ) ) &
3491  ) &
3492  + ssice * ( 1.0_rp - zerosw )
3493  !
3494  !----- change of SDF
3495  do myu = 2, nspc
3496  do mm = 1, iflg( myu,ijk )
3497  !--- advection speed
3498  do n = 1, nbin+1
3499  uadv( n,myu,ijk ) = cbnd( n,myu )*rexpxbnd( n )*gtice(ijk)*sicetnd ! U of eq. (A.18) of Suzuki (2004)
3500  enddo
3501  uadv( 0, myu,ijk ) = 0.0_rp
3502  uadv( nbin+2,myu,ijk ) = 0.0_rp
3503  enddo
3504  enddo
3505 
3506  enddo
3507  enddo
3508 
3509  call prof_rapstart('_SBM_AdvIce', 3)
3510 
3511  !$omp parallel do private(acoef,aip,aim,ai,flq,cplus,cmins)
3512  do ijk = 1, ijkmax
3513  do nn = 1, loopflg(ijk)
3514  do myu = 2, nspc
3515  do mm = 1, iflg( myu,ijk )
3516  do n = 0, nbin+2
3517  crn( n,myu,ijk ) = uadv( n,myu,ijk )*dtcnd(ijk)/dxmic
3518  enddo
3519 
3520  do n = 0, nbin+1
3521  acoef(n,0) = - ( gcn( n+1,myu,ijk )-26.0_rp*gcn( n,myu,ijk )+gcn( n-1,myu,ijk ) ) / 48.0_rp
3522  acoef(n,1) = ( gcn( n+1,myu,ijk ) -gcn( n-1,myu,ijk ) ) / 16.0_rp
3523  acoef(n,2) = ( gcn( n+1,myu,ijk )- 2.0_rp*gcn( n,myu,ijk )+gcn( n-1,myu,ijk ) ) / 48.0_rp
3524  enddo
3525 
3526  do n = 0, nbin+1
3527  cplus = 1.0_rp - ( crn(n+1,myu,ijk) + abs(crn(n+1,myu,ijk)) )
3528 
3529  aip(n) = acoef(n,0) * ( 1.0_rp-cplus**1 ) &
3530  + acoef(n,1) * ( 1.0_rp-cplus**2 ) &
3531  + acoef(n,2) * ( 1.0_rp-cplus**3 )
3532  aip(n) = max( aip(n), 0.0_rp )
3533  enddo
3534 
3535  do n = 0, nbin+1
3536  cmins = 1.0_rp - ( abs(crn(n,myu,ijk)) - crn(n,myu,ijk) )
3537 
3538  aim(n) = acoef(n,0) * ( 1.0_rp-cmins**1 ) &
3539  - acoef(n,1) * ( 1.0_rp-cmins**2 ) &
3540  + acoef(n,2) * ( 1.0_rp-cmins**3 )
3541  aim(n) = max( aim(n), 0.0_rp )
3542  enddo
3543 
3544  do n = 0, nbin+1
3545  ai(n) = acoef(n,0) * 2.0_rp &
3546  + acoef(n,2) * 2.0_rp
3547  ai(n) = max( ai(n), aip(n)+aim(n)+cldmin )
3548  enddo
3549 
3550  do n = 1, nbin+1
3551  flq(n) = ( aip(n-1)/ai(n-1)*gcn( n-1,myu,ijk ) &
3552  - aim(n )/ai(n )*gcn( n ,myu,ijk ) )*dxmic/dtcnd(ijk)
3553  enddo
3554 
3555  do n = 1, nbin
3556  gcn( n,myu,ijk ) = gcn( n,myu,ijk ) - ( flq(n+1)-flq(n) )*dtcnd(ijk)/dxmic
3557  enddo
3558 
3559 !!$ dumm_regene(ijk) = dumm_regene(ijk)+( -flq(1)*dtcnd(ijk)/dxmic ) &
3560 !!$ * min( uadv(1,myu,ijk),0.0_RP )/( uadv(1,myu,ijk) + EPS )
3561  enddo
3562  enddo
3563  enddo
3564  enddo
3565 
3566 
3567  if ( flg_lt ) then
3568 
3569  !$omp parallel do private(acoef_c,aip_c,aim_c,ai_c,flq_c,cplus,cmins)
3570  do ijk = 1, ijkmax
3571  do nn = 1, loopflg(ijk)
3572  do myu = 2, nspc
3573  do mm = 1, iflg( myu,ijk )
3574  do n = 0, nbin+1
3575  acoef_c(n,0) = - ( gcrgn( n+1,myu,ijk )-26.0_rp*gcrgn( n,myu,ijk )+gcrgn( n-1,myu,ijk ) ) / 48.0_rp
3576  acoef_c(n,1) = ( gcrgn( n+1,myu,ijk ) -gcrgn( n-1,myu,ijk ) ) / 16.0_rp
3577  acoef_c(n,2) = ( gcrgn( n+1,myu,ijk )- 2.0_rp*gcrgn( n,myu,ijk )+gcrgn( n-1,myu,ijk ) ) / 48.0_rp
3578  enddo
3579 
3580  do n = 0, nbin+1
3581  cplus = 1.0_rp - ( crn(n+1,myu,ijk) + abs(crn(n+1,myu,ijk)) )
3582 
3583  aip_c(n) = acoef_c(n,0) * ( 1.0_rp-cplus**1 ) &
3584  + acoef_c(n,1) * ( 1.0_rp-cplus**2 ) &
3585  + acoef_c(n,2) * ( 1.0_rp-cplus**3 )
3586  aip_c(n) = max( aip_c(n), 0.0_rp )
3587  enddo
3588 
3589  do n = 0, nbin+1
3590  cmins = 1.0_rp - ( abs(crn(n,myu,ijk)) - crn(n,myu,ijk) )
3591 
3592  aim_c(n) = acoef_c(n,0) * ( 1.0_rp-cmins**1 ) &
3593  - acoef_c(n,1) * ( 1.0_rp-cmins**2 ) &
3594  + acoef_c(n,2) * ( 1.0_rp-cmins**3 )
3595  aim_c(n) = max( aim_c(n), 0.0_rp )
3596  enddo
3597 
3598  do n = 0, nbin+1
3599  ai_c(n) = acoef_c(n,0) * 2.0_rp &
3600  + acoef_c(n,2) * 2.0_rp
3601  ai_c(n) = max( ai_c(n), aip_c(n)+aim_c(n)+cldmin )
3602  enddo
3603 
3604  do n = 1, nbin+1
3605  flq_c(n) = ( aip_c(n-1)/ai_c(n-1)*gcrgn( n-1,myu,ijk ) &
3606  - aim_c(n )/ai_c(n )*gcrgn( n ,myu,ijk ) )/dtcnd(ijk)
3607  enddo
3608 
3609  do n = 1, nbin
3610  gcrgn( n,myu,ijk ) = gcrgn( n,myu,ijk ) - ( flq_c(n+1)-flq_c(n) )*dtcnd(ijk)
3611  enddo
3612  enddo
3613  enddo
3614  enddo
3615  enddo
3616 
3617  end if
3618 
3619  call prof_rapend ('_SBM_AdvIce', 3)
3620 
3621 !OCL LOOP_NOFUSION
3622  !$omp parallel do private(sblmss,dqv)
3623  do ijk = 1, ijkmax
3624  do nn = 1, loopflg(ijk)
3625 
3626  !----- new mass
3627  gcinew(ijk) = 0.0_rp
3628  do myu = 2, nspc
3629  do n = 1, nbin
3630  gcinew(ijk) = gcinew(ijk) + gcn( n,myu,ijk )*expxctr( n )*dxmic
3631  enddo
3632  enddo
3633  !
3634  !----- change of humidity and temperature
3635  sblmss = gcinew(ijk) - gciold(ijk)
3636  dqv = sblmss/dens(ijk)
3637  qvap(ijk) = qvap(ijk) - dqv
3638  temp(ijk) = temp(ijk) + dqv*qlsbl(ijk)/cv(ijk)
3639  cp(ijk) = cp(ijk) + ( cp_ice - cp_vapor ) * dqv
3640  cv(ijk) = cv(ijk) + ( cv_ice - cv_vapor ) * dqv
3641 
3642  gciold(ijk) = gcinew(ijk)
3643  !
3644  !----- continue / end
3645  enddo
3646  enddo
3647 
3648  enddo ! nloop
3649 
3650 
3651  !
3652  !$omp parallel do
3653  do ijk = 1, ijkmax
3654  !------- number -> mass
3655  do myu = 2, nspc
3656  do n = 1, nbin
3657  gc( n,myu,ijk ) = gcn( n,myu,ijk )*expxctr( n )
3658  enddo
3659  enddo
3660  enddo
3661  if ( flg_lt ) then
3662  !$omp parallel do
3663  do ijk = 1, ijkmax
3664  !------- number -> mass
3665  do myu = 2, nspc
3666  do n = 1, nbin
3667  gcrg( n,myu,ijk ) = gcrgn( n,myu,ijk )
3668  enddo
3669  enddo
3670  enddo
3671  end if
3672 
3673  call prof_rapend ('_SBM_Icephase', 3)
3674 
3675  return
3676  end subroutine icephase
3677 
3678  !-----------------------------------------------------------------------------
3679  subroutine mixphase( &
3680  ijkmax, &
3681  flg_lt, &
3682  dens, &
3683  pres, &
3684  qdry, &
3685  temp, &
3686  qvap, &
3687  gc, &
3688  gcrg, &
3689  cp, &
3690  cv, &
3691  dtime )
3692  use scale_const, only: &
3693  pi => const_pi, &
3694  eps => const_eps, &
3695  esat0 => const_psat0, &
3696  tem00 => const_tem00, &
3697  rvap => const_rvap, &
3698  tmlt => const_tmelt
3699  use scale_atmos_hydrometeor, only: &
3700  atmos_hydrometeor_lhv_para, &
3701  atmos_hydrometeor_lhs_para, &
3702  cp_vapor, &
3703  cp_water, &
3704  cp_ice, &
3705  cv_vapor, &
3706  cv_water, &
3707  cv_ice
3708  use scale_atmos_saturation, only: &
3709  atmos_saturation_pres2qsat_liq_para, &
3710  atmos_saturation_pres2qsat_ice_para
3711  implicit none
3712 
3713  integer, intent(in) :: ijkmax
3714  logical, intent(in) :: flg_lt
3715  real(rp), intent(in) :: dens(ijkmax) ! Density [kg/m3]
3716  real(rp), intent(in) :: pres(ijkmax) ! Pressure [Pa]
3717  real(rp), intent(in) :: qdry(ijkmax) ! dry air mass ratio [kg/kg]
3718  real(rp), intent(inout) :: temp(ijkmax) ! Temperature [K]
3719  real(rp), intent(inout) :: qvap(ijkmax) ! Specific humidity [kg/kg]
3720  real(rp), intent(inout) :: gc (nbin,nspc,ijkmax) ! Mass size distribution function of hydrometeor
3721  real(rp), intent(inout) :: gcrg(nbin,nspc,ijkmax)
3722  real(rp), intent(inout) :: cp(ijkmax) ! specific heat
3723  real(rp), intent(inout) :: cv(ijkmax) ! specific heat
3724  real(dp), intent(in) :: dtime ! Time step interval
3725 
3726  integer :: n, myu, mm, ncount
3727  integer :: nloop(ijkmax)
3728  real(rp) :: gclold(ijkmax), gclnew(ijkmax)
3729  real(rp) :: gciold(ijkmax), gcinew(ijkmax)
3730  real(rp) :: cndmss, sblmss
3731  real(rp) :: gtliq(ijkmax), gtice(ijkmax)
3732  real(rp) :: umax(ijkmax), uval, dtcnd(ijkmax)
3733  real(rp) :: qlevp(ijkmax), qlsbl(ijkmax)
3734  real(rp) :: cef1, cef2, cef3, cef4, a, b, c, d
3735  real(rp) :: rmdplus, rmdmins, ssplus, ssmins, tplus, tmins
3736  real(rp) :: sliqtnd, sicetnd
3737  real(rp) :: ssliq, ssice
3738  real(rp) :: sumliq(ijkmax), sumice(ijkmax)
3739  real(rp) :: gcn( -1:nbin+2,nspc,ijkmax )
3740  real(rp), parameter :: cflfct = 0.50_rp
3741 ! real(RP) :: dumm_regene(ijkmax)
3742  real(rp) :: csum
3743  integer :: iflg( nspc,ijkmax )
3744  real(rp) :: f1, f2, emu, cefd, cefk, festl, festi
3745  real(rp) :: qsatl(ijkmax), qsati(ijkmax)
3746  real(rp), parameter :: afmyu = 1.72e-05_rp, bfmyu = 3.93e+2_rp
3747  real(rp), parameter :: cfmyu = 1.2e+02_rp, fct = 1.4e+3_rp
3748  real(rp) :: qvtmp, zerosw
3749  integer :: ijk, nn, loopflg(ijkmax)
3750 
3751  !--- local for advection
3752 ! real(RP) :: qadv( -1:nbin+2 ), uadv( 0:nbin+2 )
3753  real(rp) :: uadv ( 0:nbin+2,nspc,ijkmax )
3754  real(rp) :: crn ( 0:nbin+2,nspc,ijkmax )
3755  real(rp) :: flq ( 1:nbin+1 )
3756  real(rp) :: acoef( 0:nbin+1,0:2 )
3757  real(rp) :: aip ( 0:nbin+1 )
3758  real(rp) :: aim ( 0:nbin+1 )
3759  real(rp) :: ai ( 0:nbin+1 )
3760  real(rp) :: cmins, cplus
3761  integer :: nloopmax
3762 
3763  !--- for lithgning component (if not use lightning component, arrays shown below work dummy array)
3764  real(rp) :: gcrgn( -1:nbin+2,nspc,ijkmax )
3765  real(rp) :: flq_c( 1:nbin+1 )
3766  real(rp) :: acoef_c( 0:nbin+1,0:2 )
3767  real(rp) :: aip_c ( 0:nbin+1 )
3768  real(rp) :: aim_c ( 0:nbin+1 )
3769  real(rp) :: ai_c ( 0:nbin+1 )
3770 
3771  call prof_rapstart('_SBM_Mixphase', 3)
3772 
3773  !$omp workshare
3774 ! dumm_regene( : ) = 0.0_RP
3775  iflg( :,: ) = 0
3776  !$omp end workshare
3777  !$omp parallel do private(csum)
3778  do ijk = 1, ijkmax
3779 
3780  do myu = 1, nspc
3781  csum = 0.0_rp
3782  do n = 1, nbin
3783  csum = csum + gc( n,myu,ijk )*dxmic
3784  enddo
3785  if ( csum > cldmin ) iflg( myu,ijk ) = 1
3786  enddo
3787 
3788  enddo
3789 
3790  ! lhv
3791  call atmos_hydrometeor_lhv_para( ijkmax, 1, ijkmax, temp(:), qlevp(:) )
3792  ! lhs
3793  call atmos_hydrometeor_lhs_para( ijkmax, 1, ijkmax, temp(:), qlsbl(:) )
3794 
3795  call atmos_saturation_pres2qsat_liq_para( ijkmax, 1, ijkmax, &
3796  temp(:), pres(:), qdry(:), & ! [IN]
3797  qsatl(:) ) ! [OUT]
3798  call atmos_saturation_pres2qsat_ice_para( ijkmax, 1, ijkmax, &
3799  temp(:), pres(:), qdry(:), & ! [IN]
3800  qsati(:) ) ! [OUT]
3801 
3802  !$omp workshare
3803  gclold(:) = 0.0_rp
3804  gciold(:) = 0.0_rp
3805  nloop(:) = 0
3806  !$omp end workshare
3807  !$omp parallel do private(ssliq,ssice,emu,cefd,cefk,festl,f1,f2,festi,uval)
3808  do ijk = 1, ijkmax
3809  !----- old mass
3810  do n = 1, nbin
3811  gclold(ijk) = gclold(ijk) + gc( n,il,ijk )*dxmic
3812  enddo
3813 
3814  do myu = 2, nspc
3815  do n = 1, nbin
3816  gciold(ijk) = gciold(ijk) + gc( n,myu,ijk )*dxmic
3817  enddo
3818  enddo
3819 
3820  !----- mass -> number
3821  do myu = 1, nspc
3822  do n = 1, nbin
3823  gcn( n,myu,ijk ) = gc( n,myu,ijk ) * rexpxctr( n )
3824  gcrgn( n,myu,ijk ) = gcrg( n,myu,ijk )
3825  enddo
3826  gcn( -1,myu,ijk ) = 0.0_rp
3827  gcn( 0,myu,ijk ) = 0.0_rp
3828  gcn( nbin+1,myu,ijk ) = 0.0_rp
3829  gcn( nbin+2,myu,ijk ) = 0.0_rp
3830  gcrgn( -1,myu,ijk ) = 0.0_rp
3831  gcrgn( 0,myu,ijk ) = 0.0_rp
3832  gcrgn( nbin+1,myu,ijk ) = 0.0_rp
3833  gcrgn( nbin+2,myu,ijk ) = 0.0_rp
3834  enddo
3835 
3836  !-- supersaturation
3837  ssliq = qvap(ijk)/qsatl(ijk) - 1.0_rp
3838  ssice = qvap(ijk)/qsati(ijk) - 1.0_rp
3839 
3840  emu = afmyu*( bfmyu/( temp(ijk)+cfmyu ) )*( temp(ijk)/tmlt )**1.50_rp
3841  cefd = emu/dens(ijk)
3842  cefk = fct*emu
3843 
3844  festl = esat0*exp( qlevp(ijk)/rvap*( 1.0_rp/tem00 - 1.0_rp/temp(ijk) ) )
3845  f1 = rvap*temp(ijk)/festl/cefd
3846  f2 = qlevp(ijk)/cefk/temp(ijk)*( qlevp(ijk)/rvap/temp(ijk) - 1.0_rp )
3847  gtliq(ijk) = 4.0_rp*pi/( f1+f2 )
3848 
3849  festi = esat0*exp( qlsbl(ijk)/rvap*( 1.0_rp/tem00 - 1.0_rp/temp(ijk) ) )
3850  f1 = rvap*temp(ijk)/festi/cefd
3851  f2 = qlsbl(ijk)/cefk/temp(ijk)*( qlsbl(ijk)/rvap/temp(ijk) - 1.0_rp )
3852  gtice(ijk) = 4.0_rp*pi/( f1+f2 )
3853 
3854  !----- CFL condition
3855  umax(ijk) = cbnd( 1,il )*rexpxbnd( 1 )*gtliq(ijk)*abs( ssliq )
3856  do myu = 2, nspc
3857  uval = cbnd( 1,myu )*rexpxbnd( 1 )*gtice(ijk)*abs( ssice )
3858  umax(ijk) = max( umax(ijk),uval )
3859  enddo
3860 
3861  dtcnd(ijk) = cflfct*dxmic/umax(ijk)
3862  nloop(ijk) = int( dtime/dtcnd(ijk) ) + 1
3863  dtcnd(ijk) = dtime/nloop(ijk)
3864 
3865  nloop(ijk) = nloop(ijk) * iflg( il,ijk ) !--- for determing trivial loop
3866  nloop(ijk) = nloop(ijk) * maxval( iflg( 2:nspc,ijk ) ) !--- for determing trivial loop
3867 ! nloop(ijk) = nloop(ijk) * maxval( iflg( 1:nspc,ijk ) ) !--- for determing trivial loop
3868  enddo
3869  nloopmax = maxval(nloop,1)
3870 
3871 
3872 !OCL LOOP_NOFISSION
3873 !OCL LOOP_NOINTERCHANGE
3874  do ncount = 1, nloopmax
3875 
3876  !$omp parallel do
3877  do ijk = 1, ijkmax
3878  loopflg(ijk) = min( 1, int(nloop(ijk)/ncount) ) ! 0 or 1
3879  enddo
3880 
3881  ! lhv
3882  call atmos_hydrometeor_lhv_para( ijkmax, 1, ijkmax, temp(:), qlevp(:) )
3883  ! lhs
3884  call atmos_hydrometeor_lhs_para( ijkmax, 1, ijkmax, temp(:), qlsbl(:) )
3885 
3886 
3887 !OCL LOOP_NOFUSION
3888  !$omp parallel do private(emu,cefd,cefk,festl,f1,f2,festi)
3889  do ijk = 1, ijkmax
3890  do nn = 1, loopflg(ijk)
3891 
3892  !-- matrix for supersaturation tendency
3893  emu = afmyu*( bfmyu/( temp(ijk)+cfmyu ) )*( temp(ijk)/tmlt )**1.50_rp
3894  cefd = emu/dens(ijk)
3895  cefk = fct*emu
3896  festl = esat0*exp( qlevp(ijk)/rvap*( 1.0_rp/tem00 - 1.0_rp/temp(ijk) ) )
3897  f1 = rvap*temp(ijk)/festl/cefd
3898  f2 = qlevp(ijk)/cefk/temp(ijk)*( qlevp(ijk)/rvap/temp(ijk) - 1.0_rp )
3899  gtliq(ijk) = 4.0_rp*pi/( f1+f2 )
3900  festi = esat0*exp( qlsbl(ijk)/rvap*( 1.0_rp/tem00 - 1.0_rp/temp(ijk) ) )
3901  f1 = rvap*temp(ijk)/festi/cefd
3902  f2 = qlsbl(ijk)/cefk/temp(ijk)*( qlsbl(ijk)/rvap/temp(ijk) - 1.0_rp )
3903  gtice(ijk) = 4.0_rp*pi/( f1+f2 )
3904 
3905  sumliq(ijk) = 0.0_rp
3906  do n = 1, nbin
3907  sumliq(ijk) = sumliq(ijk) + gcn( n,il,ijk )*cctr( n,il )*dxmic
3908  enddo
3909 
3910  sumice(ijk) = 0.0_rp
3911  do myu = 2, nspc
3912  do n = 1, nbin
3913  sumice(ijk) = sumice(ijk) + gcn( n,myu,ijk )*cctr( n,myu )*dxmic
3914  enddo
3915  enddo
3916  enddo
3917  enddo
3918 
3919 
3920  call atmos_saturation_pres2qsat_liq_para( ijkmax, 1, ijkmax, &
3921  temp(:), pres(:), qdry(:), & ! [IN]
3922  qsatl(:) ) ! [OUT]
3923  call atmos_saturation_pres2qsat_ice_para( ijkmax, 1, ijkmax, &
3924  temp(:), pres(:), qdry(:), & ! [IN]
3925  qsati(:) ) ! [OUT]
3926 !OCL LOOP_NOFUSION
3927  !$omp parallel do private(ssliq,ssice,zerosw,qvtmp,cef1,cef2,cef3,cef4,a,b,c,d,rmdplus,rmdmins,ssplus,ssmins,tplus,tmins,sliqtnd,sicetnd)
3928  do ijk = 1, ijkmax
3929  do nn = 1, loopflg(ijk)
3930 
3931  !-- supersaturation
3932  ssliq = qvap(ijk)/qsatl(ijk) - 1.0_rp
3933  ssice = qvap(ijk)/qsati(ijk) - 1.0_rp
3934 
3935  zerosw = 0.5_rp + sign( 0.5_rp,qvap(ijk)-eps ) !--- zerosw = 1 (qv>0), zerosw=0 (qv=0)
3936  qvtmp = qvap(ijk) * zerosw + ( qvap(ijk)+eps ) * ( 1.0_rp-zerosw )
3937  cef1 = ( ssliq+1.0_rp )*( 1.0_rp/qvtmp + qlevp(ijk)/rvap/temp(ijk)/temp(ijk)*qlevp(ijk)/cv(ijk) )
3938  cef2 = ( ssliq+1.0_rp )*( 1.0_rp/qvtmp + qlevp(ijk)/rvap/temp(ijk)/temp(ijk)*qlsbl(ijk)/cv(ijk) )
3939  cef3 = ( ssice+1.0_rp )*( 1.0_rp/qvtmp + qlsbl(ijk)/rvap/temp(ijk)/temp(ijk)*qlevp(ijk)/cv(ijk) )
3940  cef4 = ( ssice+1.0_rp )*( 1.0_rp/qvtmp + qlsbl(ijk)/rvap/temp(ijk)/temp(ijk)*qlsbl(ijk)/cv(ijk) )
3941 
3942  a = - cef1*sumliq(ijk)*gtliq(ijk)/dens(ijk) ! a of (A.19) of Suzuki (2004)
3943  b = - cef2*sumice(ijk)*gtice(ijk)/dens(ijk) ! b of (A.19) of Suzuki (2004)
3944  c = - cef3*sumliq(ijk)*gtliq(ijk)/dens(ijk) ! c of (A.19) of Suzuki (2004)
3945  d = - cef4*sumice(ijk)*gtice(ijk)/dens(ijk) ! d of (A.19) of Suzuki (2004)
3946 
3947  b = b + eps * ( 1.0_rp - zerosw ) !--- avoiding division by zero when qv = 0
3948  !--- eigenvalues
3949  rmdplus = ( ( a+d ) + sqrt( ( a-d )**2 + 4.0_rp*b*c ) ) * 0.50_rp
3950  rmdmins = ( ( a+d ) - sqrt( ( a-d )**2 + 4.0_rp*b*c ) ) * 0.50_rp
3951 
3952  rmdplus = rmdplus + eps * ( 1.0_rp - zerosw ) !--- avoiding division by zero when qv = 0
3953  rmdmins = rmdmins + eps * ( 1.0_rp - zerosw ) !--- avoiding division by zero when qv = 0
3954 
3955  !--- supersaturation tendency
3956  ssplus = ( ( rmdmins-a )*ssliq - b*ssice )/b/( rmdmins-rmdplus + eps * ( 1.0_rp - zerosw ) )
3957  ssmins = ( ( a-rmdplus )*ssliq + b*ssice )/b/( rmdmins-rmdplus + eps * ( 1.0_rp - zerosw ) )
3958 
3959  tplus = ( exp( rmdplus*dtcnd(ijk) )-1.0_rp )/( rmdplus*dtcnd(ijk) ) &
3960  * ( 0.5_rp + sign( 0.5_rp,abs(rmdplus*dtcnd(ijk)-0.1_rp) ) ) &
3961  + 1.0_rp &
3962  * ( 0.5_rp - sign( 0.5_rp,abs(rmdplus*dtcnd(ijk)-0.1_rp) ) )
3963  tmins = ( exp( rmdmins*dtcnd(ijk) )-1.0_rp )/( rmdmins*dtcnd(ijk) ) &
3964  * ( 0.5_rp + sign( 0.5_rp,abs(rmdmins*dtcnd(ijk)-0.1_rp) ) ) &
3965  + 1.0_rp &
3966  * ( 0.5_rp - sign( 0.5_rp,abs(rmdmins*dtcnd(ijk)-0.1_rp) ) )
3967 
3968  sliqtnd = ssliq * ( 1.0_rp - zerosw ) &
3969  + zerosw * &
3970  ( b*tplus*ssplus + b*tmins*ssmins ) ! sliwtnd in page 116 of Suzuki (2004)
3971  sicetnd = ssice * ( 1.0_rp - zerosw ) &
3972  + zerosw * &
3973  ( ( rmdplus-a )*tplus*ssplus & ! sicetnd in page 116 of Suzuki (2004)
3974  + ( rmdmins-a )*tmins*ssmins )
3975 
3976  !--- change of SDF
3977  do myu = 1, nspc
3978  do mm = 1, iflg( myu,ijk )
3979  !--- advection speed
3980  do n = 1, nbin+1
3981  !--- myu = 1 -> ssliq, myu > 1 -> ssice
3982  uadv( n,myu,ijk ) = cbnd( n,myu )*rexpxbnd( n )*gtliq(ijk)*sliqtnd & ! U of eq. (A.18) of Suzuki (2004)
3983  * ( 0.5_rp - sign( 0.5_rp,real(myu)-1.5_rp) ) &
3984  + cbnd( n,myu )*rexpxbnd( n )*gtice(ijk)*sicetnd & ! U of eq. (A.18) of Suzuki (2004)
3985  * ( 0.5_rp + sign( 0.5_rp,real(myu)-1.5_rp) )
3986  enddo
3987  uadv( 0, myu,ijk ) = 0.0_rp
3988  uadv( nbin+2,myu,ijk ) = 0.0_rp
3989  enddo
3990  enddo
3991  enddo
3992  enddo
3993 
3994  call prof_rapstart('_SBM_AdvMix', 3)
3995 
3996  !$omp parallel do private(acoef,aip,aim,ai,flq,cplus,cmins)
3997  do ijk = 1, ijkmax
3998  do nn = 1, loopflg(ijk)
3999  do myu = 1, nspc
4000  do mm = 1, iflg( myu,ijk )
4001  do n = 0, nbin+2
4002  crn( n,myu,ijk ) = uadv( n,myu,ijk )*dtcnd(ijk)/dxmic
4003  enddo
4004 
4005  do n = 0, nbin+1
4006  acoef(n,0) = - ( gcn( n+1,myu,ijk )-26.0_rp*gcn( n,myu,ijk )+gcn( n-1,myu,ijk ) ) / 48.0_rp
4007  acoef(n,1) = ( gcn( n+1,myu,ijk ) -gcn( n-1,myu,ijk ) ) / 16.0_rp
4008  acoef(n,2) = ( gcn( n+1,myu,ijk )- 2.0_rp*gcn( n,myu,ijk )+gcn( n-1,myu,ijk ) ) / 48.0_rp
4009  enddo
4010 
4011  do n = 0, nbin+1
4012  cplus = 1.0_rp - ( crn(n+1,myu,ijk) + abs(crn(n+1,myu,ijk)) )
4013 
4014  aip(n) = acoef(n,0) * ( 1.0_rp-cplus**1 ) &
4015  + acoef(n,1) * ( 1.0_rp-cplus**2 ) &
4016  + acoef(n,2) * ( 1.0_rp-cplus**3 )
4017  aip(n) = max( aip(n), 0.0_rp )
4018  enddo
4019 
4020  do n = 0, nbin+1
4021  cmins = 1.0_rp - ( abs(crn(n,myu,ijk)) - crn(n,myu,ijk) )
4022 
4023  aim(n) = acoef(n,0) * ( 1.0_rp-cmins**1 ) &
4024  - acoef(n,1) * ( 1.0_rp-cmins**2 ) &
4025  + acoef(n,2) * ( 1.0_rp-cmins**3 )
4026  aim(n) = max( aim(n), 0.0_rp )
4027  enddo
4028 
4029  do n = 0, nbin+1
4030  ai(n) = acoef(n,0) * 2.0_rp &
4031  + acoef(n,2) * 2.0_rp
4032  ai(n) = max( ai(n), aip(n)+aim(n)+cldmin )
4033  enddo
4034 
4035  do n = 1, nbin+1
4036  flq(n) = ( aip(n-1)/ai(n-1)*gcn( n-1,myu,ijk ) &
4037  - aim(n )/ai(n )*gcn( n ,myu,ijk ) )*dxmic/dtcnd(ijk)
4038  enddo
4039 
4040  do n = 1, nbin
4041  gcn( n,myu,ijk ) = gcn( n,myu,ijk ) - ( flq(n+1)-flq(n) )*dtcnd(ijk)/dxmic
4042  enddo
4043  enddo
4044  enddo
4045  enddo
4046  enddo
4047 
4048 !!$ !$omp parallel do
4049 !!$ do ijk = 1, ijkmax
4050 !!$ do nn = 1, loopflg(ijk)
4051 !!$ do myu = 1, nspc
4052 !!$ do mm = 1, iflg( myu,ijk )
4053 !!$ dumm_regene(ijk) = dumm_regene(ijk)+( -flq(1,myu,ijk)*dtcnd(ijk)/dxmic ) &
4054 !!$ * min( uadv(1,myu,ijk),0.0_RP )/( uadv(1,myu,ijk)+EPS )
4055 !!$ enddo
4056 !!$ enddo
4057 !!$ enddo
4058 !!$ enddo
4059 
4060  if ( flg_lt ) then
4061  !$omp parallel do private(acoef_c,aip_c,aim_c,ai_c,flq_c,cplus,cmins)
4062  do ijk = 1, ijkmax
4063  do nn = 1, loopflg(ijk)
4064  do myu = 1, nspc
4065  do mm = 1, iflg( myu,ijk )
4066  do n = 0, nbin+1
4067  acoef_c(n,0) = - ( gcrgn( n+1,myu,ijk )-26.0_rp*gcrgn( n,myu,ijk )+gcrgn( n-1,myu,ijk ) ) / 48.0_rp
4068  acoef_c(n,1) = ( gcrgn( n+1,myu,ijk ) -gcrgn( n-1,myu,ijk ) ) / 16.0_rp
4069  acoef_c(n,2) = ( gcrgn( n+1,myu,ijk )- 2.0_rp*gcrgn( n,myu,ijk )+gcrgn( n-1,myu,ijk ) ) / 48.0_rp
4070  enddo
4071 
4072  do n = 0, nbin+1
4073  cplus = 1.0_rp - ( crn(n+1,myu,ijk) + abs(crn(n+1,myu,ijk)) )
4074 
4075  aip_c(n) = acoef_c(n,0) * ( 1.0_rp-cplus**1 ) &
4076  + acoef_c(n,1) * ( 1.0_rp-cplus**2 ) &
4077  + acoef_c(n,2) * ( 1.0_rp-cplus**3 )
4078  aip_c(n) = max( aip_c(n), 0.0_rp )
4079  enddo
4080 
4081  do n = 0, nbin+1
4082  cmins = 1.0_rp - ( abs(crn(n,myu,ijk)) - crn(n,myu,ijk) )
4083 
4084  aim_c(n) = acoef_c(n,0) * ( 1.0_rp-cmins**1 ) &
4085  - acoef_c(n,1) * ( 1.0_rp-cmins**2 ) &
4086  + acoef_c(n,2) * ( 1.0_rp-cmins**3 )
4087  aim_c(n) = max( aim_c(n), 0.0_rp )
4088  enddo
4089 
4090  do n = 0, nbin+1
4091  ai_c(n) = acoef_c(n,0) * 2.0_rp &
4092  + acoef_c(n,2) * 2.0_rp
4093  ai_c(n) = max( ai_c(n), aip_c(n)+aim_c(n)+cldmin )
4094  enddo
4095 
4096  do n = 1, nbin+1
4097  flq_c(n) = ( aip_c(n-1)/ai_c(n-1)*gcrgn( n-1,myu,ijk ) &
4098  - aim_c(n )/ai_c(n )*gcrgn( n ,myu,ijk ) )/dtcnd(ijk)
4099  enddo
4100 
4101  do n = 1, nbin
4102  gcrgn( n,myu,ijk ) = gcrgn( n,myu,ijk ) - ( flq_c(n+1)-flq_c(n) )*dtcnd(ijk)
4103  enddo
4104  enddo
4105  enddo
4106  enddo
4107  enddo
4108  end if
4109 
4110  call prof_rapend ('_SBM_AdvMix', 3)
4111 
4112 !OCL LOOP_NOFUSION
4113  !$omp parallel do private(cndmss,sblmss)
4114  do ijk = 1, ijkmax
4115  do nn = 1, loopflg(ijk)
4116  !--- new mass
4117  gclnew(ijk) = 0.0_rp
4118  do n = 1, nbin
4119  gclnew(ijk) = gclnew(ijk) + gcn( n,il,ijk )*expxctr( n )*dxmic
4120  enddo
4121 
4122  gcinew(ijk) = 0.0_rp
4123  do myu = 2, nspc
4124  do n = 1, nbin
4125  gcinew(ijk) = gcinew(ijk) + gcn( n,myu,ijk )*expxctr( n )*dxmic
4126  enddo
4127  enddo
4128 
4129  !--- change of humidity and temperature
4130  cndmss = gclnew(ijk) - gclold(ijk)
4131  sblmss = gcinew(ijk) - gciold(ijk)
4132 
4133  qvap(ijk) = qvap(ijk) - ( cndmss + sblmss ) / dens(ijk)
4134  temp(ijk) = temp(ijk) + ( cndmss*qlevp(ijk)+sblmss*qlsbl(ijk) ) / dens(ijk) / cv(ijk)
4135  cp(ijk) = cp(ijk) + ( ( cp_water - cp_vapor ) * cndmss + ( cp_ice - cp_vapor ) * sblmss ) / dens(ijk)
4136  cv(ijk) = cv(ijk) + ( ( cv_water - cv_vapor ) * cndmss + ( cv_ice - cv_vapor ) * sblmss ) / dens(ijk)
4137 
4138  gclold(ijk) = gclnew(ijk)
4139  gciold(ijk) = gcinew(ijk)
4140  enddo
4141  enddo
4142 
4143  enddo ! ncount
4144 
4145  !$omp parallel do
4146  do ijk = 1, ijkmax
4147  !----- number -> mass
4148  do myu = 1, nspc
4149  do n = 1, nbin
4150  gc( n,myu,ijk ) = gcn( n,myu,ijk )*expxctr( n )
4151  enddo
4152  enddo
4153  enddo
4154 
4155  if ( flg_lt ) then
4156  !$omp parallel do
4157  do ijk = 1, ijkmax
4158  !----- number -> mass
4159  do myu = 1, nspc
4160  do n = 1, nbin
4161  gcrg( n,myu,ijk ) = gcrgn( n,myu,ijk )
4162  enddo
4163  enddo
4164  enddo
4165  end if
4166 
4167 ! if ( .NOT. flg_regeneration ) then
4168 ! dumm_regene(:) = 0.0_RP
4169 ! endif
4170 
4171  call prof_rapend ('_SBM_Mixphase', 3)
4172 
4173  return
4174  end subroutine mixphase
4175 
4176  !-----------------------------------------------------------------------------
4177  subroutine ice_nucleat( &
4178  ijkmax, &
4179  num_cold, &
4180  index_cold, &
4181  dens, &
4182  pres, &
4183  qdry, &
4184  temp, &
4185  qvap, &
4186  gc, &
4187  cp, &
4188  cv, &
4189  dtime )
4190  use scale_const, only: &
4191  qlmlt => const_emelt
4192  use scale_atmos_saturation, only: &
4193  atmos_saturation_pres2qsat_ice_para
4194  use scale_atmos_hydrometeor, only: &
4195  atmos_hydrometeor_lhs_para, &
4196  cp_vapor, &
4197  cp_ice, &
4198  cv_vapor, &
4199  cv_ice
4200  implicit none
4201 
4202  integer, intent(in) :: ijkmax
4203  integer, intent(in) :: num_cold
4204  integer, intent(in) :: index_cold(ijkmax)
4205  real(rp), intent(in) :: dens(ijkmax) ! Density [kg/m3]
4206  real(rp), intent(in) :: pres(ijkmax) ! Pressure [Pa]
4207  real(rp), intent(in) :: qdry(ijkmax) ! dry air mass ratio [kg/kg]
4208  real(rp), intent(inout) :: temp(ijkmax) ! Temperature [K]
4209  real(rp), intent(inout) :: qvap(ijkmax) ! Specific humidity [kg/kg]
4210  real(rp), intent(inout) :: gc (nbin,nspc,ijkmax) ! Mass size distribution function of hydrometeor
4211  real(rp), intent(inout) :: cp(ijkmax) ! specific heat
4212  real(rp), intent(inout) :: cv(ijkmax) ! specific heat
4213  real(dp), intent(in) :: dtime ! Time step interval
4214 
4215  real(rp) :: ssice
4216  real(rp) :: numin, tdel, qdel
4217 ! real(RP), parameter :: n0 = 1.E+3_RP ! N_{IN0} in page 112 of Suzuki (2004)
4218  real(rp), parameter :: acoef = -0.639_rp, bcoef = 12.96_rp ! A and B in paeg 112 of Suzuki (2004)
4219  ! threshould to determine the type of freezed hydrometeor (the detail is described in page 113 of Suzuki(2004))
4220  real(rp), parameter :: tcolmu = 269.0_rp, tcolml = 265.0_rp! -4[degC], -8[degC]
4221  real(rp), parameter :: tdendu = 257.0_rp, tdendl = 255.0_rp! -14[degC], -18[degC]
4222  real(rp), parameter :: tplatu = 250.6_rp ! -22.4[degC]
4223  !
4224  real(rp) :: qsati(ijkmax), qlsbl(ijkmax)
4225  integer :: ijk, indirect
4226  real(rp) :: numice
4227  integer :: myu, n
4228  integer :: ispc
4229 
4230  call prof_rapstart('_SBM_IceNucleat', 3)
4231 
4232  call atmos_saturation_pres2qsat_ice_para( ijkmax, 1, ijkmax, &
4233  temp(:), pres(:), qdry(:), & ! [IN]
4234  qsati(:) ) ! [OUT]
4235  ! lhs
4236  call atmos_hydrometeor_lhs_para( ijkmax, 1, ijkmax, temp(:), qlsbl(:) )
4237 
4238  !$omp parallel do private(ijk,ssice,numice,numin,ispc,tdel,qdel)
4239  do indirect = 1, num_cold
4240  ijk = index_cold(indirect)
4241 
4242  !--- supersaturation
4243  ssice = qvap(ijk)/qsati(ijk) - 1.0_rp
4244 
4245  if( ssice <= 0.0_rp ) cycle
4246 
4247  numice = 0.0_rp
4248  do myu = ic, ih
4249  do n = 1, nbin
4250  numice = numice + gc( n,myu,ijk )*rexpxctr( n )*dxmic
4251  enddo
4252  enddo
4253 
4254  numin = n0_icenucl * exp( acoef + bcoef * ssice * 1.e+2_rp )
4255  if( numin > numice ) then
4256  !--- -4 [deg] > T >= -8 [deg] and T < -22.4 [deg] -> column
4257  if ( temp(ijk) <= tplatu .OR. ( temp(ijk) >= tcolml .AND. temp(ijk) < tcolmu ) ) then
4258  ispc = ic
4259  !--- -14 [deg] > T >= -18 [deg] -> dendrite
4260  elseif( temp(ijk) <= tdendu .AND. temp(ijk) >= tdendl ) then
4261  ispc = id
4262  !--- else -> plate
4263  else
4264  ispc = ip
4265  endif
4266 
4267  numin = (numin-numice) * expxctr( 1 )
4268  numin = min( numin,qvap(ijk)*dens(ijk) )
4269  gc( 1,ispc,ijk ) = gc( 1,ispc,ijk ) + numin / dxmic
4270 
4271  tdel = numin/dens(ijk)*qlsbl(ijk)/cv(ijk)
4272  temp(ijk) = temp(ijk) + tdel
4273  qdel = numin/dens(ijk)
4274  qvap(ijk) = qvap(ijk) - qdel
4275  cp(ijk) = cp(ijk) + ( cp_ice - cp_vapor ) * qdel
4276  cv(ijk) = cv(ijk) + ( cv_ice - cv_vapor ) * qdel
4277  endif
4278 
4279  enddo
4280 
4281 
4282  call prof_rapend ('_SBM_IceNucleat', 3)
4283 
4284  return
4285  end subroutine ice_nucleat
4286 
4287  !-----------------------------------------------------------------------------
4288  subroutine freezing( &
4289  ijkmax, &
4290  num_cold, &
4291  index_cold, &
4292  flg_lt, &
4293  dens, &
4294  temp, &
4295  gc, &
4296  gcrg, &
4297  cp, &
4298  cv, &
4299  dtime )
4300  use scale_const, only: &
4301  eps => const_eps, &
4302  pi => const_pi, &
4303  tmlt => const_tmelt, &
4304  rhow => const_dwatr
4305  use scale_atmos_hydrometeor, only: &
4306  atmos_hydrometeor_lhf_para, &
4307  cp_water, &
4308  cp_ice, &
4309  cv_water, &
4310  cv_ice
4311  implicit none
4312 
4313  integer, intent(in) :: ijkmax
4314  integer, intent(in) :: num_cold
4315  integer, intent(in) :: index_cold(ijkmax)
4316  logical, intent(in) :: flg_lt
4317  real(rp), intent(in) :: dens(ijkmax) ! Density [kg/m3]
4318  real(rp), intent(inout) :: temp(ijkmax) ! Temperature [K]
4319  real(rp), intent(inout) :: gc (nbin,nspc,ijkmax) ! Mass size distribution function of hydrometeor
4320  real(rp), intent(inout) :: gcrg(nbin,nspc,ijkmax)
4321  real(rp), intent(inout) :: cp(ijkmax) ! specific heat
4322  real(rp), intent(inout) :: cv(ijkmax) ! specific heat
4323  real(dp), intent(in) :: dtime ! Time step interval
4324 
4325  integer :: nbound, n
4326  real(rp) :: xbound, tc, rate, dmp, frz, sumfrz, tdel
4327  real(rp), parameter :: coefa = 1.0e-01_rp ! a_{fr} of eq.(3.19) of Suzuki (2004)
4328  real(rp), parameter :: coefb = 0.66_rp ! b_{fr} of eq.(3.19) of Suzuki (2004)
4329  real(rp), parameter :: rbound = 2.0e-04_rp ! 200 um
4330 ! real(RP), parameter :: tthreth = 235.0_RP ! -38 [degC] threshold for using Bigg's parameterization
4331 ! real(RP), parameter :: ncoefim = 1.0E+7_RP ! N_{im0} of eq.(3.18) of Suzuki (2004)
4332 ! real(RP), parameter :: gamm = 3.3_RP ! gamma of eq.(3.18) of Suzuki (2004)
4333  integer :: ijk, indirect
4334  real(rp) :: qlmlt(ijkmax)
4335  real(rp) :: dcrg
4336 
4337  call prof_rapstart('_SBM_Freezing', 3)
4338 
4339  call atmos_hydrometeor_lhf_para( ijkmax, 1, ijkmax, temp(:), qlmlt(:) )
4340 
4341  xbound = log( rhow * 4.0_rp*pi/3.0_rp * rbound**3 )
4342  nbound = int( ( xbound-xbnd( 1 ) )/dxmic ) + 1
4343 
4344  !$omp parallel do private(ijk,tc,rate,dcrg,sumfrz,dmp,frz,tdel)
4345  do indirect = 1, num_cold
4346  ijk = index_cold(indirect)
4347 
4348 ! if ( temp(ijk) <= tthreth ) then !--- Bigg (1975)
4349 
4350  tc = temp(ijk)-tmlt
4351  rate = coefa*exp( -coefb*tc )
4352 
4353  dcrg = 0.0_rp
4354  sumfrz = 0.0_rp
4355  do n = 1, nbound-1
4356  dmp = rate*expxctr( n )
4357  frz = gc( n,il,ijk )*( 1.0_rp-exp( -dmp*dtime ) )
4358  frz = min( frz, gc( n,il,ijk ) )
4359 
4360  if ( flg_lt ) then
4361  dcrg = frz / ( gc(n,il,ijk)+eps ) * gcrg( n,il,ijk )
4362  gcrg( n,il,ijk ) = gcrg( n,il,ijk ) - dcrg
4363  gcrg( n,ip,ijk ) = gcrg( n,ip,ijk ) + dcrg
4364  end if
4365 
4366  gc( n,il,ijk ) = gc( n,il,ijk ) - frz
4367  gc( n,ip,ijk ) = gc( n,ip,ijk ) + frz
4368 
4369  sumfrz = sumfrz + frz
4370  enddo
4371  do n = nbound, nbin
4372  dmp = rate*expxctr( n )
4373  frz = gc( n,il,ijk )*( 1.0_rp-exp( -dmp*dtime ) )
4374  frz = min( frz, gc( n,il,ijk ) )
4375 
4376  if ( flg_lt ) then
4377  dcrg = frz / ( gc(n,il,ijk)+eps ) * gcrg( n,il,ijk )
4378  gcrg( n,il,ijk ) = gcrg( n,il,ijk ) - dcrg
4379  gcrg( n,ih,ijk ) = gcrg( n,ih,ijk ) + dcrg
4380  end if
4381 
4382  gc( n,il,ijk ) = gc( n,il,ijk ) - frz
4383  gc( n,ih,ijk ) = gc( n,ih,ijk ) + frz
4384 
4385  sumfrz = sumfrz + frz
4386  enddo
4387 
4388 ! elseif( temp(ijk) > tthreth ) then !--- Vali (1975)
4389 !
4390 ! tc = temp(ijk)-tmlt
4391 ! dmp = ncoefim * ( 0.1_RP * ( -tc )**gamm )
4392 ! sumfrz = 0.0_RP
4393 ! do n = 1, nbound-1
4394 ! frz = gc( n,il,ijk )*expxctr( n )/rhow*dmp
4395 ! frz = min( frz,gc( n,il,ijk ) )
4396 ! gc( n,il,ijk ) = gc( n,il,ijk ) - frz
4397 ! gc( n,ip,ijk ) = gc( n,ip,ijk ) + frz
4398 !
4399 ! if ( flg_lt ) then
4400 ! dcrg = frz / ( gc(n,il,ijk)+EPS ) * gcrg( n,il,ijk )
4401 ! gcrg( n,il,ijk ) = gcrg( n,il,ijk ) - dcrg
4402 ! gcrg( n,ip,ijk ) = gcrg( n,ip,ijk ) + dcrg
4403 ! endif
4404 !
4405 ! sumfrz = sumfrz + frz
4406 ! enddo
4407 ! do n = nbound, nbin
4408 ! frz = gc( n,il,ijk )*dmp*expxctr( n )/rhow
4409 ! frz = min( frz,gc( n,il,ijk ) )
4410 ! gc( n,il,ijk ) = gc( n,il,ijk ) - frz
4411 ! gc( n,ih,ijk ) = gc( n,ih,ijk ) + frz
4412 !
4413 ! if ( flg_lt ) then
4414 ! dcrg = frz / ( gc(n,il,ijk)+EPS ) * gcrg( n,il,ijk )
4415 ! gcrg( n,il,ijk ) = gcrg( n,il,ijk ) - dcrg
4416 ! gcrg( n,ih,ijk ) = gcrg( n,ih,ijk ) + dcrg
4417 ! endif
4418 !
4419 ! sumfrz = sumfrz + frz
4420 ! enddo
4421 ! endif
4422  sumfrz = sumfrz*dxmic
4423 
4424  tdel = sumfrz/dens(ijk)*qlmlt(ijk)/cv(ijk)
4425  temp(ijk) = temp(ijk) + tdel
4426  cp(ijk) = cp(ijk) + ( cp_ice - cp_water ) * sumfrz/dens(ijk)
4427  cv(ijk) = cv(ijk) + ( cv_ice - cv_water ) * sumfrz/dens(ijk)
4428  enddo
4429 
4430  call prof_rapend ('_SBM_Freezing', 3)
4431 
4432  return
4433  end subroutine freezing
4434 
4435  !-----------------------------------------------------------------------------
4436  subroutine melting( &
4437  ijkmax, &
4438  num_warm, &
4439  index_warm, &
4440  flg_lt, &
4441  dens, &
4442  temp, &
4443  gc, &
4444  gcrg, &
4445  cp, &
4446  cv, &
4447  dtime )
4448  use scale_atmos_hydrometeor, only: &
4449  atmos_hydrometeor_lhf_para, &
4450  cp_water, &
4451  cp_ice, &
4452  cv_water, &
4453  cv_ice
4454  implicit none
4455 
4456  integer, intent(in) :: ijkmax
4457  integer, intent(in) :: num_warm
4458  integer, intent(in) :: index_warm(ijkmax)
4459  logical, intent(in) :: flg_lt
4460  real(rp), intent(in) :: dens(ijkmax) ! Density [kg/m3]
4461  real(rp), intent(inout) :: temp(ijkmax) ! Temperature [K]
4462  real(rp), intent(inout) :: gc (nbin,nspc,ijkmax) ! Mass size distribution function of hydrometeor
4463  real(rp), intent(inout) :: gcrg(nbin,nspc,ijkmax)
4464  real(rp), intent(inout) :: cp(ijkmax) ! specific heat
4465  real(rp), intent(inout) :: cv(ijkmax) ! specific heat
4466  real(dp), intent(in) :: dtime ! Time step interval
4467 
4468  integer :: n, m
4469  real(rp) :: summlt, sumice, tdel
4470  integer :: ijk, indirect
4471  real(rp) :: qlmlt(ijkmax)
4472  real(rp) :: dcrg
4473 
4474  call prof_rapstart('_SBM_Melting', 3)
4475 
4476  call atmos_hydrometeor_lhf_para( ijkmax, 1, ijkmax, temp(:), qlmlt(:) )
4477 
4478  !$omp parallel do private(ijk,summlt,sumice,dcrg,tdel)
4479  do indirect = 1, num_warm
4480  ijk = index_warm(indirect)
4481 
4482  summlt = 0.0_rp
4483  do n = 1, nbin
4484  sumice = 0.0_rp
4485  do m = ic, ih
4486  sumice = sumice + gc( n,m,ijk )
4487  gc( n,m,ijk ) = 0.0_rp
4488  enddo
4489  gc( n,il,ijk ) = gc( n,il,ijk ) + sumice
4490  if ( flg_lt ) then
4491  dcrg = 0.0_rp
4492  do m = ic, ih
4493  dcrg = dcrg + gcrg( n,m,ijk )
4494  gcrg( n,m,ijk ) = 0.0_rp
4495  end do
4496  gcrg( n,il,ijk ) = gcrg( n,il,ijk ) + dcrg
4497  end if
4498  summlt = summlt + sumice !--- All freezed particle melt instantaneously
4499  enddo
4500  summlt = summlt*dxmic
4501 
4502  tdel = - summlt/dens(ijk)*qlmlt(ijk)/cv(ijk)
4503  temp(ijk) = temp(ijk) + tdel
4504  cp(ijk) = cp(ijk) + ( cp_water - cp_ice ) * summlt/dens(ijk)
4505  cv(ijk) = cv(ijk) + ( cv_water - cv_ice ) * summlt/dens(ijk)
4506  !
4507  enddo
4508 
4509  call prof_rapend ('_SBM_Melting', 3)
4510 
4511  return
4512  end subroutine melting
4513 
4514  !-----------------------------------------------------------------------------
4515  subroutine collmain( &
4516  KA, IA, JA , &
4517  ijkmax, &
4518  flg_lt, &
4519  d0_crg, &
4520  v0_crg, &
4521  dq, &
4522  beta, &
4523  temp, &
4524  ghyd, &
4525  gcrg, &
4526  crg_sep,&
4527  dt )
4528  implicit none
4529 
4530  integer, intent(in) :: ka
4531  integer, intent(in) :: ia
4532  integer, intent(in) :: ja
4533 
4534  integer, intent(in) :: ijkmax
4535  logical, intent(in) :: flg_lt
4536  real(rp), intent(in) :: v0_crg, d0_crg
4537  real(rp), intent(in) :: dq(ijkmax)
4538  real(rp), intent(in) :: beta(ijkmax)
4539  real(rp), intent(in) :: temp(ijkmax) ! Temperature [K]
4540  real(rp), intent(inout) :: ghyd(nbin,nspc,ijkmax) ! Mass size distribution function of hydrometeor
4541  real(rp), intent(inout) :: gcrg(nbin,nspc,ijkmax)
4542  real(rp), intent(out) :: crg_sep(nspc,ijkmax)
4543  real(dp), intent(in) :: dt ! Time step interval
4544  !---------------------------------------------------------------------------
4545 
4546  if ( flg_rndm ) then ! stochastic method
4547  call r_collcoag( ka, ia, ja, & ! [IN]
4548  ijkmax, & ! [IN]
4549  wgtbin, & ! [IN]
4550  temp(:), & ! [IN]
4551  ghyd(:,:,:), & ! [INOUT]
4552  dt ) ! [IN]
4553  else ! default
4554  call collcoag( ijkmax, & ! [IN]
4555  flg_lt, & ! [IN]
4556  d0_crg, & ! [IN]
4557  v0_crg, & ! [IN]
4558  dq(:), & ! [IN]
4559  beta(:), & ! [IN]
4560  temp(:), & ! [IN]
4561  ghyd(:,:,:), & ! [INOUT]
4562  gcrg(:,:,:), & ! [INOUT]
4563  crg_sep(:,:),& ! [OUT]
4564  dt ) ! [IN]
4565  endif
4566 
4567  return
4568  end subroutine collmain
4569 
4570  !-----------------------------------------------------------------------------
4571  subroutine collmainf( &
4572  KA, IA, JA, &
4573  ijkmax, &
4574  flg_lt, &
4575  d0_crg, &
4576  v0_crg, &
4577  dq, &
4578  beta, &
4579  temp, &
4580  ghyd, &
4581  gcrg, &
4582  crg_sep,&
4583  dt )
4584  implicit none
4585 
4586  integer, intent(in) :: ka
4587  integer, intent(in) :: ia
4588  integer, intent(in) :: ja
4589 
4590  integer, intent(in) :: ijkmax
4591  logical, intent(in) :: flg_lt
4592  real(rp), intent(in) :: v0_crg, d0_crg
4593  real(rp), intent(in) :: dq(ijkmax)
4594  real(rp), intent(in) :: beta(ijkmax)
4595  real(rp), intent(in) :: temp(ijkmax) ! Temperature [K]
4596  real(rp), intent(inout) :: ghyd(nbin,nspc,ijkmax) ! Mass size distribution function of hydrometeor
4597  real(rp), intent(inout) :: gcrg(nbin,nspc,ijkmax)
4598  real(rp), intent(out) :: crg_sep(nspc,ijkmax)
4599  real(dp), intent(in) :: dt ! Time step interval
4600  !---------------------------------------------------------------------------
4601 
4602  if ( flg_rndm ) then ! stochastic method
4603  call r_collcoag( ka, ia, ja, & ! [IN]
4604  ijkmax, & ! [IN]
4605  wgtbin, & ! [IN]
4606  temp(:), & ! [IN]
4607  ghyd(:,:,:), & ! [INOUT]
4608  dt ) ! [IN]
4609  else ! default
4610  call collcoag( ijkmax, & ! [IN]
4611  flg_lt, & ! [IN]
4612  d0_crg, & ! [IN]
4613  v0_crg, & ! [IN]
4614  dq(:), & ! [IN]
4615  beta(:), & ! [IN]
4616  temp(:), & ! [IN]
4617  ghyd(:,:,:), & ! [INOUT]
4618  gcrg(:,:,:), & ! [INOUT]
4619  crg_sep(:,:),& ! [OUT]
4620  dt ) ! [IN]
4621  endif
4622 
4623  return
4624  end subroutine collmainf
4625 
4626  !-----------------------------------------------------------------------------
4627  !--- reference paper
4628  ! Bott et al. (1998) J. Atmos. Sci. vol.55, pp. 2284-
4629  subroutine collcoag( &
4630  ijkmax, &
4631  flg_lt, &
4632  d0_crg, &
4633  v0_crg, &
4634  dq, &
4635  beta, &
4636  temp, &
4637  gc, &
4638  gcrg, &
4639  crg_sep,&
4640  dtime )
4641  implicit none
4642 
4643  integer, intent(in) :: ijkmax
4644  logical, intent(in) :: flg_lt
4645  real(rp), intent(in) :: d0_crg, v0_crg
4646  real(rp), intent(in) :: dq(ijkmax)
4647  real(rp), intent(in) :: beta(ijkmax)
4648  real(rp), intent(in) :: temp(ijkmax) ! Temperature [K]
4649  real(rp), intent(inout) :: gc (nbin,nspc,ijkmax) ! Mass size distribution function of hydrometeor
4650  real(rp), intent(inout) :: gcrg(nbin,nspc,ijkmax)
4651  real(rp), intent(out) :: crg_sep(nspc,ijkmax)
4652  real(dp), intent(in) :: dtime ! Time step interval
4653 
4654  integer :: i, j, k
4655  real(rp) :: xi, xj, xnew, dmpi, dmpj, frci, frcj
4656  real(rp) :: gprime, gprimk, wgt, crn, sum, flux
4657  real(rp) :: acoef0, acoef1, acoef2
4658  real(rp), parameter :: dmpmin = 1.e-01_rp
4659  real(rp) :: suri, surj
4660 
4661  integer :: myu, n, irsl, ilrg, isml
4662  integer :: ibnd( ijkmax )
4663  integer :: iflg( nspc,ijkmax )
4664  integer :: iexst( nbin,nspc,ijkmax )
4665  real(rp) :: csum
4666  integer :: ijk, nn, mm, pp, qq
4667 
4668  !--- for non-inducting charge separation
4669  real(rp) :: drhoi, drhoj, drhok, alpha, dgenei, dgenej
4670  integer :: ispc
4671  !---------------------------------------------------------------------------
4672 
4673  call prof_rapstart('_SBM_CollCoag', 3)
4674 
4675  !$omp workshare
4676  iflg( :,: ) = 0
4677  iexst( :,:,: ) = 0
4678  crg_sep( :,: ) = 0.0_rp
4679  !$omp end workshare
4680  !$omp parallel do private(csum)
4681  do ijk = 1, ijkmax
4682  !--- judgement of particle existence
4683  do myu = 1, nspc
4684  csum = 0.0_rp
4685  do n = 1, nbin
4686  csum = csum + gc( n,myu,ijk )*dxmic
4687  enddo
4688  if ( csum > cldmin ) iflg( myu,ijk ) = 1
4689  enddo
4690 
4691  if ( temp(ijk) < tcrit ) then
4692  ibnd(ijk) = 1
4693  else
4694  ibnd(ijk) = 2
4695  endif
4696 
4697  do myu = 1, nspc
4698  do n = 1, nbin
4699  if ( gc( n,myu,ijk ) > cldmin ) then
4700  iexst( n,myu,ijk ) = 1
4701  endif
4702  enddo
4703  enddo
4704 
4705  enddo
4706 
4707 !OCL NORECURRENCE
4708  !$omp parallel do &
4709  !$omp private(irsl,ispc,k,xi,xj,xnew,dmpi,dmpj,frci,frcj, &
4710  !$omp gprime,suri,surj,gprimk,wgt,crn,acoef0,acoef1,acoef2,sum,flux,alpha, &
4711  !$omp drhoi,drhoj,drhok,dgenei,dgenej)
4712  do ijk = 1, ijkmax
4713  do isml = 1, nspc
4714  do nn = 1, iflg( isml,ijk )
4715 
4716  do ilrg = 1, nspc
4717  do mm = 1, iflg( ilrg,ijk )
4718  !--- rule of interaction
4719  irsl = ifrsl( ibnd(ijk),isml,ilrg )
4720 
4721  if ( isml /= ilrg ) then
4722  ispc = 0
4723  else if ( isml == ilrg ) then
4724  ispc = 1
4725  endif
4726 
4727  do i = 1, nbin-1 ! small
4728  do pp = 1, iexst( i,isml,ijk )
4729 
4730 !OCL SIMD
4731  do j = i+1, nbin ! large
4732 ! do qq = 1, iexst( j,ilrg,ijk )
4733  if ( iexst(j,ilrg,ijk) == 1 ) then
4734 
4735  k = kindx( i,j )
4736  xi = expxctr( i )
4737  xj = expxctr( j )
4738  xnew = log( xi+xj )
4739 
4740  dmpi = ck( isml,ilrg,i,j )*gc( j,ilrg,ijk )/xj*dxmic*dtime ! dg_{i}/dt*dt in page 119 of Suzuki (2004)
4741  dmpj = ck( ilrg,isml,i,j )*gc( i,isml,ijk )/xi*dxmic*dtime ! dg_{j}/dt*dt in page 119 of Suzuki (2004)
4742 
4743  if ( dmpi <= dmpmin ) then
4744  frci = gc( i,isml,ijk )*dmpi ! Dg_{i} in page 119 of Suzuki (2004)
4745  else
4746  frci = gc( i,isml,ijk )*( 1.0_rp-exp( -dmpi ) ) ! Dg_{i} in page 119 of Suzuki (2004)
4747  endif
4748 
4749  if ( dmpj <= dmpmin ) then
4750  frcj = gc( j,ilrg,ijk )*dmpj ! Dg_{j} in page 119 of Suzuki (2004)
4751  else
4752  frcj = gc( j,ilrg,ijk )*( 1.0_rp-exp( -dmpj ) ) ! Dg_{j} in page 119 of Suzuki (2004)
4753  endif
4754 
4755  gprime = frci+frcj
4756 ! if ( gprime <= 0.0_RP ) cycle large
4757  if ( gprime > 0.0_rp .AND. k < nbin ) then
4758 
4759  suri = gc( i,isml,ijk )
4760  surj = gc( j,ilrg,ijk )
4761  gc( i,isml,ijk ) = gc( i,isml,ijk )-frci
4762  gc( j,ilrg,ijk ) = gc( j,ilrg,ijk )-frcj
4763  gc( i,isml,ijk ) = max( gc( i,isml,ijk )-frci, 0.0_rp )
4764  gc( j,ilrg,ijk ) = max( gc( j,ilrg,ijk )-frcj, 0.0_rp )
4765  frci = suri - gc( i,isml,ijk )
4766  frcj = surj - gc( j,ilrg,ijk )
4767  gprime = frci+frcj ! g' in page 119 of Suzuki (2004)
4768 
4769  if ( flg_lt ) then
4770  !--- charge transfer from two particles to generated large particle
4771  drhoi = frci/suri*gcrg( i,isml,ijk )
4772  drhoj = frcj/surj*gcrg( j,ilrg,ijk )
4773  drhok = drhoi + drhoj
4774 
4775  !--- charge density generated by non-inductive charging
4776  alpha = 5.0_rp * ( 2.0_rp*radc( i )/d0_crg )**2 &
4777  * abs( vt(ilrg,j)-vt(isml,i) )/v0_crg ! alpha in eq. (12) of Mansell et al. (2005)
4778  alpha = min( 10.0_rp, alpha )
4779 
4780  dgenei = frci / xi * rcoll( isml,ilrg,i,j ) * ( -dq( ijk ) ) &
4781  * alpha * beta( ijk ) * flg_noninduct( isml,ilrg ) ! eq. (8) of Mansell et al. (2005)
4782  dgenej = - dgenei
4783  gcrg( i,isml,ijk ) = gcrg( i,isml,ijk ) + ( dgenei-drhoi )
4784  gcrg( j,ilrg,ijk ) = gcrg( j,ilrg,ijk ) + ( dgenej-drhoj )
4785  crg_sep( isml,ijk ) = crg_sep( isml,ijk ) + dgenei
4786  crg_sep( ilrg,ijk ) = crg_sep( ilrg,ijk ) + dgenej
4787  end if
4788 
4789  gprimk = gc( k,irsl,ijk ) + gprime ! g'_{k} in page 119 of Suzuki (2004)
4790  wgt = gprime / gprimk ! w in page 119 of Suzuki (2004)
4791  crn = ( xnew-xctr( k ) )/( xctr( k+1 )-xctr( k ) ) ! c_{k} in page 119 of Suzuki (2004)
4792 
4793  acoef0 = -( gc( k+1,irsl,ijk )-26.0_rp*gprimk+gc( k-1,irsl,ijk ) )/24.0_rp ! a_{k,0} in page 119 of Suzuki (2004)
4794  acoef1 = ( gc( k+1,irsl,ijk )-gc( k-1,irsl,ijk ) ) *0.50_rp ! a_{k,1} in page 119 of Suzuki (2004)
4795  acoef2 = ( gc( k+1,irsl,ijk )-2.0_rp*gprimk+gc( k-1,irsl,ijk ) ) *0.50_rp ! a_{k,2} in page 119 of Suzuki (2004)
4796 
4797 !!$ sum = 0.0_RP
4798 !!$ do l = 0, ldeg
4799 !!$ sum = sum + acoef( l )/( l+1 )/2.0_RP**( l+1 ) &
4800 !!$ *( 1.0_RP-( 1.0_RP-2.0_RP*crn )**( l+1 ) )
4801 !!$ enddo
4802 
4803  sum = acoef0 / 2.0_rp * ( 1.0_rp - ( 1.0_rp-2.0_rp*crn ) ) &
4804  + acoef1 / 8.0_rp * ( 1.0_rp - ( 1.0_rp-2.0_rp*crn )**2 ) &
4805  + acoef2 / 24.0_rp * ( 1.0_rp - ( 1.0_rp-2.0_rp*crn )**3 )
4806 
4807  flux = wgt*sum ! f_{k+1/2} in page 119 of Suzuki (2004)
4808  flux = min( max( flux,0.0_rp ),gprime )
4809 
4810  gc( k,irsl,ijk ) = gprimk - flux ! tilda{g_{k}} in page 119 of Suzuki (2004)
4811  gc( k+1,irsl,ijk ) = gc( k+1,irsl,ijk ) + flux ! tilda{g_{k+1}} in page 119 of Suzuki (2004)
4812 
4813  if ( flg_lt ) then
4814  !--- charge transfer from two particles to generated large particle
4815  if( gprime /= 0.0_rp ) then
4816  gcrg( k,irsl,ijk ) = gcrg( k,irsl,ijk ) + drhok * ( gprime-flux )/gprime
4817  gcrg( k+1,irsl,ijk ) = gcrg( k+1,irsl,ijk ) + drhok * flux/gprime
4818  endif
4819  end if
4820 
4821  endif
4822 
4823  endif
4824 ! enddo
4825  enddo !large
4826  enddo
4827  enddo !small
4828  !
4829  enddo
4830  enddo
4831 
4832  enddo
4833  enddo
4834 
4835  enddo
4836 
4837  call prof_rapend ('_SBM_CollCoag', 3)
4838 
4839  return
4840  end subroutine collcoag
4841 
4842  !-------------------------------------------------------------------
4843  subroutine getrule( ifrsl, indx ) !--- out
4844  ! subroutine for creating lookup table (Table A.2 of Suzuki 2004)
4845  integer, intent(out) :: indx( nbin,nbin )
4846  integer, intent(out) :: ifrsl( 2,nspc_mk,nspc_mk )
4847  integer :: i, j, k
4848  real(rp) :: xnew
4849  !
4850  !$omp parallel do private(xnew,k)
4851  do i = 1, nbin
4852  do j = 1, nbin
4853  xnew = log( expxctr( i )+expxctr( j ) )
4854  k = int( ( xnew-xctr( 1 ) )/dxmic ) + 1
4855  k = max( max( k,j ),i )
4856  indx( i,j ) = k
4857  enddo
4858  enddo
4859  !
4860  !
4861  !--- liquid + liquid -> liquid
4862  ifrsl( 1:2,il,il ) = il
4863  !
4864  !--- liquid + column -> ( graupel, hail ) + column
4865  ifrsl( 1:2,il,ic ) = ic
4866  ifrsl( 1,ic,il ) = ig
4867  ifrsl( 2,ic,il ) = ih
4868  !
4869  !--- liquid + plate -> ( graupel, hail ) + plate
4870  ifrsl( 1:2,il,ip ) = ip
4871  ifrsl( 1,ip,il ) = ig
4872  ifrsl( 2,ip,il ) = ih
4873  !
4874  !--- liquid + dendrite -> ( graupel, hail ) + dendrite
4875  ifrsl( 1:2,il,id ) = id
4876  ifrsl( 1,id,il ) = ig
4877  ifrsl( 2,id,il ) = ih
4878  !
4879  !--- liquid + snowflake -> ( graupel, hail ) + snowflake
4880  ifrsl( 1:2,il,iss ) = iss
4881  ifrsl( 1,iss,il ) = ig
4882  ifrsl( 2,iss,il ) = ih
4883  !
4884  !--- liquid + graupel -> ( graupel, hail )
4885  ifrsl( 1:2,il,ig ) = ig
4886  ifrsl( 1,ig,il ) = ig
4887  ifrsl( 2,ig,il ) = ih
4888  !
4889  !--- liquid + hail -> ( graupel, hail )
4890  ifrsl( 1:2,il,ih ) = ih
4891  ifrsl( 1,ih,il ) = ig
4892  ifrsl( 2,ih,il ) = ih
4893  !
4894  !
4895  !--- column + column -> snowflake
4896  ifrsl( 1:2,ic,ic ) = iss
4897  !
4898  !--- column + plate -> snowflake
4899  ifrsl( 1:2,ic,ip ) = iss
4900  ifrsl( 1:2,ip,ic ) = iss
4901  !
4902  !--- column + dendrite -> snowflake
4903  ifrsl( 1:2,ic,id ) = iss
4904  ifrsl( 1:2,id,ic ) = iss
4905  !
4906  !--- column + snowflake -> snowflake
4907  ifrsl( 1:2,ic,iss ) = iss
4908  ifrsl( 1:2,iss,ic ) = iss
4909  !
4910  !--- column + graupel -> column + graupel
4911  ifrsl( 1:2,ic,ig ) = ig
4912  ifrsl( 1:2,ig,ic ) = ic
4913  !
4914  !--- column + hail -> column + ( graupel, hail )
4915  ifrsl( 1:2,ih,ic ) = ic
4916  ifrsl( 1,ic,ih ) = ig
4917  ifrsl( 2,ic,ih ) = ih
4918  !
4919  !
4920  !--- plate + plate -> snowflake
4921  ifrsl( 1:2,ip,ip ) = iss
4922  !
4923  !--- plate + dendrite -> snowflake
4924  ifrsl( 1:2,ip,id ) = iss
4925  ifrsl( 1:2,id,ip ) = iss
4926  !
4927  !--- plate + snowflake -> snowflake
4928  ifrsl( 1:2,ip,iss ) = iss
4929  ifrsl( 1:2,iss,ip ) = iss
4930  !
4931  !--- plate + graupel -> plate + graupel
4932  ifrsl( 1:2,ip,ig ) = ig
4933  ifrsl( 1:2,ig,ip ) = ip
4934  !
4935  !--- plate + hail -> plate + ( graupel, hail )
4936  ifrsl( 1:2,ih,ip ) = ip
4937  ifrsl( 1,ip,ih ) = ig
4938  ifrsl( 2,ip,ih ) = ih
4939  !
4940  !
4941  !--- dendrite + dendrite -> snowflake
4942  ifrsl( 1:2,id,id ) = iss
4943  !
4944  !--- dendrite + snowflake -> snowflake
4945  ifrsl( 1:2,id,iss ) = iss
4946  ifrsl( 1:2,iss,id ) = iss
4947  !
4948  !--- dendrite + graupel -> dendrite + graupel
4949  ifrsl( 1:2,id,ig ) = ig
4950  ifrsl( 1:2,ig,id ) = id
4951  !
4952  !--- dendrite + hail -> dendrite + ( graupel, hail )
4953  ifrsl( 1:2,ih,id ) = id
4954  ifrsl( 1,id,ih ) = ig
4955  ifrsl( 2,id,ih ) = ih
4956  !
4957  !
4958  !--- snowflake + snowflake -> snowflake
4959  ifrsl( 1:2,iss,iss ) = iss
4960  !
4961  !--- snowflake + graupel -> snowflake + graupel
4962  ifrsl( 1:2,iss,ig ) = ig
4963  ifrsl( 1:2,ig,iss ) = iss
4964  !
4965  !--- snowflake + hail -> snowflake + ( graupel, hail )
4966  ifrsl( 1:2,ih,iss ) = iss
4967  ifrsl( 1,iss,ih ) = ig
4968  ifrsl( 2,iss,ih ) = ih
4969  !
4970  !
4971  !--- graupel + graupel -> graupel
4972  ifrsl( 1:2,ig,ig ) = ig
4973  !
4974  !--- graupel + hail -> ( graupel, hail )
4975  ifrsl( 1,ig,ih ) = ig
4976  ifrsl( 1,ih,ig ) = ig
4977  ifrsl( 2,ig,ih ) = ih
4978  ifrsl( 2,ih,ig ) = ih
4979  !
4980  !--- hail + hail -> hail
4981  ifrsl( 1:2,ih,ih ) = ih
4982 
4983  return
4984  end subroutine getrule
4985 
4986  !-----------------------------------------------------------------------------
4987  subroutine faero( &
4988  ijkmax, &
4989  f0, &
4990  ga )
4991  implicit none
4992 
4993  integer , intent(in) :: ijkmax
4994  real(rp), intent(in) :: f0(ijkmax)
4995  real(rp), intent(inout) :: ga(nccn,ijkmax)
4996 
4997 ! real(RP), parameter :: alpha = 3.0_RP
4998 
4999  integer :: ijk, n
5000  !---------------------------------------------------------------------------
5001 
5002  call prof_rapstart('_SBM_FAero', 3)
5003 
5004  !$omp parallel do
5005  do ijk = 1, ijkmax
5006  do n = 1, nccn
5007  ga(n,ijk) = ga(n,ijk) + f0(ijk) * marate(n) * expxactr(n) / dxaer ! [note] marate is never set.
5008  enddo
5009  enddo
5010 
5011  call prof_rapend ('_SBM_FAero', 3)
5012 
5013  return
5014  end subroutine faero
5015 
5016  !-----------------------------------------------------------------------------
5017  ! + Y. Sato added for stochastic method
5018  ! + Reference Sato et al. (2009) JGR, doi:10.1029/2008JD011247
5019  subroutine random_setup( mset ) !--- in
5020  use scale_random, only: &
5021  random_uniform
5022  use scale_prc, only: &
5023  prc_abort
5024 
5025  integer, intent(in) :: mset
5026 
5027  !--- local ----
5028  integer :: n
5029  real(rp) :: nbinr, tmp1
5030  real(rp) :: rans( mbin ), ranl( mbin )
5031  integer :: pq
5032  real(rp), allocatable :: ranstmp( : )
5033  real(rp), allocatable :: ranltmp( : )
5034  integer :: p, q
5035  integer :: k, temp
5036  integer, allocatable :: orderb( : )
5037  real(rp) :: abq1
5038  real(rp) :: a
5039  real(rp), allocatable :: randnum(:,:,:)
5040  !-------------------------------------------------------
5041  pq = nbin*(nbin-1)/2
5042  allocate( blrg( mset, mbin ) )
5043  allocate( bsml( mset, mbin ) )
5044  allocate( ranstmp( pq ) )
5045  allocate( ranltmp( pq ) )
5046  allocate( orderb( pq ) )
5047  allocate( randnum(1,1,pq) )
5048 
5049  a = real( nbin )*real( nbin-1 )*0.50_rp
5050  if ( a < mbin ) then
5051  log_error("ATMOS_PHY_MP_SUZUKI10_random_setup",*) "mbin should be smaller than {nbin}_C_{2}"
5052  call prc_abort
5053  endif
5054 
5055  wgtbin = a/real( mbin )
5056  nbinr = real( nbin )
5057 
5058  do p = 1, pq
5059  orderb( p ) = p
5060  enddo
5061 
5062  do p = 1, nbin-1
5063  ranstmp( (p-1)*nbin-(p*(p-1))/2+1 : p*nbin-(p*(p+1))/2 ) = p
5064  do q = 1, nbin-p
5065  ranltmp( (p-1)*nbin-(p*(p-1))/2+q ) = p+q
5066  enddo
5067  enddo
5068 
5069  do n = 1, mset
5070  call random_uniform( randnum )
5071  do p = 1, pq
5072  abq1 = randnum( 1,1,p )
5073  k = int( abq1*( pq-p-1 ) ) + p
5074  temp = orderb( p )
5075  orderb( p ) = orderb( k )
5076  orderb( k ) = temp
5077  enddo
5078 
5079  do p = 1, mbin
5080  if ( p <= pq ) then
5081  rans( p ) = ranstmp( orderb( p ) )
5082  ranl( p ) = ranltmp( orderb( p ) )
5083  else
5084  rans( p ) = ranstmp( orderb( p-pq ) )
5085  ranl( p ) = ranltmp( orderb( p-pq ) )
5086  endif
5087  if ( rans( p ) >= ranl( p ) ) then
5088  tmp1 = rans( p )
5089  rans( p ) = ranl( p )
5090  ranl( p ) = tmp1
5091  endif
5092  enddo
5093  blrg( n,1:mbin ) = int( ranl( 1:mbin ) )
5094  bsml( n,1:mbin ) = int( rans( 1:mbin ) )
5095  enddo
5096 
5097  deallocate( ranstmp )
5098  deallocate( ranltmp )
5099  deallocate( orderb )
5100  deallocate( randnum )
5101 
5102  end subroutine random_setup
5103 
5104  !-----------------------------------------------------------------------------
5105  !--- reference paper
5106  ! Bott et al. (1998) J. Atmos. Sci. vol.55, pp. 2284-
5107  ! Bott et al. (2000) J. Atmos. Sci. Vol.57, pp. 284-
5108  subroutine r_collcoag( &
5109  KA, IA, JA, &
5110  ijkmax, &
5111  swgt, &
5112  temp, &
5113  gc, &
5114  dtime )
5115  use scale_random, only: &
5116  random_uniform
5117  implicit none
5118 
5119  integer, intent(in) :: ka
5120  integer, intent(in) :: ia
5121  integer, intent(in) :: ja
5122 
5123  integer, intent(in) :: ijkmax
5124  real(rp), intent(in) :: swgt
5125  real(rp), intent(in) :: temp(ijkmax) ! Temperature [K]
5126  real(rp), intent(inout) :: gc (nbin,nspc,ijkmax) ! Mass size distribution function of hydrometeor
5127  real(dp), intent(in) :: dtime ! Time step interval
5128 
5129  integer :: i, j, k, l
5130  real(rp) :: xi, xj, xnew, dmpi, dmpj, frci, frcj
5131  real(rp) :: gprime, gprimk, wgt, crn, sum, flux
5132  integer, parameter :: ldeg = 2
5133  real(rp), parameter :: dmpmin = 1.e-01_rp, cmin = 1.e-10_rp
5134  real(rp) :: acoef( 0:ldeg )
5135  !
5136  !--- Y.sato added to use code6
5137  integer :: nums( mbin ), numl( mbin )
5138  real(rp), parameter :: gt = 1.0_rp
5139  integer :: s, det
5140  real(rp) :: nbinr, mbinr ! use to weight
5141 ! real(RP) :: beta
5142  real(rp) :: tmpi, tmpj
5143 
5144  real(rp) :: rndm(nspc,nspc)
5145 
5146  integer :: ibnd( ijkmax )
5147  integer :: iflg( nspc,ijkmax )
5148  integer :: iexst( nbin,nspc,ijkmax )
5149  real(rp) :: csum
5150  integer :: ijk, nn, mm, pp, qq, myu, n, isml, ilrg, irsl
5151  !---------------------------------------------------------------------------
5152 
5153  call prof_rapstart('_SBM_CollCoagR', 3)
5154 
5155  !$omp workshare
5156  iflg( :,: ) = 0
5157  iexst( :,:,: ) = 0
5158  !$omp end workshare
5159  !$omp parallel do private(csum)
5160  do ijk = 1, ijkmax
5161  !--- judgement of particle existence
5162  do isml = 1, nspc
5163  csum = 0.0_rp
5164  do n = 1, nbin
5165  csum = csum + gc( n,isml,ijk )*dxmic
5166  enddo
5167  if ( csum > cldmin ) iflg( isml,ijk ) = 1
5168  enddo
5169 
5170  if ( temp(ijk) < tcrit ) then
5171  ibnd(ijk) = 1
5172  else
5173  ibnd(ijk) = 2
5174  endif
5175 
5176  do myu = 1, nspc
5177  do n = 1, nbin
5178  if ( gc( n,myu,ijk ) > cldmin ) then
5179  iexst( n,myu,ijk ) = 1
5180  endif
5181  enddo
5182  enddo
5183 
5184  enddo
5185 
5186  !$omp parallel do private(rndm,irsl,det,nbinr,mbinr,nums,numl,i,j,k,xi,xj,xnew,dmpi,dmpj,frci,frcj,tmpi,tmpj,gprime,gprimk,wgt,crn,acoef,sum,flux)
5187  do ijk = 1, ijkmax
5188 
5189  call random_uniform( rndm(:,:) )
5190 
5191  do isml = 1, nspc
5192  do nn = 1, iflg( isml,ijk )
5193 
5194  do ilrg = 1, nspc
5195  do mm = 1, iflg( ilrg,ijk )
5196  !--- rule of interaction
5197  irsl = ifrsl( ibnd(ijk),isml,ilrg )
5198 
5199  det = int( rndm(ilrg,isml)*ia*ja*ka + 1 )
5200  nbinr = real( nbin )
5201  mbinr = real( mbin )
5202  nums( 1:mbin ) = bsml( det,1:mbin )
5203  numl( 1:mbin ) = blrg( det,1:mbin )
5204 
5205  do s = 1, mbin
5206  i = nums( s )
5207  j = numl( s )
5208 
5209  do pp = 1, iexst( i,isml,ijk )
5210  do qq = 1, iexst( j,ilrg,ijk )
5211 
5212  k = kindx( i,j )
5213  xi = expxctr( i )
5214  xj = expxctr( j )
5215  xnew = log( xi+xj )
5216 
5217  dmpi = ck( isml,ilrg,i,j )*gc( j,ilrg,ijk )/xj*dxmic*dtime
5218  dmpj = ck( ilrg,isml,i,j )*gc( i,isml,ijk )/xi*dxmic*dtime
5219 
5220  if ( dmpi <= dmpmin ) then
5221  frci = gc( i,isml,ijk )*dmpi
5222  else
5223  frci = gc( i,isml,ijk )*( 1.0_rp-exp( -dmpi ) )
5224  endif
5225 
5226  if ( dmpj <= dmpmin ) then
5227  frcj = gc( j,ilrg,ijk )*dmpj
5228  else
5229  frcj = gc( j,ilrg,ijk )*( 1.0_rp-exp( -dmpj ) )
5230  endif
5231  tmpi = gc( i,isml,ijk )
5232  tmpj = gc( j,ilrg,ijk )
5233 
5234  gc( i,isml,ijk ) = gc( i,isml,ijk )-frci*swgt
5235  gc( j,ilrg,ijk ) = gc( j,ilrg,ijk )-frcj*swgt
5236 
5237  if ( j /= k ) then
5238  gc( j,ilrg,ijk ) = max( gc( j,ilrg,ijk ), 0.0_rp )
5239  endif
5240  gc( i,isml,ijk ) = max( gc( i,isml,ijk ), 0.0_rp )
5241 
5242  frci = tmpi - gc( i,isml,ijk )
5243  frcj = tmpj - gc( j,ilrg,ijk )
5244 
5245  gprime = frci+frcj
5246 
5247  !-----------------------------------------------
5248  !--- Exponential Flux Method (Bott, 2000, JAS)
5249  !-----------------------------------------------
5250  ! if ( gprime <= 0.0_RP ) cycle !large
5251  ! if ( gprime > 0.0_RP .AND. k < nbin ) then
5252  ! gprimk = gc( (irsl-1)*nbin+k ) + gprime
5253  !
5254  ! beta = log( gc( (irsl-1)*nbin+k+1 )/gprimk+1.E-60_RP )
5255  ! crn = ( xnew-xctr( k ) )/( xctr( k+1 )-xctr( k ) )
5256  !
5257  ! flux = ( gprime/beta )*( exp( beta*0.50_RP ) -exp( beta*( 0.50_RP-crn ) ) )
5258  ! flux = min( gprimk ,gprime )
5259  !
5260  ! gc( (irsl-1)*nbin+k ) = gprimk - flux
5261  ! gc( (irsl-1)*nbin+k+1 ) = gc( (irsl-1)*nbin+k+1 ) + flux
5262  ! endif
5263 
5264  !-----------------------------------------------
5265  !--- Flux Method (Bott, 1998, JAS)
5266  !-----------------------------------------------
5267 ! if ( gprime <= 0.0_RP ) cycle !large
5268  if ( gprime > 0.0_rp .AND. k < nbin ) then
5269 
5270  gprimk = gc( k,irsl,ijk ) + gprime
5271  wgt = gprime / gprimk
5272  crn = ( xnew-xctr( k ) )/( xctr( k+1 )-xctr( k ) )
5273 
5274  acoef( 0 ) = -( gc( k+1,irsl,ijk )-26.0_rp*gprimk+gc( k-1,irsl,ijk ) )/24.0_rp
5275  acoef( 1 ) = ( gc( k+1,irsl,ijk )-gc( k-1,irsl,ijk ) ) *0.5_rp
5276  acoef( 2 ) = ( gc( k+1,irsl,ijk )-2.0_rp*gprimk+gc( k-1,irsl,ijk ) ) *0.50_rp
5277 
5278  sum = 0.0_rp
5279  do l = 0, ldeg
5280  sum = sum + acoef( l )/( l+1 )/2.0_rp**( l+1 ) &
5281  *( 1.0_rp-( 1.0_rp-2.0_rp*crn )**( l+1 ) )
5282  enddo
5283 
5284  flux = wgt*sum
5285  flux = min( max( flux,0.0_rp ),gprime )
5286 
5287  gc( k,irsl,ijk ) = gprimk - flux
5288  gc( k+1,irsl,ijk ) = gc( k+1,irsl,ijk ) + flux
5289  endif
5290 
5291  enddo
5292  enddo
5293 
5294  enddo ! bin
5295  !
5296  enddo
5297  enddo
5298 
5299  enddo
5300  enddo
5301 
5302  enddo
5303 
5304  call prof_rapend ('_SBM_CollCoagR', 3)
5305 
5306  return
5307  end subroutine r_collcoag
5308 
5309  !-----------------------------------------------------------------------------
5310  !----- mkpara is module to create micpara.dat, which is parameter file of
5311  !----- micrphyisical proprties of hydrometeors (collision kernel, radius...).
5312  !----- Imported from preprocess/mk_para2 at 2013/12/26 (Y.Sato)
5313  subroutine mkpara
5314 
5315  implicit none
5316 
5317  integer :: i, j
5318  !-----------------------------------------------------------------------------
5319 
5320  allocate( radc_mk( nbin ) )
5321  allocate( xctr_mk( nbin ) )
5322  allocate( xbnd_mk( nbin+1 ) )
5323  allocate( cctr_mk( nspc_mk,nbin ) )
5324  allocate( cbnd_mk( nspc_mk,nbin+1 ) )
5325  allocate( ck_mk( nspc_mk,nspc_mk,nbin,nbin ) )
5326  allocate( vt_mk( nspc_mk,nbin ) )
5327  allocate( br_mk( nspc_mk,nbin ) )
5328 
5329  !--- file reading
5330  call rdkdat
5331 
5332  !--- grid setting
5333  call sdfgrid
5334 
5335  !--- capacity
5336  call getcp
5337 
5338  !--- collection kernel
5339  call getck
5340 
5341  !--- terminal velocity
5342  call getvt
5343 
5344  !--- bulk radius
5345  call getbr
5346 
5347  !--- output
5348  call paraout
5349 
5350  deallocate( radc_mk )
5351  deallocate( xctr_mk )
5352  deallocate( xbnd_mk )
5353  deallocate( cctr_mk )
5354  deallocate( cbnd_mk )
5355  deallocate( ck_mk )
5356  deallocate( vt_mk )
5357  deallocate( br_mk )
5358 
5359  end subroutine mkpara
5360 
5361  !---------------------------------------------------------------------------------------
5362  subroutine rdkdat
5363 
5364  implicit none
5365  integer, parameter :: il = 1, ic = 2, ip = 3, id = 4
5366  integer, parameter :: is = 5, ig = 6, ih = 7
5367  integer, parameter :: icemax = 3
5368 
5369  integer :: k, kk, i, j
5370  real(dp) :: xl( ndat ), rlec( ndat ), vrl( ndat )
5371  real(dp) :: blkradl( ndat ), blkdnsl( ndat )
5372 
5373  real(dp) :: xi( ndat,icemax ), riec( ndat,icemax ), vri( ndat,icemax )
5374  real(dp) :: blkradi( ndat,icemax ), blkdnsi( ndat,icemax )
5375 
5376  real(dp) :: xs( ndat ), rsec( ndat ), vrs( ndat )
5377  real(dp) :: blkrads( ndat ), blkdnss( ndat )
5378 
5379  real(dp) :: xg( ndat ), rgec( ndat ), vrg( ndat )
5380  real(dp) :: blkradg( ndat ), blkdnsg( ndat )
5381 
5382  real(dp) :: xh( ndat ), rhec( ndat ), vrh( ndat )
5383  real(dp) :: blkradh( ndat ), blkdnsh( ndat )
5384 
5385  data xl(1:ndat) / &
5386  0.33510e-10_dp,0.67021e-10_dp,0.13404e-09_dp,0.26808e-09_dp,0.53617e-09_dp,0.10723e-08_dp, &
5387  0.21447e-08_dp,0.42893e-08_dp,0.85786e-08_dp,0.17157e-07_dp,0.34315e-07_dp,0.68629e-07_dp, &
5388  0.13726e-06_dp,0.27452e-06_dp,0.54903e-06_dp,0.10981e-05_dp,0.21961e-05_dp,0.43923e-05_dp, &
5389  0.87845e-05_dp,0.17569e-04_dp,0.35138e-04_dp,0.70276e-04_dp,0.14055e-03_dp,0.28110e-03_dp, &
5390  0.56221e-03_dp,0.11244e-02_dp,0.22488e-02_dp,0.44977e-02_dp,0.89954e-02_dp,0.17991e-01_dp, &
5391  0.35981e-01_dp,0.71963e-01_dp,0.14393e+00_dp /
5392  data xi(1:ndat,1) / &
5393  0.33510e-10_dp,0.67021e-10_dp,0.13404e-09_dp,0.26808e-09_dp,0.53617e-09_dp,0.10723e-08_dp, &
5394  0.21447e-08_dp,0.42893e-08_dp,0.85786e-08_dp,0.17157e-07_dp,0.34315e-07_dp,0.68629e-07_dp, &
5395  0.13726e-06_dp,0.27452e-06_dp,0.54903e-06_dp,0.10981e-05_dp,0.21961e-05_dp,0.43923e-05_dp, &
5396  0.87845e-05_dp,0.17569e-04_dp,0.35138e-04_dp,0.70276e-04_dp,0.14055e-03_dp,0.28110e-03_dp, &
5397  0.56221e-03_dp,0.11244e-02_dp,0.22488e-02_dp,0.44977e-02_dp,0.89954e-02_dp,0.17991e-01_dp, &
5398  0.35981e-01_dp,0.71963e-01_dp,0.14393e+00_dp /
5399  data xi(1:ndat,2) / &
5400  0.33510e-10_dp,0.67021e-10_dp,0.13404e-09_dp,0.26808e-09_dp,0.53617e-09_dp,0.10723e-08_dp, &
5401  0.21447e-08_dp,0.42893e-08_dp,0.85786e-08_dp,0.17157e-07_dp,0.34315e-07_dp,0.68629e-07_dp, &
5402  0.13726e-06_dp,0.27452e-06_dp,0.54903e-06_dp,0.10981e-05_dp,0.21961e-05_dp,0.43923e-05_dp, &
5403  0.87845e-05_dp,0.17569e-04_dp,0.35138e-04_dp,0.70276e-04_dp,0.14055e-03_dp,0.28110e-03_dp, &
5404  0.56221e-03_dp,0.11244e-02_dp,0.22488e-02_dp,0.44977e-02_dp,0.89954e-02_dp,0.17991e-01_dp, &
5405  0.35981e-01_dp,0.71963e-01_dp,0.14393e+00 /
5406  data xi(1:ndat,3) / &
5407  0.33510e-10_dp,0.67021e-10_dp,0.13404e-09_dp,0.26808e-09_dp,0.53617e-09_dp,0.10723e-08_dp, &
5408  0.21447e-08_dp,0.42893e-08_dp,0.85786e-08_dp,0.17157e-07_dp,0.34315e-07_dp,0.68629e-07_dp, &
5409  0.13726e-06_dp,0.27452e-06_dp,0.54903e-06_dp,0.10981e-05_dp,0.21961e-05_dp,0.43923e-05_dp, &
5410  0.87845e-05_dp,0.17569e-04_dp,0.35138e-04_dp,0.70276e-04_dp,0.14055e-03_dp,0.28110e-03_dp, &
5411  0.56221e-03_dp,0.11244e-02_dp,0.22488e-02_dp,0.44977e-02_dp,0.89954e-02_dp,0.17991e-01_dp, &
5412  0.35981e-01_dp,0.71963e-01_dp,0.14393e+00_dp /
5413  data xs(1:ndat) / &
5414  0.33510e-10_dp,0.67021e-10_dp,0.13404e-09_dp,0.26808e-09_dp,0.53617e-09_dp,0.10723e-08_dp, &
5415  0.21447e-08_dp,0.42893e-08_dp,0.85786e-08_dp,0.17157e-07_dp,0.34315e-07_dp,0.68629e-07_dp, &
5416  0.13726e-06_dp,0.27452e-06_dp,0.54903e-06_dp,0.10981e-05_dp,0.21961e-05_dp,0.43923e-05_dp, &
5417  0.87845e-05_dp,0.17569e-04_dp,0.35138e-04_dp,0.70276e-04_dp,0.14055e-03_dp,0.28110e-03_dp, &
5418  0.56221e-03_dp,0.11244e-02_dp,0.22488e-02_dp,0.44977e-02_dp,0.89954e-02_dp,0.17991e-01_dp, &
5419  0.35981e-01_dp,0.71963e-01_dp,0.14393e+00_dp /
5420  data xg(1:ndat) / &
5421  0.33510e-10_dp,0.67021e-10_dp,0.13404e-09_dp,0.26808e-09_dp,0.53617e-09_dp,0.10723e-08_dp, &
5422  0.21447e-08_dp,0.42893e-08_dp,0.85786e-08_dp,0.17157e-07_dp,0.34315e-07_dp,0.68629e-07_dp, &
5423  0.13726e-06_dp,0.27452e-06_dp,0.54903e-06_dp,0.10981e-05_dp,0.21961e-05_dp,0.43923e-05_dp, &
5424  0.87845e-05_dp,0.17569e-04_dp,0.35138e-04_dp,0.70276e-04_dp,0.14055e-03_dp,0.28110e-03_dp, &
5425  0.56221e-03_dp,0.11244e-02_dp,0.22488e-02_dp,0.44977e-02_dp,0.89954e-02_dp,0.17991e-01_dp, &
5426  0.35981e-01_dp,0.71963e-01_dp,0.14393e+00_dp /
5427  data xh(1:ndat) / &
5428  0.33510e-10_dp,0.67021e-10_dp,0.13404e-09_dp,0.26808e-09_dp,0.53617e-09_dp,0.10723e-08_dp, &
5429  0.21447e-08_dp,0.42893e-08_dp,0.85786e-08_dp,0.17157e-07_dp,0.34315e-07_dp,0.68629e-07_dp, &
5430  0.13726e-06_dp,0.27452e-06_dp,0.54903e-06_dp,0.10981e-05_dp,0.21961e-05_dp,0.43923e-05_dp, &
5431  0.87845e-05_dp,0.17569e-04_dp,0.35138e-04_dp,0.70276e-04_dp,0.14055e-03_dp,0.28110e-03_dp, &
5432  0.56221e-03_dp,0.11244e-02_dp,0.22488e-02_dp,0.44977e-02_dp,0.89954e-02_dp,0.17991e-01_dp, &
5433  0.35981e-01_dp,0.71963e-01_dp,0.14393e+00_dp /
5434  data rlec(1:ndat) / &
5435  0.20000e-03_dp,0.25198e-03_dp,0.31748e-03_dp,0.40000e-03_dp,0.50397e-03_dp,0.63496e-03_dp, &
5436  0.80000e-03_dp,0.10079e-02_dp,0.12699e-02_dp,0.16000e-02_dp,0.20159e-02_dp,0.25398e-02_dp, &
5437  0.32000e-02_dp,0.40317e-02_dp,0.50797e-02_dp,0.64000e-02_dp,0.80635e-02_dp,0.10159e-01_dp, &
5438  0.12800e-01_dp,0.16127e-01_dp,0.20319e-01_dp,0.25600e-01_dp,0.32254e-01_dp,0.40637e-01_dp, &
5439  0.51200e-01_dp,0.64508e-01_dp,0.81275e-01_dp,0.10240e+00_dp,0.12902e+00_dp,0.16255e+00_dp, &
5440  0.20480e+00_dp,0.25803e+00_dp,0.32510e+00_dp /
5441  data riec(1:ndat,1) / &
5442  0.31936e-03_dp,0.40397e-03_dp,0.51099e-03_dp,0.64638e-03_dp,0.81764e-03_dp,0.10343e-02_dp, &
5443  0.13084e-02_dp,0.16551e-02_dp,0.20937e-02_dp,0.26486e-02_dp,0.33506e-02_dp,0.42387e-02_dp, &
5444  0.64360e-02_dp,0.81426e-02_dp,0.10302e-01_dp,0.13035e-01_dp,0.16494e-01_dp,0.20872e-01_dp, &
5445  0.26412e-01_dp,0.33426e-01_dp,0.42304e-01_dp,0.53543e-01_dp,0.67770e-01_dp,0.85783e-01_dp, &
5446  0.10859e+00_dp,0.13746e+00_dp,0.17403e+00_dp,0.22032e+00_dp,0.27895e+00_dp,0.35319e+00_dp, &
5447  0.44722e+00_dp,0.56630e+00_dp,0.71712e+00_dp /
5448  data riec(1:ndat,2) / &
5449  0.13188e-03_dp,0.16615e-03_dp,0.20953e-03_dp,0.27728e-03_dp,0.36694e-03_dp,0.48559e-03_dp, &
5450  0.64261e-03_dp,0.85040e-03_dp,0.11254e-02_dp,0.14893e-02_dp,0.19709e-02_dp,0.26082e-02_dp, &
5451  0.34515e-02_dp,0.45676e-02_dp,0.60446e-02_dp,0.79991e-02_dp,0.10586e-01_dp,0.14009e-01_dp, &
5452  0.18539e-01_dp,0.24533e-01_dp,0.32466e-01_dp,0.42964e-01_dp,0.56857e-01_dp,0.75242e-01_dp, &
5453  0.99573e-01_dp,0.13177e+00_dp,0.17438e+00_dp,0.23077e+00_dp,0.30539e+00_dp,0.40414e+00_dp, &
5454  0.53482e+00_dp,0.70775e+00_dp,0.93661e+00_dp /
5455  data riec(1:ndat,3) / 0.14998e-03_dp,0.18896e-03_dp,0.23808e-03_dp, &
5456  0.29996e-03_dp,0.37793e-03_dp,0.47616e-03_dp,0.61048e-03_dp,0.81343e-03_dp,0.10839e-02_dp, &
5457  0.14442e-02_dp,0.19243e-02_dp,0.25640e-02_dp,0.34164e-02_dp,0.45522e-02_dp,0.60656e-02_dp, &
5458  0.80820e-02_dp,0.10769e-01_dp,0.14349e-01_dp,0.19119e-01_dp,0.25475e-01_dp,0.44576e-01_dp, &
5459  0.62633e-01_dp,0.88006e-01_dp,0.12366e+00_dp,0.17375e+00_dp,0.24414e+00_dp,0.34304e+00_dp, &
5460  0.48201e+00_dp,0.67728e+00_dp,0.95164e+00_dp,0.13372e+01_dp,0.18788e+01_dp,0.26400e+01_dp /
5461  data rsec(1:ndat) / &
5462  0.92832e-03_dp,0.11696e-02_dp,0.14736e-02_dp,0.18566e-02_dp,0.23392e-02_dp,0.29472e-02_dp, &
5463  0.37133e-02_dp,0.46784e-02_dp,0.58944e-02_dp,0.74265e-02_dp,0.93569e-02_dp,0.11789e-01_dp, &
5464  0.14853e-01_dp,0.18714e-01_dp,0.23578e-01_dp,0.29706e-01_dp,0.37427e-01_dp,0.47156e-01_dp, &
5465  0.59412e-01_dp,0.74855e-01_dp,0.94311e-01_dp,0.11882e+00_dp,0.14971e+00_dp,0.18862e+00_dp, &
5466  0.23765e+00_dp,0.29942e+00_dp,0.37724e+00_dp,0.47530e+00_dp,0.59884e+00_dp,0.75449e+00_dp, &
5467  0.95060e+00_dp,0.11977e+01_dp,0.15090e+01_dp /
5468  data rgec(1:ndat) / 0.27144e-03_dp,0.34200e-03_dp,0.43089e-03_dp, &
5469  0.54288e-03_dp,0.68399e-03_dp,0.86177e-03_dp,0.10858e-02_dp,0.13680e-02_dp,0.17235e-02_dp, &
5470  0.21715e-02_dp,0.27360e-02_dp,0.34471e-02_dp,0.43431e-02_dp,0.54719e-02_dp,0.68942e-02_dp, &
5471  0.86861e-02_dp,0.10944e-01_dp,0.13788e-01_dp,0.17372e-01_dp,0.21888e-01_dp,0.27577e-01_dp, &
5472  0.34745e-01_dp,0.43775e-01_dp,0.55154e-01_dp,0.69489e-01_dp,0.87551e-01_dp,0.11031e+00_dp, &
5473  0.13898e+00_dp,0.17510e+00_dp,0.22061e+00_dp,0.27796e+00_dp,0.35020e+00_dp,0.44123e+00_dp /
5474  data rhec(1:ndat) / &
5475  0.20715e-03_dp,0.26099e-03_dp,0.32883e-03_dp,0.41430e-03_dp,0.52198e-03_dp,0.65766e-03_dp, &
5476  0.82860e-03_dp,0.10440e-02_dp,0.13153e-02_dp,0.16572e-02_dp,0.20879e-02_dp,0.26306e-02_dp, &
5477  0.33144e-02_dp,0.41759e-02_dp,0.52613e-02_dp,0.66288e-02_dp,0.83517e-02_dp,0.10523e-01_dp, &
5478  0.13258e-01_dp,0.16703e-01_dp,0.21045e-01_dp,0.26515e-01_dp,0.33407e-01_dp,0.42090e-01_dp, &
5479  0.53030e-01_dp,0.66814e-01_dp,0.84180e-01_dp,0.10606e+00_dp,0.13363e+00_dp,0.16836e+00_dp, &
5480  0.21212e+00_dp,0.26725e+00_dp,0.33672e+00_dp /
5481  data vrl(1:ndat) / &
5482  0.50000e-01_dp,0.78000e-01_dp,0.12000e+00_dp,0.19000e+00_dp,0.31000e+00_dp,0.49000e+00_dp, &
5483  0.77000e+00_dp,0.12000e+01_dp,0.19000e+01_dp,0.30000e+01_dp,0.48000e+01_dp,0.74000e+01_dp, &
5484  0.11000e+02_dp,0.17000e+02_dp,0.26000e+02_dp,0.37000e+02_dp,0.52000e+02_dp,0.71000e+02_dp, &
5485  0.94000e+02_dp,0.12000e+03_dp,0.16000e+03_dp,0.21000e+03_dp,0.26000e+03_dp,0.33000e+03_dp, &
5486  0.41000e+03_dp,0.48000e+03_dp,0.57000e+03_dp,0.66000e+03_dp,0.75000e+03_dp,0.82000e+03_dp, &
5487  0.88000e+03_dp,0.90000e+03_dp,0.90000e+03_dp /
5488  data vri(1:ndat,1) / 0.30000e-01_dp,0.40000e-01_dp,0.60000e-01_dp, &
5489  0.80000e-01_dp,0.11000e+00_dp,0.15000e+00_dp,0.17000e+00_dp,0.18000e+00_dp,0.20000e+00_dp, &
5490  0.25000e+00_dp,0.40000e+00_dp,0.60000e+01_dp,0.10000e+02_dp,0.15000e+02_dp,0.20000e+02_dp, &
5491  0.25000e+02_dp,0.31000e+02_dp,0.37000e+02_dp,0.41000e+02_dp,0.46000e+02_dp,0.51000e+02_dp, &
5492  0.55000e+02_dp,0.59000e+02_dp,0.62000e+02_dp,0.64000e+02_dp,0.67000e+02_dp,0.68000e+02_dp, &
5493  0.69000e+02_dp,0.70000e+02_dp,0.71000e+02_dp,0.71500e+02_dp,0.71750e+02_dp,0.72000e+02_dp /
5494  data vri(1:ndat,2) / &
5495  0.30000e-01_dp,0.40000e-01_dp,0.50000e-01_dp,0.70000e-01_dp,0.90000e-01_dp,0.12000e+00_dp, &
5496  0.50000e+00_dp,0.80000e+00_dp,0.16000e+01_dp,0.18000e+01_dp,0.20000e+01_dp,0.30000e+01_dp, &
5497  0.40000e+01_dp,0.50000e+01_dp,0.80000e+01_dp,0.13000e+02_dp,0.19000e+02_dp,0.26000e+02_dp, &
5498  0.32000e+02_dp,0.38000e+02_dp,0.47000e+02_dp,0.55000e+02_dp,0.65000e+02_dp,0.73000e+02_dp, &
5499  0.77000e+02_dp,0.79000e+02_dp,0.80000e+02_dp,0.81000e+02_dp,0.81000e+02_dp,0.82000e+02_dp, &
5500  0.82000e+02_dp,0.82000e+02_dp,0.82000e+02_dp /
5501  data vri(1:ndat,3) / 0.35000e-01_dp,0.45000e-01_dp,0.55000e-01_dp, &
5502  0.75000e-01_dp,0.95000e-01_dp,0.13000e+00_dp,0.60000e+00_dp,0.90000e+00_dp,0.17000e+01_dp, &
5503  0.20000e+01_dp,0.25000e+01_dp,0.38000e+01_dp,0.50000e+01_dp,0.70000e+01_dp,0.90000e+01_dp, &
5504  0.11000e+02_dp,0.14000e+02_dp,0.17000e+02_dp,0.21000e+02_dp,0.25000e+02_dp,0.32000e+02_dp, &
5505  0.38000e+02_dp,0.44000e+02_dp,0.49000e+02_dp,0.53000e+02_dp,0.55000e+02_dp,0.58000e+02_dp, &
5506  0.59000e+02_dp,0.61000e+02_dp,0.62000e+02_dp,0.63000e+02_dp,0.64000e+02_dp,0.65000e+02_dp /
5507  data vrs(1:ndat) / &
5508  0.20000e-01_dp,0.31000e-01_dp,0.49000e-01_dp,0.77000e-01_dp,0.12000e+00_dp,0.19000e+00_dp, &
5509  0.30000e+00_dp,0.48000e+00_dp,0.76000e+00_dp,0.12000e+01_dp,0.19000e+01_dp,0.30000e+01_dp, &
5510  0.48000e+01_dp,0.75000e+01_dp,0.11000e+02_dp,0.16000e+02_dp,0.21000e+02_dp,0.26000e+02_dp, &
5511  0.34000e+02_dp,0.41000e+02_dp,0.49000e+02_dp,0.57000e+02_dp,0.65000e+02_dp,0.73000e+02_dp, &
5512  0.81000e+02_dp,0.87000e+02_dp,0.93000e+02_dp,0.99000e+02_dp,0.10750e+03_dp,0.11500e+03_dp, &
5513  0.12500e+03_dp,0.13500e+03_dp,0.14500e+03_dp /
5514  data vrg(1:ndat) / 0.39000e-01_dp,0.62000e-01_dp,0.97000e-01_dp, &
5515  0.15000e+00_dp,0.24000e+00_dp,0.38000e+00_dp,0.61000e+00_dp,0.96000e+00_dp,0.15000e+01_dp, &
5516  0.24000e+01_dp,0.38000e+01_dp,0.61000e+01_dp,0.96000e+01_dp,0.15000e+02_dp,0.23000e+02_dp, &
5517  0.31000e+02_dp,0.39000e+02_dp,0.49000e+02_dp,0.59000e+02_dp,0.68000e+02_dp,0.79000e+02_dp, &
5518  0.88000e+02_dp,0.10000e+03_dp,0.11000e+03_dp,0.13000e+03_dp,0.15000e+03_dp,0.17000e+03_dp, &
5519  0.20000e+03_dp,0.23000e+03_dp,0.26000e+03_dp,0.30000e+03_dp,0.35000e+03_dp,0.40000e+03_dp /
5520  data vrh(1:ndat) / &
5521  0.53000e-01_dp,0.84000e-01_dp,0.13000e+00_dp,0.21000e+00_dp,0.33000e+00_dp,0.52000e+00_dp, &
5522  0.82000e+00_dp,0.13000e+01_dp,0.21000e+01_dp,0.33000e+01_dp,0.52000e+01_dp,0.82000e+01_dp, &
5523  0.13000e+02_dp,0.20000e+02_dp,0.28000e+02_dp,0.36000e+02_dp,0.46000e+02_dp,0.56000e+02_dp, &
5524  0.67000e+02_dp,0.80000e+02_dp,0.97000e+02_dp,0.12000e+03_dp,0.14000e+03_dp,0.17000e+03_dp, &
5525  0.20000e+03_dp,0.24000e+03_dp,0.29000e+03_dp,0.35000e+03_dp,0.42000e+03_dp,0.51000e+03_dp, &
5526  0.61000e+03_dp,0.74000e+03_dp,0.89000e+03_dp /
5527  data blkradl(1:ndat) / &
5528  0.20000e-03_dp,0.25198e-03_dp,0.31748e-03_dp,0.40000e-03_dp,0.50397e-03_dp,0.63496e-03_dp, &
5529  0.80000e-03_dp,0.10079e-02_dp,0.12699e-02_dp,0.16000e-02_dp,0.20159e-02_dp,0.25398e-02_dp, &
5530  0.32000e-02_dp,0.40317e-02_dp,0.50797e-02_dp,0.64000e-02_dp,0.80635e-02_dp,0.10159e-01_dp, &
5531  0.12800e-01_dp,0.16127e-01_dp,0.20319e-01_dp,0.25600e-01_dp,0.32254e-01_dp,0.40637e-01_dp, &
5532  0.51200e-01_dp,0.64508e-01_dp,0.81275e-01_dp,0.10240e+00_dp,0.12902e+00_dp,0.16255e+00_dp, &
5533  0.20480e+00_dp,0.25803e+00_dp,0.32510e+00_dp /
5534  data blkradi(1:ndat,1) / 0.57452e-03_dp,0.72384e-03_dp,0.91199e-03_dp, &
5535  0.11490e-02_dp,0.14477e-02_dp,0.18240e-02_dp,0.22981e-02_dp,0.28954e-02_dp,0.36479e-02_dp, &
5536  0.45961e-02_dp,0.57908e-02_dp,0.72959e-02_dp,0.11572e-01_dp,0.14770e-01_dp,0.18851e-01_dp, &
5537  0.24060e-01_dp,0.30709e-01_dp,0.39194e-01_dp,0.50025e-01_dp,0.63848e-01_dp,0.81491e-01_dp, &
5538  0.10401e+00_dp,0.13275e+00_dp,0.16943e+00_dp,0.21625e+00_dp,0.27601e+00_dp,0.35228e+00_dp, &
5539  0.44962e+00_dp,0.57387e+00_dp,0.73244e+00_dp,0.93484e+00_dp,0.11932e+01_dp,0.15229e+01_dp /
5540  data blkradi(1:ndat,2) / &
5541  0.20715e-03_dp,0.26099e-03_dp,0.32912e-03_dp,0.43555e-03_dp,0.57638e-03_dp,0.76276e-03_dp, &
5542  0.10094e-02_dp,0.13358e-02_dp,0.17677e-02_dp,0.23394e-02_dp,0.30958e-02_dp,0.40969e-02_dp, &
5543  0.54216e-02_dp,0.71748e-02_dp,0.94948e-02_dp,0.12565e-01_dp,0.16628e-01_dp,0.22005e-01_dp, &
5544  0.29120e-01_dp,0.38537e-01_dp,0.50998e-01_dp,0.67488e-01_dp,0.89311e-01_dp,0.11819e+00_dp, &
5545  0.15641e+00_dp,0.20698e+00_dp,0.27391e+00_dp,0.36249e+00_dp,0.47970e+00_dp,0.63481e+00_dp, &
5546  0.84009e+00_dp,0.11117e+01_dp,0.14712e+01_dp /
5547  data blkradi(1:ndat,3) / 0.23559e-03_dp,0.29682e-03_dp,0.37397e-03_dp, &
5548  0.47118e-03_dp,0.59365e-03_dp,0.74795e-03_dp,0.95894e-03_dp,0.12777e-02_dp,0.17025e-02_dp, &
5549  0.22685e-02_dp,0.30227e-02_dp,0.40275e-02_dp,0.53665e-02_dp,0.71506e-02_dp,0.95278e-02_dp, &
5550  0.12695e-01_dp,0.16916e-01_dp,0.22539e-01_dp,0.30032e-01_dp,0.40017e-01_dp,0.70019e-01_dp, &
5551  0.98384e-01_dp,0.13824e+00_dp,0.19424e+00_dp,0.27293e+00_dp,0.38350e+00_dp,0.53885e+00_dp, &
5552  0.75714e+00_dp,0.10639e+01_dp,0.14948e+01_dp,0.21004e+01_dp,0.29513e+01_dp,0.41469e+01_dp /
5553  data blkrads(1:ndat) / &
5554  0.20715e-03_dp,0.26148e-03_dp,0.33067e-03_dp,0.41710e-03_dp,0.52691e-03_dp,0.66640e-03_dp, &
5555  0.84289e-03_dp,0.10674e-02_dp,0.13513e-02_dp,0.17129e-02_dp,0.21670e-02_dp,0.27521e-02_dp, &
5556  0.34989e-02_dp,0.44777e-02_dp,0.57347e-02_dp,0.75389e-02_dp,0.99020e-02_dp,0.13161e-01_dp, &
5557  0.17372e-01_dp,0.23337e-01_dp,0.31058e-01_dp,0.41194e-01_dp,0.55153e-01_dp,0.74854e-01_dp, &
5558  0.99806e-01_dp,0.13463e+00_dp,0.18136e+00_dp,0.24282e+00_dp,0.32955e+00_dp,0.44123e+00_dp, &
5559  0.59884e+00_dp,0.77090e+00_dp,0.99387e+00_dp /
5560  data blkradg(1:ndat) / 0.27144e-03_dp,0.34200e-03_dp,0.43089e-03_dp, &
5561  0.54288e-03_dp,0.68399e-03_dp,0.86177e-03_dp,0.10858e-02_dp,0.13680e-02_dp,0.17235e-02_dp, &
5562  0.21715e-02_dp,0.27360e-02_dp,0.34471e-02_dp,0.43431e-02_dp,0.54719e-02_dp,0.68942e-02_dp, &
5563  0.86861e-02_dp,0.10944e-01_dp,0.13788e-01_dp,0.17372e-01_dp,0.21888e-01_dp,0.27577e-01_dp, &
5564  0.34745e-01_dp,0.43775e-01_dp,0.55154e-01_dp,0.69489e-01_dp,0.87551e-01_dp,0.11031e+00_dp, &
5565  0.13898e+00_dp,0.17510e+00_dp,0.22061e+00_dp,0.27796e+00_dp,0.35020e+00_dp,0.44123e+00_dp /
5566  data blkradh(1:ndat) / &
5567  0.20715e-03_dp,0.26099e-03_dp,0.32883e-03_dp,0.41430e-03_dp,0.52198e-03_dp,0.65766e-03_dp, &
5568  0.82860e-03_dp,0.10440e-02_dp,0.13153e-02_dp,0.16572e-02_dp,0.20879e-02_dp,0.26306e-02_dp, &
5569  0.33144e-02_dp,0.41759e-02_dp,0.52613e-02_dp,0.66288e-02_dp,0.83517e-02_dp,0.10523e-01_dp, &
5570  0.13258e-01_dp,0.16703e-01_dp,0.21045e-01_dp,0.26515e-01_dp,0.33407e-01_dp,0.42090e-01_dp, &
5571  0.53030e-01_dp,0.66814e-01_dp,0.84180e-01_dp,0.10606e+00_dp,0.13363e+00_dp,0.16836e+00_dp, &
5572  0.21212e+00_dp,0.26725e+00_dp,0.33672e+00_dp /
5573  data blkdnsl(1:ndat) / &
5574  0.10000e+01_dp,0.10000e+01_dp,0.10000e+01_dp,0.10000e+01_dp,0.10000e+01_dp,0.10000e+01_dp, &
5575  0.10000e+01_dp,0.10000e+01_dp,0.10000e+01_dp,0.10000e+01_dp,0.10000e+01_dp,0.10000e+01_dp, &
5576  0.10000e+01_dp,0.10000e+01_dp,0.10000e+01_dp,0.10000e+01_dp,0.10000e+01_dp,0.10000e+01_dp, &
5577  0.10000e+01_dp,0.10000e+01_dp,0.10000e+01_dp,0.10000e+01_dp,0.10000e+01_dp,0.10000e+01_dp, &
5578  0.10000e+01_dp,0.10000e+01_dp,0.10000e+01_dp,0.10000e+01_dp,0.10000e+01_dp,0.10000e+01_dp, &
5579  0.10000e+01_dp,0.10000e+01_dp,0.10000e+01_dp /
5580  data blkdnsi(1:ndat,1) / 0.90000e+00_dp,0.90000e+00_dp,0.90000e+00_dp, &
5581  0.90000e+00_dp,0.90000e+00_dp,0.90000e+00_dp,0.90000e+00_dp,0.90000e+00_dp,0.90000e+00_dp, &
5582  0.90000e+00_dp,0.90000e+00_dp,0.90000e+00_dp,0.87368e+00_dp,0.87072e+00_dp,0.86777e+00_dp, &
5583  0.86483e+00_dp,0.86189e+00_dp,0.85897e+00_dp,0.85606e+00_dp,0.85316e+00_dp,0.85026e+00_dp, &
5584  0.84738e+00_dp,0.84451e+00_dp,0.84164e+00_dp,0.83879e+00_dp,0.83595e+00_dp,0.83311e+00_dp, &
5585  0.83029e+00_dp,0.82747e+00_dp,0.82467e+00_dp,0.82187e+00_dp,0.81908e+00_dp,0.81631e+00_dp /
5586  data blkdnsi(1:ndat,2) / &
5587  0.90000e+00_dp,0.90000e+00_dp,0.90000e+00_dp,0.90000e+00_dp,0.90000e+00_dp,0.90000e+00_dp, &
5588  0.90000e+00_dp,0.90000e+00_dp,0.90000e+00_dp,0.90000e+00_dp,0.90000e+00_dp,0.90000e+00_dp, &
5589  0.90000e+00_dp,0.90000e+00_dp,0.90000e+00_dp,0.90000e+00_dp,0.90000e+00_dp,0.90000e+00_dp, &
5590  0.90000e+00_dp,0.90000e+00_dp,0.90000e+00_dp,0.90000e+00_dp,0.90000e+00_dp,0.90000e+00_dp, &
5591  0.90000e+00_dp,0.90000e+00_dp,0.90000e+00_dp,0.90000e+00_dp,0.90000e+00_dp,0.90000e+00_dp, &
5592  0.90000e+00_dp,0.90000e+00_dp,0.90000e+00_dp /
5593  data blkdnsi(1:ndat,3) / 0.61183e+00_dp,0.61183e+00_dp,0.61183e+00_dp, &
5594  0.61183e+00_dp,0.61183e+00_dp,0.61183e+00_dp,0.61183e+00_dp,0.61183e+00_dp,0.61183e+00_dp, &
5595  0.61183e+00_dp,0.61183e+00_dp,0.61183e+00_dp,0.61183e+00_dp,0.61183e+00_dp,0.61183e+00_dp, &
5596  0.61183e+00_dp,0.61183e+00_dp,0.61183e+00_dp,0.61183e+00_dp,0.61183e+00_dp,0.51790e+00_dp, &
5597  0.45557e+00_dp,0.40075e+00_dp,0.35252e+00_dp,0.31010e+00_dp,0.27278e+00_dp,0.23995e+00_dp, &
5598  0.21108e+00_dp,0.18567e+00_dp,0.16333e+00_dp,0.14367e+00_dp,0.12638e+00_dp,0.11118e+00_dp /
5599  data blkdnss(1:ndat) / &
5600  0.90000e+00_dp,0.89500e+00_dp,0.88500e+00_dp,0.88200e+00_dp,0.87500e+00_dp,0.86500e+00_dp, &
5601  0.85500e+00_dp,0.84200e+00_dp,0.83000e+00_dp,0.81500e+00_dp,0.80500e+00_dp,0.78600e+00_dp, &
5602  0.76500e+00_dp,0.73000e+00_dp,0.69500e+00_dp,0.61183e+00_dp,0.54000e+00_dp,0.46000e+00_dp, &
5603  0.40000e+00_dp,0.33000e+00_dp,0.28000e+00_dp,0.24000e+00_dp,0.20000e+00_dp,0.16000e+00_dp, &
5604  0.13500e+00_dp,0.11000e+00_dp,0.90000e-01_dp,0.75000e-01_dp,0.60000e-01_dp,0.50000e-01_dp, &
5605  0.40000e-01_dp,0.37500e-01_dp,0.35000e-01_dp /
5606  data blkdnsg(1:ndat) / &
5607  0.40000e+00_dp,0.40000e+00_dp,0.40000e+00_dp,0.40000e+00_dp,0.40000e+00_dp,0.40000e+00_dp, &
5608  0.40000e+00_dp,0.40000e+00_dp,0.40000e+00_dp,0.40000e+00_dp,0.40000e+00_dp,0.40000e+00_dp, &
5609  0.40000e+00_dp,0.40000e+00_dp,0.40000e+00_dp,0.40000e+00_dp,0.40000e+00_dp,0.40000e+00_dp, &
5610  0.40000e+00_dp,0.40000e+00_dp,0.40000e+00_dp,0.40000e+00_dp,0.40000e+00_dp,0.40000e+00_dp, &
5611  0.40000e+00_dp,0.40000e+00_dp,0.40000e+00_dp,0.40000e+00_dp,0.40000e+00_dp,0.40000e+00_dp, &
5612  0.40000e+00_dp,0.40000e+00_dp,0.40000e+00_dp /
5613  data blkdnsh(1:ndat) / &
5614  0.90000e+00_dp,0.90000e+00_dp,0.90000e+00_dp,0.90000e+00_dp,0.90000e+00_dp,0.90000e+00_dp, &
5615  0.90000e+00_dp,0.90000e+00_dp,0.90000e+00_dp,0.90000e+00_dp,0.90000e+00_dp,0.90000e+00_dp, &
5616  0.90000e+00_dp,0.90000e+00_dp,0.90000e+00_dp,0.90000e+00_dp,0.90000e+00_dp,0.90000e+00_dp, &
5617  0.90000e+00_dp,0.90000e+00_dp,0.90000e+00_dp,0.90000e+00_dp,0.90000e+00_dp,0.90000e+00_dp, &
5618  0.90000e+00_dp,0.90000e+00_dp,0.90000e+00_dp,0.90000e+00_dp,0.90000e+00_dp,0.90000e+00_dp, &
5619  0.90000e+00_dp,0.90000e+00_dp,0.90000e+00_dp /
5620 
5621  !--- mass
5622  do k = 1, ndat
5623  xmss( il,k ) = log( dble( xl( k ) )*1.e-03_dp )
5624  xmss( ic,k ) = log( dble( xi( k,1 ) )*1.e-03_dp )
5625  xmss( ip,k ) = log( dble( xi( k,2 ) )*1.e-03_dp )
5626  xmss( id,k ) = log( dble( xi( k,3 ) )*1.e-03_dp )
5627  xmss( is,k ) = log( dble( xs( k ) )*1.e-03_dp )
5628  xmss( ig,k ) = log( dble( xg( k ) )*1.e-03_dp )
5629  xmss( ih,k ) = log( dble( xh( k ) )*1.e-03_dp )
5630  end do
5631 
5632  !--- capacity
5633  do k = 1, ndat ! cm -> m
5634  zcap( il,k ) = dble(rlec( k ))*1.e-02_dp
5635  zcap( ic,k ) = dble(riec( k,1 ))*1.e-02_dp
5636  zcap( ip,k ) = dble(riec( k,2 ))*1.e-02_dp
5637  zcap( id,k ) = dble(riec( k,3 ))*1.e-02_dp
5638  zcap( is,k ) = dble(rsec( k ))*1.e-02_dp
5639  zcap( ig,k ) = dble(rgec( k ))*1.e-02_dp
5640  zcap( ih,k ) = dble(rhec( k ))*1.e-02_dp
5641  end do
5642 
5643  !--- terminal velocity
5644  do k = 1, ndat ! cm/s -> m/s
5645  vtrm( il,k ) = dble(vrl( k ))*1.e-02_dp
5646  vtrm( ic,k ) = dble(vri( k,1 ))*1.e-02_dp
5647  vtrm( ip,k ) = dble(vri( k,2 ))*1.e-02_dp
5648  vtrm( id,k ) = dble(vri( k,3 ))*1.e-02_dp
5649  vtrm( is,k ) = dble(vrs( k ))*1.e-02_dp
5650  vtrm( ig,k ) = dble(vrg( k ))*1.e-02_dp
5651  vtrm( ih,k ) = dble(vrh( k ))*1.e-02_dp
5652  end do
5653 
5654  !--- bulk radii
5655  do k = 1, ndat ! cm -> mm
5656  blkr( il,k ) = dble(blkradl( k ))*10.0_dp
5657  blkr( ic,k ) = dble(blkradi( k,1 ))*10.0_dp
5658  blkr( ip,k ) = dble(blkradi( k,2 ))*10.0_dp
5659  blkr( id,k ) = dble(blkradi( k,3 ))*10.0_dp
5660  blkr( is,k ) = dble(blkrads( k ))*10.0_dp
5661  blkr( ig,k ) = dble(blkradg( k ))*10.0_dp
5662  blkr( ih,k ) = dble(blkradh( k ))*10.0_dp
5663  end do
5664 
5665  !--- bulk density
5666  do k = 1, ndat ! g/cm^3 -> kg/m^3
5667  blkd( il,k ) = dble(blkdnsl( k ))*1000._dp
5668  blkd( ic,k ) = dble(blkdnsi( k,1 ))*1000._dp
5669  blkd( ip,k ) = dble(blkdnsi( k,2 ))*1000._dp
5670  blkd( id,k ) = dble(blkdnsi( k,3 ))*1000._dp
5671  blkd( is,k ) = dble(blkdnss( k ))*1000._dp
5672  blkd( ig,k ) = dble(blkdnsg( k ))*1000._dp
5673  blkd( ih,k ) = dble(blkdnsh( k ))*1000._dp
5674  end do
5675 
5676  !--- collection kernel
5677  ! cm**3/s -> m**3/s
5678  do i = 1, ndat
5679  do j = 1, ndat
5680  do k = 1, 7
5681  do kk = 1, 7
5682  ykrn( k,kk,i,j ) = kernels( k,kk,i,j )*1.e-06_dp
5683  enddo
5684  enddo
5685 
5686  end do
5687  end do
5688 
5689  end subroutine rdkdat
5690  !---------------------------------------------------------------------------------------
5691  subroutine sdfgrid
5692 
5693  real(dp) :: xsta, xend
5694  integer :: n
5695  real(dp):: pi_dp = 3.1415920_dp
5696  real(dp):: rhow_dp = 1.0e+03_dp
5697 
5698  xsta = log( rhow_dp * 4.0_dp*pi_dp/3.0_dp * ( 3.e-06_dp )**3.0_dp )
5699  xend = log( rhow_dp * 4.0_dp*pi_dp/3.0_dp * ( 3.e-03_dp )**3.0_dp )
5700 
5701  dxmic_mk = ( xend-xsta )/nbin
5702  dxmic = dxmic_mk
5703  do n = 1, nbin+1
5704  xbnd_mk( n ) = xsta + dxmic_mk*( n-1 )
5705  end do
5706  do n = 1, nbin
5707  xctr_mk( n ) = ( xbnd_mk( n )+xbnd_mk( n+1 ) )*0.50_dp
5708  radc_mk( n ) = ( exp( xctr_mk( n ) )*3.0_dp/4.0_dp/pi_dp/rhow_dp )**( 1.0_dp/3.0_dp )
5709  end do
5710 
5711  end subroutine sdfgrid
5712  !---------------------------------------------------------------------------------------
5713  subroutine getcp
5714 
5715  integer :: myu, n
5716 
5717  !$omp parallel do
5718  do n = 1, nbin
5719  do myu = 1, nspc_mk
5720  cctr_mk( myu,n ) = fcpc( myu,xctr_mk( n ) )
5721  end do
5722  end do
5723  !$omp parallel do
5724  do n = 1, nbin+1
5725  do myu = 1, nspc_mk
5726  cbnd_mk( myu,n ) = fcpc( myu,xbnd_mk( n ) )
5727  end do
5728  end do
5729 
5730  end subroutine getcp
5731  !---------------------------------------------------------------------------------------
5732  function fcpc( myu,x )
5733 
5734  integer, intent(in) :: myu
5735  real(dp), intent(in) :: x
5736  real(dp) :: fcpc
5737 
5738  real(dp) :: qknt( ndat+kdeg ), elm( ndat,ndat ), coef( ndat )
5739 
5740  call getknot &
5741  ( ndat, kdeg, myu, xmss(:,:), & !--- in
5742  qknt ) !--- out
5743 
5744  call getmatrx &
5745  ( ndat, kdeg, qknt, myu, xmss(:,:), & !--- in
5746  elm ) !--- out
5747 
5748  call getcoef &
5749  ( ndat, kdeg, elm, myu, zcap(:,:), & !--- in
5750  coef ) !--- out
5751 
5752  fcpc = fspline( ndat, kdeg, coef, qknt, myu, xmss(:,:), x )
5753 
5754  end function fcpc
5755  !---------------------------------------------------------------------------------------
5756  subroutine getck
5757 
5758  integer :: myu, nyu, i, j
5759 
5760  log_info("ATMOS_PHY_MP_SUZUKI10_getck",*) 'Create micpara.dat'
5761  if( kphase == 0 ) then
5762  log_info("ATMOS_PHY_MP_SUZUKI10_getck",*) 'Hydro-dynamic kernel'
5763  else if( kphase == 1 ) then
5764  log_info("ATMOS_PHY_MP_SUZUKI10_getck",*) 'Long Kernel'
5765  else if( kphase == 2 ) then
5766  log_info("ATMOS_PHY_MP_SUZUKI10_getck",*) 'Golovin Kernel'
5767  endif
5768 
5769  !$omp parallel do
5770  do j = 1, nbin
5771  do i = 1, nbin
5772  do nyu = 1, nspc_mk
5773  do myu = 1, nspc_mk
5774  ck_mk( myu,nyu,i,j ) = fckrn( myu,nyu,xctr_mk( i ),xctr_mk( j ) )
5775  end do
5776  end do
5777  end do
5778  end do
5779 
5780  return
5781 
5782  end subroutine getck
5783  !---------------------------------------------------------------------------------------
5784 !OCL SERIAL
5785  function fckrn( myu,nyu,x,y )
5786 
5787  integer, intent(in) :: myu, nyu
5788  real(dp), intent(in) :: x, y
5789  real(dp) :: fckrn
5790 
5791  real(dp) :: qknt( ndat+kdeg ), rknt( ndat+kdeg )
5792  real(dp) :: coef( ndat,ndat )
5793 
5794  real(dp) :: xlrg, xsml, vlrg, vsml, rlrg
5795  real(dp):: pi_dp = 3.1415920_dp
5796  real(dp):: rhow_dp = 1.0e+03_dp
5797 
5798  if( kphase == 0 ) then
5799  call getknot &
5800  ( ndat, kdeg, myu, xmss(:,:), & !--- in
5801  qknt ) !--- out
5802 
5803  rknt( : ) = qknt( : )
5804 
5805  call getcoef2 &
5806  ( ndat, ndat, kdeg, kdeg, &
5807  myu, myu, & !--- in
5808  xmss(:,:), xmss(:,:), & !--- in
5809  qknt, rknt, & !--- in
5810  ykrn(:,:,:,:), & !--- in
5811  coef ) !--- out
5812 
5813  fckrn = fspline2 &
5814  ( ndat, ndat, kdeg, kdeg, & !--- in
5815  coef, qknt, rknt, myu, & !--- in
5816  xmss(:,:), xmss(:,:), & !--- in
5817  x, y ) !--- in
5818  else if( kphase == 1 ) then
5819  xlrg = max( x, y )
5820  xsml = min( x, y )
5821 
5822  vlrg = (exp( xlrg ) / rhow_dp )*1.e+06_dp
5823  vsml = (exp( xsml ) / rhow_dp )*1.e+06_dp
5824 
5825  rlrg = ( exp( xlrg )/( 4.0_dp*pi_dp*rhow_dp )*3.0_dp )**(1.0_dp/3.0_dp )*1.e+06_dp
5826 
5827  if( rlrg <=50.0_dp ) then
5828  fckrn = 9.44e+03_dp*( vlrg*vlrg + vsml*vsml )
5829  else
5830  fckrn = 5.78e-03_dp*( vlrg+vsml )
5831  end if
5832  else if( kphase == 2 ) then
5833  fckrn = 1.5_dp*( exp(x) +exp(y) )
5834  end if
5835 
5836  return
5837 
5838  end function fckrn
5839  !---------------------------------------------------------------------------------------
5840  subroutine getvt
5841 
5842  integer :: myu, n
5843 
5844  !$omp parallel do
5845  do n = 1, nbin
5846  do myu = 1, nspc_mk
5847  vt_mk( myu,n ) = max( fvterm( myu,xctr_mk( n ) ), 0.0_dp )
5848  end do
5849  end do
5850 
5851  end subroutine getvt
5852  !---------------------------------------------------------------------------------------
5853  function fvterm( myu,x )
5854 
5855  integer, intent(in) :: myu
5856  real(dp), intent(in) :: x
5857  real(dp) :: fvterm
5858 
5859  real(dp) :: qknt( ndat+kdeg ), elm( ndat,ndat ), coef( ndat )
5860 
5861  call getknot &
5862  ( ndat, kdeg, myu, xmss(:,:), & !--- in
5863  qknt ) !--- out
5864 
5865  call getmatrx &
5866  ( ndat, kdeg, qknt, myu, xmss(:,:), & !--- in
5867  elm ) !--- out
5868 
5869  call getcoef &
5870  ( ndat, kdeg, elm, myu, vtrm(:,:), & !--- in
5871  coef ) !--- out
5872 
5873  fvterm = fspline( ndat, kdeg, coef, qknt, myu, xmss(:,:), x )
5874 
5875  end function fvterm
5876  !---------------------------------------------------------------------------------------
5877  subroutine getbr
5878 
5879  integer :: myu, n
5880 
5881  !$omp parallel do
5882  do n = 1, nbin
5883  do myu = 1, nspc_mk
5884  br_mk( myu,n ) = fbulkrad( myu, xctr_mk( n ) )
5885  end do
5886  end do
5887 
5888  end subroutine getbr
5889  !---------------------------------------------------------------------------------------
5890  function fbulkrad( myu,x )
5891 
5892  integer, intent(in) :: myu
5893  real(dp), intent(in) :: x
5894  real(dp) :: fbulkrad
5895 
5896  real(dp) :: qknt( ndat+kdeg ), elm( ndat,ndat ), coef( ndat )
5897 
5898  call getknot &
5899  ( ndat, kdeg, myu, xmss(:,:), & !--- in
5900  qknt ) !--- out
5901 
5902  call getmatrx &
5903  ( ndat, kdeg, qknt, myu, xmss(:,:), & !--- in
5904  elm ) !--- out
5905 
5906  call getcoef &
5907  ( ndat, kdeg, elm, myu, blkr(:,:), & !--- in
5908  coef ) !--- out
5909 
5910  fbulkrad = fspline( ndat, kdeg, coef, qknt, myu, xmss(:,:), x )
5911 
5912  end function fbulkrad
5913  !---------------------------------------------------------------------------------------
5914  subroutine paraout
5915 
5916  integer :: myu, nyu, i, j, n
5917 
5918  open ( fid_micpara, file = fname_micpara, form = 'formatted', status='new' )
5919 
5920  write( fid_micpara,* ) nspc_mk, nbin
5921 
5922  ! grid parameter
5923  do n = 1, nbin
5924  xctr( n ) = xctr_mk( n )
5925  radc( n ) = radc_mk( n )
5926  write( fid_micpara,* ) n, xctr( n ), radc( n )
5927  end do
5928  do n = 1, nbin+1
5929  xbnd( n ) = xbnd_mk( n )
5930  write( fid_micpara,* ) n, xbnd( n )
5931  end do
5932  write( fid_micpara,* ) dxmic_mk
5933 
5934  ! capacity
5935  do myu = 1, nspc_mk
5936  do n = 1, nbin
5937  cctr( n,myu ) = cctr_mk( myu,n )
5938  write( fid_micpara,* ) myu, n, cctr( n,myu )
5939  end do
5940  do n = 1, nbin+1
5941  cbnd( n,myu ) = cbnd_mk( myu,n )
5942  write( fid_micpara,* ) myu, n, cbnd( n,myu )
5943  end do
5944  end do
5945 
5946  ! collection kernel
5947  do myu = 1, nspc_mk
5948  do nyu = 1, nspc_mk
5949  do i = 1, nbin
5950  do j = 1, nbin
5951  ck( myu,nyu,i,j ) = ck_mk( myu,nyu,i,j )
5952  write( fid_micpara,* ) myu, nyu, i, j, ck( myu,nyu,i,j )
5953  end do
5954  end do
5955  end do
5956  end do
5957 
5958  ! falling velocity
5959  do myu = 1, nspc_mk
5960  do n = 1, nbin
5961  vt( myu,n ) = vt_mk( myu,n )
5962  write( fid_micpara,* ) myu, n, vt( myu,n )
5963  end do
5964  end do
5965 
5966  ! bulk radius
5967  do myu = 1, nspc_mk
5968  do n = 1, nbin
5969  br( myu,n ) = br_mk( myu,n )
5970  write( fid_micpara,* ) myu, n, br( myu,n )
5971  end do
5972  end do
5973 
5974  close ( fid_micpara )
5975 
5976  end subroutine paraout
5977  !---------------------------------------------------------------------------------------
5978  !---- unify from other files
5979  !---------------------------------------------------------------------------------------
5980 !OCL SERIAL
5981  subroutine tinvss(n,a,dt,e,nn,iw,inder)
5982 
5983  implicit none
5984 
5985  integer, intent(in) :: n, nn
5986  integer, intent(inout) :: inder
5987  real(dp), intent(inout) :: a(nn,n)
5988  real(dp), intent(inout) :: dt
5989  real(dp), intent(in) :: e
5990  integer, intent(inout) :: iw( 2*n )
5991  integer :: i, j, k, kk, ij, nnn
5992  real(dp) :: work, aa, az, eps
5993  integer :: ipiv, jpiv
5994  real(dp) :: piv
5995  !-----------------------------------------
5996 
5997  inder = 0
5998  if( n < 1 ) then
5999  goto 910
6000  elseif( n == 1 ) then
6001  goto 930
6002  elseif( n > 1 ) then
6003  goto 101
6004  endif
6005 
6006 101 continue
6007  if( n > nn ) then
6008  inder = -1
6009  write(6,690) "n= ", n, "nn= ", nn, &
6010  "n should be less than or equal to nn in TINVSS"
6011 690 format(a8,i5,5x,a4,i5,a55)
6012  return
6013  endif
6014 
6015  eps = 0.0_dp
6016  dt = 1.0_dp
6017  do k = 1, n
6018  piv = 0.0_dp
6019  do i = k, n
6020  do j = k, n
6021  if( abs(a(i,j)) <= abs(piv) ) goto 110 !
6022  ipiv = i
6023  jpiv = j
6024  piv = a(i,j)
6025 110 continue
6026  enddo
6027  enddo
6028  dt = dt * piv
6029  if( abs(piv) <= eps ) goto 920
6030  if( k == 1 ) eps = abs(piv)*e
6031  if( ipiv == k ) goto 130
6032  dt = -dt
6033  do j = 1, n
6034  work = a(ipiv,j)
6035  a(ipiv,j) = a(k,j)
6036  a(k,j) = work
6037  enddo
6038 130 continue
6039  if( jpiv == k ) goto 150
6040  dt = -dt
6041  do i = 1, n
6042  work = a(i,jpiv)
6043  a(i,jpiv) = a(i,k)
6044  a(i,k) = work
6045  enddo
6046 150 continue
6047  iw(2*k-1) = ipiv
6048  aa=1.0_dp/piv
6049  iw(2*k) = jpiv
6050  do j = 1, n
6051  a(k,j) = a(k,j)*aa
6052  enddo
6053  do i = 1, n
6054  if( i == k ) goto 220
6055  az = a(i,k)
6056  if( az == 0.0_dp ) goto 220
6057  do j = 1, n
6058  a(i,j) = a(i,j)-a(k,j)*az
6059  enddo
6060  a(i,k) = -aa*az
6061 220 continue
6062  enddo
6063  a(k,k) = aa
6064  enddo
6065  do kk = 2, n
6066  k=n+1-kk
6067  ij=iw(2*k)
6068  if( ij == k ) goto 420
6069  do j = 1, n
6070  work=a(ij,j)
6071  a(ij,j) = a(k,j)
6072  a(k,j) = work
6073  enddo
6074 420 continue
6075  ij = iw(2*k-1)
6076  if( ij == k ) goto 400
6077  do i = 1, n
6078  work=a(i,ij)
6079  a(i,ij)=a(i,k)
6080  a(i,k)=work
6081  enddo
6082 400 continue
6083  enddo
6084 
6085  return
6086 
6087 910 continue
6088  inder = -1
6089  write(*,691) "n= ", n, "should be positive in TINVSS"
6090 691 format(a8,i5,5x,a30)
6091  return
6092 
6093 
6094 920 continue
6095  dt = 0.0_dp
6096  inder = n-k+1
6097  nnn = k-1
6098  write(*,692) 'given matrix A to TINVSS is ill conditioned, or sigular withrank =', nnn, &
6099  'return with no further calculation'
6100 692 format(a,1x,i4,1x,a)
6101  return
6102 
6103 930 continue
6104  dt=a(1,1)
6105  k=1
6106  if( dt == 0.0_dp ) goto 920
6107  a(1,1) = 1.0_dp/a(1,1)
6108  return
6109 
6110  end subroutine tinvss
6111  !---------------------------------------------------------------
6112 !OCL SERIAL
6113  subroutine getknot &
6114  ( ndat, kdeg, myu, xdat, & !--- in
6115  qknt ) !--- out
6116 
6117  integer, intent(in) :: ndat ! number of data
6118  integer, intent(in) :: kdeg ! degree of Spline + 1
6119  integer, intent(in) :: myu
6120  real(dp), intent(in) :: xdat(nspc_mk,ndat) ! data of independent var.
6121 
6122  real(dp), intent(out) :: qknt( ndat+kdeg ) ! knots for B-Spline
6123 
6124  !--- local
6125  integer :: i
6126 
6127  do i = 1, kdeg
6128  qknt( i ) = xdat(myu,1)
6129  end do
6130 
6131  do i = 1, ndat-kdeg
6132  qknt( i+kdeg ) = ( xdat(myu,i)+xdat(myu,i+kdeg) )*0.50_dp
6133  end do
6134 
6135  do i = 1, kdeg
6136  qknt( ndat+i ) = xdat(myu,ndat)
6137  end do
6138 
6139  return
6140 
6141  end subroutine getknot
6142  !---------------------------------------------------------------
6143  recursive function fbspl ( ndat, inum, kdeg, qknt, myu, xdat, x ) &
6144  result(bspl)
6145 
6146  real(dp) :: bspl
6147 
6148  integer, intent(in) :: ndat ! number of data
6149  integer, intent(in) :: inum ! index of B-Spline
6150  integer, intent(in) :: kdeg ! degree of B-Spline + 1
6151  real(dp), intent(in) :: qknt(ndat+kdeg) ! knot of B-Spline
6152  integer, intent(in) :: myu
6153  real(dp), intent(in) :: xdat(nspc_mk,ndat) ! data of independent variable
6154  real(dp), intent(in) :: x ! interpolation point
6155 
6156  !--- local
6157  real(dp) :: bsp1, bsp2
6158 
6159  if ( ( inum == 1 .AND. x == xdat(myu,1) ) .OR. &
6160  ( inum == ndat .AND. x == xdat(myu,ndat) ) ) then
6161  bspl = 1.
6162  return
6163  end if
6164 
6165  if ( kdeg == 1 ) then
6166  if ( x >= qknt( inum ) .AND. x < qknt( inum+1 ) ) then
6167  bspl = 1.0_dp
6168  else
6169  bspl = 0.0_dp
6170  end if
6171  else
6172  if ( qknt( inum+kdeg-1 ) /= qknt( inum ) ) then
6173  bsp1 = ( x-qknt( inum ) ) &
6174  /( qknt( inum+kdeg-1 )-qknt( inum ) ) &
6175  * fbspl( ndat, inum, kdeg-1, qknt, myu, xdat, x )
6176  else
6177  bsp1 = 0.0_dp
6178  end if
6179  if ( qknt( inum+kdeg ) /= qknt( inum+1 ) ) then
6180  bsp2 = ( qknt( inum+kdeg )-x ) &
6181  /( qknt( inum+kdeg )-qknt( inum+1 ) ) &
6182  * fbspl( ndat, inum+1, kdeg-1, qknt, myu, xdat, x )
6183  else
6184  bsp2 = 0.0_dp
6185  end if
6186  bspl = bsp1 + bsp2
6187  end if
6188 
6189  end function fbspl
6190 !!$ !---------------------------------------------------------------
6191 !!$ function fpb( ndat, i, kdeg, qknt, myu, xdat, elm, x )
6192 !!$
6193 !!$ real :: fpb
6194 !!$ integer :: ndat, i, kdeg
6195 !!$ integer :: myu
6196 !!$ real(DP) :: qknt(ndat+kdeg), xdat(nspc_mk,ndat), elm(ndat,ndat)
6197 !!$ real(DP) :: x
6198 !!$
6199 !!$ integer :: l
6200 !!$ real(DP) :: sum
6201 !!$
6202 !!$ sum = 0.0_DP
6203 !!$ do l = 1, ndat
6204 !!$ sum = sum + elm( l,i )*fbspl( ndat, l, kdeg, qknt, myu, xdat, x )
6205 !!$ end do
6206 !!$
6207 !!$ fpb = sum
6208 !!$
6209 !!$ end function fpb
6210  !---------------------------------------------------------------
6211 !OCL SERIAL
6212  subroutine getmatrx &
6213  ( ndat, kdeg, qknt, myu, xdat, & !--- in
6214  elm ) !--- out
6215 
6216 ! use scale_tinvss, only: TINVSS
6217 
6218  integer, intent(in) :: ndat
6219  integer, intent(in) :: kdeg
6220  real(dp), intent(in) :: qknt(ndat+kdeg)
6221  integer, intent(in) :: myu
6222  real(dp), intent(in) :: xdat(nspc_mk,ndat)
6223 
6224  real(dp), intent(out) :: elm(ndat,ndat)
6225 
6226  !--- local
6227  real(dp) :: dt
6228  integer :: iw( 2*ndat ), i, j, inder
6229  real(dp), parameter :: eps = 0.
6230 
6231  do i = 1, ndat
6232  do j = 1, ndat
6233  elm( i,j ) = fbspl( ndat, j, kdeg, qknt, myu, xdat, xdat(myu,i) )
6234  end do
6235  end do
6236 
6237  call tinvss( ndat, elm, dt, eps, ndat, iw, inder )
6238 
6239  return
6240 
6241  end subroutine getmatrx
6242  !---------------------------------------------------------------
6243 !OCL SERIAL
6244  subroutine getcoef &
6245  ( ndat, kdeg, elm, myu, ydat, & !--- in
6246  coef ) !--- out
6247 
6248  integer, intent(in) :: ndat ! number of data
6249  integer, intent(in) :: kdeg ! degree of Spline + 1
6250  real(dp), intent(in) :: elm(ndat,ndat) ! matrix ( inverse )
6251  integer, intent(in) :: myu
6252  real(dp), intent(in) :: ydat(nspc_mk,ndat) ! data of dependent var.
6253 
6254  real(dp), intent(out) :: coef(ndat) ! expansion coefficient
6255 
6256  !--- local
6257  integer :: i, j
6258  real(dp) :: sum
6259 
6260  do i = 1, ndat
6261  sum = 0.0_dp
6262  do j = 1, ndat
6263  sum = sum + elm( i,j )*ydat(myu,j)
6264  end do
6265  coef( i ) = sum
6266  end do
6267 
6268  return
6269 
6270  end subroutine getcoef
6271  !---------------------------------------------------------------
6272 !OCL SERIAL
6273  function fspline ( ndat, kdeg, coef, qknt, myu, xdat, x )
6274 
6275  integer, intent(in) :: ndat
6276  integer, intent(in) :: kdeg
6277  real(dp), intent(in) :: coef( ndat )
6278  real(dp), intent(in) :: qknt( ndat+kdeg )
6279  integer, intent(in) :: myu
6280  real(dp), intent(in) :: xdat( ndat )
6281  real(dp), intent(in) :: x
6282 
6283  real(dp) :: fspline
6284 
6285  !--- local
6286  real(dp) :: sum
6287  integer :: i
6288 
6289  sum = 0.0_dp
6290  do i = 1, ndat
6291  sum = sum + coef( i )*fbspl( ndat, i, kdeg, qknt, myu, xdat, x )
6292  end do
6293 
6294  fspline = sum
6295 
6296  return
6297 
6298  end function fspline
6299  !---------------------------------------------------------------
6300 !OCL SERIAL
6301  subroutine getcoef2 &
6302  ( mdat, ndat, kdeg, ldeg, myu, nyu, & !--- in
6303  xdat, ydat, qknt, rknt, zdat, & !--- in
6304  coef ) !--- out
6305 
6306 ! use scale_tinvss, only: TINVSS
6307 
6308  integer, intent(in) :: mdat ! number of data (x-direction)
6309  integer, intent(in) :: ndat ! number of data (y-direction)
6310  integer, intent(in) :: kdeg ! degree of Spline + 1 (x)
6311  integer, intent(in) :: ldeg ! degree of Spline + 1 (y)
6312  integer, intent(in) :: myu
6313  integer, intent(in) :: nyu
6314  real(dp), intent(in) :: xdat(nspc_mk,mdat) ! data of independent var. (x)
6315  real(dp), intent(in) :: ydat(nspc_mk,ndat) ! data of independent var. (y)
6316  real(dp), intent(in) :: qknt(mdat+kdeg) ! knots of B-Spline (x)
6317  real(dp), intent(in) :: rknt(ndat+ldeg) ! knots of B-Spline (y)
6318  real(dp), intent(in) :: zdat(nspc_mk,nspc_mk,mdat,ndat) ! data of dependent var.
6319 
6320  real(dp), intent(out) :: coef( mdat,ndat ) ! expansion coefficient
6321 
6322  !--- local
6323  real(dp) :: elmx( mdat,mdat ), elmy( ndat,ndat )
6324  integer :: iw1( 2*mdat ), iw2( 2*ndat )
6325  real(dp) :: beta( mdat,ndat ), sum, dt
6326  real(dp), parameter :: eps = 0.0_dp
6327  integer :: i, j, k, l, inder
6328 
6329  do i = 1, mdat
6330  do j = 1, mdat
6331  elmx( i,j ) = fbspl( mdat, j, kdeg, qknt, myu, xdat, xdat(myu,i) )
6332  end do
6333  end do
6334  call tinvss( mdat, elmx, dt, eps, mdat, iw1, inder )
6335 
6336  do l = 1, ndat
6337  do i = 1, mdat
6338  sum = 0.0_dp
6339  do j = 1, mdat
6340  sum = sum + elmx( i,j )*zdat(myu,nyu,j,l)
6341  end do
6342  beta( i,l ) = sum
6343  end do
6344  end do
6345 
6346  do i = 1, ndat
6347  do j = 1, ndat
6348  elmy( i,j ) = fbspl( ndat, j, ldeg, rknt, myu, ydat, ydat(myu,i) )
6349  end do
6350  end do
6351  call tinvss( ndat, elmy, dt, eps, ndat, iw2, inder )
6352 
6353  do k = 1, mdat
6354  do i = 1, ndat
6355  sum = 0.0_dp
6356  do j = 1, ndat
6357  sum = sum + elmy( i,j )*beta( k,j )
6358  end do
6359  coef( k,i ) = sum
6360  end do
6361  end do
6362 
6363  return
6364 
6365  end subroutine getcoef2
6366  !---------------------------------------------------------------
6367 !OCL SERIAL
6368  function fspline2 &
6369  ( mdat, ndat, kdeg, ldeg, & !--- in
6370  coef, qknt, rknt, myu, &
6371  xdat, ydat, & !--- in
6372  x, y ) !--- in
6373 
6374  integer, intent(in) :: mdat
6375  integer, intent(in) :: ndat
6376  integer, intent(in) :: kdeg
6377  integer, intent(in) :: ldeg
6378  real(dp), intent(in) :: coef( mdat,ndat )
6379  real(dp), intent(in) :: qknt( mdat+kdeg )
6380  real(dp), intent(in) :: rknt( ndat+ldeg )
6381  integer, intent(in) :: myu
6382  real(dp), intent(in) :: xdat(nspc_mk,mdat), ydat(nspc_mk,ndat)
6383  real(dp), intent(in) :: x, y
6384 
6385  real(dp) :: fspline2
6386 
6387  !--- local
6388  real(dp) :: sum, add
6389  integer :: i, j
6390 
6391  sum = 0.0_dp
6392  do i = 1, mdat
6393  do j = 1, ndat
6394  add = coef( i,j )*fbspl( mdat, i, kdeg, qknt, myu, xdat, x ) &
6395  *fbspl( ndat, j, ldeg, rknt, myu, ydat, y )
6396  sum = sum + add
6397  end do
6398  end do
6399 
6400  fspline2 = sum
6401 
6402  return
6403 
6404  end function fspline2
6405 
6406 end module scale_atmos_phy_mp_suzuki10
scale_prc::prc_abort
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:350
scale_atmos_hydrometeor::i_hr
integer, parameter, public i_hr
liquid water rain
Definition: scale_atmos_hydrometeor.F90:98
scale_atmos_hydrometeor::cp_water
real(rp), public cp_water
CP for water [J/kg/K].
Definition: scale_atmos_hydrometeor.F90:152
scale_atmos_hydrometeor::i_hs
integer, parameter, public i_hs
snow
Definition: scale_atmos_hydrometeor.F90:100
scale_precision
module PRECISION
Definition: scale_precision.F90:14
scale_const::const_rvap
real(rp), parameter, public const_rvap
specific gas constant (water vapor) [J/kg/K]
Definition: scale_const.F90:68
scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_crg_qtrc2qhyd
subroutine, public atmos_phy_mp_suzuki10_crg_qtrc2qhyd(KA, KS, KE, IA, IS, IE, JA, JS, JE, QTRC0, Qecrg)
get charge density ratio of each category
Definition: scale_atmos_phy_mp_suzuki10.F90:1991
scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_qtrc2qhyd
subroutine, public atmos_phy_mp_suzuki10_qtrc2qhyd(KA, KS, KE, IA, IS, IE, JA, JS, JE, QTRC0, Qe)
Calculate mass ratio of each category.
Definition: scale_atmos_phy_mp_suzuki10.F90:1646
scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_tendency
subroutine, public atmos_phy_mp_suzuki10_tendency(KA, KS, KE, IA, IS, IE, JA, JS, JE, KIJMAX, dt, DENS, PRES, TEMP, QTRC, QDRY, CPtot, CVtot, CCN, RHOQ_t, RHOE_t, CPtot_t, CVtot_t, EVAPORATE, flg_lt, d0_crg, v0_crg, dqcrg, beta_crg, QTRC_crg, QSPLT_in, Sarea, RHOC_t_mp)
Cloud Microphysics.
Definition: scale_atmos_phy_mp_suzuki10.F90:954
scale_const::const_emelt
real(rp), parameter, public const_emelt
Definition: scale_const.F90:79
scale_const::const_eps
real(rp), public const_eps
small number
Definition: scale_const.F90:35
scale_prof::prof_rapstart
subroutine, public prof_rapstart(rapname_base, level, disable_barrier)
Start raptime.
Definition: scale_prof.F90:174
scale_atmos_hydrometeor
module atmosphere / hydrometeor
Definition: scale_atmos_hydrometeor.F90:12
scale_atmos_hydrometeor::i_hh
integer, parameter, public i_hh
hail
Definition: scale_atmos_hydrometeor.F90:102
scale_random
module RANDOM
Definition: scale_random.F90:11
scale_atmos_phy_mp_suzuki10::nbin
integer, public nbin
Definition: scale_atmos_phy_mp_suzuki10.F90:138
scale_io::io_get_available_fid
integer function, public io_get_available_fid()
search & get available file ID
Definition: scale_io.F90:373
scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_terminal_velocity
subroutine, public atmos_phy_mp_suzuki10_terminal_velocity(KA, vterm_o)
get terminal velocity
Definition: scale_atmos_phy_mp_suzuki10.F90:1423
scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_nbnd
integer, public atmos_phy_mp_suzuki10_nbnd
Definition: scale_atmos_phy_mp_suzuki10.F90:60
scale_prc
module PROCESS
Definition: scale_prc.F90:11
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:149
scale_atmos_hydrometeor::i_hi
integer, parameter, public i_hi
ice water cloud
Definition: scale_atmos_hydrometeor.F90:99
scale_io
module STDIO
Definition: scale_io.F90:10
scale_prc::prc_masterrank
integer, parameter, public prc_masterrank
master process in each communicator
Definition: scale_prc.F90:67
scale_const
module CONSTANT
Definition: scale_const.F90:11
scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_nwaters
integer, public atmos_phy_mp_suzuki10_nwaters
Definition: scale_atmos_phy_mp_suzuki10.F90:57
scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_qtrc2nhyd
subroutine, public atmos_phy_mp_suzuki10_qtrc2nhyd(KA, KS, KE, IA, IS, IE, JA, JS, JE, DENS, QTRC0, Ne)
Calculate number concentration of each category.
Definition: scale_atmos_phy_mp_suzuki10.F90:1717
scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_qhyd2qtrc
subroutine, public atmos_phy_mp_suzuki10_qhyd2qtrc(KA, KS, KE, IA, IS, IE, JA, JS, JE, Qe, QTRC, QNUM)
get mass ratio of each category
Definition: scale_atmos_phy_mp_suzuki10.F90:1788
scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_tracer_units
character(len=h_short), dimension(:), allocatable, public atmos_phy_mp_suzuki10_tracer_units
Definition: scale_atmos_phy_mp_suzuki10.F90:65
scale_prof
module profiler
Definition: scale_prof.F90:11
scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_effective_radius
subroutine, public atmos_phy_mp_suzuki10_effective_radius(KA, KS, KE, IA, IS, IE, JA, JS, JE, DENS0, TEMP0, QTRC0, Re)
Calculate Effective Radius.
Definition: scale_atmos_phy_mp_suzuki10.F90:1508
scale_atmos_phy_mp_suzuki10
module Spectran Bin Microphysics
Definition: scale_atmos_phy_mp_suzuki10.F90:23
scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_cloud_fraction
subroutine, public atmos_phy_mp_suzuki10_cloud_fraction(KA, KS, KE, IA, IS, IE, JA, JS, JE, QTRC0, mask_criterion, cldfrac)
Calculate Cloud Fraction.
Definition: scale_atmos_phy_mp_suzuki10.F90:1449
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:88
scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_tracer_descriptions
character(len=h_mid), dimension(:), allocatable, public atmos_phy_mp_suzuki10_tracer_descriptions
Definition: scale_atmos_phy_mp_suzuki10.F90:64
scale_const::const_pi
real(rp), parameter, public const_pi
pi
Definition: scale_const.F90:32
scale_atmos_hydrometeor::i_hc
integer, parameter, public i_hc
liquid water cloud
Definition: scale_atmos_hydrometeor.F90:97
scale_atmos_hydrometeor::i_qv
integer, public i_qv
Definition: scale_atmos_hydrometeor.F90:93
scale_const::const_dwatr
real(rp), parameter, public const_dwatr
density of water [kg/m3]
Definition: scale_const.F90:89
scale_const::const_tem00
real(rp), parameter, public const_tem00
temperature reference (0C) [K]
Definition: scale_const.F90:99
scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_tracer_setup
subroutine, public atmos_phy_mp_suzuki10_tracer_setup
Config.
Definition: scale_atmos_phy_mp_suzuki10.F90:300
scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_ntracers
integer, public atmos_phy_mp_suzuki10_ntracers
Definition: scale_atmos_phy_mp_suzuki10.F90:56
scale_comm_cartesc
module COMMUNICATION
Definition: scale_comm_cartesC.F90:11
scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_tracer_names
character(len=h_short), dimension(:), allocatable, public atmos_phy_mp_suzuki10_tracer_names
Definition: scale_atmos_phy_mp_suzuki10.F90:63
scale_const::const_tmelt
real(rp), parameter, public const_tmelt
Definition: scale_const.F90:80
scale_atmos_phy_mp_suzuki10::nccn
integer, public nccn
Definition: scale_atmos_phy_mp_suzuki10.F90:140
scale_prof::prof_rapend
subroutine, public prof_rapend(rapname_base, level, disable_barrier)
Save raptime.
Definition: scale_prof.F90:246
scale_const::const_dice
real(rp), parameter, public const_dice
density of ice [kg/m3]
Definition: scale_const.F90:90
scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_finalize
subroutine, public atmos_phy_mp_suzuki10_finalize
finalize
Definition: scale_atmos_phy_mp_suzuki10.F90:895
scale_specfunc::sf_gamma
real(rp) function, public sf_gamma(x)
Gamma function.
Definition: scale_specfunc.F90:50
scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_setup
subroutine, public atmos_phy_mp_suzuki10_setup(KA, IA, JA, flg_lt)
Setup.
Definition: scale_atmos_phy_mp_suzuki10.F90:394
scale_atmos_saturation
module atmosphere / saturation
Definition: scale_atmos_saturation.F90:12
scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_nices
integer, public atmos_phy_mp_suzuki10_nices
Definition: scale_atmos_phy_mp_suzuki10.F90:58
scale_io::io_fid_conf
integer, public io_fid_conf
Config file ID.
Definition: scale_io.F90:57
scale_atmos_hydrometeor::n_hyd
integer, parameter, public n_hyd
Definition: scale_atmos_hydrometeor.F90:95
scale_atmos_hydrometeor::cp_vapor
real(rp), public cp_vapor
CP for vapor [J/kg/K].
Definition: scale_atmos_hydrometeor.F90:150
scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_nccn
integer, public atmos_phy_mp_suzuki10_nccn
Definition: scale_atmos_phy_mp_suzuki10.F90:59
scale_atmos_hydrometeor::cv_water
real(rp), public cv_water
CV for water [J/kg/K].
Definition: scale_atmos_hydrometeor.F90:151
scale_prc::prc_ismaster
logical, public prc_ismaster
master process in local communicator?
Definition: scale_prc.F90:92
scale_atmos_hydrometeor::cv_ice
real(rp), public cv_ice
CV for ice [J/kg/K].
Definition: scale_atmos_hydrometeor.F90:153
scale_atmos_hydrometeor::cp_ice
real(rp), public cp_ice
CP for ice [J/kg/K].
Definition: scale_atmos_hydrometeor.F90:154
scale_atmos_hydrometeor::i_hg
integer, parameter, public i_hg
graupel
Definition: scale_atmos_hydrometeor.F90:101