SCALE-RM
Functions/Subroutines | Variables
scale_atmos_phy_ae_kajino13 Module Reference

module ATMOSPHERE / Physics Aerosol Microphysics More...

Functions/Subroutines

subroutine, public atmos_phy_ae_kajino13_config (AE_TYPE, QA, QS)
 Config. More...
 
subroutine, public atmos_phy_ae_kajino13_setup
 Setup. More...
 
subroutine, public atmos_phy_ae_kajino13 (QQA, DENS, MOMZ, MOMX, MOMY, RHOT, EMIT, NREG, QTRC, CN, CCN, RHOQ_t_AE)
 Aerosol Microphysics. More...
 
subroutine aerosol_zerochem (deltt, temp_k, pres_pa, super, flag_npf, flag_cond, flag_coag, aerosol_procs, conc_gas, emis_procs, emis_gas, aerosol_activ)
 
subroutine aerosol_nucleation (conc_h2so4, J_1nm)
 
subroutine, public atmos_phy_ae_kajino13_effectiveradius (Re, QTRC, RH)
 Calculate Effective Radius. More...
 
subroutine, public atmos_phy_ae_kajino13_mkinit (QTRC, CCN, DENS, RHOT, m0_init, dg_init, sg_init, d_min_inp, d_max_inp, k_min_inp, k_max_inp, n_kap_inp)
 

Variables

character(len=h_short), dimension(:), allocatable, target, public atmos_phy_ae_kajino13_name
 
character(len=h_mid), dimension(:), allocatable, target, public atmos_phy_ae_kajino13_desc
 
character(len=h_short), dimension(:), allocatable, target, public atmos_phy_ae_kajino13_unit
 
real(rp), dimension(n_ae), target, public atmos_phy_ae_kajino13_dens
 

Detailed Description

module ATMOSPHERE / Physics Aerosol Microphysics

Description
kajino13 code
Author
Team SCALE
History
  • 2015-03-25 (Y.Sato) [new]
NAMELIST
  • PARAM_TRACER_KAJINO13
    nametypedefault valuecomment
    AE_CTG integer 1
    NASIZ integer, dimension(3)
    NAKAP integer, dimension(3)

  • PARAM_ATMOS_PHY_AE_KAJINO13
    nametypedefault valuecomment
    H2SO4DT real(RP) 5.E-6_RP h2so4 production rate (Temporal) [ug/m3/s]
    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
trim(ofilename) Total number mixing ratio of emitted aerosol num/kg total_emit_aerosol_number

Function/Subroutine Documentation

◆ atmos_phy_ae_kajino13_config()

subroutine, public scale_atmos_phy_ae_kajino13::atmos_phy_ae_kajino13_config ( character(len=*), intent(in)  AE_TYPE,
integer, intent(out)  QA,
integer, intent(out)  QS 
)

Config.

Definition at line 188 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_stdio::io_fid_conf, scale_stdio::io_fid_log, scale_stdio::io_fid_nml, scale_stdio::io_l, scale_stdio::io_nml, scale_process::prc_mpistop(), and scale_tracer::tracer_regist().

Referenced by scale_atmos_phy_ae::atmos_phy_ae_config().

188  use scale_process, only: &
190  use scale_tracer, only: &
192  implicit none
193 
194  character(len=*), intent(in) :: AE_TYPE
195  integer, intent(out) :: QA
196  integer, intent(out) :: QS
197 
198  integer, allocatable :: aero_idx(:,:,:,:)
199  integer :: n_kap_max, n_siz_max, ncat_max
200  integer :: NASIZ(3), NAKAP(3)
201  character(len=H_SHORT) :: attribute, catego, aunit
202 
203  namelist / param_tracer_kajino13 / &
204  ae_ctg, &
205  nasiz, &
206  nakap
207 
208  integer :: m, ierr, ik, ic, ia0, is0
209 
210  !---------------------------------------------------------------------------
211 
212  if( io_l ) write(io_fid_log,*)
213  if( io_l ) write(io_fid_log,*) '++++++ Module[Aerosol Tracer] / Categ[ATMOS PHYSICS] / Origin[SCALElib]'
214  if( io_l ) write(io_fid_log,*) '*** Tracers for Kajino(2013) scheme'
215 
216  if ( ae_type /= 'KAJINO13' .AND. ae_type /= 'NONE' ) then
217  write(*,*) 'xxx ATMOS_PHY_AE_TYPE is not KAJINO13. Check!'
218  call prc_mpistop
219  endif
220 
221  ncat_max = max( ic_mix, ic_sea, ic_dus )
222 
223  nasiz(:) = 64
224  nakap(:) = 1
225 
226  rewind(io_fid_conf)
227  read(io_fid_conf,nml=param_tracer_kajino13,iostat=ierr)
228 
229  if( ierr < 0 ) then !--- missing
230  if( io_l ) write(io_fid_log,*) '*** Not found namelist. Default used.'
231  elseif( ierr > 0 ) then !--- fatal error
232  write(*,*) 'xxx Not appropriate names in namelist PARAM_TRACER_KAJINO13, Check!'
233  call prc_mpistop
234  end if
235 
236  if( io_nml ) write(io_fid_nml,nml=param_tracer_kajino13)
237 
238  if( ae_ctg > ncat_max ) then
239  write(*,*) 'xxx AE_CTG should be smaller than', ncat_max+1, 'stop'
240  call prc_mpistop
241  endif
242 
243  allocate( nsiz(ae_ctg) )
244  allocate( nkap(ae_ctg) )
245 
246  nkap(1:ae_ctg) = nakap(1:ae_ctg)
247  nsiz(1:ae_ctg) = nasiz(1:ae_ctg)
248 
249  if( maxval( nkap ) /= 1 .OR. minval( nkap ) /= 1 ) then
250  write(*,*) 'xxx NKAP(:) /= 1 is not supported now, Stop!'
251  call prc_mpistop
252  end if
253 
254 ! do ia0 = 1, N_ATR
255  do ic = 1, ae_ctg
256  do ik = 1, nkap(ic)
257  do is0 = 1, nsiz(ic)
258  qa_ae = qa_ae + n_atr
259  enddo
260  enddo
261  enddo
262  qa_ae = qa_ae + gas_ctg
263 
264  allocate( atmos_phy_ae_kajino13_name(qa_ae) )
265  allocate( atmos_phy_ae_kajino13_desc(qa_ae) )
266  allocate( atmos_phy_ae_kajino13_unit(qa_ae) )
267 
268  n_siz_max = 0
269  n_kap_max = 0
270  do ic = 1, ae_ctg
271  n_siz_max = max(n_siz_max, nsiz(ic))
272  n_kap_max = max(n_kap_max, nkap(ic))
273  enddo
274 
275  allocate( aero_idx(n_atr,ae_ctg,n_kap_max,n_siz_max) )
276  m = 0
277 ! do ia0 = 1, N_ATR
278  do ic = 1, ae_ctg
279  do ik = 1, nkap(ic)
280  do is0 = 1, nsiz(ic)
281  do ia0 = 1, n_atr
282  m = m+1
283  aero_idx(ia0,ic,ik,is0) = m
284  enddo
285  enddo
286  enddo
287  enddo
288 
289  !-----------------------------------------------------------------------------
290  !
291  !++ calculate each category and aerosol
292  !
293  !-----------------------------------------------------------------------------
294  ic = qa_ae-gas_ctg+ig_h2so4
295  write(atmos_phy_ae_kajino13_unit(ic),'(a)') 'kg/kg'
296  ic = qa_ae-gas_ctg+ig_cgas
297  write(atmos_phy_ae_kajino13_unit(ic),'(a)') 'kg/kg'
298 
299 ! do ia0 = 1, N_ATR
300 ! if( ia0 == 1 ) then
301 ! write(attribute,'(a)') "Number"
302 ! write(aunit,'(a)') "num/kg"
303 ! elseif( ia0 == 2 ) then
304 ! write(attribute,'(a)') "Section"
305 ! write(aunit,'(a)') "m2/kg"
306 ! elseif( ia0 == 3 ) then
307 ! write(attribute,'(a)') "Volume"
308 ! write(aunit,'(a)') "m3/kg"
309 ! elseif( ia0 == 4 ) then
310 ! write(attribute,'(a)') "Mass"
311 ! write(aunit,'(a)') "kg/kg"
312 ! elseif( ia0 == 5 ) then
313 ! write(attribute,'(a)') "kpXmass"
314 ! write(aunit,'(a)') "kg/kg"
315 ! endif
316  do ic = 1, ae_ctg !aerosol category
317  do ik = 1, nkap(ic) !kappa bin
318  do is0 = 1, nsiz(ic)
319  do ia0 = 1, n_atr
320  if( ia0 == 1 ) then
321  write(attribute,'(a)') "Number"
322  write(aunit,'(a)') "num/kg"
323  elseif( ia0 == 2 ) then
324  write(attribute,'(a)') "Section"
325  write(aunit,'(a)') "m2/kg"
326  elseif( ia0 == 3 ) then
327  write(attribute,'(a)') "Volume"
328  write(aunit,'(a)') "m3/kg"
329  elseif( ia0 == 4 ) then
330  write(attribute,'(a)') "Mass"
331  write(aunit,'(a)') "kg/kg"
332  elseif( ia0 == 5 ) then
333  write(attribute,'(a)') "kXm"
334  write(aunit,'(a)') "kg/kg"
335  endif
336  if( ic == ic_mix ) then
337  write(catego,'(a)') "Sulf_"
338  elseif( ic == ic_sea ) then
339  write(catego,'(a)') "Salt_"
340  elseif( ic == ic_dus ) then
341  write(catego,'(a)') "Dust_"
342  endif
343  write(atmos_phy_ae_kajino13_unit(aero_idx(ia0,ic,ik,is0)),'(a)') trim(aunit)
344  write(atmos_phy_ae_kajino13_name(aero_idx(ia0,ic,ik,is0)),'(a,a,i0)') trim(catego), trim(attribute), is0
345  write(atmos_phy_ae_kajino13_desc(aero_idx(ia0,ic,ik,is0)),'(a,a,a,i0)') trim(attribute), ' mixing radio of ', trim(catego), is0
346  enddo
347  enddo
348  enddo
349  enddo
350  ic = qa_ae-gas_ctg+ig_h2so4
351  write(atmos_phy_ae_kajino13_name(ic),'(a)') 'H2SO4_Gas'
352  ic = qa_ae-gas_ctg+ig_cgas
353  write(atmos_phy_ae_kajino13_name(ic),'(a)') 'Condensable_GAS'
354 
355  ic = qa_ae-gas_ctg+ig_h2so4
356  write(atmos_phy_ae_kajino13_desc(ic),'(a)') 'Mixing ratio of H2SO4 Gas'
357  ic = qa_ae-gas_ctg+ig_cgas
358  write(atmos_phy_ae_kajino13_desc(ic),'(a)') 'Mixing ratio of Condensable GAS'
359 
360  deallocate(aero_idx)
361 
362  call tracer_regist( qs, & ! [OUT]
363  qa_ae, & ! [IN]
367 
368  qa = qa_ae
369  qaes = qs
370  qaee = qs + qa_ae - 1
371 
372  return
subroutine, public prc_mpistop
Abort MPI.
character(len=h_mid), dimension(:), allocatable, target, public atmos_phy_ae_kajino13_desc
character(len=h_short), dimension(:), allocatable, target, public atmos_phy_ae_kajino13_name
module TRACER
character(len=h_short), dimension(:), allocatable, target, public atmos_phy_ae_kajino13_unit
module PROCESS
subroutine, public tracer_regist(QS, NQ, NAME, DESC, UNIT, CV, CP, R, ADVC, MASS)
Regist tracer.
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 378 of file scale_atmos_phy_ae_kajino13.F90.

References atmos_phy_ae_kajino13_dens, float(), scale_stdio::io_fid_conf, scale_stdio::io_fid_log, scale_stdio::io_fid_nml, scale_stdio::io_l, scale_stdio::io_nml, scale_process::prc_mpistop(), and scale_time::time_dtsec_atmos_phy_ae.

Referenced by scale_atmos_phy_ae::atmos_phy_ae_config().

378  use scale_process, only: &
380  use scale_time, only: &
382  implicit none
383 
384  real(RP), allocatable :: d_min_inp(:)
385  real(RP), allocatable :: d_max_inp(:)
386  real(RP), allocatable :: k_min_inp(:)
387  real(RP), allocatable :: k_max_inp(:)
388  integer , allocatable :: n_kap_inp(:)
389 
390  real(RP), parameter :: d_min_def = 1.e-9_rp ! default lower bound of 1st size bin
391  real(RP), parameter :: d_max_def = 1.e-5_rp ! upper bound of last size bin
392  integer , parameter :: n_kap_def = 1 ! number of kappa bins
393  real(RP), parameter :: k_min_def = 0.e0_rp ! lower bound of 1st kappa bin
394  real(RP), parameter :: k_max_def = 1.e0_rp ! upper bound of last kappa bin
395 
396  namelist / param_atmos_phy_ae_kajino13 / &
397  h2so4dt, &
398  ocgasdt, &
399 ! c_ratio, &
400  c_kappa, &
401  t_npf, &
402 ! d_min_inp, &
403 ! d_max_inp, &
404 ! k_min_inp, &
405 ! k_max_inp, &
406 ! n_kap_inp, &
407  flag_npf, &
408  flag_cond, &
409  flag_coag, &
410  flag_ccn_interactive, &
411  flag_regeneration, &
412  dg_reg, &
413  sg_reg, &
414  logk_aenucl, &
415  nbins_out
416 
417  integer :: it, ierr
418  !---------------------------------------------------------------------------
419 
420  if( io_l ) write(io_fid_log,*)
421  if( io_l ) write(io_fid_log,*) '++++++ Module[Aerosol] / Categ[ATMOS PHYSICS] / Origin[SCALElib]'
422  if( io_l ) write(io_fid_log,*) '*** Kajino(2013) scheme'
423 
424  !--- setup parameter
425  pi6 = pi / 6._rp ! pi/6
426  sixpi = 6._rp / pi ! 6/pi
427  forpi = 4._rp / pi ! 4/pi
428 
429 
430 ! deltt = TIME_DTSEC_ATMOS_PHY_AE
431  n_ctg = ae_ctg
432  allocate( rnum_out(nbins_out) )
433  allocate( n_siz(n_ctg) )
434  allocate( d_min(n_ctg) )
435  allocate( d_max(n_ctg) )
436  allocate( n_kap(n_ctg) )
437  allocate( k_min(n_ctg) )
438  allocate( k_max(n_ctg) )
439  allocate( d_min_inp(n_ctg) )
440  allocate( d_max_inp(n_ctg) )
441  allocate( n_kap_inp(n_ctg) )
442  allocate( k_min_inp(n_ctg) )
443  allocate( k_max_inp(n_ctg) )
444  allocate( ctg_name(n_ctg) )
445 
446  n_siz(1:n_ctg) = nsiz(1:n_ctg) ! number of size bins
447  d_min(1:n_ctg) = d_min_def ! lower bound of 1st size bin
448  d_max(1:n_ctg) = d_max_def ! upper bound of last size bin
449  n_kap(1:n_ctg) = n_kap_def ! number of kappa bins
450  k_min(1:n_ctg) = k_min_def ! lower bound of 1st kappa bin
451  k_max(1:n_ctg) = k_max_def ! upper bound of last kappa bin
452 
453  do it = 1, n_ctg
454  if( n_ctg == 1 ) then
455  write(ctg_name(it),'(a)') "Sulfate"
456  elseif( n_ctg == 2 ) then
457  write(ctg_name(it),'(a)') "Seasalt"
458  elseif( n_ctg == 3 ) then
459  write(ctg_name(it),'(a)') "Dust"
460  endif
461  enddo
462 
463  !--- read namelist
464  rewind(io_fid_conf)
465  read(io_fid_conf,nml=param_atmos_phy_ae_kajino13,iostat=ierr)
466  if( ierr < 0 ) then !--- missing
467  if( io_l ) write(io_fid_log,*) '*** Not found namelist. Default used.'
468  elseif( ierr > 0 ) then !--- fatal error
469  write(*,*) 'xxx Not appropriate names in namelist PARAM_ATMOS_PHY_AE_KAJINO13. Check!'
470  call prc_mpistop
471  endif
472  if( io_nml ) write(io_fid_nml,nml=param_atmos_phy_ae_kajino13)
473 
474  !--- now only the default setting is supported
475 ! n_siz(1:n_ctg) = NSIZ(1:n_ctg) ! number of size bins
476 ! d_min(1:n_ctg) = d_min_inp(1:n_ctg) ! lower bound of 1st size bin
477 ! d_max(1:n_ctg) = d_max_inp(1:n_ctg) ! upper bound of last size bin
478 ! n_kap(1:n_ctg) = n_kap_inp(1:n_ctg) ! number of kappa bins
479 ! k_min(1:n_ctg) = k_min_inp(1:n_ctg) ! lower bound of 1st kappa bin
480 ! k_max(1:n_ctg) = k_max_inp(1:n_ctg) ! upper bound of last kappa bin
481 
482  !--- diagnose parameters (n_trans, n_siz_max, n_kap_max)
483  n_trans = 0
484  n_siz_max = 0
485  n_kap_max = 0
486  do ic = 1, n_ctg
487  n_trans = n_trans + n_siz(ic) * n_kap(ic) * n_atr
488  n_siz_max = max(n_siz_max, n_siz(ic))
489  n_kap_max = max(n_kap_max, n_kap(ic))
490  enddo
491 
492  !--- bin settings
493  allocate(d_lw(n_siz_max,n_ctg))
494  allocate(d_ct(n_siz_max,n_ctg))
495  allocate(d_up(n_siz_max,n_ctg))
496  allocate(k_lw(n_kap_max,n_ctg))
497  allocate(k_ct(n_kap_max,n_ctg))
498  allocate(k_up(n_kap_max,n_ctg))
499  d_lw(:,:) = 0.0_rp
500  d_ct(:,:) = 0.0_rp
501  d_up(:,:) = 0.0_rp
502  k_lw(:,:) = 0.0_rp
503  k_ct(:,:) = 0.0_rp
504  k_up(:,:) = 0.0_rp
505 
506  do ic = 1, n_ctg
507  dlogd = (log(d_max(ic)) - log(d_min(ic)))/float(n_siz(ic))
508  do is0 = 1, n_siz(ic) !size bin
509  d_lw(is0,ic) = exp(log(d_min(ic))+dlogd* float(is0-1) )
510  d_ct(is0,ic) = exp(log(d_min(ic))+dlogd*(float(is0)-0.5_rp))
511  d_up(is0,ic) = exp(log(d_min(ic))+dlogd* float(is0) )
512  enddo !is (1:n_siz(ic))
513 
514  dk = (k_max(ic) - k_min(ic))/float(n_kap(ic))
515  do ik = 1, n_kap(ic) !size bin
516  k_lw(ik,ic) = k_min(ic) + dk * float(ik-1)
517  k_ct(ik,ic) = k_min(ic) + dk *(float(ik)-0.5_rp)
518  k_up(ik,ic) = k_min(ic) + dk * float(ik)
519  enddo !ik (1:n_kap(ic))
520 
521  enddo !ic (1:n_ctg)
522 
523 !find size bin of regenerated aerosols
524  do is0 = 1, n_siz(ic_mix)
525  if (dg_reg >= d_lw(is0,ic_mix) .AND. &
526  dg_reg < d_up(is0,ic_mix) ) then
527  is0_reg = is0
528  endif !d_lw < dg_reg < d_up
529  enddo
530 
531  !--- coagulation rule
532  ! [ NOTE: current version has one category and single
533  ! hygroscopicity bins and thus does not consider inter-category
534  ! nor inter-hygroscopicity-section coagulation ]
535  mcomb = 0
536  do ic = 1, n_ctg !=1
537  do ik = 1, n_kap(ic) !=1
538  do is2 = 1, n_siz(ic)
539  do is1 = 1, n_siz(ic)
540  if ( d_ct(is2,ic) >= d_ct(is1,ic) ) then
541  mcomb = mcomb + 1
542  endif !d_ct(is2) >= d_ct(is1)
543  enddo !is1 (1:n_siz(ic) )
544  enddo !is2 (1:n_siz(ic) )
545  enddo !ik(1:n_kap(ic))
546  enddo !ic(1:n_ctg)
547 
548  allocate(is_i(mcomb))
549  allocate(is_j(mcomb))
550  allocate(is_k(mcomb))
551  allocate(ik_i(mcomb))
552  allocate(ik_j(mcomb))
553  allocate(ik_k(mcomb))
554  allocate(ic_i(mcomb))
555  allocate(ic_j(mcomb))
556  allocate(ic_k(mcomb))
557 
558  mc = 0
559  do ic = 1, n_ctg !=1
560  do ik = 1, n_kap(ic) !=1
561  do is2 = 1, n_siz(ic)
562  do is1 = 1, n_siz(ic)
563  if ( d_ct(is2,ic) >= d_ct(is1,ic) ) then
564  mc = mc + 1
565  is_i(mc) = is1
566  ik_i(mc) = ik
567  ic_i(mc) = ic
568  is_j(mc) = is2
569  ik_j(mc) = ik
570  ic_j(mc) = ic
571  is_k(mc) = is2
572  ik_k(mc) = ik
573  ic_k(mc) = ic
574  endif !d_ct(is2) >= d_ct(is1)
575  enddo !is1 (1:n_siz(ic) )
576  enddo !is2 (1:n_siz(ic) )
577  enddo !ik(1:n_kap(ic))
578  enddo !ic(1:n_ctg)
579 
580  !--- gas concentration
581 ! conc_h2so4 = 0.0_RP
582 ! conc_cgas = 0.0_RP
583 
584  allocate( it_procs2trans(n_atr,n_siz_max,n_kap_max,n_ctg) )
585  allocate( ia_trans2procs(n_trans) )
586  allocate( is_trans2procs(n_trans) )
587  allocate( ik_trans2procs(n_trans) )
588  allocate( ic_trans2procs(n_trans) )
589 
590  it_procs2trans(:,:,:,:)= -999
591  ia_trans2procs(:) = 0
592  is_trans2procs(:) = 0
593  ik_trans2procs(:) = 0
594  ic_trans2procs(:) = 0
595 
596  !--- get pointer for trans2procs, procs2trans
597  it = 0
598  do ic = 1, n_ctg !aerosol category
599  do ik = 1, n_kap(ic) !kappa bin
600  do is0 = 1, n_siz(ic) !size bin
601  do ia0 = 1, n_atr !attributes
602  it = it + 1
603  it_procs2trans(ia0,is0,ik,ic)= it
604  ia_trans2procs(it) = ia0
605  is_trans2procs(it) = is0
606  ik_trans2procs(it) = ik
607  ic_trans2procs(it) = ic
608  enddo !ia (1:n_atr_prog )
609  enddo !is (1:n_siz(ic) )
610  enddo !ik (1:n_kap(ic) )
611  enddo !ic (1:n_ctg )
612 
613  atmos_phy_ae_kajino13_dens(:) = rhod_ae
614 
615  return
subroutine, public prc_mpistop
Abort MPI.
real(rp), dimension(n_ae), target, public atmos_phy_ae_kajino13_dens
typedef float(real32_t)
module TIME
Definition: scale_time.F90:15
module PROCESS
real(dp), public time_dtsec_atmos_phy_ae
time interval of physics(aerosol ) [sec]
Definition: scale_time.F90:46
Here is the call graph for this function:
Here is the caller graph for this function:

◆ atmos_phy_ae_kajino13()

subroutine, public scale_atmos_phy_ae_kajino13::atmos_phy_ae_kajino13 ( integer, intent(in)  QQA,
real(rp), dimension(ka,ia,ja), intent(inout)  DENS,
real(rp), dimension(ka,ia,ja), intent(inout)  MOMZ,
real(rp), dimension(ka,ia,ja), intent(inout)  MOMX,
real(rp), dimension(ka,ia,ja), intent(inout)  MOMY,
real(rp), dimension(ka,ia,ja), intent(inout)  RHOT,
real(rp), dimension(ka,ia,ja,qqa), intent(inout)  EMIT,
real(rp), dimension(ka,ia,ja), intent(in)  NREG,
real(rp), dimension(ka,ia,ja,qa), intent(inout)  QTRC,
real(rp), dimension(ka,ia,ja), intent(out)  CN,
real(rp), dimension(ka,ia,ja), intent(out)  CCN,
real(rp), dimension(ka,ia,ja,qa), intent(inout)  RHOQ_t_AE 
)

Aerosol Microphysics.

Definition at line 633 of file scale_atmos_phy_ae_kajino13.F90.

References aerosol_zerochem(), scale_const::const_cpdry, scale_const::const_cvdry, scale_const::const_pre00, scale_const::const_rdry, scale_const::const_rvap, scale_atmos_hydrometeor::i_qv, scale_grid_index::ie, scale_stdio::io_fid_log, scale_stdio::io_l, scale_grid_index::is, scale_grid_index::je, scale_grid_index::js, scale_grid_index::ke, scale_grid_index::ks, scale_tracer::qa, scale_time::time_dtsec_atmos_phy_ae, scale_tracer::tracer_cv, scale_tracer::tracer_mass, and scale_tracer::tracer_r.

Referenced by scale_atmos_phy_ae::atmos_phy_ae_config().

633  use scale_grid_index
634  use scale_tracer
635  use scale_const, only: &
636  const_cpdry, &
637  const_cvdry, &
638  const_rvap, &
639  const_pre00, &
640  const_rdry
641  use scale_atmos_hydrometeor, only: &
642  i_qv
643  use scale_atmos_saturation, only: &
644  pres2qsat_liq => atmos_saturation_pres2qsat_liq
645  use scale_time, only: &
647  use scale_history, only: &
648  hist_in
649  implicit none
650  integer, intent(in) :: QQA
651  real(RP), intent(inout) :: DENS(KA,IA,JA)
652  real(RP), intent(inout) :: MOMZ(KA,IA,JA)
653  real(RP), intent(inout) :: MOMX(KA,IA,JA)
654  real(RP), intent(inout) :: MOMY(KA,IA,JA)
655  real(RP), intent(inout) :: RHOT(KA,IA,JA)
656  real(RP), intent(inout) :: EMIT(KA,IA,JA,QQA)
657  real(RP), intent(in) :: NREG(KA,IA,JA)
658  real(RP), intent(inout) :: QTRC(KA,IA,JA,QA)
659  real(RP), intent(out) :: CN(KA,IA,JA)
660  real(RP), intent(out) :: CCN(KA,IA,JA)
661  real(RP), intent(inout) :: RHOQ_t_AE(KA,IA,JA,QA)
662 
663  real(RP) :: QTRC0(KA,IA,JA,QA)
664  real(RP) :: QTRC1(KA,IA,JA,QA)
665 
666  !--- local
667  real(RP) :: pres_ae(KA,IA,JA)
668  real(RP) :: temp_ae(KA,IA,JA)
669  real(RP) :: qv_ae(KA,IA,JA)
670  real(RP) :: ssliq_ae(KA,IA,JA)
671  real(RP) :: th(KA,IA,JA)
672  real(RP) :: q(KA,IA,JA,QA)
673  real(RP) :: qdry(KA,IA,JA)
674  real(RP) :: rrhog(KA,IA,JA)
675  real(RP) :: cva(KA,IA,JA)
676  real(RP) :: cpa(KA,IA,JA)
677  real(RP) :: t_ccn, t_cn
678  real(RP) :: Rmoist
679  real(RP) :: qsat_tmp
680  real(RP) :: m0_reg, m2_reg, m3_reg !regenerated aerosols [m^k/m3]
681  real(RP) :: ms_reg !regenerated aerosol mass [ug/m3]
682  real(RP) :: reg_factor_m2,reg_factor_m3 !to save cpu time for moment conversion
683  !--- aerosol variables
684  real(RP),allocatable :: aerosol_procs(:,:,:,:) !(n_atr,n_siz_max,n_kap_max,n_ctg)
685  real(RP),allocatable :: aerosol_activ(:,:,:,:) !(n_atr,n_siz_max,n_kap_max,n_ctg)
686  real(RP),allocatable :: emis_procs(:,:,:,:) !(n_atr,n_siz_max,n_kap_max,n_ctg)
687  real(RP),allocatable :: emis_gas(:) !emission of gas
688  real(RP) :: total_aerosol_mass(KA,IA,JA,n_ctg)
689  real(RP) :: total_aerosol_number(KA,IA,JA,n_ctg)
690  real(RP) :: total_emit_aerosol_mass(KA,IA,JA,n_ctg)
691  real(RP) :: total_emit_aerosol_number(KA,IA,JA,n_ctg)
692  !--- gas
693 ! real(RP),allocatable :: conc_h2so4(:,:,:,:) !concentration [ug/m3]
694 ! real(RP) :: conc_gas(KA,IA,JA,GAS_CTG) !concentration [ug/m3]
695  real(RP) :: conc_gas(GAS_CTG) !concentration [ug/m3]
696  character(len=H_LONG) :: ofilename
697  integer :: i, j, k, iq, it
698 
699  if( io_l ) write(io_fid_log,*) '*** Atmos physics step: Aerosol(kajino13)'
700 
701  !--- Negative fixer
702  do j = js, je
703  do i = is, ie
704  do k = ks, ke
705  do ic = 1, n_ctg !aerosol category
706  do ik = 1, n_kap(ic) !kappa bin
707  do is0 = 1, n_siz(ic) !size bin
708  if (qtrc(k,i,j,qaes-1+it_procs2trans(ia_m0,is0,ik,ic)) < 0.0_rp .or. &
709  qtrc(k,i,j,qaes-1+it_procs2trans(ia_m2,is0,ik,ic)) < 0.0_rp .or. &
710  qtrc(k,i,j,qaes-1+it_procs2trans(ia_m3,is0,ik,ic)) < 0.0_rp .or. &
711  qtrc(k,i,j,qaes-1+it_procs2trans(ia_ms,is0,ik,ic)) < 0.0_rp .or. &
712  qtrc(k,i,j,qaes-1+it_procs2trans(ia_kp,is0,ik,ic)) < 0.0_rp ) then
713  qtrc(k,i,j,qaes-1+it_procs2trans(1:n_atr,is0,ik,ic)) = 0.0_rp
714  endif
715  enddo
716  enddo
717  enddo
718  enddo
719  enddo
720  enddo
721 
722  do iq = 1, qa
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  allocate( aerosol_procs(n_atr,n_siz_max,n_kap_max,n_ctg) )
733  allocate( aerosol_activ(n_atr,n_siz_max,n_kap_max,n_ctg) )
734  allocate( emis_procs(n_atr,n_siz_max,n_kap_max,n_ctg) )
735  allocate( emis_gas(gas_ctg) )
736 
737  reg_factor_m2 = dg_reg**2._rp * exp( 2.0_rp *(log(sg_reg)**2._rp)) !m0_reg to m2_reg
738  reg_factor_m3 = dg_reg**3._rp * exp( 4.5_rp *(log(sg_reg)**2._rp)) !m0_reg to m3_reg
739 
740  !--- convert SCALE variable to zerochem variable
741 
742  aerosol_procs(:,:,:,:) = 0.0_rp
743  aerosol_activ(:,:,:,:) = 0.0_rp
744  emis_procs(:,:,:,:) = 0.0_rp
745  emis_gas(:) = 0.0_rp
746  pres_ae(:,:,:) = 0.0_rp
747  temp_ae(:,:,:) = 0.0_rp
748  qv_ae(:,:,:) = 0.0_rp
749 
750  do j = js, je
751  do i = is, ie
752 
753  do k = ks, ke
754  rrhog(k,i,j) = 1.0_rp / dens(k,i,j)
755  th(k,i,j) = rhot(k,i,j) * rrhog(k,i,j)
756  enddo
757  do k = ks, ke
758  calc_qdry( qdry(k,i,j), qtrc, tracer_mass, k, i, j, iq )
759  enddo
760  do k = ks, ke
761  calc_cv( cva(k,i,j), qdry(k,i,j), qtrc, k, i, j, iq, const_cvdry, tracer_cv )
762  enddo
763  do k = ks, ke
764  calc_r( rmoist, qdry(k,i,j), qtrc, k, i, j, iq, const_rdry, tracer_r )
765  cpa(k,i,j) = cva(k,i,j) + rmoist
766  calc_pre( pres_ae(k,i,j), dens(k,i,j), th(k,i,j), rmoist, cpa(k,i,j), const_pre00 )
767  temp_ae(k,i,j) = pres_ae(k,i,j) / ( dens(k,i,j) * rmoist )
768  if ( i_qv > 0 ) then
769  qv_ae(k,i,j) = qtrc(k,i,j,i_qv)
770  !--- calculate super saturation of water
771  call pres2qsat_liq( qsat_tmp,temp_ae(k,i,j),pres_ae(k,i,j) )
772  ssliq_ae(k,i,j) = qv_ae(k,i,j)/qsat_tmp - 1.0_rp
773  else
774  ssliq_ae(k,i,j) = - 1.0_rp
775  end if
776  enddo
777 
778  enddo
779  enddo
780 
781  ! tiny number, tiny mass
782  do j = js, je
783  do i = is, ie
784  do k = ks, ke
785  do ic = 1, n_ctg !aerosol category
786  do ik = 1, n_kap(ic) !kappa bin
787  do is0 = 1, n_siz(ic) !size bin
788  if (qtrc0(k,i,j,qaes-1+it_procs2trans(ia_m0,is0,ik,ic))*dens(k,i,j) < cleannumber) then
789  do ia0 = 1, n_atr
790  qtrc0(k,i,j,qaes-1+it_procs2trans(ia0,is0,ik,ic)) = 0._rp !to save cpu time and avoid underflow
791  enddo !ia0 (1:n_atr )
792  endif
793  enddo !is (1:n_siz(ic) )
794  enddo !ik (1:n_kap(ic) )
795  enddo !ic (1:n_ctg )
796  enddo
797  enddo
798  enddo
799 
800  do iq = 1, qa
801  do j = js, je
802  do i = is, ie
803  do k = ks, ke
804  qtrc1(k,i,j,iq) = qtrc0(k,i,j,iq) ! save
805  enddo
806  enddo
807  enddo
808  enddo
809 
810  !---- Calculate aerosol processs
811  cn(:,:,:) = 0.0_rp
812  ccn(:,:,:) = 0.0_rp
813  do k = ks, ke
814  do j = js, je
815  do i = is, ie
816 
817  !--- aerosol_trans at initial time
818  ! [xx/kg] -> [xx/m3]
819  do it = 1, n_trans
820  aerosol_procs(ia_trans2procs(it), &
821  is_trans2procs(it), &
822  ik_trans2procs(it), &
823  ic_trans2procs(it)) = qtrc1(k,i,j,qaes-1+it)*dens(k,i,j)
824  emis_procs(ia_trans2procs(it), &
825  is_trans2procs(it), &
826  ik_trans2procs(it), &
827  ic_trans2procs(it)) = emit(k,i,j,it)*dens(k,i,j)
828  enddo !it(1:n_trans)
829  ! mixing ratio [kg/kg] -> concentration [ug/m3]
830  conc_gas(1:gas_ctg) &
831  = qtrc1(k,i,j,qaee-gas_ctg+1:qaee-gas_ctg+gas_ctg)*dens(k,i,j)*1.e+9_rp
832 
833  emis_gas(1:gas_ctg) = emit(k,i,j,qa_ae-gas_ctg+ig_h2so4:qa_ae-gas_ctg+ig_cgas)
834 
835  call aerosol_zerochem( &
836  time_dtsec_atmos_phy_ae, & !--- in
837  temp_ae(k,i,j), & !--- in
838  pres_ae(k,i,j), & !--- in
839  ssliq_ae(k,i,j), & !--- in
840  flag_npf, flag_cond, flag_coag,& !--- in
841  aerosol_procs, & !--- inout
842  conc_gas, & !--- inout
843  emis_procs, & !--- out
844  emis_gas, & !--- out
845  aerosol_activ ) !--- out
846 
847 ! aerosol loss due to activation to cloud droplets
848  if (flag_ccn_interactive) then
849  do is0 = 1, n_siz(ic_mix)
850  do ia0 = 1, n_atr !attributes
851  aerosol_activ(ia0,is0,ik_out,ic_mix) = min(max(0._rp, aerosol_activ(ia0,is0,ik_out,ic_mix)), &
852  aerosol_procs(ia0,is0,ik_out,ic_mix) )
853  aerosol_procs(ia0,is0,ik_out,ic_mix) = &
854  aerosol_procs(ia0,is0,ik_out,ic_mix) - aerosol_activ(ia0,is0,ik_out,ic_mix)
855  enddo
856  enddo
857  endif !flag_ccn_interactive
858 
859 ! aerosol regeneration due to evaporation of cloud droplets
860 ! using prescribed size parameters and to internal mixture category (ic_mix)
861  if (flag_regeneration) then
862  m0_reg = nreg(k,i,j) !#/m3
863 ! m2_reg = m0_reg * dg_reg**2._RP * exp( 2.0_RP *(log(sg_reg)**2._RP)) !m2/m3
864 ! m3_reg = m0_reg * dg_reg**3._RP * exp( 4.5_RP *(log(sg_reg)**2._RP)) !m3/m3
865  m2_reg = m0_reg * reg_factor_m2 !m2/m3
866  m3_reg = m0_reg * reg_factor_m3 !m3/m3
867  ms_reg = m3_reg * pi6 * conv_vl_ms !ug/m3
868  aerosol_procs(ia_m0,is0_reg,ik_out,ic_mix) = &
869  aerosol_procs(ia_m0,is0_reg,ik_out,ic_mix) + m0_reg !#/m3
870  aerosol_procs(ia_m2,is0_reg,ik_out,ic_mix) = &
871  aerosol_procs(ia_m2,is0_reg,ik_out,ic_mix) + m2_reg !m2/m3
872  aerosol_procs(ia_m3,is0_reg,ik_out,ic_mix) = &
873  aerosol_procs(ia_m3,is0_reg,ik_out,ic_mix) + m3_reg !m3/m3
874  aerosol_procs(ia_ms,is0_reg,ik_out,ic_mix) = &
875  aerosol_procs(ia_ms,is0_reg,ik_out,ic_mix) + ms_reg !ug/m3
876 ! additional attirbute to be added (ia_kp)
877  endif !flag_regeneration
878 
879 ! diagnosed variables
880  do is0 = 1, n_siz(ic_mix)
881  ccn(k,i,j) = ccn(k,i,j) + aerosol_activ(ia_m0,is0,ik_out,ic_mix)
882  cn(k,i,j) = cn(k,i,j) + aerosol_procs(ia_m0,is0,ik_out,ic_mix)
883  enddo
884 
885 ! call trans_ccn(aerosol_procs, aerosol_activ, t_ccn, t_cn, &
886 ! n_ctg, n_kap_max, n_siz_max, N_ATR, &
887 ! ic_mix, ia_m0, ia_m2, ia_m3, ik_out, n_siz, &
888 ! rnum_out, nbins_out)
889 
890 ! CN(k,i,j) = t_cn
891 ! CCN(k,i,j) = t_ccn
892 
893  ! [xx/m3] -> [xx/kg]
894  do ic = 1, n_ctg !category
895  do ik = 1, n_kap(ic) !kappa bin
896  do is0 = 1, n_siz(ic) !size bin
897  do ia0 = 1, n_atr !attributes
898  qtrc1(k,i,j,qaes-1+it_procs2trans(ia0,is0,ik,ic)) = aerosol_procs(ia0,is0,ik,ic) / dens(k,i,j)
899  enddo !ia (1:N_ATR )
900  enddo !is (1:n_siz(ic) )
901  enddo !ik (1:n_kap(ic) )
902  enddo !ic (1:n_ctg )
903  ! [ug/m3] -> mixing ratio [kg/kg]
904  qtrc1(k,i,j,qaee-gas_ctg+1:qaee-gas_ctg+gas_ctg) &
905  = conc_gas(1:gas_ctg) / dens(k,i,j)*1.e-9_rp
906 
907  ! tiny number, tiny mass
908  do ic = 1, n_ctg !aerosol category
909  do ik = 1, n_kap(ic) !kappa bin
910  do is0 = 1, n_siz(ic) !size bin
911  if (qtrc1(k,i,j,qaes-1+it_procs2trans(ia_m0,is0,ik,ic))*dens(k,i,j) < cleannumber) then
912  do ia0 = 1, n_atr
913  qtrc1(k,i,j,qaes-1+it_procs2trans(ia0,is0,ik,ic)) = 0._rp !to save cpu time and avoid underflow
914  enddo !ia0 (1:n_atr )
915  endif
916  enddo !is (1:n_siz(ic) )
917  enddo !ik (1:n_kap(ic) )
918  enddo !ic (1:n_ctg )
919 
920  !--- Negative fixer
921  do ic = 1, n_ctg !aerosol category
922  do ik = 1, n_kap(ic) !kappa bin
923  do is0 = 1, n_siz(ic) !size bin
924  if (qtrc(k,i,j,qaes-1+it_procs2trans(ia_m0,is0,ik,ic)) < 0.0_rp .or. &
925  qtrc(k,i,j,qaes-1+it_procs2trans(ia_m2,is0,ik,ic)) < 0.0_rp .or. &
926  qtrc(k,i,j,qaes-1+it_procs2trans(ia_m3,is0,ik,ic)) < 0.0_rp .or. &
927  qtrc(k,i,j,qaes-1+it_procs2trans(ia_ms,is0,ik,ic)) < 0.0_rp .or. &
928  qtrc(k,i,j,qaes-1+it_procs2trans(ia_kp,is0,ik,ic)) < 0.0_rp ) then
929  qtrc(k,i,j,qaes-1+it_procs2trans(1:n_atr,is0,ik,ic)) = 0.0_rp
930  endif
931  enddo
932  enddo
933  enddo
934 
935  ! for history
936  total_aerosol_mass(k,i,j,:) = 0.0_rp
937  total_aerosol_number(k,i,j,:) = 0.0_rp
938  total_emit_aerosol_mass(k,i,j,:) = 0.0_rp
939  total_emit_aerosol_number(k,i,j,:) = 0.0_rp
940  do ic = 1, n_ctg
941  do ik = 1, n_kap(ic)
942  do is0 = 1, n_siz(ic)
943  total_aerosol_mass(k,i,j,ic) = total_aerosol_mass(k,i,j,ic) &
944  + qtrc1(k,i,j,qaes-1+it_procs2trans(ia_ms,is0,ik,ic))
945  total_aerosol_number(k,i,j,ic) = total_aerosol_number(k,i,j,ic) &
946  + qtrc1(k,i,j,qaes-1+it_procs2trans(ia_m0,is0,ik,ic))
947  total_emit_aerosol_mass(k,i,j,ic) = total_emit_aerosol_mass(k,i,j,ic) &
948  + emit(k,i,j,it_procs2trans(ia_ms,is0,ik,ic))
949  total_emit_aerosol_number(k,i,j,ic) = total_emit_aerosol_number(k,i,j,ic) &
950  + emit(k,i,j,it_procs2trans(ia_m0,is0,ik,ic))
951  enddo
952  enddo
953  enddo
954 
955  enddo
956  enddo
957  enddo
958 
959  do ic = 1, n_ctg
960  write(ofilename,'(a,a)') trim(ctg_name(ic)), 'mass'
961  call hist_in( total_aerosol_mass(:,:,:,ic), trim(ofilename), 'Total mass mixing ratio of aerosol', 'kg/kg' )
962  write(ofilename,'(a,a)') trim(ctg_name(ic)), 'number'
963  call hist_in( total_aerosol_number(:,:,:,ic), trim(ofilename), 'Total number mixing ratio of aerosol', 'num/kg' )
964  write(ofilename,'(a,a)') trim(ctg_name(ic)), 'mass_emit'
965  call hist_in( total_emit_aerosol_mass(:,:,:,ic), trim(ofilename), 'Total mass mixing ratio of emitted aerosol', 'kg/kg' )
966  write(ofilename,'(a,a)') trim(ctg_name(ic)), 'number_emit'
967  call hist_in( total_emit_aerosol_number(:,:,:,ic), trim(ofilename), 'Total number mixing ratio of emitted aerosol', 'num/kg' )
968  enddo
969 
970  do ic = 1, n_ctg
971  write(ofilename,'(a,a)') trim(ctg_name(ic)), 'mass'
972  call hist_in( total_aerosol_mass(:,:,:,ic), trim(ofilename), 'Total mass mixing ratio of aerosol', 'kg/kg' )
973  write(ofilename,'(a,a)') trim(ctg_name(ic)), 'number'
974  call hist_in( total_aerosol_number(:,:,:,ic), trim(ofilename), 'Total number mixing ratio of aerosol', 'num/kg' )
975  write(ofilename,'(a,a)') trim(ctg_name(ic)), 'mass_emit'
976  call hist_in( total_emit_aerosol_mass(:,:,:,ic), trim(ofilename), 'Total mass mixing ratio of emitted aerosol', 'kg/kg' )
977  write(ofilename,'(a,a)') trim(ctg_name(ic)), 'number_emit'
978  call hist_in( total_emit_aerosol_number(:,:,:,ic), trim(ofilename), 'Total number mixing ratio of emitted aerosol', 'num/kg' )
979  enddo
980 
981  call hist_in( emit(:,:,:,qa_ae-gas_ctg+ig_h2so4), 'H2SO4_emit', 'Emission ratio of H2SO4 gas', 'ug/m3/s' )
982  call hist_in( emit(:,:,:,qa_ae-gas_ctg+ig_cgas), 'CGAS_emit', 'Emission ratio of Condensabule gas', 'ug/m3/s' )
983 
984  deallocate( aerosol_procs )
985  deallocate( aerosol_activ )
986  deallocate( emis_procs )
987  deallocate( emis_gas )
988 
989  do iq = 1, qa
990  do j = js, je
991  do i = is, ie
992  do k = ks, ke
993  rhoq_t_ae(k,i,j,iq) = ( qtrc1(k,i,j,iq) - qtrc0(k,i,j,iq) ) * dens(k,i,j) / time_dtsec_atmos_phy_ae
994  enddo
995  enddo
996  enddo
997  enddo
998 
999  return
integer, public is
start point of inner domain: x, local
real(rp), public const_cvdry
specific heat (dry air,constant volume) [J/kg/K]
Definition: scale_const.F90:59
integer, public je
end point of inner domain: y, local
real(rp), public const_cpdry
specific heat (dry air,constant pressure) [J/kg/K]
Definition: scale_const.F90:58
module ATMOSPHERE / Saturation adjustment
real(rp), dimension(qa_max), public tracer_r
integer, public ke
end point of inner domain: z, local
real(rp), dimension(qa_max), public tracer_cv
real(rp), public const_rdry
specific gas constant (dry air) [J/kg/K]
Definition: scale_const.F90:57
module grid index
module TRACER
real(rp), public const_pre00
pressure reference [Pa]
Definition: scale_const.F90:90
integer, public js
start point of inner domain: y, local
module TIME
Definition: scale_time.F90:15
real(dp), public time_dtsec_atmos_phy_ae
time interval of physics(aerosol ) [sec]
Definition: scale_time.F90:46
real(rp), parameter, public const_rvap
specific gas constant (water vapor) [J/kg/K]
Definition: scale_const.F90:65
module CONSTANT
Definition: scale_const.F90:14
integer, public ks
start point of inner domain: z, local
integer, public ie
end point of inner domain: x, local
module HISTORY
real(rp), dimension(qa_max), public tracer_mass
Here is the call graph for this function:
Here is the caller graph for this function:

◆ aerosol_zerochem()

subroutine scale_atmos_phy_ae_kajino13::aerosol_zerochem ( real(dp), intent(in)  deltt,
real(rp), intent(in)  temp_k,
real(rp), intent(in)  pres_pa,
real(rp), intent(in)  super,
logical, intent(in)  flag_npf,
logical, intent(in)  flag_cond,
logical, intent(in)  flag_coag,
real(rp), dimension(n_atr,n_siz_max,n_kap_max,n_ctg), intent(inout)  aerosol_procs,
real(rp), dimension(gas_ctg), intent(inout)  conc_gas,
real(rp), dimension(n_atr,n_siz_max,n_kap_max,n_ctg), intent(in)  emis_procs,
real(rp), dimension(gas_ctg), intent(in)  emis_gas,
real(rp), dimension(n_atr,n_siz_max,n_kap_max,n_ctg), intent(out)  aerosol_activ 
)

Definition at line 1014 of file scale_atmos_phy_ae_kajino13.F90.

References aerosol_nucleation(), scale_prof::prof_rapend(), scale_prof::prof_rapstart(), scale_time::time_nowsec, and scale_time::time_startdaysec.

Referenced by atmos_phy_ae_kajino13().

1014 
1015  use scale_time, only: &
1016  time_nowsec, &
1018 
1019  implicit none
1020  ! i/o variables
1021  real(DP),intent(in) :: deltt ! delta t [sec]
1022  real(RP),intent(in) :: temp_k ! temperature [K]
1023  real(RP),intent(in) :: pres_pa ! pressure [Pa]
1024  real(RP),intent(in) :: super ! supersaturation [-]
1025 ! real(RP),intent(in) :: h2so4dt ! h2so4 production rate [ug/m3/s]
1026 ! real(RP),intent(in) :: c_ratio ! ratio of condensable mass to h2so4 (after NPF) [-]
1027 ! real(RP),intent(in) :: c_kappa ! hygroscopicity of condensable mass [-]
1028  logical, intent(in) :: flag_npf ! (on/off) new particle formation
1029  logical, intent(in) :: flag_cond ! (on/off) condensation
1030  logical, intent(in) :: flag_coag ! (on/off) coagulation
1031  real(RP),intent(inout) :: aerosol_procs(N_ATR,n_siz_max,n_kap_max,n_ctg)
1032  real(RP),intent(inout) :: conc_gas(GAS_CTG)
1033  real(RP),intent(out) :: aerosol_activ(N_ATR,n_siz_max,n_kap_max,n_ctg)
1034  real(RP),intent(in) :: emis_procs(N_ATR,n_siz_max,n_kap_max,n_ctg)
1035  real(RP),intent(in) :: emis_gas(GAS_CTG)
1036  ! local variables
1037  real(RP) :: J_1nm ! nucleation rate of 1nm particles [#/cm3/s]
1038  integer :: ic_nuc ! category for 1nm new particles
1039  integer :: ik_nuc ! kappa bin for 1nm new particles
1040  integer :: is_nuc ! size bin for 1nm new particles
1041  real(RP) :: c_ratio, t_elaps
1042  real(RP) :: chem_gas(GAS_CTG)
1043 
1044  chem_gas(ig_h2so4) = h2so4dt
1045  chem_gas(ig_cgas) = ocgasdt
1046 
1047  !--- convert unit of aerosol mass [ia=ia_ms]
1048  do ic = 1, n_ctg !aerosol category
1049  do ik = 1, n_kap(ic) !kappa bin
1050  do is0 = 1, n_siz(ic) !size bin
1051  aerosol_procs(ia_ms,is0,ik,ic) = aerosol_procs(ia_ms,is0,ik,ic) * 1.e+9_rp ! [kg/m3] -> [ug/m3]
1052  enddo !is (1:n_siz(ic) )
1053  enddo !ik (1:n_kap(ic) )
1054  enddo !ic (1:n_ctg )
1055 
1056  ! emission
1057  do ic = 1, n_ctg !aerosol category
1058  do ik = 1, n_kap(ic) !kappa bin
1059  do is0 = 1, n_siz(ic) !size bin
1060  do ia0 = 1, n_atr !attributes (prognostic)
1061  aerosol_procs(ia0,is0,ik,ic) = aerosol_procs(ia0,is0,ik,ic) &
1062  + emis_procs(ia0,is0,ik,ic) * deltt
1063  enddo !ia0 (1:n_atr )
1064  enddo !is (1:n_siz(ic) )
1065  enddo !ik (1:n_kap(ic) )
1066  enddo !ic (1:n_ctg )
1067 
1068  ! tiny number, tiny mass
1069  do ic = 1, n_ctg !aerosol category
1070  do ik = 1, n_kap(ic) !kappa bin
1071  do is0 = 1, n_siz(ic) !size bin
1072  if (aerosol_procs(ia_m0,is0,ik,ic) < cleannumber) then
1073  do ia0 = 1, n_atr
1074  aerosol_procs(ia0,is0,ik,ic) = 0._rp !to save cpu time and avoid underflow
1075  enddo !ia0 (1:n_atr )
1076  endif
1077  enddo !is (1:n_siz(ic) )
1078  enddo !ik (1:n_kap(ic) )
1079  enddo !ic (1:n_ctg )
1080 
1081  ! update conc_h2so4
1082 ! conc_h2so4 = conc_h2so4 + h2so4dt * deltt ! [ug/m3]
1083  conc_gas(ig_h2so4) = conc_gas(ig_h2so4) + chem_gas(ig_h2so4) * deltt ! [ug/m3]
1084  conc_gas(ig_cgas) = conc_gas(ig_cgas) + chem_gas(ig_cgas) * deltt ! [ug/m3]
1085 
1086  conc_gas(ig_h2so4) = conc_gas(ig_h2so4) + emis_gas(ig_h2so4) * deltt ! [ug/m3]
1087  conc_gas(ig_cgas) = conc_gas(ig_cgas) + emis_gas(ig_cgas) * deltt ! [ug/m3]
1088 
1089  if( chem_gas(ig_h2so4)+emis_gas(ig_h2so4) /= 0.0_rp ) then
1090  c_ratio = ( chem_gas(ig_h2so4)+emis_gas(ig_h2so4)+chem_gas(ig_cgas)+emis_gas(ig_cgas) ) &
1091  / ( chem_gas(ig_h2so4)+emis_gas(ig_h2so4) )
1092  else
1093  c_ratio = 0.0_rp
1094  endif
1095 
1096  ! new particle formation
1097  t_elaps = time_nowsec - time_startdaysec
1098  j_1nm = 0._rp
1099  ic_nuc = ic_mix
1100  ik_nuc = 1
1101  is_nuc = 1
1102  if (flag_npf) then
1103  call prof_rapstart('ATM_Aerosol_NPF',1)
1104  if( t_elaps <= t_npf ) then ! call aerosol_nucleation(conc_h2so4, J_1nm) !(i) conc_h2so4 (o) J_1nm
1105  call aerosol_nucleation(conc_gas(ig_h2so4), j_1nm) !(i) conc_h2so4 (o) J_1nm
1106  endif
1107  call prof_rapend('ATM_Aerosol_NPF',1)
1108  else
1109  j_1nm = 0._rp !( new particle formation does not occur )
1110  endif !if flag_npf=.true.
1111 
1112  ! condensation
1113  if (flag_cond) then
1114  call prof_rapstart('ATM_Aerosol_cond',1)
1115  call aerosol_condensation(j_1nm, temp_k, pres_pa, deltt, &
1116  ic_nuc, ik_nuc, is_nuc, n_ctg, n_kap, n_siz, n_atr, &
1117  n_siz_max, n_kap_max, ia_m0, ia_m2, ia_m3, ia_ms, &
1118  ! d_lw, d_ct, d_up, k_lw, k_ct, k_up, &
1119  d_lw, d_ct, d_up, &
1120  conc_gas(ig_h2so4), c_ratio, aerosol_procs )
1121  call prof_rapend('ATM_Aerosol_cond',1)
1122  endif !if flag_cond=.true.
1123 
1124  ! coagulation
1125  if (flag_coag) then
1126  call prof_rapstart('ATM_Aerosol_coag',1)
1127  call aerosol_coagulation(deltt, temp_k, pres_pa, &
1128  mcomb,is_i,is_j,is_k,ik_i,ik_j,ik_k,ic_i,ic_j,ic_k, &
1129  n_atr,n_siz_max,n_kap_max,n_ctg,ia_m0,ia_m2,ia_m3,ia_ms, &
1130  n_siz,n_kap,d_lw,d_ct,d_up,aerosol_procs)
1131  call prof_rapend('ATM_Aerosol_coag',1)
1132  endif !if flag_coag=.true.
1133 
1134  call aerosol_activation(c_kappa, super, temp_k, ia_m0, ia_m2, ia_m3, &
1135  n_atr,n_siz_max,n_kap_max,n_ctg,n_siz,n_kap, &
1136  d_ct,aerosol_procs, aerosol_activ)
1137 
1138 ! conc_gas(IG_CGAS) = conc_gas(IG_H2SO4)*( c_ratio-1.0_RP )
1139 
1140  !--- convert unit of aerosol mass [ia=ia_ms]
1141  do ic = 1, n_ctg !aerosol category
1142  do ik = 1, n_kap(ic) !kappa bin
1143  do is0 = 1, n_siz(ic) !size bin
1144  aerosol_procs(ia_ms,is0,ik,ic) = aerosol_procs(ia_ms,is0,ik,ic) * 1.e-9_rp ! [ug/m3] -> [kg/m3]
1145  enddo !is (1:n_siz(ic) )
1146  enddo !ik (1:n_kap(ic) )
1147  enddo !ic (1:n_ctg )
1148 
1149  return
real(dp), public time_nowsec
subday part of current time [sec]
Definition: scale_time.F90:68
real(dp), public time_startdaysec
second of start time [sec]
Definition: scale_time.F90:74
module TIME
Definition: scale_time.F90:15
Here is the call graph for this function:
Here is the caller graph for this function:

◆ aerosol_nucleation()

subroutine scale_atmos_phy_ae_kajino13::aerosol_nucleation ( real(rp), intent(in)  conc_h2so4,
real(rp), intent(inout)  J_1nm 
)

Definition at line 1155 of file scale_atmos_phy_ae_kajino13.F90.

References float().

Referenced by aerosol_zerochem().

1155  implicit none
1156  !i/o variables
1157  real(RP), intent(in) :: conc_h2so4 !H2SO4 concentration [ug/m3]
1158  real(RP), intent(inout) :: J_1nm !nucleation rate of 1nm particles [#/cm3/s]
1159  !local variables
1160  real(RP) :: conc_num_h2so4 !H2SO4 number concentration [#/cm3]
1161  ! real(RP), parameter :: logk = -12.4_RP !constant coefficient for kinetic nucleation [-]
1162  ! real(RP), parameter :: logk = -11.4_RP !constant coefficient for kinetic nucleation [-]
1163 
1164  conc_num_h2so4 = conc_h2so4 * conv_m2n ![ug/m3] -> [#/cm3]
1165 
1166  !emperical formula of new particle formation rate [Kuang et al., 2008]
1167  j_1nm = (10._rp**(logk_aenucl)) * conc_num_h2so4 ** 2._rp
1168 
1169  return
Here is the call graph for this function:
Here is the caller graph for this function:

◆ atmos_phy_ae_kajino13_effectiveradius()

subroutine, public scale_atmos_phy_ae_kajino13::atmos_phy_ae_kajino13_effectiveradius ( real(rp), dimension (ka,ia,ja,n_ae), intent(out)  Re,
real(rp), dimension(ka,ia,ja,qa), intent(in)  QTRC,
real(rp), dimension (ka,ia,ja), intent(in)  RH 
)

Calculate Effective Radius.

Definition at line 2243 of file scale_atmos_phy_ae_kajino13.F90.

References scale_const::const_undef, and scale_atmos_aerosol::n_ae.

Referenced by scale_atmos_phy_ae::atmos_phy_ae_config().

2243  use scale_grid_index
2244  use scale_tracer
2245  use scale_const, only: &
2246  undef => const_undef
2247  use scale_atmos_aerosol, only: &
2248  n_ae
2249  implicit none
2250  real(RP), intent(out) :: Re (KA,IA,JA,N_AE) ! effective radius
2251  real(RP), intent(in) :: QTRC(KA,IA,JA,QA) ! tracer mass concentration [kg/kg]
2252  real(RP), intent(in) :: RH (KA,IA,JA) ! relative humidity (0-1)
2253  !---------------------------------------------------------------------------
2254 
2255  re(:,:,:,:) = 0.0_rp
2256 
2257 ! Re(:,:,:,I_ae_seasalt) = 2.E-4_RP
2258 ! Re(:,:,:,I_ae_dust ) = 4.E-6_RP
2259 ! Re(:,:,:,I_ae_bc ) = 4.E-8_RP
2260 ! Re(:,:,:,I_ae_oc ) = RH(:,:,:)
2261 ! Re(:,:,:,I_ae_sulfate) = RH(:,:,:)
2262 
2263  return
real(rp), public const_undef
Definition: scale_const.F90:43
module grid index
module TRACER
module CONSTANT
Definition: scale_const.F90:14
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 ( real(rp), dimension(ka,ia,ja,qa), intent(inout)  QTRC,
real(rp), dimension (ka,ia,ja), intent(out)  CCN,
real(rp), dimension(ka,ia,ja), intent(in)  DENS,
real(rp), dimension(ka,ia,ja), intent(in)  RHOT,
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 
)

Definition at line 2274 of file scale_atmos_phy_ae_kajino13.F90.

References scale_const::const_cpdry, scale_const::const_cvdry, scale_const::const_pi, scale_const::const_pre00, scale_const::const_rdry, scale_atmos_hydrometeor::i_qv, scale_grid_index::ia, scale_grid_index::ie, scale_stdio::io_fid_log, scale_stdio::io_l, scale_grid_index::is, scale_grid_index::ja, scale_grid_index::je, scale_grid_index::js, scale_grid_index::ka, scale_grid_index::ke, scale_grid_index::ks, scale_tracer::tracer_cp, scale_tracer::tracer_cv, scale_tracer::tracer_mass, and scale_tracer::tracer_r.

Referenced by mod_mkinit::rect_setup().

2274  use scale_const, only: &
2275  pi => const_pi, &
2276  const_cvdry, &
2277  const_cpdry, &
2278  const_rdry, &
2279  const_pre00
2280  use scale_atmos_hydrometeor, only: &
2281  i_qv
2282  use scale_atmos_saturation, only: &
2283  saturation_pres2qsat_liq => atmos_saturation_pres2qsat_liq
2284  implicit none
2285 
2286  real(RP), intent(inout) :: QTRC(KA,IA,JA,QA)
2287  real(RP), intent(out) :: CCN (KA,IA,JA)
2288  real(RP), intent(in) :: DENS(KA,IA,JA)
2289  real(RP), intent(in) :: RHOT(KA,IA,JA)
2290  real(RP), intent(in) :: m0_init ! initial total num. conc. of modes (Atk,Acm,Cor) [#/m3]
2291  real(RP), intent(in) :: dg_init ! initial number equivalen diameters of modes [m]
2292  real(RP), intent(in) :: sg_init ! initial standard deviation [-]
2293  real(RP), intent(in) :: d_min_inp(3)
2294  real(RP), intent(in) :: d_max_inp(3)
2295  real(RP), intent(in) :: k_min_inp(3)
2296  real(RP), intent(in) :: k_max_inp(3)
2297  integer, intent(in) :: n_kap_inp(3)
2298 
2299  integer, parameter :: ia_m0 = 1 ! 1. number conc [#/m3]
2300  integer, parameter :: ia_m2 = 2 ! 2. 2nd mom conc [m2/m3]
2301  integer, parameter :: ia_m3 = 3 ! 3. 3rd mom conc [m3/m3]
2302  integer, parameter :: ia_ms = 4 ! 4. mass conc [ug/m3]
2303  integer, parameter :: ia_kp = 5
2304  integer, parameter :: ik_out = 1
2305 
2306  real(RP) :: c_kappa = 0.3_rp ! hygroscopicity of condensable mass [-]
2307  real(RP), parameter :: cleannumber = 1.e-3_rp
2308  ! local variables
2309  real(RP), parameter :: rhod_ae = 1.83_rp ! particle density [g/cm3] sulfate assumed
2310  real(RP), parameter :: conv_vl_ms = rhod_ae/1.e-12_rp ! M3(volume)[m3/m3] to mass[m3/m3]
2311 
2312  integer :: n_trans
2313  real(RP) :: m0t, dgt, sgt, m2t, m3t, mst
2314  real(RP), allocatable :: aerosol_procs(:,:,:,:) ! (n_atr,n_siz_max,n_kap_max,n_ctg)
2315  real(RP), allocatable :: aerosol_activ(:,:,:,:) ! (n_atr,n_siz_max,n_kap_max,n_ctg)
2316  real(RP), allocatable :: emis_procs (:,:,:,:) ! (n_atr,n_siz_max,n_kap_max,n_ctg)
2317  !--- gas
2318  real(RP) :: conc_gas(GAS_CTG) !concentration [ug/m3]
2319  integer :: n_siz_max, n_kap_max, n_ctg
2320 ! integer, allocatable :: it_procs2trans(:,:,:,:) !procs to trans conversion
2321 ! integer, allocatable :: ia_trans2procs(:) !trans to procs conversion
2322 ! integer, allocatable :: is_trans2procs(:) !trans to procs conversion
2323 ! integer, allocatable :: ik_trans2procs(:) !trans to procs conversion
2324 ! integer, allocatable :: ic_trans2procs(:)
2325  !--- bin settings (lower bound, center, upper bound)
2326  real(RP), allocatable :: d_lw(:,:), d_ct(:,:), d_up(:,:) !diameter [m]
2327  real(RP), allocatable :: k_lw(:,:), k_ct(:,:), k_up(:,:) !kappa [-]
2328  real(RP) :: dlogd, dk !delta log(D), delta K
2329  real(RP), allocatable :: d_min(:) !lower bound of 1st size bin (n_ctg)
2330  real(RP), allocatable :: d_max(:) !upper bound of last size bin (n_ctg)
2331  integer, allocatable :: n_kap(:) !number of kappa bins (n_ctg)
2332  real(RP), allocatable :: k_min(:) !lower bound of 1st kappa bin (n_ctg)
2333  real(RP), allocatable :: k_max(:)
2334 
2335  real(RP) :: pott, qdry, pres
2336  real(RP) :: temp, cpa, cva, qsat_tmp, ssliq, Rmoist
2337  real(RP) :: pi6
2338 
2339  integer :: ia0, ik, is0, ic, k, i, j, it, iq
2340  !---------------------------------------------------------------------------
2341 
2342 
2343  pi6 = pi / 6._rp
2344  n_ctg = ae_ctg
2345 
2346  allocate( d_min(n_ctg) )
2347  allocate( d_max(n_ctg) )
2348  allocate( n_kap(n_ctg) )
2349  allocate( k_min(n_ctg) )
2350  allocate( k_max(n_ctg) )
2351 
2352 
2353  d_min(1:n_ctg) = d_min_inp(1:n_ctg) ! lower bound of 1st size bin
2354  d_max(1:n_ctg) = d_max_inp(1:n_ctg) ! upper bound of last size bin
2355  n_kap(1:n_ctg) = n_kap_inp(1:n_ctg) ! number of kappa bins
2356  k_min(1:n_ctg) = k_min_inp(1:n_ctg) ! lower bound of 1st kappa bin
2357  k_max(1:n_ctg) = k_max_inp(1:n_ctg) ! upper bound of last kappa bin
2358 
2359  !--- diagnose parameters (n_trans, n_siz_max, n_kap_max)
2360  n_trans = 0
2361  n_siz_max = 0
2362  n_kap_max = 0
2363  do ic = 1, n_ctg
2364  n_trans = n_trans + nsiz(ic) * nkap(ic) * n_atr
2365  n_siz_max = max(n_siz_max, nsiz(ic))
2366  n_kap_max = max(n_kap_max, nkap(ic))
2367  enddo
2368 
2369  allocate( aerosol_procs(n_atr,n_siz_max,n_kap_max,n_ctg) )
2370  allocate( aerosol_activ(n_atr,n_siz_max,n_kap_max,n_ctg) )
2371  allocate( emis_procs(n_atr,n_siz_max,n_kap_max,n_ctg) )
2372  aerosol_procs(:,:,:,:) = 0._rp
2373  aerosol_activ(:,:,:,:) = 0._rp
2374  emis_procs(:,:,:,:) = 0._rp
2375 
2376 ! allocate( it_procs2trans(N_ATR,n_siz_max,n_kap_max,n_ctg) )
2377 ! allocate( ia_trans2procs(n_trans) )
2378 ! allocate( is_trans2procs(n_trans) )
2379 ! allocate( ik_trans2procs(n_trans) )
2380 ! allocate( ic_trans2procs(n_trans) )
2381 
2382  !bin setting
2383  allocate(d_lw(n_siz_max,n_ctg))
2384  allocate(d_ct(n_siz_max,n_ctg))
2385  allocate(d_up(n_siz_max,n_ctg))
2386  allocate(k_lw(n_kap_max,n_ctg))
2387  allocate(k_ct(n_kap_max,n_ctg))
2388  allocate(k_up(n_kap_max,n_ctg))
2389  d_lw(:,:) = 0._rp
2390  d_ct(:,:) = 0._rp
2391  d_up(:,:) = 0._rp
2392  k_lw(:,:) = 0._rp
2393  k_ct(:,:) = 0._rp
2394  k_up(:,:) = 0._rp
2395 
2396  do ic = 1, n_ctg
2397 
2398  dlogd = (log(d_max(ic)) - log(d_min(ic)))/real(NSIZ(ic),kind=rp)
2399 
2400  do is0 = 1, nsiz(ic) !size bin
2401  d_lw(is0,ic) = exp(log(d_min(ic))+dlogd* real(is0-1,kind=RP) )
2402  d_ct(is0,ic) = exp(log(d_min(ic))+dlogd*(real(is0 ,kind=rp)-0.5_rp))
2403  d_up(is0,ic) = exp(log(d_min(ic))+dlogd* real(is0 ,kind=RP) )
2404  enddo !is (1:n_siz(ic))
2405  dk = (k_max(ic) - k_min(ic))/real(n_kap(ic),kind=rp)
2406  do ik = 1, n_kap(ic) !size bin
2407  k_lw(ik,ic) = k_min(ic) + dk * real(ik-1,kind=rp)
2408  k_ct(ik,ic) = k_min(ic) + dk *(real(ik ,kind=rp)-0.5_rp)
2409  k_up(ik,ic) = k_min(ic) + dk * real(ik ,kind=rp)
2410  enddo !ik (1:n_kap(ic))
2411 
2412  enddo !ic (1:n_ctg)
2413 ! ik = 1 !only one kappa bin
2414 
2415  m0t = m0_init !total M0 [#/m3]
2416  dgt = dg_init ![m]
2417  sgt = sg_init ![-]
2418 
2419  if ( m0t <= cleannumber ) then
2420  m0t = cleannumber
2421  dgt = 0.1e-6_rp
2422  sgt = 1.3_rp
2423  if( io_l ) write(io_fid_log,*) '*** WARNING! Initial aerosol number is set as ', cleannumber, '[#/m3]'
2424  endif
2425 
2426  m2t = m0t*dgt**(2.d0) *dexp(2.0d0 *(dlog(real(sgt,kind=dp))**2.d0)) !total M2 [m2/m3]
2427  m3t = m0t*dgt**(3.d0) *dexp(4.5d0 *(dlog(real(sgt,kind=dp))**2.d0)) !total M3 [m3/m3]
2428  mst = m3t*pi6*conv_vl_ms !total Ms [ug/m3]
2429 
2430  do ic = 1, n_ctg
2431  !aerosol_procs initial condition
2432  do ik = 1, n_kap(ic) !kappa bin
2433  do is0 = 1, nsiz(ic)
2434  if (dgt >= d_lw(is0,ic) .and. dgt < d_up(is0,ic)) then
2435  aerosol_procs(ia_m0,is0,ik,ic) = aerosol_procs(ia_m0,is0,ik,ic) + m0t ![#/m3]
2436  aerosol_procs(ia_m2,is0,ik,ic) = aerosol_procs(ia_m2,is0,ik,ic) + m2t ![m2/m3]
2437  aerosol_procs(ia_m3,is0,ik,ic) = aerosol_procs(ia_m3,is0,ik,ic) + m3t ![m3/m3]
2438  aerosol_procs(ia_ms,is0,ik,ic) = aerosol_procs(ia_ms,is0,ik,ic) + mst*1.e-9_rp ![kg/m3]
2439  elseif (dgt < d_lw(1,ic)) then
2440  aerosol_procs(ia_m0,1 ,ik,ic) = aerosol_procs(ia_m0,1 ,ik,ic) + m0t ![#/m3]
2441  aerosol_procs(ia_m2,1 ,ik,ic) = aerosol_procs(ia_m2,1 ,ik,ic) + m2t ![m2/m3]
2442  aerosol_procs(ia_m3,1 ,ik,ic) = aerosol_procs(ia_m3,1 ,ik,ic) + m3t ![m3/m3]
2443  aerosol_procs(ia_ms,1 ,ik,ic) = aerosol_procs(ia_ms,1 ,ik,ic) + mst*1.e-9_rp ![kg/m3]
2444  elseif (dgt >= d_up(nsiz(ic),ic)) then
2445  aerosol_procs(ia_m0,nsiz(ic),ik,ic) = aerosol_procs(ia_m0,nsiz(ic),ik,ic) + m0t ![#/m3]
2446  aerosol_procs(ia_m2,nsiz(ic),ik,ic) = aerosol_procs(ia_m2,nsiz(ic),ik,ic) + m2t ![m2/m3]
2447  aerosol_procs(ia_m3,nsiz(ic),ik,ic) = aerosol_procs(ia_m3,nsiz(ic),ik,ic) + m3t ![m3/m3]
2448  aerosol_procs(ia_ms,nsiz(ic),ik,ic) = aerosol_procs(ia_ms,nsiz(ic),ik,ic) + mst*1.e-9_rp ![kg/m3]
2449  endif
2450  enddo
2451  enddo
2452  enddo
2453 
2454  conc_gas(:) = 0.0_rp
2455  do j = 1, ja
2456  do i = 1, ia
2457  do k = 1, ka
2458  it = 0
2459  do ic = 1, n_ctg ! category
2460  do ik = 1, nkap(ic) ! kappa bin
2461  do is0 = 1, nsiz(ic) ! size bin
2462  do ia0 = 1, n_atr ! attributes
2463  qtrc(k,i,j,qaes+it) = aerosol_procs(ia0,is0,ik,ic) / dens(k,i,j) !#,m2,m3,kg/m3 -> #,m2,m3,kg/kg
2464  it = it + 1
2465  enddo !ia0 (1:N_ATR )
2466  enddo !is (1:n_siz(ic) )
2467  enddo !ik (1:n_kap(ic) )
2468  enddo !ic (1:n_ctg )
2469 
2470  do ic = 1, gas_ctg ! GAS category
2471  qtrc(k,i,j,qaes+it) = conc_gas(ic) / dens(k,i,j) * 1.e-9_rp !mixing ratio [kg/kg]
2472  it = it + 1
2473  enddo
2474  enddo
2475  enddo
2476  enddo
2477 
2478  ccn(:,:,:) = 0.0_rp
2479  do j = js, je
2480  do i = is, ie
2481  do k = ks, ke
2482  !--- calculate super saturation of water
2483  if ( i_qv > 0 ) then
2484  calc_qdry(qdry, qtrc, tracer_mass, k, i, j, iq)
2485  calc_r(rmoist, qdry, qtrc, k, i, j, iq, const_rdry, tracer_r )
2486  calc_cv(cva, qdry, qtrc, k, i, j, iq, const_cvdry, tracer_cv)
2487  calc_cp(cpa, qdry, qtrc, k, i, j, iq, const_cpdry, tracer_cp)
2488 
2489  pres = const_pre00 * ( rhot(k,i,j) * rmoist / const_pre00 )**(cpa/cva)
2490  temp = pres / ( dens(k,i,j) * rmoist )
2491 
2492  call saturation_pres2qsat_liq( qsat_tmp, temp, pres )
2493 
2494  ssliq = qtrc(k,i,j,i_qv) / qsat_tmp - 1.0_rp
2495  else
2496  ssliq = -1.0_rp
2497  endif
2498 
2499  call aerosol_activation( c_kappa, ssliq, temp, ia_m0, ia_m2, ia_m3, &
2500  n_atr, n_siz_max, n_kap_max, n_ctg, nsiz, n_kap, &
2501  d_ct, aerosol_procs, aerosol_activ )
2502 
2503  do is0 = 1, nsiz(ic_mix)
2504  ccn(k,i,j) = ccn(k,i,j) + aerosol_activ(ia_m0,is0,ik_out,ic_mix)
2505  enddo
2506  enddo
2507  enddo
2508  enddo
2509 
2510  return
integer, public is
start point of inner domain: x, local
real(rp), public const_cvdry
specific heat (dry air,constant volume) [J/kg/K]
Definition: scale_const.F90:59
integer, public je
end point of inner domain: y, local
real(rp), public const_cpdry
specific heat (dry air,constant pressure) [J/kg/K]
Definition: scale_const.F90:58
module ATMOSPHERE / Saturation adjustment
real(rp), dimension(qa_max), public tracer_r
integer, public ke
end point of inner domain: z, local
real(rp), dimension(qa_max), public tracer_cv
real(rp), public const_rdry
specific gas constant (dry air) [J/kg/K]
Definition: scale_const.F90:57
real(rp), dimension(qa_max), public tracer_cp
real(rp), public const_pre00
pressure reference [Pa]
Definition: scale_const.F90:90
integer, public js
start point of inner domain: y, local
module CONSTANT
Definition: scale_const.F90:14
integer, public ks
start point of inner domain: z, local
integer, public ie
end point of inner domain: x, local
real(rp), public const_pi
pi
Definition: scale_const.F90:34
integer, parameter, public rp
integer, parameter, public dp
real(rp), dimension(qa_max), public tracer_mass
Here is the caller graph for this function:

Variable Documentation

◆ atmos_phy_ae_kajino13_name

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

Definition at line 57 of file scale_atmos_phy_ae_kajino13.F90.

Referenced by scale_atmos_phy_ae::atmos_phy_ae_config(), and atmos_phy_ae_kajino13_config().

57  character(len=H_SHORT), public, target, allocatable :: ATMOS_PHY_AE_kajino13_NAME(:)

◆ atmos_phy_ae_kajino13_desc

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

Definition at line 58 of file scale_atmos_phy_ae_kajino13.F90.

Referenced by scale_atmos_phy_ae::atmos_phy_ae_config(), and atmos_phy_ae_kajino13_config().

58  character(len=H_MID) , public, target, allocatable :: ATMOS_PHY_AE_kajino13_DESC(:)

◆ atmos_phy_ae_kajino13_unit

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

Definition at line 59 of file scale_atmos_phy_ae_kajino13.F90.

Referenced by scale_atmos_phy_ae::atmos_phy_ae_config(), and atmos_phy_ae_kajino13_config().

59  character(len=H_SHORT), public, target, allocatable :: ATMOS_PHY_AE_kajino13_UNIT(:)

◆ atmos_phy_ae_kajino13_dens

real(rp), dimension(n_ae), target, public scale_atmos_phy_ae_kajino13::atmos_phy_ae_kajino13_dens

Definition at line 61 of file scale_atmos_phy_ae_kajino13.F90.

Referenced by scale_atmos_phy_ae::atmos_phy_ae_config(), and atmos_phy_ae_kajino13_setup().

61  real(RP), public, target :: ATMOS_PHY_AE_kajino13_DENS(N_AE) ! hydrometeor density [kg/m3]=[g/L]