SCALE-RM
Functions/Subroutines | Variables
scale_atmos_phy_ae_kajino13 Module Reference

module atmosphere / physics / aerosol / Kajino13 More...

Functions/Subroutines

subroutine, public atmos_phy_ae_kajino13_tracer_setup (QA_AE)
 Tracer setup. More...
 
subroutine, public atmos_phy_ae_kajino13_setup
 Setup. More...
 
subroutine, public atmos_phy_ae_kajino13_tendency (KA, KS, KE, IA, IS, IE, JA, JS, JE, QA_AE, TEMP, PRES, QDRY, NREG, DENS, QV, QTRC, EMIT, dt, RHOQ_t_AE, CN, CCN)
 Aerosol Microphysics. More...
 
subroutine, public atmos_phy_ae_kajino13_effective_radius (KA, IA, JA, QA_AE, QTRC, RH, Re)
 Calculate Effective Radius. More...
 
subroutine, public atmos_phy_ae_kajino13_mkinit (KA, KS, KE, IA, IS, IE, JA, JS, JE, QA_AE, DENS, TEMP, PRES, QDRY, QV, m0_init, dg_init, sg_init, d_min_inp, d_max_inp, k_min_inp, k_max_inp, n_kap_inp, QTRC, CCN)
 
subroutine, public atmos_phy_ae_kajino13_negative_fixer (KA, KS, KE, IA, IS, IE, JA, JS, JE, QA_AE, QTRC)
 

Variables

character(len=h_short), dimension(:), allocatable, public atmos_phy_ae_kajino13_name
 
character(len=h_mid), dimension(:), allocatable, public atmos_phy_ae_kajino13_desc
 
character(len=h_short), dimension(:), allocatable, public atmos_phy_ae_kajino13_unit
 
real(rp), public atmos_phy_ae_kajino13_h2so4dt = 5.E-6_RP
 
real(rp), public atmos_phy_ae_kajino13_ocgasdt = 8.E-5_RP
 
real(rp), public atmos_phy_ae_kajino13_c_kappa = 0.3_RP
 
logical, public atmos_phy_ae_kajino13_flag_npf = .false.
 
logical, public atmos_phy_ae_kajino13_flag_cond = .true.
 
logical, public atmos_phy_ae_kajino13_flag_coag = .true.
 
logical, public atmos_phy_ae_kajino13_flag_ccn_interactive = .true.
 
logical, public atmos_phy_ae_kajino13_flag_regeneration = .true.
 
real(rp), public atmos_phy_ae_kajino13_dg_reg = 5.E-7_RP
 
real(rp), public atmos_phy_ae_kajino13_sg_reg = 1.6_RP
 
real(rp), public atmos_phy_ae_kajino13_logk_aenucl = -12.4_RP
 

Detailed Description

module atmosphere / physics / aerosol / Kajino13

Description
kajino13 aerosol microphysics scheme
Author
Team SCALE
NAMELIST
  • PARAM_ATMOS_PHY_AE_KAJINO13_TRACER
    nametypedefault valuecomment
    AE_CTG integer 1
    NASIZ integer, dimension(3)
    NAKAP integer, dimension(3)

  • PARAM_ATMOS_PHY_AE_KAJINO13
    nametypedefault valuecomment
    ATMOS_PHY_AE_KAJINO13_H2SO4DT real(RP) 5.E-6_RP h2so4 production rate (Temporal) [ug/m3/s]
    ATMOS_PHY_AE_KAJINO13_OCGASDT real(RP) 8.E-5_RP other condensational bas production rate 16*h2so4dt (see Kajino et al. 2013)

History Output
namedescriptionunitvariable
CGAS_emit Emission ratio of Condensabule gas ug/m3/s EMIT
H2SO4_emit Emission ratio of H2SO4 gas ug/m3/s EMIT
{ctg_name}_mass Total mass mixing ratio of aerosol kg/kg total_aerosol_mass
{ctg_name}_mass_emit Total mass mixing ratio of emitted aerosol kg/kg total_emit_aerosol_mass
{ctg_name}_number Total number mixing ratio of aerosol num/kg total_aerosol_number
{ctg_name}_number_emit Total number mixing ratio of emitted aerosol num/kg total_emit_aerosol_number

Function/Subroutine Documentation

◆ atmos_phy_ae_kajino13_tracer_setup()

subroutine, public scale_atmos_phy_ae_kajino13::atmos_phy_ae_kajino13_tracer_setup ( integer, intent(out)  QA_AE)

Tracer setup.

Definition at line 183 of file scale_atmos_phy_ae_kajino13.F90.

References atmos_phy_ae_kajino13_desc, atmos_phy_ae_kajino13_name, atmos_phy_ae_kajino13_unit, scale_io::io_fid_conf, and scale_prc::prc_abort().

Referenced by mod_atmos_phy_ae_driver::atmos_phy_ae_driver_tracer_setup().

183  use scale_prc, only: &
184  prc_abort
185 
186  integer, intent(out) :: qa_ae
187 
188  integer, allocatable :: aero_idx(:,:,:,:)
189  integer :: ncat_max
190  integer :: nasiz(3), nakap(3)
191  character(len=H_SHORT) :: attribute, catego, aunit
192 
193  namelist / param_atmos_phy_ae_kajino13_tracer / &
194  ae_ctg, &
195  nasiz, &
196  nakap
197 
198  integer :: m, ierr, ik, ic, ia0, is0
199 
200  !---------------------------------------------------------------------------
201 
202  log_newline
203  log_info("ATMOS_PHY_AE_kajino13_tracer_setup",*) 'Setup'
204 
205  ! note: tentatively, aerosol module should be called at all time. we need dummy subprogram.
206  !if ( ATMOS_sw_phy_ae ) then
207 
208  log_warn("ATMOS_PHY_AE_kajino13_tracer_setup",*) '### Kajino13 scheme is still experimental'
209 
210  ncat_max = max( ic_mix, ic_sea, ic_dus )
211 
212  nasiz(:) = 64
213  nakap(:) = 1
214 
215  rewind(io_fid_conf)
216  read(io_fid_conf,nml=param_atmos_phy_ae_kajino13_tracer,iostat=ierr)
217 
218  if( ierr < 0 ) then !--- missing
219  log_info("ATMOS_PHY_AE_kajino13_tracer_setup",*) 'Not found namelist. Default used.'
220  elseif( ierr > 0 ) then !--- fatal error
221  log_error("ATMOS_PHY_AE_kajino13_tracer_setup",*) 'Not appropriate names in namelist PARAM_ATMOS_PHY_AE_KAJINO13_TRACER, Check!'
222  call prc_abort
223  end if
224 
225  log_nml(param_atmos_phy_ae_kajino13_tracer)
226 
227  if( ae_ctg > ncat_max ) then
228  log_error("ATMOS_PHY_AE_kajino13_tracer_setup",*) 'AE_CTG should be smaller than', ncat_max+1, 'stop'
229  call prc_abort
230  endif
231 
232  allocate( nsiz(ae_ctg) )
233  allocate( nkap(ae_ctg) )
234 
235  nkap(1:ae_ctg) = nakap(1:ae_ctg)
236  nsiz(1:ae_ctg) = nasiz(1:ae_ctg)
237 
238  if( maxval( nkap ) /= 1 .OR. minval( nkap ) /= 1 ) then
239  log_error("ATMOS_PHY_AE_kajino13_tracer_setup",*) 'NKAP(:) /= 1 is not supported now, Stop!'
240  call prc_abort
241  end if
242 
243  ! do ia0 = 1, N_ATR
244  do ic = 1, ae_ctg
245  do ik = 1, nkap(ic)
246  do is0 = 1, nsiz(ic)
247  qa_ae = qa_ae + n_atr
248  enddo
249  enddo
250  enddo
251  qa_ae = qa_ae + gas_ctg
252 
253  allocate( atmos_phy_ae_kajino13_name(qa_ae) )
254  allocate( atmos_phy_ae_kajino13_desc(qa_ae) )
255  allocate( atmos_phy_ae_kajino13_unit(qa_ae) )
256 
257  n_siz_max = 0
258  n_kap_max = 0
259  do ic = 1, ae_ctg
260  n_siz_max = max(n_siz_max, nsiz(ic))
261  n_kap_max = max(n_kap_max, nkap(ic))
262  enddo
263 
264  allocate( aero_idx(n_atr,ae_ctg,n_kap_max,n_siz_max) )
265  m = 0
266  ! do ia0 = 1, N_ATR
267  do ic = 1, ae_ctg
268  do ik = 1, nkap(ic)
269  do is0 = 1, nsiz(ic)
270  do ia0 = 1, n_atr
271  m = m+1
272  aero_idx(ia0,ic,ik,is0) = m
273  enddo
274  enddo
275  enddo
276  enddo
277 
278  !-----------------------------------------------------------------------------
279  !
280  !++ calculate each category and aerosol
281  !
282  !-----------------------------------------------------------------------------
283  ic = qa_ae-gas_ctg+ig_h2so4
284  write(atmos_phy_ae_kajino13_unit(ic),'(a)') 'kg/kg'
285  ic = qa_ae-gas_ctg+ig_cgas
286  write(atmos_phy_ae_kajino13_unit(ic),'(a)') 'kg/kg'
287 
288  !do ia0 = 1, N_ATR
289  ! if( ia0 == 1 ) then
290  ! write(attribute,'(a)') "Number"
291  ! write(aunit,'(a)') "num/kg"
292  ! elseif( ia0 == 2 ) then
293  ! write(attribute,'(a)') "Section"
294  ! write(aunit,'(a)') "m2/kg"
295  ! elseif( ia0 == 3 ) then
296  ! write(attribute,'(a)') "Volume"
297  ! write(aunit,'(a)') "m3/kg"
298  ! elseif( ia0 == 4 ) then
299  ! write(attribute,'(a)') "Mass"
300  ! write(aunit,'(a)') "kg/kg"
301  ! elseif( ia0 == 5 ) then
302  ! write(attribute,'(a)') "kpXmass"
303  ! write(aunit,'(a)') "kg/kg"
304  ! endif
305  do ic = 1, ae_ctg !aerosol category
306  do ik = 1, nkap(ic) !kappa bin
307  do is0 = 1, nsiz(ic)
308  do ia0 = 1, n_atr
309  if( ia0 == 1 ) then
310  write(attribute,'(a)') "Number"
311  write(aunit,'(a)') "num/kg"
312  elseif( ia0 == 2 ) then
313  write(attribute,'(a)') "Section"
314  write(aunit,'(a)') "m2/kg"
315  elseif( ia0 == 3 ) then
316  write(attribute,'(a)') "Volume"
317  write(aunit,'(a)') "m3/kg"
318  elseif( ia0 == 4 ) then
319  write(attribute,'(a)') "Mass"
320  write(aunit,'(a)') "kg/kg"
321  elseif( ia0 == 5 ) then
322  write(attribute,'(a)') "kXm"
323  write(aunit,'(a)') "kg/kg"
324  endif
325  if( ic == ic_mix ) then
326  write(catego,'(a)') "Sulf_"
327  elseif( ic == ic_sea ) then
328  write(catego,'(a)') "Salt_"
329  elseif( ic == ic_dus ) then
330  write(catego,'(a)') "Dust_"
331  endif
332  write(atmos_phy_ae_kajino13_unit(aero_idx(ia0,ic,ik,is0)),'(a)') trim(aunit)
333  write(atmos_phy_ae_kajino13_name(aero_idx(ia0,ic,ik,is0)),'(a,a,i0)') trim(catego), trim(attribute), is0
334  write(atmos_phy_ae_kajino13_desc(aero_idx(ia0,ic,ik,is0)),'(a,a,a,i0)') trim(attribute), ' mixing radio of ', trim(catego), is0
335  enddo
336  enddo
337  enddo
338  enddo
339  ic = qa_ae-gas_ctg+ig_h2so4
340  write(atmos_phy_ae_kajino13_name(ic),'(a)') 'H2SO4_Gas'
341  ic = qa_ae-gas_ctg+ig_cgas
342  write(atmos_phy_ae_kajino13_name(ic),'(a)') 'Condensable_GAS'
343 
344  ic = qa_ae-gas_ctg+ig_h2so4
345  write(atmos_phy_ae_kajino13_desc(ic),'(a)') 'Mixing ratio of H2SO4 Gas'
346  ic = qa_ae-gas_ctg+ig_cgas
347  write(atmos_phy_ae_kajino13_desc(ic),'(a)') 'Mixing ratio of Condensable GAS'
348 
349  deallocate(aero_idx)
350 
351  return
integer, public io_fid_conf
Config file ID.
Definition: scale_io.F90:55
module PROCESS
Definition: scale_prc.F90:11
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:338
Here is the call graph for this function:
Here is the caller graph for this function:

◆ atmos_phy_ae_kajino13_setup()

subroutine, public scale_atmos_phy_ae_kajino13::atmos_phy_ae_kajino13_setup ( )

Setup.

Definition at line 357 of file scale_atmos_phy_ae_kajino13.F90.

References atmos_phy_ae_kajino13_c_kappa, atmos_phy_ae_kajino13_dg_reg, atmos_phy_ae_kajino13_flag_ccn_interactive, atmos_phy_ae_kajino13_flag_coag, atmos_phy_ae_kajino13_flag_cond, atmos_phy_ae_kajino13_flag_npf, atmos_phy_ae_kajino13_flag_regeneration, atmos_phy_ae_kajino13_h2so4dt, atmos_phy_ae_kajino13_logk_aenucl, atmos_phy_ae_kajino13_ocgasdt, atmos_phy_ae_kajino13_sg_reg, float(), scale_io::io_fid_conf, and scale_prc::prc_abort().

Referenced by mod_atmos_phy_ae_driver::atmos_phy_ae_driver_setup(), and atmos_phy_ae_kajino13_mkinit().

357  use scale_prc, only: &
358  prc_abort
359  implicit none
360 
361  real(RP), allocatable :: d_min_inp(:)
362  real(RP), allocatable :: d_max_inp(:)
363  real(RP), allocatable :: k_min_inp(:)
364  real(RP), allocatable :: k_max_inp(:)
365  integer , allocatable :: n_kap_inp(:)
366 
367  real(RP), parameter :: d_min_def = 1.e-9_rp ! default lower bound of 1st size bin
368  real(RP), parameter :: d_max_def = 1.e-5_rp ! upper bound of last size bin
369  integer , parameter :: n_kap_def = 1 ! number of kappa bins
370  real(RP), parameter :: k_min_def = 0.e0_rp ! lower bound of 1st kappa bin
371  real(RP), parameter :: k_max_def = 1.e0_rp ! upper bound of last kappa bin
372 
373  namelist / param_atmos_phy_ae_kajino13 / &
374  atmos_phy_ae_kajino13_h2so4dt, &
375  atmos_phy_ae_kajino13_ocgasdt, &
376 ! ATMOS_PHY_AE_KAJINO13_c_ratio, &
377  atmos_phy_ae_kajino13_c_kappa, &
378 ! d_min_inp, &
379 ! d_max_inp, &
380 ! k_min_inp, &
381 ! k_max_inp, &
382 ! n_kap_inp, &
383  atmos_phy_ae_kajino13_flag_npf, &
384  atmos_phy_ae_kajino13_flag_cond, &
385  atmos_phy_ae_kajino13_flag_coag, &
386  atmos_phy_ae_kajino13_flag_ccn_interactive, &
387  atmos_phy_ae_kajino13_flag_regeneration, &
388  atmos_phy_ae_kajino13_dg_reg, &
389  atmos_phy_ae_kajino13_sg_reg, &
390  atmos_phy_ae_kajino13_logk_aenucl, &
391  atmos_phy_ae_kajino13_nbins_out
392 
393  integer :: it, ierr
394  !---------------------------------------------------------------------------
395 
396  log_newline
397  log_info("ATMOS_PHY_AE_kajino13_setup",*) 'Setup'
398  log_info("ATMOS_PHY_AE_kajino13_setup",*) 'Kajino(2013) scheme'
399 
400  !--- setup parameter
401  pi6 = pi / 6._rp ! pi/6
402  sixpi = 6._rp / pi ! 6/pi
403  forpi = 4._rp / pi ! 4/pi
404 
405 
406  n_ctg = ae_ctg
407  allocate( rnum_out(atmos_phy_ae_kajino13_nbins_out) )
408  allocate( n_siz(n_ctg) )
409  allocate( d_min(n_ctg) )
410  allocate( d_max(n_ctg) )
411  allocate( n_kap(n_ctg) )
412  allocate( k_min(n_ctg) )
413  allocate( k_max(n_ctg) )
414  allocate( d_min_inp(n_ctg) )
415  allocate( d_max_inp(n_ctg) )
416  allocate( n_kap_inp(n_ctg) )
417  allocate( k_min_inp(n_ctg) )
418  allocate( k_max_inp(n_ctg) )
419  allocate( ctg_name(n_ctg) )
420 
421  n_siz(1:n_ctg) = nsiz(1:n_ctg) ! number of size bins
422  d_min(1:n_ctg) = d_min_def ! lower bound of 1st size bin
423  d_max(1:n_ctg) = d_max_def ! upper bound of last size bin
424  n_kap(1:n_ctg) = n_kap_def ! number of kappa bins
425  k_min(1:n_ctg) = k_min_def ! lower bound of 1st kappa bin
426  k_max(1:n_ctg) = k_max_def ! upper bound of last kappa bin
427 
428  do it = 1, n_ctg
429  if( n_ctg == 1 ) then
430  write(ctg_name(it),'(a)') "Sulfate"
431  elseif( n_ctg == 2 ) then
432  write(ctg_name(it),'(a)') "Seasalt"
433  elseif( n_ctg == 3 ) then
434  write(ctg_name(it),'(a)') "Dust"
435  endif
436  enddo
437 
438  !--- read namelist
439  rewind(io_fid_conf)
440  read(io_fid_conf,nml=param_atmos_phy_ae_kajino13,iostat=ierr)
441  if( ierr < 0 ) then !--- missing
442  log_info("ATMOS_PHY_AE_kajino13_setup",*) 'Not found namelist. Default used.'
443  elseif( ierr > 0 ) then !--- fatal error
444  log_error("ATMOS_PHY_AE_kajino13_setup",*) 'Not appropriate names in namelist PARAM_ATMOS_PHY_AE_KAJINO13. Check!'
445  call prc_abort
446  endif
447  log_nml(param_atmos_phy_ae_kajino13)
448 
449  !--- now only the default setting is supported
450 ! n_siz(1:n_ctg) = NSIZ(1:n_ctg) ! number of size bins
451 ! d_min(1:n_ctg) = d_min_inp(1:n_ctg) ! lower bound of 1st size bin
452 ! d_max(1:n_ctg) = d_max_inp(1:n_ctg) ! upper bound of last size bin
453 ! n_kap(1:n_ctg) = n_kap_inp(1:n_ctg) ! number of kappa bins
454 ! k_min(1:n_ctg) = k_min_inp(1:n_ctg) ! lower bound of 1st kappa bin
455 ! k_max(1:n_ctg) = k_max_inp(1:n_ctg) ! upper bound of last kappa bin
456 
457  !--- diagnose parameters (n_trans, n_siz_max, n_kap_max)
458  n_trans = 0
459  n_siz_max = 0
460  n_kap_max = 0
461  do ic = 1, n_ctg
462  n_trans = n_trans + n_siz(ic) * n_kap(ic) * n_atr
463  n_siz_max = max(n_siz_max, n_siz(ic))
464  n_kap_max = max(n_kap_max, n_kap(ic))
465  enddo
466 
467  !--- bin settings
468  allocate(d_lw(n_siz_max,n_ctg))
469  allocate(d_ct(n_siz_max,n_ctg))
470  allocate(d_up(n_siz_max,n_ctg))
471  allocate(k_lw(n_kap_max,n_ctg))
472  allocate(k_ct(n_kap_max,n_ctg))
473  allocate(k_up(n_kap_max,n_ctg))
474  d_lw(:,:) = 0.0_rp
475  d_ct(:,:) = 0.0_rp
476  d_up(:,:) = 0.0_rp
477  k_lw(:,:) = 0.0_rp
478  k_ct(:,:) = 0.0_rp
479  k_up(:,:) = 0.0_rp
480 
481  do ic = 1, n_ctg
482  dlogd = (log(d_max(ic)) - log(d_min(ic)))/float(n_siz(ic))
483  do is0 = 1, n_siz(ic) !size bin
484  d_lw(is0,ic) = exp(log(d_min(ic))+dlogd* float(is0-1) )
485  d_ct(is0,ic) = exp(log(d_min(ic))+dlogd*(float(is0)-0.5_rp))
486  d_up(is0,ic) = exp(log(d_min(ic))+dlogd* float(is0) )
487  enddo !is (1:n_siz(ic))
488 
489  dk = (k_max(ic) - k_min(ic))/float(n_kap(ic))
490  do ik = 1, n_kap(ic) !size bin
491  k_lw(ik,ic) = k_min(ic) + dk * float(ik-1)
492  k_ct(ik,ic) = k_min(ic) + dk *(float(ik)-0.5_rp)
493  k_up(ik,ic) = k_min(ic) + dk * float(ik)
494  enddo !ik (1:n_kap(ic))
495 
496  enddo !ic (1:n_ctg)
497 
498 !find size bin of regenerated aerosols
499  do is0 = 1, n_siz(ic_mix)
500  if ( atmos_phy_ae_kajino13_dg_reg >= d_lw(is0,ic_mix) .AND. &
501  atmos_phy_ae_kajino13_dg_reg < d_up(is0,ic_mix) ) then
502  is0_reg = is0
503  endif !d_lw < ATMOS_PHY_AE_KAJINO13_dg_reg < d_up
504  enddo
505 
506  !--- coagulation rule
507  ! [ NOTE: current version has one category and single
508  ! hygroscopicity bins and thus does not consider inter-category
509  ! nor inter-hygroscopicity-section coagulation ]
510  mcomb = 0
511  do ic = 1, n_ctg !=1
512  do ik = 1, n_kap(ic) !=1
513  do is2 = 1, n_siz(ic)
514  do is1 = 1, n_siz(ic)
515  if ( d_ct(is2,ic) >= d_ct(is1,ic) ) then
516  mcomb = mcomb + 1
517  endif !d_ct(is2) >= d_ct(is1)
518  enddo !is1 (1:n_siz(ic) )
519  enddo !is2 (1:n_siz(ic) )
520  enddo !ik(1:n_kap(ic))
521  enddo !ic(1:n_ctg)
522 
523  allocate(is_i(mcomb))
524  allocate(is_j(mcomb))
525  allocate(is_k(mcomb))
526  allocate(ik_i(mcomb))
527  allocate(ik_j(mcomb))
528  allocate(ik_k(mcomb))
529  allocate(ic_i(mcomb))
530  allocate(ic_j(mcomb))
531  allocate(ic_k(mcomb))
532 
533  mc = 0
534  do ic = 1, n_ctg !=1
535  do ik = 1, n_kap(ic) !=1
536  do is2 = 1, n_siz(ic)
537  do is1 = 1, n_siz(ic)
538  if ( d_ct(is2,ic) >= d_ct(is1,ic) ) then
539  mc = mc + 1
540  is_i(mc) = is1
541  ik_i(mc) = ik
542  ic_i(mc) = ic
543  is_j(mc) = is2
544  ik_j(mc) = ik
545  ic_j(mc) = ic
546  is_k(mc) = is2
547  ik_k(mc) = ik
548  ic_k(mc) = ic
549  endif !d_ct(is2) >= d_ct(is1)
550  enddo !is1 (1:n_siz(ic) )
551  enddo !is2 (1:n_siz(ic) )
552  enddo !ik(1:n_kap(ic))
553  enddo !ic(1:n_ctg)
554 
555  !--- gas concentration
556 ! conc_h2so4 = 0.0_RP
557 ! conc_cgas = 0.0_RP
558 
559  allocate( it_procs2trans(n_atr,n_siz_max,n_kap_max,n_ctg) )
560  allocate( ia_trans2procs(n_trans) )
561  allocate( is_trans2procs(n_trans) )
562  allocate( ik_trans2procs(n_trans) )
563  allocate( ic_trans2procs(n_trans) )
564 
565  it_procs2trans(:,:,:,:)= -999
566  ia_trans2procs(:) = 0
567  is_trans2procs(:) = 0
568  ik_trans2procs(:) = 0
569  ic_trans2procs(:) = 0
570 
571  !--- get pointer for trans2procs, procs2trans
572  it = 0
573  do ic = 1, n_ctg !aerosol category
574  do ik = 1, n_kap(ic) !kappa bin
575  do is0 = 1, n_siz(ic) !size bin
576  do ia0 = 1, n_atr !attributes
577  it = it + 1
578  it_procs2trans(ia0,is0,ik,ic)= it
579  ia_trans2procs(it) = ia0
580  is_trans2procs(it) = is0
581  ik_trans2procs(it) = ik
582  ic_trans2procs(it) = ic
583  enddo !ia (1:n_atr_prog )
584  enddo !is (1:n_siz(ic) )
585  enddo !ik (1:n_kap(ic) )
586  enddo !ic (1:n_ctg )
587 
588  allocate( aerosol_procs(n_atr,n_siz_max,n_kap_max,n_ctg) )
589  allocate( aerosol_activ(n_atr,n_siz_max,n_kap_max,n_ctg) )
590  allocate( emis_procs(n_atr,n_siz_max,n_kap_max,n_ctg) )
591  allocate( emis_gas(gas_ctg) )
592  aerosol_procs(:,:,:,:) = 0.0_rp
593  aerosol_activ(:,:,:,:) = 0.0_rp
594  emis_procs(:,:,:,:) = 0.0_rp
595  emis_gas(:) = 0.0_rp
596 
597  return
integer, public io_fid_conf
Config file ID.
Definition: scale_io.F90:55
module PROCESS
Definition: scale_prc.F90:11
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:338
typedef float(real32_t)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ atmos_phy_ae_kajino13_tendency()

subroutine, public scale_atmos_phy_ae_kajino13::atmos_phy_ae_kajino13_tendency ( integer, intent(in)  KA,
integer, intent(in)  KS,
integer, intent(in)  KE,
integer, intent(in)  IA,
integer, intent(in)  IS,
integer, intent(in)  IE,
integer, intent(in)  JA,
integer, intent(in)  JS,
integer, intent(in)  JE,
integer, intent(in)  QA_AE,
real(rp), dimension(ka,ia,ja), intent(in)  TEMP,
real(rp), dimension(ka,ia,ja), intent(in)  PRES,
real(rp), dimension(ka,ia,ja), intent(in)  QDRY,
real(rp), dimension(ka,ia,ja), intent(in)  NREG,
real(rp), dimension(ka,ia,ja), intent(in)  DENS,
real(rp), dimension (ka,ia,ja), intent(in)  QV,
real(rp), dimension(ka,ia,ja,qa_ae), intent(in)  QTRC,
real(rp), dimension(ka,ia,ja,qa_ae), intent(in)  EMIT,
real(dp), intent(in)  dt,
real(rp), dimension(ka,ia,ja,qa_ae), intent(out)  RHOQ_t_AE,
real(rp), dimension(ka,ia,ja), intent(out)  CN,
real(rp), dimension(ka,ia,ja), intent(out)  CCN 
)

Aerosol Microphysics.

Definition at line 619 of file scale_atmos_phy_ae_kajino13.F90.

References atmos_phy_ae_kajino13_dg_reg, atmos_phy_ae_kajino13_flag_ccn_interactive, atmos_phy_ae_kajino13_flag_coag, atmos_phy_ae_kajino13_flag_cond, atmos_phy_ae_kajino13_flag_npf, atmos_phy_ae_kajino13_flag_regeneration, and atmos_phy_ae_kajino13_sg_reg.

Referenced by mod_atmos_phy_ae_driver::atmos_phy_ae_driver_calc_tendency().

619  use scale_atmos_saturation, only: &
620  pres2qsat_liq => atmos_saturation_pres2qsat_liq
621  use scale_file_history, only: &
622  file_history_in
623  implicit none
624  integer, intent(in) :: ka, ks, ke
625  integer, intent(in) :: ia, is, ie
626  integer, intent(in) :: ja, js, je
627  integer, intent(in) :: qa_ae
628  real(RP), intent(in) :: temp(ka,ia,ja)
629  real(RP), intent(in) :: pres(ka,ia,ja)
630  real(RP), intent(in) :: qdry(ka,ia,ja)
631  real(RP), intent(in) :: nreg(ka,ia,ja)
632  real(RP), intent(in) :: dens(ka,ia,ja)
633  real(RP), intent(in) :: qv (ka,ia,ja)
634  real(RP), intent(in) :: qtrc(ka,ia,ja,qa_ae)
635  real(RP), intent(in) :: emit(ka,ia,ja,qa_ae)
636  real(DP), intent(in) :: dt
637 
638  real(RP), intent(out) :: rhoq_t_ae(ka,ia,ja,qa_ae)
639  real(RP), intent(out) :: cn(ka,ia,ja)
640  real(RP), intent(out) :: ccn(ka,ia,ja)
641 
642  real(RP) :: qtrc0(ka,ia,ja,qa_ae)
643  real(RP) :: qtrc1(ka,ia,ja,qa_ae)
644 
645  !--- local
646  real(RP) :: ssliq_ae(ka,ia,ja)
647  real(RP) :: rrhog(ka,ia,ja)
648 ! real(RP) :: t_ccn, t_cn
649  real(RP) :: qsat_tmp
650  real(RP) :: m0_reg, m2_reg, m3_reg !regenerated aerosols [m^k/m3]
651  real(RP) :: ms_reg !regenerated aerosol mass [ug/m3]
652  real(RP) :: reg_factor_m2,reg_factor_m3 !to save cpu time for moment conversion
653  !--- aerosol variables
654  real(RP) :: total_aerosol_mass(ka,ia,ja,n_ctg)
655  real(RP) :: total_aerosol_number(ka,ia,ja,n_ctg)
656  real(RP) :: total_emit_aerosol_mass(ka,ia,ja,n_ctg)
657  real(RP) :: total_emit_aerosol_number(ka,ia,ja,n_ctg)
658  !--- gas
659 ! real(RP),allocatable :: conc_h2so4(:,:,:,:) !concentration [ug/m3]
660 ! real(RP) :: conc_gas(KA,IA,JA,GAS_CTG) !concentration [ug/m3]
661  real(RP) :: conc_gas(gas_ctg) !concentration [ug/m3]
662  integer :: i, j, k, iq, it
663 
664  log_progress(*) 'atmosphere / physics / aerosol / kajino13'
665 
666  aerosol_procs(:,:,:,:) = 0.0_rp
667  aerosol_activ(:,:,:,:) = 0.0_rp
668  emis_procs(:,:,:,:) = 0.0_rp
669  emis_gas(:) = 0.0_rp
670 
671  do iq = 1, qa_ae
672  do j = js, je
673  do i = is, ie
674  do k = ks, ke
675  qtrc0(k,i,j,iq) = qtrc(k,i,j,iq) ! save
676  enddo
677  enddo
678  enddo
679  enddo
680 
681  reg_factor_m2 = atmos_phy_ae_kajino13_dg_reg**2._rp * exp( 2.0_rp *(log(atmos_phy_ae_kajino13_sg_reg)**2._rp)) !m0_reg to m2_reg
682  reg_factor_m3 = atmos_phy_ae_kajino13_dg_reg**3._rp * exp( 4.5_rp *(log(atmos_phy_ae_kajino13_sg_reg)**2._rp)) !m0_reg to m3_reg
683 
684  !--- convert SCALE variable to zerochem variable
685 
686  do j = js, je
687  do i = is, ie
688 
689  do k = ks, ke
690  rrhog(k,i,j) = 1.0_rp / dens(k,i,j)
691  enddo
692  do k = ks, ke
693  !--- calculate super saturation of water
694  call pres2qsat_liq( temp(k,i,j), pres(k,i,j), qdry(k,i,j), & ! [IN]
695  qsat_tmp ) ! [OUT]
696  ssliq_ae(k,i,j) = qv(k,i,j)/qsat_tmp - 1.0_rp
697  enddo
698 
699  enddo
700  enddo
701 
702  ! tiny number, tiny mass
703  do j = js, je
704  do i = is, ie
705  do k = ks, ke
706  do ic = 1, n_ctg !aerosol category
707  do ik = 1, n_kap(ic) !kappa bin
708  do is0 = 1, n_siz(ic) !size bin
709  if (qtrc0(k,i,j,it_procs2trans(ia_m0,is0,ik,ic))*dens(k,i,j) < cleannumber) then
710  do ia0 = 1, n_atr
711  qtrc0(k,i,j,it_procs2trans(ia0,is0,ik,ic)) = 0._rp !to save cpu time and avoid underflow
712  enddo !ia0 (1:n_atr )
713  endif
714  enddo !is (1:n_siz(ic) )
715  enddo !ik (1:n_kap(ic) )
716  enddo !ic (1:n_ctg )
717  enddo
718  enddo
719  enddo
720 
721  do iq = 1, qa_ae
722  do j = js, je
723  do i = is, ie
724  do k = ks, ke
725  qtrc1(k,i,j,iq) = qtrc0(k,i,j,iq) ! save
726  enddo
727  enddo
728  enddo
729  enddo
730 
731  !---- Calculate aerosol processs
732  cn(:,:,:) = 0.0_rp
733  ccn(:,:,:) = 0.0_rp
734  do k = ks, ke
735  do j = js, je
736  do i = is, ie
737 
738  !--- aerosol_trans at initial time
739  ! [xx/kg] -> [xx/m3]
740  do it = 1, n_trans
741  aerosol_procs(ia_trans2procs(it), &
742  is_trans2procs(it), &
743  ik_trans2procs(it), &
744  ic_trans2procs(it)) = qtrc1(k,i,j,it)*dens(k,i,j)
745  emis_procs(ia_trans2procs(it), &
746  is_trans2procs(it), &
747  ik_trans2procs(it), &
748  ic_trans2procs(it)) = emit(k,i,j,it)*dens(k,i,j)
749  enddo !it(1:n_trans)
750  ! mixing ratio [kg/kg] -> concentration [ug/m3]
751  conc_gas(1:gas_ctg) &
752  = qtrc1(k,i,j,qa_ae-gas_ctg+1:qa_ae-gas_ctg+gas_ctg)*dens(k,i,j)*1.e+9_rp
753 
754  emis_gas(1:gas_ctg) = emit(k,i,j,qa_ae-gas_ctg+ig_h2so4:qa_ae-gas_ctg+ig_cgas)
755 
756  call aerosol_zerochem( &
757  dt, & !--- in
758  temp(k,i,j), & !--- in
759  pres(k,i,j), & !--- in
760  ssliq_ae(k,i,j), & !--- in
761  atmos_phy_ae_kajino13_flag_npf, & !--- in
762  atmos_phy_ae_kajino13_flag_cond, & !--- in
763  atmos_phy_ae_kajino13_flag_coag, & !--- in
764  aerosol_procs, & !--- inout
765  conc_gas, & !--- inout
766  emis_procs, & !--- out
767  emis_gas, & !--- out
768  aerosol_activ ) !--- out
769 
770 ! aerosol loss due to activation to cloud droplets
771  if (atmos_phy_ae_kajino13_flag_ccn_interactive) then
772  do is0 = 1, n_siz(ic_mix)
773  do ia0 = 1, n_atr !attributes
774  aerosol_activ(ia0,is0,ik_out,ic_mix) = min(max(0._rp, aerosol_activ(ia0,is0,ik_out,ic_mix)), &
775  aerosol_procs(ia0,is0,ik_out,ic_mix) )
776 
777  aerosol_procs(ia0,is0,ik_out,ic_mix) = &
778  aerosol_procs(ia0,is0,ik_out,ic_mix) - aerosol_activ(ia0,is0,ik_out,ic_mix)
779  enddo
780  enddo
781  endif !flag_ccn_interactive
782 
783 ! aerosol regeneration due to evaporation of cloud droplets
784 ! using prescribed size parameters and to internal mixture category (ic_mix)
785  if (atmos_phy_ae_kajino13_flag_regeneration) then
786  m0_reg = nreg(k,i,j) !#/m3
787 ! m2_reg = m0_reg * ATMOS_PHY_AE_KAJINO13_dg_reg**2._RP * exp( 2.0_RP *(log(ATMOS_PHY_AE_KAJINO13_sg_reg)**2._RP)) !m2/m3
788 ! m3_reg = m0_reg * ATMOS_PHY_AE_KAJINO13_dg_reg**3._RP * exp( 4.5_RP *(log(ATMOS_PHY_AE_KAJINO13_sg_reg)**2._RP)) !m3/m3
789  m2_reg = m0_reg * reg_factor_m2 !m2/m3
790  m3_reg = m0_reg * reg_factor_m3 !m3/m3
791  ms_reg = m3_reg * pi6 * conv_vl_ms !ug/m3
792  aerosol_procs(ia_m0,is0_reg,ik_out,ic_mix) = &
793  aerosol_procs(ia_m0,is0_reg,ik_out,ic_mix) + m0_reg !#/m3
794  aerosol_procs(ia_m2,is0_reg,ik_out,ic_mix) = &
795  aerosol_procs(ia_m2,is0_reg,ik_out,ic_mix) + m2_reg !m2/m3
796  aerosol_procs(ia_m3,is0_reg,ik_out,ic_mix) = &
797  aerosol_procs(ia_m3,is0_reg,ik_out,ic_mix) + m3_reg !m3/m3
798  aerosol_procs(ia_ms,is0_reg,ik_out,ic_mix) = &
799  aerosol_procs(ia_ms,is0_reg,ik_out,ic_mix) + ms_reg !ug/m3
800 ! additional attirbute to be added (ia_kp)
801  endif !flag_regeneration
802 
803 ! diagnosed variables
804  do is0 = 1, n_siz(ic_mix)
805  ccn(k,i,j) = ccn(k,i,j) + aerosol_activ(ia_m0,is0,ik_out,ic_mix)
806  cn(k,i,j) = cn(k,i,j) + aerosol_procs(ia_m0,is0,ik_out,ic_mix)
807  enddo
808 
809 ! call trans_ccn(aerosol_procs, aerosol_activ, t_ccn, t_cn, &
810 ! n_ctg, n_kap_max, n_siz_max, N_ATR, &
811 ! ic_mix, ia_m0, ia_m2, ia_m3, ik_out, n_siz, &
812 ! rnum_out, ATMOS_PHY_AE_KAJINO13_nbins_out)
813 
814 ! CN(k,i,j) = t_cn
815 ! CCN(k,i,j) = t_ccn
816 
817  ! [xx/m3] -> [xx/kg]
818  do ic = 1, n_ctg !category
819  do ik = 1, n_kap(ic) !kappa bin
820  do is0 = 1, n_siz(ic) !size bin
821  do ia0 = 1, n_atr !attributes
822  qtrc1(k,i,j,it_procs2trans(ia0,is0,ik,ic)) = aerosol_procs(ia0,is0,ik,ic) / dens(k,i,j)
823  enddo !ia (1:N_ATR )
824  enddo !is (1:n_siz(ic) )
825  enddo !ik (1:n_kap(ic) )
826  enddo !ic (1:n_ctg )
827  ! [ug/m3] -> mixing ratio [kg/kg]
828  qtrc1(k,i,j,qa_ae-gas_ctg+1:qa_ae-gas_ctg+gas_ctg) = conc_gas(1:gas_ctg) / dens(k,i,j)*1.e-9_rp
829 
830  ! tiny number, tiny mass
831  do ic = 1, n_ctg !aerosol category
832  do ik = 1, n_kap(ic) !kappa bin
833  do is0 = 1, n_siz(ic) !size bin
834  if (qtrc1(k,i,j,it_procs2trans(ia_m0,is0,ik,ic))*dens(k,i,j) < cleannumber) then
835  do ia0 = 1, n_atr
836  qtrc1(k,i,j,it_procs2trans(ia0,is0,ik,ic)) = 0._rp !to save cpu time and avoid underflow
837  enddo !ia0 (1:n_atr )
838  endif
839  enddo !is (1:n_siz(ic) )
840  enddo !ik (1:n_kap(ic) )
841  enddo !ic (1:n_ctg )
842 
843  ! for history
844  total_aerosol_mass(k,i,j,:) = 0.0_rp
845  total_aerosol_number(k,i,j,:) = 0.0_rp
846  total_emit_aerosol_mass(k,i,j,:) = 0.0_rp
847  total_emit_aerosol_number(k,i,j,:) = 0.0_rp
848  do ic = 1, n_ctg
849  do ik = 1, n_kap(ic)
850  do is0 = 1, n_siz(ic)
851  total_aerosol_mass(k,i,j,ic) = total_aerosol_mass(k,i,j,ic) &
852  + qtrc1(k,i,j,it_procs2trans(ia_ms,is0,ik,ic))
853  total_aerosol_number(k,i,j,ic) = total_aerosol_number(k,i,j,ic) &
854  + qtrc1(k,i,j,it_procs2trans(ia_m0,is0,ik,ic))
855  total_emit_aerosol_mass(k,i,j,ic) = total_emit_aerosol_mass(k,i,j,ic) &
856  + emit(k,i,j,it_procs2trans(ia_ms,is0,ik,ic))
857  total_emit_aerosol_number(k,i,j,ic) = total_emit_aerosol_number(k,i,j,ic) &
858  + emit(k,i,j,it_procs2trans(ia_m0,is0,ik,ic))
859  enddo
860  enddo
861  enddo
862 
863  enddo
864  enddo
865  enddo
866 
867  do ic = 1, n_ctg
868  call file_history_in( total_aerosol_mass(:,:,:,ic), trim(ctg_name(ic))//'_mass', 'Total mass mixing ratio of aerosol', 'kg/kg' )
869  call file_history_in( total_aerosol_number(:,:,:,ic), trim(ctg_name(ic))//'_number', 'Total number mixing ratio of aerosol', 'num/kg' )
870  call file_history_in( total_emit_aerosol_mass(:,:,:,ic), trim(ctg_name(ic))//'_mass_emit', 'Total mass mixing ratio of emitted aerosol', 'kg/kg' )
871  call file_history_in( total_emit_aerosol_number(:,:,:,ic), trim(ctg_name(ic))//'_number_emit', 'Total number mixing ratio of emitted aerosol', 'num/kg' )
872  enddo
873 
874  call file_history_in( emit(:,:,:,qa_ae-gas_ctg+ig_h2so4), 'H2SO4_emit', 'Emission ratio of H2SO4 gas', 'ug/m3/s' )
875  call file_history_in( emit(:,:,:,qa_ae-gas_ctg+ig_cgas), 'CGAS_emit', 'Emission ratio of Condensabule gas', 'ug/m3/s' )
876 
877 
878  do iq = 1, qa_ae
879  do j = js, je
880  do i = is, ie
881  do k = ks, ke
882  rhoq_t_ae(k,i,j,iq) = ( qtrc1(k,i,j,iq) - qtrc0(k,i,j,iq) ) * dens(k,i,j) / dt
883  enddo
884  enddo
885  enddo
886  enddo
887 
888  return
module atmosphere / saturation
integer, public ia
of whole cells: x, local, with HALO
integer, public ja
of whole cells: y, local, with HALO
integer, public is
start point of inner domain: x, local
integer, public ie
end point of inner domain: x, local
integer, public ke
end point of inner domain: z, local
integer, public je
end point of inner domain: y, local
integer, public ks
start point of inner domain: z, local
integer, public js
start point of inner domain: y, local
integer, public ka
of whole cells: z, local, with HALO
module file_history
Here is the caller graph for this function:

◆ atmos_phy_ae_kajino13_effective_radius()

subroutine, public scale_atmos_phy_ae_kajino13::atmos_phy_ae_kajino13_effective_radius ( integer, intent(in)  KA,
integer, intent(in)  IA,
integer, intent(in)  JA,
integer, intent(in)  QA_AE,
real(rp), dimension(ka,ia,ja,qa_ae), intent(in)  QTRC,
real(rp), dimension (ka,ia,ja), intent(in)  RH,
real(rp), dimension (ka,ia,ja,n_ae), intent(out)  Re 
)

Calculate Effective Radius.

Definition at line 897 of file scale_atmos_phy_ae_kajino13.F90.

References scale_atmos_aerosol::n_ae.

Referenced by mod_atmos_phy_ae_vars::atmos_phy_ae_vars_get_diagnostic().

897  use scale_atmos_aerosol, only: &
898  n_ae
899  implicit none
900  integer, intent(in) :: ka, ia, ja, qa_ae
901 
902  real(RP), intent(in) :: qtrc(ka,ia,ja,qa_ae) ! tracer mass concentration [kg/kg]
903  real(RP), intent(in) :: rh (ka,ia,ja) ! relative humidity (0-1)
904 
905  real(RP), intent(out) :: re (ka,ia,ja,n_ae) ! effective radius
906  !---------------------------------------------------------------------------
907 
908  re(:,:,:,:) = 0.0_rp
909 
910 ! Re(:,:,:,I_ae_seasalt) = 2.E-4_RP
911 ! Re(:,:,:,I_ae_dust ) = 4.E-6_RP
912 ! Re(:,:,:,I_ae_bc ) = 4.E-8_RP
913 ! Re(:,:,:,I_ae_oc ) = RH(:,:,:)
914 ! Re(:,:,:,I_ae_sulfate) = RH(:,:,:)
915 
916  return
integer, parameter, public n_ae
integer, public ia
of whole cells: x, local, with HALO
integer, public ja
of whole cells: y, local, with HALO
module atmosphere / aerosol
integer, public ka
of whole cells: z, local, with HALO
Here is the caller graph for this function:

◆ atmos_phy_ae_kajino13_mkinit()

subroutine, public scale_atmos_phy_ae_kajino13::atmos_phy_ae_kajino13_mkinit ( integer, intent(in)  KA,
integer, intent(in)  KS,
integer, intent(in)  KE,
integer, intent(in)  IA,
integer, intent(in)  IS,
integer, intent(in)  IE,
integer, intent(in)  JA,
integer, intent(in)  JS,
integer, intent(in)  JE,
integer, intent(in)  QA_AE,
real(rp), dimension(ka,ia,ja), intent(in)  DENS,
real(rp), dimension(ka,ia,ja), intent(in)  TEMP,
real(rp), dimension(ka,ia,ja), intent(in)  PRES,
real(rp), dimension(ka,ia,ja), intent(in)  QDRY,
real(rp), dimension (ka,ia,ja), intent(in)  QV,
real(rp), intent(in)  m0_init,
real(rp), intent(in)  dg_init,
real(rp), intent(in)  sg_init,
real(rp), dimension(3), intent(in)  d_min_inp,
real(rp), dimension(3), intent(in)  d_max_inp,
real(rp), dimension(3), intent(in)  k_min_inp,
real(rp), dimension(3), intent(in)  k_max_inp,
integer, dimension(3), intent(in)  n_kap_inp,
real(rp), dimension(ka,ia,ja,qa_ae), intent(out)  QTRC,
real(rp), dimension (ka,ia,ja), intent(out)  CCN 
)

Definition at line 940 of file scale_atmos_phy_ae_kajino13.F90.

References atmos_phy_ae_kajino13_setup(), and scale_const::const_pi.

Referenced by mod_mkinit::rect_setup().

940  use scale_const, only: &
941  pi => const_pi
942  use scale_atmos_saturation, only: &
943  saturation_pres2qsat_liq => atmos_saturation_pres2qsat_liq
944  implicit none
945 
946  integer, intent(in) :: ka, ks, ke
947  integer, intent(in) :: ia, is, ie
948  integer, intent(in) :: ja, js, je
949  integer, intent(in) :: qa_ae
950 
951  real(RP), intent(in) :: dens(ka,ia,ja)
952  real(RP), intent(in) :: temp(ka,ia,ja)
953  real(RP), intent(in) :: pres(ka,ia,ja)
954  real(RP), intent(in) :: qdry(ka,ia,ja)
955  real(RP), intent(in) :: qv (ka,ia,ja)
956  real(RP), intent(in) :: m0_init ! initial total num. conc. of modes (Atk,Acm,Cor) [#/m3]
957  real(RP), intent(in) :: dg_init ! initial number equivalen diameters of modes [m]
958  real(RP), intent(in) :: sg_init ! initial standard deviation [-]
959  real(RP), intent(in) :: d_min_inp(3)
960  real(RP), intent(in) :: d_max_inp(3)
961  real(RP), intent(in) :: k_min_inp(3)
962  real(RP), intent(in) :: k_max_inp(3)
963  integer, intent(in) :: n_kap_inp(3)
964 
965  real(RP), intent(out) :: qtrc(ka,ia,ja,qa_ae)
966  real(RP), intent(out) :: ccn (ka,ia,ja)
967 
968  integer, parameter :: ia_m0 = 1 ! 1. number conc [#/m3]
969  integer, parameter :: ia_m2 = 2 ! 2. 2nd mom conc [m2/m3]
970  integer, parameter :: ia_m3 = 3 ! 3. 3rd mom conc [m3/m3]
971  integer, parameter :: ia_ms = 4 ! 4. mass conc [ug/m3]
972  integer, parameter :: ia_kp = 5
973  integer, parameter :: ik_out = 1
974 
975  real(RP) :: c_kappa = 0.3_rp ! hygroscopicity of condensable mass [-]
976  real(RP), parameter :: cleannumber = 1.e-3_rp
977  ! local variables
978  real(RP), parameter :: rhod_ae = 1.83_rp ! particle density [g/cm3] sulfate assumed
979  real(RP), parameter :: conv_vl_ms = rhod_ae/1.e-12_rp ! M3(volume)[m3/m3] to mass[m3/m3]
980 
981  real(RP) :: m0t, dgt, sgt, m2t, m3t, mst
982  !--- gas
983  real(RP) :: conc_gas(gas_ctg) !concentration [ug/m3]
984  !--- bin settings (lower bound, center, upper bound)
985  real(RP), allocatable :: d_lw(:,:), d_ct(:,:), d_up(:,:) !diameter [m]
986  real(RP), allocatable :: k_lw(:,:), k_ct(:,:), k_up(:,:) !kappa [-]
987  real(RP) :: dlogd, dk !delta log(D), delta K
988  real(RP), allocatable :: d_min(:) !lower bound of 1st size bin (n_ctg)
989  real(RP), allocatable :: d_max(:) !upper bound of last size bin (n_ctg)
990  integer, allocatable :: n_kap(:) !number of kappa bins (n_ctg)
991  real(RP), allocatable :: k_min(:) !lower bound of 1st kappa bin (n_ctg)
992  real(RP), allocatable :: k_max(:)
993 
994  real(RP) :: qsat_tmp, ssliq
995  real(RP) :: pi6
996 
997  integer :: ia0, ik, is0, ic, k, i, j, it
998  !---------------------------------------------------------------------------
999 
1000 
1001  call atmos_phy_ae_kajino13_setup
1002 
1003 
1004  pi6 = pi / 6._rp
1005  n_ctg = ae_ctg
1006 
1007  allocate( d_min(n_ctg) )
1008  allocate( d_max(n_ctg) )
1009  allocate( n_kap(n_ctg) )
1010  allocate( k_min(n_ctg) )
1011  allocate( k_max(n_ctg) )
1012 
1013 
1014  d_min(1:n_ctg) = d_min_inp(1:n_ctg) ! lower bound of 1st size bin
1015  d_max(1:n_ctg) = d_max_inp(1:n_ctg) ! upper bound of last size bin
1016  n_kap(1:n_ctg) = n_kap_inp(1:n_ctg) ! number of kappa bins
1017  k_min(1:n_ctg) = k_min_inp(1:n_ctg) ! lower bound of 1st kappa bin
1018  k_max(1:n_ctg) = k_max_inp(1:n_ctg) ! upper bound of last kappa bin
1019 
1020 
1021  !bin setting
1022  allocate(d_lw(n_siz_max,n_ctg))
1023  allocate(d_ct(n_siz_max,n_ctg))
1024  allocate(d_up(n_siz_max,n_ctg))
1025  allocate(k_lw(n_kap_max,n_ctg))
1026  allocate(k_ct(n_kap_max,n_ctg))
1027  allocate(k_up(n_kap_max,n_ctg))
1028  d_lw(:,:) = 0._rp
1029  d_ct(:,:) = 0._rp
1030  d_up(:,:) = 0._rp
1031  k_lw(:,:) = 0._rp
1032  k_ct(:,:) = 0._rp
1033  k_up(:,:) = 0._rp
1034 
1035  do ic = 1, n_ctg
1036 
1037  dlogd = (log(d_max(ic)) - log(d_min(ic)))/real(NSIZ(ic),kind=rp)
1038 
1039  do is0 = 1, nsiz(ic) !size bin
1040  d_lw(is0,ic) = exp(log(d_min(ic))+dlogd* real(is0-1,kind=RP) )
1041  d_ct(is0,ic) = exp(log(d_min(ic))+dlogd*(real(is0 ,kind=rp)-0.5_rp))
1042  d_up(is0,ic) = exp(log(d_min(ic))+dlogd* real(is0 ,kind=RP) )
1043  enddo !is (1:n_siz(ic))
1044  dk = (k_max(ic) - k_min(ic))/real(n_kap(ic),kind=rp)
1045  do ik = 1, n_kap(ic) !size bin
1046  k_lw(ik,ic) = k_min(ic) + dk * real(ik-1,kind=rp)
1047  k_ct(ik,ic) = k_min(ic) + dk *(real(ik ,kind=rp)-0.5_rp)
1048  k_up(ik,ic) = k_min(ic) + dk * real(ik ,kind=rp)
1049  enddo !ik (1:n_kap(ic))
1050 
1051  enddo !ic (1:n_ctg)
1052 ! ik = 1 !only one kappa bin
1053 
1054  m0t = m0_init !total M0 [#/m3]
1055  dgt = dg_init ![m]
1056  sgt = sg_init ![-]
1057 
1058  if ( m0t <= cleannumber ) then
1059  m0t = cleannumber
1060  dgt = 0.1e-6_rp
1061  sgt = 1.3_rp
1062  log_warn("ATMOS_PHY_AE_kajino13_mkinit",*) 'Initial aerosol number is set as ', cleannumber, '[#/m3]'
1063  endif
1064 
1065  m2t = m0t*dgt**(2.d0) *dexp(2.0d0 *(dlog(real(sgt,kind=dp))**2.d0)) !total M2 [m2/m3]
1066  m3t = m0t*dgt**(3.d0) *dexp(4.5d0 *(dlog(real(sgt,kind=dp))**2.d0)) !total M3 [m3/m3]
1067  mst = m3t*pi6*conv_vl_ms !total Ms [ug/m3]
1068 
1069  do ic = 1, n_ctg
1070  !aerosol_procs initial condition
1071  do ik = 1, n_kap(ic) !kappa bin
1072  do is0 = 1, nsiz(ic)
1073  if (dgt >= d_lw(is0,ic) .and. dgt < d_up(is0,ic)) then
1074  aerosol_procs(ia_m0,is0,ik,ic) = aerosol_procs(ia_m0,is0,ik,ic) + m0t ![#/m3]
1075  aerosol_procs(ia_m2,is0,ik,ic) = aerosol_procs(ia_m2,is0,ik,ic) + m2t ![m2/m3]
1076  aerosol_procs(ia_m3,is0,ik,ic) = aerosol_procs(ia_m3,is0,ik,ic) + m3t ![m3/m3]
1077  aerosol_procs(ia_ms,is0,ik,ic) = aerosol_procs(ia_ms,is0,ik,ic) + mst*1.e-9_rp ![kg/m3]
1078  elseif (dgt < d_lw(1,ic)) then
1079  aerosol_procs(ia_m0,1 ,ik,ic) = aerosol_procs(ia_m0,1 ,ik,ic) + m0t ![#/m3]
1080  aerosol_procs(ia_m2,1 ,ik,ic) = aerosol_procs(ia_m2,1 ,ik,ic) + m2t ![m2/m3]
1081  aerosol_procs(ia_m3,1 ,ik,ic) = aerosol_procs(ia_m3,1 ,ik,ic) + m3t ![m3/m3]
1082  aerosol_procs(ia_ms,1 ,ik,ic) = aerosol_procs(ia_ms,1 ,ik,ic) + mst*1.e-9_rp ![kg/m3]
1083  elseif (dgt >= d_up(nsiz(ic),ic)) then
1084  aerosol_procs(ia_m0,nsiz(ic),ik,ic) = aerosol_procs(ia_m0,nsiz(ic),ik,ic) + m0t ![#/m3]
1085  aerosol_procs(ia_m2,nsiz(ic),ik,ic) = aerosol_procs(ia_m2,nsiz(ic),ik,ic) + m2t ![m2/m3]
1086  aerosol_procs(ia_m3,nsiz(ic),ik,ic) = aerosol_procs(ia_m3,nsiz(ic),ik,ic) + m3t ![m3/m3]
1087  aerosol_procs(ia_ms,nsiz(ic),ik,ic) = aerosol_procs(ia_ms,nsiz(ic),ik,ic) + mst*1.e-9_rp ![kg/m3]
1088  endif
1089  enddo
1090  enddo
1091  enddo
1092 
1093  conc_gas(:) = 0.0_rp
1094  do j = 1, ja
1095  do i = 1, ia
1096  do k = 1, ka
1097  it = 1
1098  do ic = 1, n_ctg ! category
1099  do ik = 1, nkap(ic) ! kappa bin
1100  do is0 = 1, nsiz(ic) ! size bin
1101  do ia0 = 1, n_atr ! attributes
1102  qtrc(k,i,j,it) = aerosol_procs(ia0,is0,ik,ic) / dens(k,i,j) !#,m2,m3,kg/m3 -> #,m2,m3,kg/kg
1103  it = it + 1
1104  enddo !ia0 (1:N_ATR )
1105  enddo !is (1:n_siz(ic) )
1106  enddo !ik (1:n_kap(ic) )
1107  enddo !ic (1:n_ctg )
1108 
1109  do ic = 1, gas_ctg ! GAS category
1110  qtrc(k,i,j,it) = conc_gas(ic) / dens(k,i,j) * 1.e-9_rp !mixing ratio [kg/kg]
1111  it = it + 1
1112  enddo
1113  enddo
1114  enddo
1115  enddo
1116 
1117  ccn(:,:,:) = 0.0_rp
1118  do j = js, je
1119  do i = is, ie
1120  do k = ks, ke
1121  !--- calculate super saturation of water
1122  call saturation_pres2qsat_liq( temp(k,i,j), pres(k,i,j), qdry(k,i,j), & ! [IN]
1123  qsat_tmp ) ! [OUT]
1124  ssliq = qv(k,i,j) / qsat_tmp - 1.0_rp
1125 
1126  call aerosol_activation( c_kappa, ssliq, temp(k,i,j), ia_m0, ia_m2, ia_m3, &
1127  n_atr, n_siz_max, n_kap_max, n_ctg, nsiz, n_kap, &
1128  d_ct, aerosol_procs, aerosol_activ )
1129 
1130  do is0 = 1, nsiz(ic_mix)
1131  ccn(k,i,j) = ccn(k,i,j) + aerosol_activ(ia_m0,is0,ik_out,ic_mix)
1132  enddo
1133  enddo
1134  enddo
1135  enddo
1136 
1137  return
module atmosphere / saturation
integer, public ia
of whole cells: x, local, with HALO
integer, public ja
of whole cells: y, local, with HALO
integer, public is
start point of inner domain: x, local
integer, public ie
end point of inner domain: x, local
integer, public ke
end point of inner domain: z, local
integer, public je
end point of inner domain: y, local
integer, public ks
start point of inner domain: z, local
module CONSTANT
Definition: scale_const.F90:11
integer, public js
start point of inner domain: y, local
integer, public ka
of whole cells: z, local, with HALO
real(rp), public const_pi
pi
Definition: scale_const.F90:31
integer, parameter, public rp
integer, parameter, public dp
Here is the call graph for this function:
Here is the caller graph for this function:

◆ atmos_phy_ae_kajino13_negative_fixer()

subroutine, public scale_atmos_phy_ae_kajino13::atmos_phy_ae_kajino13_negative_fixer ( integer, intent(in)  KA,
integer, intent(in)  KS,
integer, intent(in)  KE,
integer, intent(in)  IA,
integer, intent(in)  IS,
integer, intent(in)  IE,
integer, intent(in)  JA,
integer, intent(in)  JS,
integer, intent(in)  JE,
integer, intent(in)  QA_AE,
real(rp), dimension(ka,ia,ja,qa_ae), intent(inout)  QTRC 
)

Definition at line 1143 of file scale_atmos_phy_ae_kajino13.F90.

References atmos_phy_ae_kajino13_c_kappa, atmos_phy_ae_kajino13_h2so4dt, atmos_phy_ae_kajino13_logk_aenucl, atmos_phy_ae_kajino13_ocgasdt, float(), scale_atmos_grid_cartesc_index::ke, scale_atmos_grid_cartesc_index::ks, scale_prof::prof_rapend(), and scale_prof::prof_rapstart().

Referenced by mod_atmos_phy_ae_driver::atmos_phy_ae_driver_adjustment().

1143  integer, intent(in) :: ka, ks, ke
1144  integer, intent(in) :: ia, is, ie
1145  integer, intent(in) :: ja, js, je
1146  integer, intent(in) :: qa_ae
1147 
1148  real(RP), intent(inout) :: qtrc(ka,ia,ja,qa_ae)
1149 
1150  integer :: k ,i, j
1151  integer :: ic, ik, is0
1152 
1153  !--- Negative fixer
1154  do j = js, je
1155  do i = is, ie
1156  do k = ks, ke
1157  do ic = 1, n_ctg !aerosol category
1158  do ik = 1, n_kap(ic) !kappa bin
1159  do is0 = 1, n_siz(ic) !size bin
1160  if (qtrc(k,i,j,it_procs2trans(ia_m0,is0,ik,ic)) < 0.0_rp .or. &
1161  qtrc(k,i,j,it_procs2trans(ia_m2,is0,ik,ic)) < 0.0_rp .or. &
1162  qtrc(k,i,j,it_procs2trans(ia_m3,is0,ik,ic)) < 0.0_rp .or. &
1163  qtrc(k,i,j,it_procs2trans(ia_ms,is0,ik,ic)) < 0.0_rp .or. &
1164  qtrc(k,i,j,it_procs2trans(ia_kp,is0,ik,ic)) < 0.0_rp ) then
1165  qtrc(k,i,j,it_procs2trans(1:n_atr,is0,ik,ic)) = 0.0_rp
1166  endif
1167  enddo
1168  enddo
1169  enddo
1170  enddo
1171  enddo
1172  enddo
1173 
1174  return
integer, public ia
of whole cells: x, local, with HALO
integer, public ja
of whole cells: y, local, with HALO
integer, public is
start point of inner domain: x, local
integer, public ie
end point of inner domain: x, local
integer, public ke
end point of inner domain: z, local
integer, public je
end point of inner domain: y, local
integer, public ks
start point of inner domain: z, local
integer, public js
start point of inner domain: y, local
integer, public ka
of whole cells: z, local, with HALO
Here is the call graph for this function:
Here is the caller graph for this function:

Variable Documentation

◆ atmos_phy_ae_kajino13_name

character(len=h_short), dimension(:), allocatable, public scale_atmos_phy_ae_kajino13::atmos_phy_ae_kajino13_name

Definition at line 51 of file scale_atmos_phy_ae_kajino13.F90.

Referenced by mod_atmos_phy_ae_driver::atmos_phy_ae_driver_tracer_setup(), and atmos_phy_ae_kajino13_tracer_setup().

51  character(len=H_SHORT), public, allocatable :: atmos_phy_ae_kajino13_name(:)

◆ atmos_phy_ae_kajino13_desc

character(len=h_mid), dimension(:), allocatable, public scale_atmos_phy_ae_kajino13::atmos_phy_ae_kajino13_desc

Definition at line 52 of file scale_atmos_phy_ae_kajino13.F90.

Referenced by mod_atmos_phy_ae_driver::atmos_phy_ae_driver_tracer_setup(), and atmos_phy_ae_kajino13_tracer_setup().

52  character(len=H_MID) , public, allocatable :: atmos_phy_ae_kajino13_desc(:)

◆ atmos_phy_ae_kajino13_unit

character(len=h_short), dimension(:), allocatable, public scale_atmos_phy_ae_kajino13::atmos_phy_ae_kajino13_unit

Definition at line 53 of file scale_atmos_phy_ae_kajino13.F90.

Referenced by mod_atmos_phy_ae_driver::atmos_phy_ae_driver_tracer_setup(), and atmos_phy_ae_kajino13_tracer_setup().

53  character(len=H_SHORT), public, allocatable :: atmos_phy_ae_kajino13_unit(:)

◆ atmos_phy_ae_kajino13_h2so4dt

real(rp), public scale_atmos_phy_ae_kajino13::atmos_phy_ae_kajino13_h2so4dt = 5.E-6_RP

Definition at line 55 of file scale_atmos_phy_ae_kajino13.F90.

Referenced by atmos_phy_ae_kajino13_negative_fixer(), and atmos_phy_ae_kajino13_setup().

55  real(RP), public :: atmos_phy_ae_kajino13_h2so4dt = 5.e-6_rp ! h2so4 production rate (Temporal) [ug/m3/s]

◆ atmos_phy_ae_kajino13_ocgasdt

real(rp), public scale_atmos_phy_ae_kajino13::atmos_phy_ae_kajino13_ocgasdt = 8.E-5_RP

Definition at line 56 of file scale_atmos_phy_ae_kajino13.F90.

Referenced by atmos_phy_ae_kajino13_negative_fixer(), and atmos_phy_ae_kajino13_setup().

56  real(RP), public :: atmos_phy_ae_kajino13_ocgasdt = 8.e-5_rp ! other condensational bas production rate 16*h2so4dt (see Kajino et al. 2013)

◆ atmos_phy_ae_kajino13_c_kappa

real(rp), public scale_atmos_phy_ae_kajino13::atmos_phy_ae_kajino13_c_kappa = 0.3_RP

Definition at line 58 of file scale_atmos_phy_ae_kajino13.F90.

Referenced by atmos_phy_ae_kajino13_negative_fixer(), and atmos_phy_ae_kajino13_setup().

58  real(RP), public :: atmos_phy_ae_kajino13_c_kappa = 0.3_rp ! hygroscopicity of condensable mass [-]

◆ atmos_phy_ae_kajino13_flag_npf

logical, public scale_atmos_phy_ae_kajino13::atmos_phy_ae_kajino13_flag_npf = .false.

Definition at line 59 of file scale_atmos_phy_ae_kajino13.F90.

Referenced by atmos_phy_ae_kajino13_setup(), and atmos_phy_ae_kajino13_tendency().

59  logical, public :: atmos_phy_ae_kajino13_flag_npf = .false.

◆ atmos_phy_ae_kajino13_flag_cond

logical, public scale_atmos_phy_ae_kajino13::atmos_phy_ae_kajino13_flag_cond = .true.

Definition at line 60 of file scale_atmos_phy_ae_kajino13.F90.

Referenced by atmos_phy_ae_kajino13_setup(), and atmos_phy_ae_kajino13_tendency().

60  logical, public :: atmos_phy_ae_kajino13_flag_cond = .true.

◆ atmos_phy_ae_kajino13_flag_coag

logical, public scale_atmos_phy_ae_kajino13::atmos_phy_ae_kajino13_flag_coag = .true.

Definition at line 61 of file scale_atmos_phy_ae_kajino13.F90.

Referenced by atmos_phy_ae_kajino13_setup(), and atmos_phy_ae_kajino13_tendency().

61  logical, public :: atmos_phy_ae_kajino13_flag_coag = .true.

◆ atmos_phy_ae_kajino13_flag_ccn_interactive

logical, public scale_atmos_phy_ae_kajino13::atmos_phy_ae_kajino13_flag_ccn_interactive = .true.

Definition at line 62 of file scale_atmos_phy_ae_kajino13.F90.

Referenced by atmos_phy_ae_kajino13_setup(), and atmos_phy_ae_kajino13_tendency().

62  logical, public :: atmos_phy_ae_kajino13_flag_ccn_interactive = .true.

◆ atmos_phy_ae_kajino13_flag_regeneration

logical, public scale_atmos_phy_ae_kajino13::atmos_phy_ae_kajino13_flag_regeneration = .true.

Definition at line 63 of file scale_atmos_phy_ae_kajino13.F90.

Referenced by atmos_phy_ae_kajino13_setup(), and atmos_phy_ae_kajino13_tendency().

63  logical, public :: atmos_phy_ae_kajino13_flag_regeneration = .true.

◆ atmos_phy_ae_kajino13_dg_reg

real(rp), public scale_atmos_phy_ae_kajino13::atmos_phy_ae_kajino13_dg_reg = 5.E-7_RP

Definition at line 64 of file scale_atmos_phy_ae_kajino13.F90.

Referenced by atmos_phy_ae_kajino13_setup(), and atmos_phy_ae_kajino13_tendency().

64  real(RP), public :: atmos_phy_ae_kajino13_dg_reg = 5.e-7_rp ! dg of regenerated aerosol (droplet mode 500 nm) [m]

◆ atmos_phy_ae_kajino13_sg_reg

real(rp), public scale_atmos_phy_ae_kajino13::atmos_phy_ae_kajino13_sg_reg = 1.6_RP

Definition at line 65 of file scale_atmos_phy_ae_kajino13.F90.

Referenced by atmos_phy_ae_kajino13_setup(), and atmos_phy_ae_kajino13_tendency().

65  real(RP), public :: atmos_phy_ae_kajino13_sg_reg = 1.6_rp ! sg of regenerated aerosol [-]

◆ atmos_phy_ae_kajino13_logk_aenucl

real(rp), public scale_atmos_phy_ae_kajino13::atmos_phy_ae_kajino13_logk_aenucl = -12.4_RP

Definition at line 66 of file scale_atmos_phy_ae_kajino13.F90.

Referenced by atmos_phy_ae_kajino13_negative_fixer(), and atmos_phy_ae_kajino13_setup().

66  real(RP), public :: atmos_phy_ae_kajino13_logk_aenucl = -12.4_rp !constant coefficient for kinetic nucleation [-]