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_finalize
 finalize 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 184 of file scale_atmos_phy_ae_kajino13.F90.

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

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().

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 358 of file scale_atmos_phy_ae_kajino13.F90.

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

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().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ atmos_phy_ae_kajino13_finalize()

subroutine, public scale_atmos_phy_ae_kajino13::atmos_phy_ae_kajino13_finalize

finalize

Definition at line 604 of file scale_atmos_phy_ae_kajino13.F90.

604 
605  deallocate( nsiz )
606  deallocate( nkap )
607 
608  deallocate( atmos_phy_ae_kajino13_name )
609  deallocate( atmos_phy_ae_kajino13_desc )
610  deallocate( atmos_phy_ae_kajino13_unit )
611 
612  deallocate( rnum_out )
613  deallocate( n_siz )
614  deallocate( d_min )
615  deallocate( d_max )
616  deallocate( n_kap )
617  deallocate( k_min )
618  deallocate( k_max )
619  deallocate( ctg_name )
620 
621  deallocate( d_lw )
622  deallocate( d_ct )
623  deallocate( d_up )
624  deallocate( k_lw )
625  deallocate( k_ct )
626  deallocate( k_up )
627 
628  deallocate( is_i )
629  deallocate( is_j )
630  deallocate( is_k )
631  deallocate( ik_i )
632  deallocate( ik_j )
633  deallocate( ik_k )
634  deallocate( ic_i )
635  deallocate( ic_j )
636  deallocate( ic_k )
637 
638  deallocate( it_procs2trans )
639  deallocate( ia_trans2procs )
640  deallocate( is_trans2procs )
641  deallocate( ik_trans2procs )
642  deallocate( ic_trans2procs )
643 
644  deallocate( aerosol_procs )
645  deallocate( aerosol_activ )
646  deallocate( emis_procs )
647  deallocate( emis_gas )
648 
649  return

References atmos_phy_ae_kajino13_desc, atmos_phy_ae_kajino13_name, and atmos_phy_ae_kajino13_unit.

Referenced by mod_atmos_phy_ae_driver::atmos_phy_ae_driver_finalize().

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 670 of file scale_atmos_phy_ae_kajino13.F90.

670  use scale_atmos_saturation, only: &
671  pres2qsat_liq => atmos_saturation_pres2qsat_liq
672  use scale_file_history, only: &
673  file_history_in
674  implicit none
675  integer, intent(in) :: KA, KS, KE
676  integer, intent(in) :: IA, IS, IE
677  integer, intent(in) :: JA, JS, JE
678  integer, intent(in) :: QA_AE
679  real(RP), intent(in) :: TEMP(KA,IA,JA)
680  real(RP), intent(in) :: PRES(KA,IA,JA)
681  real(RP), intent(in) :: QDRY(KA,IA,JA)
682  real(RP), intent(in) :: NREG(KA,IA,JA)
683  real(RP), intent(in) :: DENS(KA,IA,JA)
684  real(RP), intent(in) :: QV (KA,IA,JA)
685  real(RP), intent(in) :: QTRC(KA,IA,JA,QA_AE)
686  real(RP), intent(in) :: EMIT(KA,IA,JA,QA_AE)
687  real(DP), intent(in) :: dt
688 
689  real(RP), intent(out) :: RHOQ_t_AE(KA,IA,JA,QA_AE)
690  real(RP), intent(out) :: CN(KA,IA,JA)
691  real(RP), intent(out) :: CCN(KA,IA,JA)
692 
693  real(RP) :: QTRC0(KA,IA,JA,QA_AE)
694  real(RP) :: QTRC1(KA,IA,JA,QA_AE)
695 
696  !--- local
697  real(RP) :: ssliq_ae(KA,IA,JA)
698  real(RP) :: rrhog(KA,IA,JA)
699 ! real(RP) :: t_ccn, t_cn
700  real(RP) :: qsat_tmp
701  real(RP) :: m0_reg, m2_reg, m3_reg !regenerated aerosols [m^k/m3]
702  real(RP) :: ms_reg !regenerated aerosol mass [ug/m3]
703  real(RP) :: reg_factor_m2,reg_factor_m3 !to save cpu time for moment conversion
704  !--- aerosol variables
705  real(RP) :: total_aerosol_mass(KA,IA,JA,n_ctg)
706  real(RP) :: total_aerosol_number(KA,IA,JA,n_ctg)
707  real(RP) :: total_emit_aerosol_mass(KA,IA,JA,n_ctg)
708  real(RP) :: total_emit_aerosol_number(KA,IA,JA,n_ctg)
709  !--- gas
710 ! real(RP),allocatable :: conc_h2so4(:,:,:,:) !concentration [ug/m3]
711 ! real(RP) :: conc_gas(KA,IA,JA,GAS_CTG) !concentration [ug/m3]
712  real(RP) :: conc_gas(GAS_CTG) !concentration [ug/m3]
713  integer :: i, j, k, iq, it
714 
715  log_progress(*) 'atmosphere / physics / aerosol / kajino13'
716 
717  aerosol_procs(:,:,:,:) = 0.0_rp
718  aerosol_activ(:,:,:,:) = 0.0_rp
719  emis_procs(:,:,:,:) = 0.0_rp
720  emis_gas(:) = 0.0_rp
721 
722  do iq = 1, qa_ae
723  do j = js, je
724  do i = is, ie
725  do k = ks, ke
726  qtrc0(k,i,j,iq) = qtrc(k,i,j,iq) ! save
727  enddo
728  enddo
729  enddo
730  enddo
731 
732  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
733  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
734 
735  !--- convert SCALE variable to zerochem variable
736 
737  do j = js, je
738  do i = is, ie
739 
740  do k = ks, ke
741  rrhog(k,i,j) = 1.0_rp / dens(k,i,j)
742  enddo
743  do k = ks, ke
744  !--- calculate super saturation of water
745  call pres2qsat_liq( temp(k,i,j), pres(k,i,j), qdry(k,i,j), & ! [IN]
746  qsat_tmp ) ! [OUT]
747  ssliq_ae(k,i,j) = qv(k,i,j)/qsat_tmp - 1.0_rp
748  enddo
749 
750  enddo
751  enddo
752 
753  ! tiny number, tiny mass
754  do j = js, je
755  do i = is, ie
756  do k = ks, ke
757  do ic = 1, n_ctg !aerosol category
758  do ik = 1, n_kap(ic) !kappa bin
759  do is0 = 1, n_siz(ic) !size bin
760  if (qtrc0(k,i,j,it_procs2trans(ia_m0,is0,ik,ic))*dens(k,i,j) < cleannumber) then
761  do ia0 = 1, n_atr
762  qtrc0(k,i,j,it_procs2trans(ia0,is0,ik,ic)) = 0._rp !to save cpu time and avoid underflow
763  enddo !ia0 (1:n_atr )
764  endif
765  enddo !is (1:n_siz(ic) )
766  enddo !ik (1:n_kap(ic) )
767  enddo !ic (1:n_ctg )
768  enddo
769  enddo
770  enddo
771 
772  do iq = 1, qa_ae
773  do j = js, je
774  do i = is, ie
775  do k = ks, ke
776  qtrc1(k,i,j,iq) = qtrc0(k,i,j,iq) ! save
777  enddo
778  enddo
779  enddo
780  enddo
781 
782  !---- Calculate aerosol processs
783  cn(:,:,:) = 0.0_rp
784  ccn(:,:,:) = 0.0_rp
785  do k = ks, ke
786  do j = js, je
787  do i = is, ie
788 
789  !--- aerosol_trans at initial time
790  ! [xx/kg] -> [xx/m3]
791  do it = 1, n_trans
792  aerosol_procs(ia_trans2procs(it), &
793  is_trans2procs(it), &
794  ik_trans2procs(it), &
795  ic_trans2procs(it)) = qtrc1(k,i,j,it)*dens(k,i,j)
796  emis_procs(ia_trans2procs(it), &
797  is_trans2procs(it), &
798  ik_trans2procs(it), &
799  ic_trans2procs(it)) = emit(k,i,j,it)*dens(k,i,j)
800  enddo !it(1:n_trans)
801  ! mixing ratio [kg/kg] -> concentration [ug/m3]
802  conc_gas(1:gas_ctg) &
803  = qtrc1(k,i,j,qa_ae-gas_ctg+1:qa_ae-gas_ctg+gas_ctg)*dens(k,i,j)*1.e+9_rp
804 
805  emis_gas(1:gas_ctg) = emit(k,i,j,qa_ae-gas_ctg+ig_h2so4:qa_ae-gas_ctg+ig_cgas)
806 
807  call aerosol_zerochem( &
808  dt, & !--- in
809  temp(k,i,j), & !--- in
810  pres(k,i,j), & !--- in
811  ssliq_ae(k,i,j), & !--- in
812  atmos_phy_ae_kajino13_flag_npf, & !--- in
813  atmos_phy_ae_kajino13_flag_cond, & !--- in
814  atmos_phy_ae_kajino13_flag_coag, & !--- in
815  aerosol_procs, & !--- inout
816  conc_gas, & !--- inout
817  emis_procs, & !--- out
818  emis_gas, & !--- out
819  aerosol_activ ) !--- out
820 
821 ! aerosol loss due to activation to cloud droplets
822  if (atmos_phy_ae_kajino13_flag_ccn_interactive) then
823  do is0 = 1, n_siz(ic_mix)
824  do ia0 = 1, n_atr !attributes
825  aerosol_activ(ia0,is0,ik_out,ic_mix) = min(max(0._rp, aerosol_activ(ia0,is0,ik_out,ic_mix)), &
826  aerosol_procs(ia0,is0,ik_out,ic_mix) )
827 
828  aerosol_procs(ia0,is0,ik_out,ic_mix) = &
829  aerosol_procs(ia0,is0,ik_out,ic_mix) - aerosol_activ(ia0,is0,ik_out,ic_mix)
830  enddo
831  enddo
832  endif !flag_ccn_interactive
833 
834 ! aerosol regeneration due to evaporation of cloud droplets
835 ! using prescribed size parameters and to internal mixture category (ic_mix)
836  if (atmos_phy_ae_kajino13_flag_regeneration) then
837  m0_reg = nreg(k,i,j) !#/m3
838 ! 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
839 ! 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
840  m2_reg = m0_reg * reg_factor_m2 !m2/m3
841  m3_reg = m0_reg * reg_factor_m3 !m3/m3
842  ms_reg = m3_reg * pi6 * conv_vl_ms !ug/m3
843  aerosol_procs(ia_m0,is0_reg,ik_out,ic_mix) = &
844  aerosol_procs(ia_m0,is0_reg,ik_out,ic_mix) + m0_reg !#/m3
845  aerosol_procs(ia_m2,is0_reg,ik_out,ic_mix) = &
846  aerosol_procs(ia_m2,is0_reg,ik_out,ic_mix) + m2_reg !m2/m3
847  aerosol_procs(ia_m3,is0_reg,ik_out,ic_mix) = &
848  aerosol_procs(ia_m3,is0_reg,ik_out,ic_mix) + m3_reg !m3/m3
849  aerosol_procs(ia_ms,is0_reg,ik_out,ic_mix) = &
850  aerosol_procs(ia_ms,is0_reg,ik_out,ic_mix) + ms_reg !ug/m3
851 ! additional attirbute to be added (ia_kp)
852  endif !flag_regeneration
853 
854 ! diagnosed variables
855  do is0 = 1, n_siz(ic_mix)
856  ccn(k,i,j) = ccn(k,i,j) + aerosol_activ(ia_m0,is0,ik_out,ic_mix)
857  cn(k,i,j) = cn(k,i,j) + aerosol_procs(ia_m0,is0,ik_out,ic_mix)
858  enddo
859 
860 ! call trans_ccn(aerosol_procs, aerosol_activ, t_ccn, t_cn, &
861 ! n_ctg, n_kap_max, n_siz_max, N_ATR, &
862 ! ic_mix, ia_m0, ia_m2, ia_m3, ik_out, n_siz, &
863 ! rnum_out, ATMOS_PHY_AE_KAJINO13_nbins_out)
864 
865 ! CN(k,i,j) = t_cn
866 ! CCN(k,i,j) = t_ccn
867 
868  ! [xx/m3] -> [xx/kg]
869  do ic = 1, n_ctg !category
870  do ik = 1, n_kap(ic) !kappa bin
871  do is0 = 1, n_siz(ic) !size bin
872  do ia0 = 1, n_atr !attributes
873  qtrc1(k,i,j,it_procs2trans(ia0,is0,ik,ic)) = aerosol_procs(ia0,is0,ik,ic) / dens(k,i,j)
874  enddo !ia (1:N_ATR )
875  enddo !is (1:n_siz(ic) )
876  enddo !ik (1:n_kap(ic) )
877  enddo !ic (1:n_ctg )
878  ! [ug/m3] -> mixing ratio [kg/kg]
879  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
880 
881  ! tiny number, tiny mass
882  do ic = 1, n_ctg !aerosol category
883  do ik = 1, n_kap(ic) !kappa bin
884  do is0 = 1, n_siz(ic) !size bin
885  if (qtrc1(k,i,j,it_procs2trans(ia_m0,is0,ik,ic))*dens(k,i,j) < cleannumber) then
886  do ia0 = 1, n_atr
887  qtrc1(k,i,j,it_procs2trans(ia0,is0,ik,ic)) = 0._rp !to save cpu time and avoid underflow
888  enddo !ia0 (1:n_atr )
889  endif
890  enddo !is (1:n_siz(ic) )
891  enddo !ik (1:n_kap(ic) )
892  enddo !ic (1:n_ctg )
893 
894  ! for history
895  total_aerosol_mass(k,i,j,:) = 0.0_rp
896  total_aerosol_number(k,i,j,:) = 0.0_rp
897  total_emit_aerosol_mass(k,i,j,:) = 0.0_rp
898  total_emit_aerosol_number(k,i,j,:) = 0.0_rp
899  do ic = 1, n_ctg
900  do ik = 1, n_kap(ic)
901  do is0 = 1, n_siz(ic)
902  total_aerosol_mass(k,i,j,ic) = total_aerosol_mass(k,i,j,ic) &
903  + qtrc1(k,i,j,it_procs2trans(ia_ms,is0,ik,ic))
904  total_aerosol_number(k,i,j,ic) = total_aerosol_number(k,i,j,ic) &
905  + qtrc1(k,i,j,it_procs2trans(ia_m0,is0,ik,ic))
906  total_emit_aerosol_mass(k,i,j,ic) = total_emit_aerosol_mass(k,i,j,ic) &
907  + emit(k,i,j,it_procs2trans(ia_ms,is0,ik,ic))
908  total_emit_aerosol_number(k,i,j,ic) = total_emit_aerosol_number(k,i,j,ic) &
909  + emit(k,i,j,it_procs2trans(ia_m0,is0,ik,ic))
910  enddo
911  enddo
912  enddo
913 
914  enddo
915  enddo
916  enddo
917 
918  do ic = 1, n_ctg
919  call file_history_in( total_aerosol_mass(:,:,:,ic), trim(ctg_name(ic))//'_mass', 'Total mass mixing ratio of aerosol', 'kg/kg' )
920  call file_history_in( total_aerosol_number(:,:,:,ic), trim(ctg_name(ic))//'_number', 'Total number mixing ratio of aerosol', 'num/kg' )
921  call file_history_in( total_emit_aerosol_mass(:,:,:,ic), trim(ctg_name(ic))//'_mass_emit', 'Total mass mixing ratio of emitted aerosol', 'kg/kg' )
922  call file_history_in( total_emit_aerosol_number(:,:,:,ic), trim(ctg_name(ic))//'_number_emit', 'Total number mixing ratio of emitted aerosol', 'num/kg' )
923  enddo
924 
925  call file_history_in( emit(:,:,:,qa_ae-gas_ctg+ig_h2so4), 'H2SO4_emit', 'Emission ratio of H2SO4 gas', 'ug/m3/s' )
926  call file_history_in( emit(:,:,:,qa_ae-gas_ctg+ig_cgas), 'CGAS_emit', 'Emission ratio of Condensabule gas', 'ug/m3/s' )
927 
928 
929  do iq = 1, qa_ae
930  do j = js, je
931  do i = is, ie
932  do k = ks, ke
933  rhoq_t_ae(k,i,j,iq) = ( qtrc1(k,i,j,iq) - qtrc0(k,i,j,iq) ) * dens(k,i,j) / dt
934  enddo
935  enddo
936  enddo
937  enddo
938 
939  return

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, atmos_phy_ae_kajino13_sg_reg, scale_atmos_grid_cartesc_index::ie, scale_atmos_grid_cartesc_index::is, scale_atmos_grid_cartesc_index::je, scale_atmos_grid_cartesc_index::js, scale_tracer::k, scale_atmos_grid_cartesc_index::ke, and scale_atmos_grid_cartesc_index::ks.

Referenced by mod_atmos_phy_ae_driver::atmos_phy_ae_driver_calc_tendency().

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 948 of file scale_atmos_phy_ae_kajino13.F90.

948  use scale_atmos_aerosol, only: &
949  n_ae
950  implicit none
951  integer, intent(in) :: KA, IA, JA, QA_AE
952 
953  real(RP), intent(in) :: QTRC(KA,IA,JA,QA_AE) ! tracer mass concentration [kg/kg]
954  real(RP), intent(in) :: RH (KA,IA,JA) ! relative humidity (0-1)
955 
956  real(RP), intent(out) :: Re (KA,IA,JA,N_AE) ! effective radius
957  !---------------------------------------------------------------------------
958 
959  re(:,:,:,:) = 0.0_rp
960 
961 ! Re(:,:,:,I_ae_seasalt) = 2.E-4_RP
962 ! Re(:,:,:,I_ae_dust ) = 4.E-6_RP
963 ! Re(:,:,:,I_ae_bc ) = 4.E-8_RP
964 ! Re(:,:,:,I_ae_oc ) = RH(:,:,:)
965 ! Re(:,:,:,I_ae_sulfate) = RH(:,:,:)
966 
967  return

References scale_atmos_aerosol::n_ae.

Referenced by mod_atmos_phy_ae_vars::atmos_phy_ae_vars_get_diagnostic().

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 991 of file scale_atmos_phy_ae_kajino13.F90.

991  use scale_const, only: &
992  pi => const_pi
993  use scale_atmos_saturation, only: &
994  saturation_pres2qsat_liq => atmos_saturation_pres2qsat_liq
995  implicit none
996 
997  integer, intent(in) :: KA, KS, KE
998  integer, intent(in) :: IA, IS, IE
999  integer, intent(in) :: JA, JS, JE
1000  integer, intent(in) :: QA_AE
1001 
1002  real(RP), intent(in) :: DENS(KA,IA,JA)
1003  real(RP), intent(in) :: TEMP(KA,IA,JA)
1004  real(RP), intent(in) :: PRES(KA,IA,JA)
1005  real(RP), intent(in) :: QDRY(KA,IA,JA)
1006  real(RP), intent(in) :: QV (KA,IA,JA)
1007  real(RP), intent(in) :: m0_init ! initial total num. conc. of modes (Atk,Acm,Cor) [#/m3]
1008  real(RP), intent(in) :: dg_init ! initial number equivalen diameters of modes [m]
1009  real(RP), intent(in) :: sg_init ! initial standard deviation [-]
1010  real(RP), intent(in) :: d_min_inp(3)
1011  real(RP), intent(in) :: d_max_inp(3)
1012  real(RP), intent(in) :: k_min_inp(3)
1013  real(RP), intent(in) :: k_max_inp(3)
1014  integer, intent(in) :: n_kap_inp(3)
1015 
1016  real(RP), intent(out) :: QTRC(KA,IA,JA,QA_AE)
1017  real(RP), intent(out) :: CCN (KA,IA,JA)
1018 
1019  integer, parameter :: ia_m0 = 1 ! 1. number conc [#/m3]
1020  integer, parameter :: ia_m2 = 2 ! 2. 2nd mom conc [m2/m3]
1021  integer, parameter :: ia_m3 = 3 ! 3. 3rd mom conc [m3/m3]
1022  integer, parameter :: ia_ms = 4 ! 4. mass conc [ug/m3]
1023  integer, parameter :: ia_kp = 5
1024  integer, parameter :: ik_out = 1
1025 
1026  real(RP) :: c_kappa = 0.3_rp ! hygroscopicity of condensable mass [-]
1027  real(RP), parameter :: cleannumber = 1.e-3_rp
1028  ! local variables
1029  real(RP), parameter :: rhod_ae = 1.83_rp ! particle density [g/cm3] sulfate assumed
1030  real(RP), parameter :: conv_vl_ms = rhod_ae/1.e-12_rp ! M3(volume)[m3/m3] to mass[m3/m3]
1031 
1032  real(RP) :: m0t, dgt, sgt, m2t, m3t, mst
1033  !--- gas
1034  real(RP) :: conc_gas(GAS_CTG) !concentration [ug/m3]
1035  !--- bin settings (lower bound, center, upper bound)
1036  real(RP), allocatable :: d_lw(:,:), d_ct(:,:), d_up(:,:) !diameter [m]
1037  real(RP), allocatable :: k_lw(:,:), k_ct(:,:), k_up(:,:) !kappa [-]
1038  real(RP) :: dlogd, dk !delta log(D), delta K
1039  real(RP), allocatable :: d_min(:) !lower bound of 1st size bin (n_ctg)
1040  real(RP), allocatable :: d_max(:) !upper bound of last size bin (n_ctg)
1041  integer, allocatable :: n_kap(:) !number of kappa bins (n_ctg)
1042  real(RP), allocatable :: k_min(:) !lower bound of 1st kappa bin (n_ctg)
1043  real(RP), allocatable :: k_max(:)
1044 
1045  real(RP) :: qsat_tmp, ssliq
1046  real(RP) :: pi6
1047 
1048  integer :: ia0, ik, is0, ic, k, i, j, it
1049  !---------------------------------------------------------------------------
1050 
1051 
1052  call atmos_phy_ae_kajino13_setup
1053 
1054 
1055  pi6 = pi / 6._rp
1056  n_ctg = ae_ctg
1057 
1058  allocate( d_min(n_ctg) )
1059  allocate( d_max(n_ctg) )
1060  allocate( n_kap(n_ctg) )
1061  allocate( k_min(n_ctg) )
1062  allocate( k_max(n_ctg) )
1063 
1064 
1065  d_min(1:n_ctg) = d_min_inp(1:n_ctg) ! lower bound of 1st size bin
1066  d_max(1:n_ctg) = d_max_inp(1:n_ctg) ! upper bound of last size bin
1067  n_kap(1:n_ctg) = n_kap_inp(1:n_ctg) ! number of kappa bins
1068  k_min(1:n_ctg) = k_min_inp(1:n_ctg) ! lower bound of 1st kappa bin
1069  k_max(1:n_ctg) = k_max_inp(1:n_ctg) ! upper bound of last kappa bin
1070 
1071 
1072  !bin setting
1073  allocate(d_lw(n_siz_max,n_ctg))
1074  allocate(d_ct(n_siz_max,n_ctg))
1075  allocate(d_up(n_siz_max,n_ctg))
1076  allocate(k_lw(n_kap_max,n_ctg))
1077  allocate(k_ct(n_kap_max,n_ctg))
1078  allocate(k_up(n_kap_max,n_ctg))
1079  d_lw(:,:) = 0._rp
1080  d_ct(:,:) = 0._rp
1081  d_up(:,:) = 0._rp
1082  k_lw(:,:) = 0._rp
1083  k_ct(:,:) = 0._rp
1084  k_up(:,:) = 0._rp
1085 
1086  do ic = 1, n_ctg
1087 
1088  dlogd = (log(d_max(ic)) - log(d_min(ic)))/real(nsiz(ic),kind=rp)
1089 
1090  do is0 = 1, nsiz(ic) !size bin
1091  d_lw(is0,ic) = exp(log(d_min(ic))+dlogd* real(is0-1,kind=rp) )
1092  d_ct(is0,ic) = exp(log(d_min(ic))+dlogd*(real(is0 ,kind=rp)-0.5_rp))
1093  d_up(is0,ic) = exp(log(d_min(ic))+dlogd* real(is0 ,kind=rp) )
1094  enddo !is (1:n_siz(ic))
1095  dk = (k_max(ic) - k_min(ic))/real(n_kap(ic),kind=rp)
1096  do ik = 1, n_kap(ic) !size bin
1097  k_lw(ik,ic) = k_min(ic) + dk * real(ik-1,kind=rp)
1098  k_ct(ik,ic) = k_min(ic) + dk *(real(ik ,kind=rp)-0.5_rp)
1099  k_up(ik,ic) = k_min(ic) + dk * real(ik ,kind=rp)
1100  enddo !ik (1:n_kap(ic))
1101 
1102  enddo !ic (1:n_ctg)
1103 ! ik = 1 !only one kappa bin
1104 
1105  m0t = m0_init !total M0 [#/m3]
1106  dgt = dg_init ![m]
1107  sgt = sg_init ![-]
1108 
1109  if ( m0t <= cleannumber ) then
1110  m0t = cleannumber
1111  dgt = 0.1e-6_rp
1112  sgt = 1.3_rp
1113  log_warn("ATMOS_PHY_AE_kajino13_mkinit",*) 'Initial aerosol number is set as ', cleannumber, '[#/m3]'
1114  endif
1115 
1116  m2t = m0t*dgt**(2.d0) *dexp(2.0d0 *(dlog(real(sgt,kind=dp))**2.d0)) !total M2 [m2/m3]
1117  m3t = m0t*dgt**(3.d0) *dexp(4.5d0 *(dlog(real(sgt,kind=dp))**2.d0)) !total M3 [m3/m3]
1118  mst = m3t*pi6*conv_vl_ms !total Ms [ug/m3]
1119 
1120  do ic = 1, n_ctg
1121  !aerosol_procs initial condition
1122  do ik = 1, n_kap(ic) !kappa bin
1123  do is0 = 1, nsiz(ic)
1124  if (dgt >= d_lw(is0,ic) .and. dgt < d_up(is0,ic)) then
1125  aerosol_procs(ia_m0,is0,ik,ic) = aerosol_procs(ia_m0,is0,ik,ic) + m0t ![#/m3]
1126  aerosol_procs(ia_m2,is0,ik,ic) = aerosol_procs(ia_m2,is0,ik,ic) + m2t ![m2/m3]
1127  aerosol_procs(ia_m3,is0,ik,ic) = aerosol_procs(ia_m3,is0,ik,ic) + m3t ![m3/m3]
1128  aerosol_procs(ia_ms,is0,ik,ic) = aerosol_procs(ia_ms,is0,ik,ic) + mst*1.e-9_rp ![kg/m3]
1129  elseif (dgt < d_lw(1,ic)) then
1130  aerosol_procs(ia_m0,1 ,ik,ic) = aerosol_procs(ia_m0,1 ,ik,ic) + m0t ![#/m3]
1131  aerosol_procs(ia_m2,1 ,ik,ic) = aerosol_procs(ia_m2,1 ,ik,ic) + m2t ![m2/m3]
1132  aerosol_procs(ia_m3,1 ,ik,ic) = aerosol_procs(ia_m3,1 ,ik,ic) + m3t ![m3/m3]
1133  aerosol_procs(ia_ms,1 ,ik,ic) = aerosol_procs(ia_ms,1 ,ik,ic) + mst*1.e-9_rp ![kg/m3]
1134  elseif (dgt >= d_up(nsiz(ic),ic)) then
1135  aerosol_procs(ia_m0,nsiz(ic),ik,ic) = aerosol_procs(ia_m0,nsiz(ic),ik,ic) + m0t ![#/m3]
1136  aerosol_procs(ia_m2,nsiz(ic),ik,ic) = aerosol_procs(ia_m2,nsiz(ic),ik,ic) + m2t ![m2/m3]
1137  aerosol_procs(ia_m3,nsiz(ic),ik,ic) = aerosol_procs(ia_m3,nsiz(ic),ik,ic) + m3t ![m3/m3]
1138  aerosol_procs(ia_ms,nsiz(ic),ik,ic) = aerosol_procs(ia_ms,nsiz(ic),ik,ic) + mst*1.e-9_rp ![kg/m3]
1139  endif
1140  enddo
1141  enddo
1142  enddo
1143 
1144  conc_gas(:) = 0.0_rp
1145  do j = 1, ja
1146  do i = 1, ia
1147  do k = 1, ka
1148  it = 1
1149  do ic = 1, n_ctg ! category
1150  do ik = 1, nkap(ic) ! kappa bin
1151  do is0 = 1, nsiz(ic) ! size bin
1152  do ia0 = 1, n_atr ! attributes
1153  qtrc(k,i,j,it) = aerosol_procs(ia0,is0,ik,ic) / dens(k,i,j) !#,m2,m3,kg/m3 -> #,m2,m3,kg/kg
1154  it = it + 1
1155  enddo !ia0 (1:N_ATR )
1156  enddo !is (1:n_siz(ic) )
1157  enddo !ik (1:n_kap(ic) )
1158  enddo !ic (1:n_ctg )
1159 
1160  do ic = 1, gas_ctg ! GAS category
1161  qtrc(k,i,j,it) = conc_gas(ic) / dens(k,i,j) * 1.e-9_rp !mixing ratio [kg/kg]
1162  it = it + 1
1163  enddo
1164  enddo
1165  enddo
1166  enddo
1167 
1168  ccn(:,:,:) = 0.0_rp
1169  do j = js, je
1170  do i = is, ie
1171  do k = ks, ke
1172  !--- calculate super saturation of water
1173  call saturation_pres2qsat_liq( temp(k,i,j), pres(k,i,j), qdry(k,i,j), & ! [IN]
1174  qsat_tmp ) ! [OUT]
1175  ssliq = qv(k,i,j) / qsat_tmp - 1.0_rp
1176 
1177  call aerosol_activation( c_kappa, ssliq, temp(k,i,j), ia_m0, ia_m2, ia_m3, &
1178  n_atr, n_siz_max, n_kap_max, n_ctg, nsiz, n_kap, &
1179  d_ct, aerosol_procs, aerosol_activ )
1180 
1181  do is0 = 1, nsiz(ic_mix)
1182  ccn(k,i,j) = ccn(k,i,j) + aerosol_activ(ia_m0,is0,ik_out,ic_mix)
1183  enddo
1184  enddo
1185  enddo
1186  enddo
1187 
1188  return

References atmos_phy_ae_kajino13_setup(), scale_const::const_pi, scale_precision::dp, scale_atmos_grid_cartesc_index::ia, scale_atmos_grid_cartesc_index::ie, scale_atmos_grid_cartesc_index::is, scale_atmos_grid_cartesc_index::ja, scale_atmos_grid_cartesc_index::je, scale_atmos_grid_cartesc_index::js, scale_tracer::k, scale_atmos_grid_cartesc_index::ka, scale_atmos_grid_cartesc_index::ke, scale_atmos_grid_cartesc_index::ks, and scale_precision::rp.

Referenced by mod_mkinit::rect_setup().

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 1194 of file scale_atmos_phy_ae_kajino13.F90.

1194  integer, intent(in) :: KA, KS, KE
1195  integer, intent(in) :: IA, IS, IE
1196  integer, intent(in) :: JA, JS, JE
1197  integer, intent(in) :: QA_AE
1198 
1199  real(RP), intent(inout) :: QTRC(KA,IA,JA,QA_AE)
1200 
1201  integer :: k ,i, j
1202  integer :: ic, ik, is0
1203 
1204  !--- Negative fixer
1205  do j = js, je
1206  do i = is, ie
1207  do k = ks, ke
1208  do ic = 1, n_ctg !aerosol category
1209  do ik = 1, n_kap(ic) !kappa bin
1210  do is0 = 1, n_siz(ic) !size bin
1211  if (qtrc(k,i,j,it_procs2trans(ia_m0,is0,ik,ic)) < 0.0_rp .or. &
1212  qtrc(k,i,j,it_procs2trans(ia_m2,is0,ik,ic)) < 0.0_rp .or. &
1213  qtrc(k,i,j,it_procs2trans(ia_m3,is0,ik,ic)) < 0.0_rp .or. &
1214  qtrc(k,i,j,it_procs2trans(ia_ms,is0,ik,ic)) < 0.0_rp .or. &
1215  qtrc(k,i,j,it_procs2trans(ia_kp,is0,ik,ic)) < 0.0_rp ) then
1216  qtrc(k,i,j,it_procs2trans(1:n_atr,is0,ik,ic)) = 0.0_rp
1217  endif
1218  enddo
1219  enddo
1220  enddo
1221  enddo
1222  enddo
1223  enddo
1224 
1225  return

References atmos_phy_ae_kajino13_c_kappa, atmos_phy_ae_kajino13_h2so4dt, atmos_phy_ae_kajino13_logk_aenucl, atmos_phy_ae_kajino13_ocgasdt, scale_precision::dp, float(), scale_atmos_grid_cartesc_index::ie, scale_atmos_grid_cartesc_index::is, scale_atmos_grid_cartesc_index::je, scale_atmos_grid_cartesc_index::js, scale_tracer::k, 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().

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 52 of file scale_atmos_phy_ae_kajino13.F90.

52  character(len=H_SHORT), public, allocatable :: ATMOS_PHY_AE_kajino13_NAME(:)

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

◆ 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 53 of file scale_atmos_phy_ae_kajino13.F90.

53  character(len=H_MID) , public, allocatable :: ATMOS_PHY_AE_kajino13_DESC(:)

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

◆ 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 54 of file scale_atmos_phy_ae_kajino13.F90.

54  character(len=H_SHORT), public, allocatable :: ATMOS_PHY_AE_kajino13_UNIT(:)

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

◆ 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 56 of file scale_atmos_phy_ae_kajino13.F90.

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

Referenced by atmos_phy_ae_kajino13_negative_fixer(), and atmos_phy_ae_kajino13_setup().

◆ 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 57 of file scale_atmos_phy_ae_kajino13.F90.

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

Referenced by atmos_phy_ae_kajino13_negative_fixer(), and atmos_phy_ae_kajino13_setup().

◆ 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 59 of file scale_atmos_phy_ae_kajino13.F90.

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

Referenced by atmos_phy_ae_kajino13_negative_fixer(), and atmos_phy_ae_kajino13_setup().

◆ atmos_phy_ae_kajino13_flag_npf

logical, public scale_atmos_phy_ae_kajino13::atmos_phy_ae_kajino13_flag_npf = .false.

Definition at line 60 of file scale_atmos_phy_ae_kajino13.F90.

60  logical, public :: ATMOS_PHY_AE_KAJINO13_flag_npf = .false.

Referenced by atmos_phy_ae_kajino13_setup(), and atmos_phy_ae_kajino13_tendency().

◆ atmos_phy_ae_kajino13_flag_cond

logical, public scale_atmos_phy_ae_kajino13::atmos_phy_ae_kajino13_flag_cond = .true.

Definition at line 61 of file scale_atmos_phy_ae_kajino13.F90.

61  logical, public :: ATMOS_PHY_AE_KAJINO13_flag_cond = .true.

Referenced by atmos_phy_ae_kajino13_setup(), and atmos_phy_ae_kajino13_tendency().

◆ atmos_phy_ae_kajino13_flag_coag

logical, public scale_atmos_phy_ae_kajino13::atmos_phy_ae_kajino13_flag_coag = .true.

Definition at line 62 of file scale_atmos_phy_ae_kajino13.F90.

62  logical, public :: ATMOS_PHY_AE_KAJINO13_flag_coag = .true.

Referenced by atmos_phy_ae_kajino13_setup(), and atmos_phy_ae_kajino13_tendency().

◆ 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 63 of file scale_atmos_phy_ae_kajino13.F90.

63  logical, public :: ATMOS_PHY_AE_KAJINO13_flag_ccn_interactive = .true.

Referenced by atmos_phy_ae_kajino13_setup(), and atmos_phy_ae_kajino13_tendency().

◆ atmos_phy_ae_kajino13_flag_regeneration

logical, public scale_atmos_phy_ae_kajino13::atmos_phy_ae_kajino13_flag_regeneration = .true.

Definition at line 64 of file scale_atmos_phy_ae_kajino13.F90.

64  logical, public :: ATMOS_PHY_AE_KAJINO13_flag_regeneration = .true.

Referenced by atmos_phy_ae_kajino13_setup(), and atmos_phy_ae_kajino13_tendency().

◆ 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 65 of file scale_atmos_phy_ae_kajino13.F90.

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

Referenced by atmos_phy_ae_kajino13_setup(), and atmos_phy_ae_kajino13_tendency().

◆ 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 66 of file scale_atmos_phy_ae_kajino13.F90.

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

Referenced by atmos_phy_ae_kajino13_setup(), and atmos_phy_ae_kajino13_tendency().

◆ 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 67 of file scale_atmos_phy_ae_kajino13.F90.

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

Referenced by atmos_phy_ae_kajino13_negative_fixer(), and atmos_phy_ae_kajino13_setup().

scale_prc::prc_abort
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:350
scale_atmos_aerosol::n_ae
integer, parameter, public n_ae
Definition: scale_atmos_aerosol.F90:33
scale_file_history
module file_history
Definition: scale_file_history.F90:15
float
typedef float(real32_t)
scale_prc
module PROCESS
Definition: scale_prc.F90:11
scale_const
module CONSTANT
Definition: scale_const.F90:11
scale_const::const_pi
real(rp), parameter, public const_pi
pi
Definition: scale_const.F90:32
scale_atmos_aerosol
module atmosphere / aerosol
Definition: scale_atmos_aerosol.F90:12
scale_atmos_saturation
module atmosphere / saturation
Definition: scale_atmos_saturation.F90:12