22       mwwat => const_mvap, &              
    23       dnwat => const_dwatr, &             
    25       stdatmpa =>  const_pstd, &          
    26       stdtemp  =>  const_tem00, &         
    49   integer, 
private :: qa_ae  = 0
    76   integer, 
parameter :: gas_ctg = 2
    77   integer, 
parameter :: n_atr = 5
    79   integer, 
parameter :: ic_mix   =  1
    80   integer, 
parameter :: ic_sea   =  2
    81   integer, 
parameter :: ic_dus   =  3
    83   integer, 
parameter :: ig_h2so4 =  1
    84   integer, 
parameter :: ig_cgas  =  2
    87   integer, 
allocatable :: nkap(:)
    88   integer, 
allocatable :: nsiz(:)
    90   integer  :: atmos_phy_ae_kajino13_nbins_out = 1024
   100   real(RP), 
parameter  :: avo     = 6.0221367e23_rp       
   101   real(RP), 
parameter  :: boltz   = 1.38048e-23_rp        
   102   real(RP), 
parameter  :: conv_n2m = 1.e6_rp/avo*1.e6_rp &
   104   real(RP), 
parameter  :: conv_m2n = 1._rp/conv_n2m       
   105   real(RP), 
parameter  :: rhod_ae  = 1.83_rp              
   106   real(RP), 
parameter  :: rho_kg   = rhod_ae*1.e3_rp      
   107   real(RP), 
parameter  :: grav     = 9.80622_rp           
   108   real(RP), 
parameter  :: conv_ms_vl = 1.e-12_rp/rhod_ae  
   109   real(RP), 
parameter  :: conv_vl_ms = rhod_ae/1.e-12_rp  
   110   real(RP), 
parameter  :: mwrat_s6   = 96._rp/98._rp      
   111   real(RP), 
parameter  :: mwrat_s6_i = 1._rp/mwrat_s6     
   112   real(RP), 
parameter  :: diffsulf =  9.36e-6_rp          
   113   real(RP), 
parameter  :: mwh2so4  =  98._rp              
   118   real(RP), 
parameter  :: c1=1.425e-6_rp, c2=0.5039_rp, c3=108.3_rp   
   123   integer, 
parameter :: ia_m0  = 1             
   124   integer, 
parameter :: ia_m2  = 2             
   125   integer, 
parameter :: ia_m3  = 3             
   126   integer, 
parameter :: ia_ms  = 4             
   127   integer, 
parameter :: ia_kp  = 5             
   128   integer, 
parameter :: ik_out = 1
   131   integer, 
allocatable :: n_siz(:)              
   132   real(RP),
allocatable :: d_min(:)              
   133   real(RP),
allocatable :: d_max(:)              
   134   integer, 
allocatable :: n_kap(:)              
   135   real(RP),
allocatable :: k_min(:)              
   136   real(RP),
allocatable :: k_max(:)              
   140   real(RP),
allocatable :: d_lw(:,:), d_ct(:,:), d_up(:,:)  
   141   real(RP),
allocatable :: k_lw(:,:), k_ct(:,:), k_up(:,:)  
   142   real(RP) :: dlogd, dk                                    
   145   integer, 
allocatable :: is_i(:), is_j(:), is_k(:)    
   146   integer, 
allocatable :: ik_i(:), ik_j(:), ik_k(:)    
   147   integer, 
allocatable :: ic_i(:), ic_j(:), ic_k(:)    
   149   integer :: is1, is2, mc, ik, is0, ic, ia0
   161   integer, 
allocatable :: it_procs2trans(:,:,:,:) 
   162   integer, 
allocatable :: ia_trans2procs(:) 
   163   integer, 
allocatable :: is_trans2procs(:) 
   164   integer, 
allocatable :: ik_trans2procs(:) 
   165   integer, 
allocatable :: ic_trans2procs(:) 
   166   real(RP), 
allocatable :: rnum_out(:)
   168   character(len=H_SHORT),
allocatable :: ctg_name(:)
   169   real(RP), 
parameter :: cleannumber = 1.e-3 
   173   real(RP), 
allocatable :: aerosol_procs(:,:,:,:) 
   174   real(RP), 
allocatable :: aerosol_activ(:,:,:,:) 
   175   real(RP), 
allocatable :: emis_procs   (:,:,:,:) 
   176   real(RP), 
allocatable :: emis_gas     (:)       
   186     integer, 
intent(out) :: QA_AE
   188     integer, 
allocatable :: aero_idx(:,:,:,:)
   190     integer :: NASIZ(3), NAKAP(3)
   191     character(len=H_SHORT) :: attribute, catego, aunit
   193     namelist / param_atmos_phy_ae_kajino13_tracer / &
   198     integer :: m, ierr, ik, ic, ia0, is0
   203     log_info(
"ATMOS_PHY_AE_kajino13_tracer_setup",*) 
'Setup'   208     log_warn(
"ATMOS_PHY_AE_kajino13_tracer_setup",*) 
'### Kajino13 scheme is still experimental'   210     ncat_max = max( ic_mix, ic_sea, ic_dus )
   216     read(
io_fid_conf,nml=param_atmos_phy_ae_kajino13_tracer,iostat=ierr)
   219        log_info(
"ATMOS_PHY_AE_kajino13_tracer_setup",*)  
'Not found namelist. Default used.'   220     elseif( ierr > 0 ) 
then    221        log_error(
"ATMOS_PHY_AE_kajino13_tracer_setup",*) 
'Not appropriate names in namelist PARAM_ATMOS_PHY_AE_KAJINO13_TRACER, Check!'   225     log_nml(param_atmos_phy_ae_kajino13_tracer)
   227     if( ae_ctg > ncat_max ) 
then   228        log_error(
"ATMOS_PHY_AE_kajino13_tracer_setup",*) 
'AE_CTG should be smaller than', ncat_max+1, 
'stop'   232     allocate( nsiz(ae_ctg) )
   233     allocate( nkap(ae_ctg) )
   235     nkap(1:ae_ctg) = nakap(1:ae_ctg)
   236     nsiz(1:ae_ctg) = nasiz(1:ae_ctg)
   238     if( maxval( nkap ) /= 1 .OR. minval( nkap ) /= 1 ) 
then   239        log_error(
"ATMOS_PHY_AE_kajino13_tracer_setup",*) 
'NKAP(:) /= 1 is not supported now, Stop!'   247        qa_ae   = qa_ae + n_atr
   251     qa_ae = qa_ae + gas_ctg
   260        n_siz_max = max(n_siz_max, nsiz(ic))
   261        n_kap_max = max(n_kap_max, nkap(ic))
   264     allocate( aero_idx(n_atr,ae_ctg,n_kap_max,n_siz_max) )
   272        aero_idx(ia0,ic,ik,is0) = m
   283     ic = qa_ae-gas_ctg+ig_h2so4
   285     ic = qa_ae-gas_ctg+ig_cgas
   310           write(attribute,
'(a)') 
"Number"   311           write(aunit,
'(a)') 
"num/kg"   312        elseif( ia0 == 2 ) 
then   313           write(attribute,
'(a)') 
"Section"   314           write(aunit,
'(a)') 
"m2/kg"   315        elseif( ia0 == 3 ) 
then   316           write(attribute,
'(a)') 
"Volume"   317           write(aunit,
'(a)') 
"m3/kg"   318        elseif( ia0 == 4 ) 
then   319           write(attribute,
'(a)') 
"Mass"   320           write(aunit,
'(a)') 
"kg/kg"   321        elseif( ia0 == 5 ) 
then   322           write(attribute,
'(a)') 
"kXm"   323           write(aunit,
'(a)') 
"kg/kg"   325        if( ic == ic_mix ) 
then   326           write(catego,
'(a)') 
"Sulf_"   327        elseif( ic == ic_sea ) 
then   328           write(catego,
'(a)') 
"Salt_"   329        elseif( ic == ic_dus ) 
then   330           write(catego,
'(a)') 
"Dust_"   334        write(
atmos_phy_ae_kajino13_desc(aero_idx(ia0,ic,ik,is0)),
'(a,a,a,i0)') trim(attribute), 
' mixing radio of ', trim(catego), is0
   339     ic = qa_ae-gas_ctg+ig_h2so4
   341     ic = qa_ae-gas_ctg+ig_cgas
   344     ic = qa_ae-gas_ctg+ig_h2so4
   346     ic = qa_ae-gas_ctg+ig_cgas
   361     real(RP), 
allocatable :: d_min_inp(:)
   362     real(RP), 
allocatable :: d_max_inp(:)
   363     real(RP), 
allocatable :: k_min_inp(:)
   364     real(RP), 
allocatable :: k_max_inp(:)
   365     integer , 
allocatable :: n_kap_inp(:)
   367     real(RP), 
parameter :: d_min_def = 1.e-9_rp 
   368     real(RP), 
parameter :: d_max_def = 1.e-5_rp 
   369     integer , 
parameter :: n_kap_def = 1        
   370     real(RP), 
parameter :: k_min_def = 0.e0_rp  
   371     real(RP), 
parameter :: k_max_def = 1.e0_rp  
   373     namelist / param_atmos_phy_ae_kajino13 / &
   391        atmos_phy_ae_kajino13_nbins_out
   397     log_info(
"ATMOS_PHY_AE_kajino13_setup",*) 
'Setup'   398     log_info(
"ATMOS_PHY_AE_kajino13_setup",*) 
'Kajino(2013) scheme'   407     allocate( rnum_out(atmos_phy_ae_kajino13_nbins_out) )
   408     allocate( n_siz(n_ctg) )
   409     allocate( d_min(n_ctg) )
   410     allocate( d_max(n_ctg) )
   411     allocate( n_kap(n_ctg) )
   412     allocate( k_min(n_ctg) )
   413     allocate( k_max(n_ctg) )
   414     allocate( d_min_inp(n_ctg) )
   415     allocate( d_max_inp(n_ctg) )
   416     allocate( n_kap_inp(n_ctg) )
   417     allocate( k_min_inp(n_ctg) )
   418     allocate( k_max_inp(n_ctg) )
   419     allocate( ctg_name(n_ctg) )
   421     n_siz(1:n_ctg) = nsiz(1:n_ctg)       
   422     d_min(1:n_ctg) = d_min_def 
   423     d_max(1:n_ctg) = d_max_def 
   424     n_kap(1:n_ctg) = n_kap_def 
   425     k_min(1:n_ctg) = k_min_def 
   426     k_max(1:n_ctg) = k_max_def 
   429      if( n_ctg == 1 ) 
then   430        write(ctg_name(it),
'(a)') 
"Sulfate"   431      elseif( n_ctg == 2 ) 
then   432        write(ctg_name(it),
'(a)') 
"Seasalt"   433      elseif( n_ctg == 3 ) 
then   434        write(ctg_name(it),
'(a)') 
"Dust"   440     read(
io_fid_conf,nml=param_atmos_phy_ae_kajino13,iostat=ierr)
   442        log_info(
"ATMOS_PHY_AE_kajino13_setup",*) 
'Not found namelist. Default used.'   443     elseif( ierr > 0 ) 
then    444        log_error(
"ATMOS_PHY_AE_kajino13_setup",*) 
'Not appropriate names in namelist PARAM_ATMOS_PHY_AE_KAJINO13. Check!'   447     log_nml(param_atmos_phy_ae_kajino13)
   462       n_trans   = n_trans + n_siz(ic) * n_kap(ic) * n_atr
   463       n_siz_max = max(n_siz_max, n_siz(ic))
   464       n_kap_max = max(n_kap_max, n_kap(ic))
   468     allocate(d_lw(n_siz_max,n_ctg))
   469     allocate(d_ct(n_siz_max,n_ctg))
   470     allocate(d_up(n_siz_max,n_ctg))
   471     allocate(k_lw(n_kap_max,n_ctg))
   472     allocate(k_ct(n_kap_max,n_ctg))
   473     allocate(k_up(n_kap_max,n_ctg))
   482       dlogd = (log(d_max(ic)) - log(d_min(ic)))/
float(n_siz(ic))
   483       do is0 = 1, n_siz(ic)  
   484         d_lw(is0,ic) = exp(log(d_min(ic))+dlogd* 
float(is0-1)      )
   485         d_ct(is0,ic) = exp(log(d_min(ic))+dlogd*(
float(is0)-0.5_rp))
   486         d_up(is0,ic) = exp(log(d_min(ic))+dlogd* 
float(is0)        )
   489       dk  = (k_max(ic) - k_min(ic))/
float(n_kap(ic))
   491         k_lw(ik,ic) = k_min(ic) + dk  * 
float(ik-1)
   492         k_ct(ik,ic) = k_min(ic) + dk  *(
float(ik)-0.5_rp)
   493         k_up(ik,ic) = k_min(ic) + dk  * 
float(ik)
   499     do is0 = 1, n_siz(ic_mix)
   513       do is2 = 1, n_siz(ic)
   514       do is1 = 1, n_siz(ic)
   515         if ( d_ct(is2,ic) >= d_ct(is1,ic) ) 
then   523     allocate(is_i(mcomb))
   524     allocate(is_j(mcomb))
   525     allocate(is_k(mcomb))
   526     allocate(ik_i(mcomb))
   527     allocate(ik_j(mcomb))
   528     allocate(ik_k(mcomb))
   529     allocate(ic_i(mcomb))
   530     allocate(ic_j(mcomb))
   531     allocate(ic_k(mcomb))
   536       do is2 = 1, n_siz(ic)
   537       do is1 = 1, n_siz(ic)
   538         if ( d_ct(is2,ic) >= d_ct(is1,ic) ) 
then   559     allocate( it_procs2trans(n_atr,n_siz_max,n_kap_max,n_ctg)  )
   560     allocate( ia_trans2procs(n_trans) )
   561     allocate( is_trans2procs(n_trans) )
   562     allocate( ik_trans2procs(n_trans) )
   563     allocate( ic_trans2procs(n_trans) )
   565     it_procs2trans(:,:,:,:)= -999
   566     ia_trans2procs(:)      = 0
   567     is_trans2procs(:)      = 0
   568     ik_trans2procs(:)      = 0
   569     ic_trans2procs(:)      = 0
   575     do is0 = 1, n_siz(ic)  
   578       it_procs2trans(ia0,is0,ik,ic)= it
   579       ia_trans2procs(it)         = ia0
   580       is_trans2procs(it)         = is0
   581       ik_trans2procs(it)         = ik
   582       ic_trans2procs(it)         = ic
   588     allocate( aerosol_procs(n_atr,n_siz_max,n_kap_max,n_ctg)  )
   589     allocate( aerosol_activ(n_atr,n_siz_max,n_kap_max,n_ctg)  )
   590     allocate( emis_procs(n_atr,n_siz_max,n_kap_max,n_ctg)  )
   591     allocate( emis_gas(gas_ctg)  )
   592     aerosol_procs(:,:,:,:) = 0.0_rp
   593     aerosol_activ(:,:,:,:) = 0.0_rp
   594     emis_procs(:,:,:,:) = 0.0_rp
   620        pres2qsat_liq => atmos_saturation_pres2qsat_liq
   624     integer,  
intent(in) :: KA, KS, KE
   625     integer,  
intent(in) :: IA, IS, IE
   626     integer,  
intent(in) :: JA, JS, JE
   627     integer,  
intent(in) :: QA_AE
   628     real(RP), 
intent(in) :: TEMP(ka,ia,ja)
   629     real(RP), 
intent(in) :: PRES(ka,ia,ja)
   630     real(RP), 
intent(in) :: QDRY(ka,ia,ja)
   631     real(RP), 
intent(in) :: NREG(ka,ia,ja)
   632     real(RP), 
intent(in) :: DENS(ka,ia,ja)
   633     real(RP), 
intent(in) :: QV  (ka,ia,ja)
   634     real(RP), 
intent(in) :: QTRC(ka,ia,ja,qa_ae)
   635     real(RP), 
intent(in) :: EMIT(ka,ia,ja,qa_ae)
   636     real(DP), 
intent(in) :: dt
   638     real(RP), 
intent(out) :: RHOQ_t_AE(ka,ia,ja,qa_ae)
   639     real(RP), 
intent(out) :: CN(ka,ia,ja)
   640     real(RP), 
intent(out) :: CCN(ka,ia,ja)
   642     real(RP) :: QTRC0(ka,ia,ja,qa_ae)
   643     real(RP) :: QTRC1(ka,ia,ja,qa_ae)
   646     real(RP) :: ssliq_ae(ka,ia,ja)
   647     real(RP) :: rrhog(ka,ia,ja)
   650     real(RP) :: m0_reg, m2_reg, m3_reg      
   652     real(RP) :: reg_factor_m2,reg_factor_m3 
   654     real(RP) :: total_aerosol_mass(ka,ia,ja,n_ctg)
   655     real(RP) :: total_aerosol_number(ka,ia,ja,n_ctg)
   656     real(RP) :: total_emit_aerosol_mass(ka,ia,ja,n_ctg)
   657     real(RP) :: total_emit_aerosol_number(ka,ia,ja,n_ctg)
   661     real(RP) :: conc_gas(gas_ctg)    
   662     integer :: i, j, k, iq, it
   664     log_progress(*) 
'atmosphere / physics / aerosol / kajino13'   666     aerosol_procs(:,:,:,:) = 0.0_rp
   667     aerosol_activ(:,:,:,:) = 0.0_rp
   668     emis_procs(:,:,:,:) = 0.0_rp
   675        qtrc0(k,i,j,iq) = qtrc(k,i,j,iq) 
   690           rrhog(k,i,j) = 1.0_rp / dens(k,i,j)
   694           call pres2qsat_liq( temp(k,i,j), pres(k,i,j), qdry(k,i,j), & 
   696           ssliq_ae(k,i,j) = qv(k,i,j)/qsat_tmp - 1.0_rp
   708       do is0 = 1, n_siz(ic)   
   709         if (qtrc0(k,i,j,it_procs2trans(ia_m0,is0,ik,ic))*dens(k,i,j) < cleannumber) 
then   711             qtrc0(k,i,j,it_procs2trans(ia0,is0,ik,ic)) = 0._rp 
   725       qtrc1(k,i,j,iq) = qtrc0(k,i,j,iq) 
   741             aerosol_procs(ia_trans2procs(it), &
   742                           is_trans2procs(it), &
   743                           ik_trans2procs(it), &
   744                           ic_trans2procs(it)) = qtrc1(k,i,j,it)*dens(k,i,j)
   745             emis_procs(ia_trans2procs(it), &
   746                           is_trans2procs(it), &
   747                           ik_trans2procs(it), &
   748                           ic_trans2procs(it)) = emit(k,i,j,it)*dens(k,i,j)
   751        conc_gas(1:gas_ctg) &
   752                = qtrc1(k,i,j,qa_ae-gas_ctg+1:qa_ae-gas_ctg+gas_ctg)*dens(k,i,j)*1.e+9_rp
   754        emis_gas(1:gas_ctg) = emit(k,i,j,qa_ae-gas_ctg+ig_h2so4:qa_ae-gas_ctg+ig_cgas)
   756        call aerosol_zerochem( &
   772          do is0 = 1, n_siz(ic_mix)
   774            aerosol_activ(ia0,is0,ik_out,ic_mix) = min(max(0._rp, aerosol_activ(ia0,is0,ik_out,ic_mix)), &
   775                                                                  aerosol_procs(ia0,is0,ik_out,ic_mix) )
   777            aerosol_procs(ia0,is0,ik_out,ic_mix) = &
   778            aerosol_procs(ia0,is0,ik_out,ic_mix) - aerosol_activ(ia0,is0,ik_out,ic_mix)
   789          m2_reg = m0_reg * reg_factor_m2     
   790          m3_reg = m0_reg * reg_factor_m3     
   791          ms_reg = m3_reg * pi6 * conv_vl_ms  
   792          aerosol_procs(ia_m0,is0_reg,ik_out,ic_mix) = &
   793          aerosol_procs(ia_m0,is0_reg,ik_out,ic_mix) + m0_reg 
   794          aerosol_procs(ia_m2,is0_reg,ik_out,ic_mix) = &
   795          aerosol_procs(ia_m2,is0_reg,ik_out,ic_mix) + m2_reg 
   796          aerosol_procs(ia_m3,is0_reg,ik_out,ic_mix) = &
   797          aerosol_procs(ia_m3,is0_reg,ik_out,ic_mix) + m3_reg 
   798          aerosol_procs(ia_ms,is0_reg,ik_out,ic_mix) = &
   799          aerosol_procs(ia_ms,is0_reg,ik_out,ic_mix) + ms_reg 
   804        do is0 = 1, n_siz(ic_mix)
   805          ccn(k,i,j) = ccn(k,i,j) + aerosol_activ(ia_m0,is0,ik_out,ic_mix)
   806          cn(k,i,j) = cn(k,i,j) + aerosol_procs(ia_m0,is0,ik_out,ic_mix)
   820        do is0 = 1, n_siz(ic)   
   822           qtrc1(k,i,j,it_procs2trans(ia0,is0,ik,ic)) = aerosol_procs(ia0,is0,ik,ic) / dens(k,i,j)
   828        qtrc1(k,i,j,qa_ae-gas_ctg+1:qa_ae-gas_ctg+gas_ctg) = conc_gas(1:gas_ctg) / dens(k,i,j)*1.e-9_rp
   833        do is0 = 1, n_siz(ic)   
   834          if (qtrc1(k,i,j,it_procs2trans(ia_m0,is0,ik,ic))*dens(k,i,j) < cleannumber) 
then   836              qtrc1(k,i,j,it_procs2trans(ia0,is0,ik,ic)) = 0._rp 
   844        total_aerosol_mass(k,i,j,:) = 0.0_rp
   845        total_aerosol_number(k,i,j,:) = 0.0_rp
   846        total_emit_aerosol_mass(k,i,j,:) = 0.0_rp
   847        total_emit_aerosol_number(k,i,j,:) = 0.0_rp
   850        do is0 = 1, n_siz(ic)
   851            total_aerosol_mass(k,i,j,ic) = total_aerosol_mass(k,i,j,ic) &
   852                                                + qtrc1(k,i,j,it_procs2trans(ia_ms,is0,ik,ic))
   853            total_aerosol_number(k,i,j,ic) = total_aerosol_number(k,i,j,ic) &
   854                                                + qtrc1(k,i,j,it_procs2trans(ia_m0,is0,ik,ic))
   855            total_emit_aerosol_mass(k,i,j,ic) = total_emit_aerosol_mass(k,i,j,ic) &
   856                                                + emit(k,i,j,it_procs2trans(ia_ms,is0,ik,ic))
   857            total_emit_aerosol_number(k,i,j,ic) = total_emit_aerosol_number(k,i,j,ic) &
   858                                                + emit(k,i,j,it_procs2trans(ia_m0,is0,ik,ic))
   868       call file_history_in( total_aerosol_mass(:,:,:,ic), trim(ctg_name(ic))//
'_mass', 
'Total mass mixing ratio of aerosol', 
'kg/kg' )
   869       call file_history_in( total_aerosol_number(:,:,:,ic), trim(ctg_name(ic))//
'_number', 
'Total number mixing ratio of aerosol', 
'num/kg' )
   870       call file_history_in( total_emit_aerosol_mass(:,:,:,ic), trim(ctg_name(ic))//
'_mass_emit', 
'Total mass mixing ratio of emitted aerosol', 
'kg/kg' )
   871       call file_history_in( total_emit_aerosol_number(:,:,:,ic), trim(ctg_name(ic))//
'_number_emit', 
'Total number mixing ratio of emitted aerosol', 
'num/kg' )
   874     call file_history_in( emit(:,:,:,qa_ae-gas_ctg+ig_h2so4), 
'H2SO4_emit', 
'Emission ratio of H2SO4 gas', 
'ug/m3/s' )
   875     call file_history_in( emit(:,:,:,qa_ae-gas_ctg+ig_cgas),  
'CGAS_emit',  
'Emission ratio of Condensabule gas', 
'ug/m3/s' )
   882       rhoq_t_ae(k,i,j,iq) = ( qtrc1(k,i,j,iq) - qtrc0(k,i,j,iq) ) * dens(k,i,j) / dt
   900     integer,  
intent(in)  :: KA, IA, JA, QA_AE
   902     real(RP), 
intent(in)  :: QTRC(ka,ia,ja,qa_ae) 
   903     real(RP), 
intent(in)  :: RH  (ka,ia,ja)       
   905     real(RP), 
intent(out) :: Re  (ka,ia,ja,
n_ae)  
   943        saturation_pres2qsat_liq => atmos_saturation_pres2qsat_liq
   946     integer,  
intent(in)  :: KA, KS, KE
   947     integer,  
intent(in)  :: IA, IS, IE
   948     integer,  
intent(in)  :: JA, JS, JE
   949     integer,  
intent(in)  :: QA_AE
   951     real(RP), 
intent(in)  :: DENS(ka,ia,ja)
   952     real(RP), 
intent(in)  :: TEMP(ka,ia,ja)
   953     real(RP), 
intent(in)  :: PRES(ka,ia,ja)
   954     real(RP), 
intent(in)  :: QDRY(ka,ia,ja)
   955     real(RP), 
intent(in)  :: QV  (ka,ia,ja)
   956     real(RP), 
intent(in)  :: m0_init             
   957     real(RP), 
intent(in)  :: dg_init             
   958     real(RP), 
intent(in)  :: sg_init             
   959     real(RP), 
intent(in)  :: d_min_inp(3)
   960     real(RP), 
intent(in)  :: d_max_inp(3)
   961     real(RP), 
intent(in)  :: k_min_inp(3)
   962     real(RP), 
intent(in)  :: k_max_inp(3)
   963     integer,  
intent(in)  :: n_kap_inp(3)
   965     real(RP), 
intent(out) :: QTRC(ka,ia,ja,qa_ae)
   966     real(RP), 
intent(out) :: CCN (ka,ia,ja)
   968     integer,  
parameter :: ia_m0  = 1             
   969     integer,  
parameter :: ia_m2  = 2             
   970     integer,  
parameter :: ia_m3  = 3             
   971     integer,  
parameter :: ia_ms  = 4             
   972     integer,  
parameter :: ia_kp  = 5
   973     integer,  
parameter :: ik_out = 1
   975     real(RP)            :: c_kappa     = 0.3_rp     
   976     real(RP), 
parameter :: cleannumber = 1.e-3_rp
   978     real(RP), 
parameter :: rhod_ae    = 1.83_rp              
   979     real(RP), 
parameter :: conv_vl_ms = rhod_ae/1.e-12_rp     
   981     real(RP) :: m0t, dgt, sgt, m2t, m3t, mst
   983     real(RP) :: conc_gas(gas_ctg)    
   985     real(RP), 
allocatable :: d_lw(:,:), d_ct(:,:), d_up(:,:)  
   986     real(RP), 
allocatable :: k_lw(:,:), k_ct(:,:), k_up(:,:)  
   987     real(RP)  :: dlogd, dk                                    
   988     real(RP), 
allocatable :: d_min(:)              
   989     real(RP), 
allocatable :: d_max(:)              
   990     integer,  
allocatable :: n_kap(:)              
   991     real(RP), 
allocatable :: k_min(:)              
   992     real(RP), 
allocatable :: k_max(:)
   994     real(RP) :: qsat_tmp, ssliq
   997     integer  :: ia0, ik, is0, ic, k, i, j, it
  1007     allocate( d_min(n_ctg) )
  1008     allocate( d_max(n_ctg) )
  1009     allocate( n_kap(n_ctg) )
  1010     allocate( k_min(n_ctg) )
  1011     allocate( k_max(n_ctg) )
  1014     d_min(1:n_ctg) = d_min_inp(1:n_ctg)  
  1015     d_max(1:n_ctg) = d_max_inp(1:n_ctg)  
  1016     n_kap(1:n_ctg) = n_kap_inp(1:n_ctg)  
  1017     k_min(1:n_ctg) = k_min_inp(1:n_ctg)  
  1018     k_max(1:n_ctg) = k_max_inp(1:n_ctg)  
  1022     allocate(d_lw(n_siz_max,n_ctg))
  1023     allocate(d_ct(n_siz_max,n_ctg))
  1024     allocate(d_up(n_siz_max,n_ctg))
  1025     allocate(k_lw(n_kap_max,n_ctg))
  1026     allocate(k_ct(n_kap_max,n_ctg))
  1027     allocate(k_up(n_kap_max,n_ctg))
  1037       dlogd = (log(d_max(ic)) - log(d_min(ic)))/
real(NSIZ(ic),kind=
rp)
  1039       do is0 = 1, nsiz(ic)  
  1040         d_lw(is0,ic) = exp(log(d_min(ic))+dlogd* 
real(is0-1,kind=RP)        )
  1041         d_ct(is0,ic) = exp(log(d_min(ic))+dlogd*(
real(is0  ,kind=
rp)-0.5_rp))
  1042         d_up(is0,ic) = exp(log(d_min(ic))+dlogd* 
real(is0  ,kind=RP)        )
  1044       dk    = (k_max(ic) - k_min(ic))/
real(n_kap(ic),kind=
rp)
  1045       do ik = 1, n_kap(ic)  
  1046         k_lw(ik,ic) = k_min(ic) + dk  * 
real(ik-1,kind=
rp)
  1047         k_ct(ik,ic) = k_min(ic) + dk  *(
real(ik  ,kind=
rp)-0.5_rp)
  1048         k_up(ik,ic) = k_min(ic) + dk  * 
real(ik  ,kind=
rp)
  1058     if ( m0t <= cleannumber ) 
then  1062        log_warn(
"ATMOS_PHY_AE_kajino13_mkinit",*) 
'Initial aerosol number is set as ', cleannumber, 
'[#/m3]'  1065     m2t = m0t*dgt**(2.d0) *dexp(2.0d0 *(dlog(
real(sgt,kind=
dp))**2.d0)) 
  1066     m3t = m0t*dgt**(3.d0) *dexp(4.5d0 *(dlog(
real(sgt,kind=
dp))**2.d0)) 
  1067     mst = m3t*pi6*conv_vl_ms                              
  1071     do ik = 1, n_kap(ic)   
  1072     do is0 = 1, nsiz(ic)
  1073        if (dgt >= d_lw(is0,ic) .and. dgt < d_up(is0,ic)) 
then  1074          aerosol_procs(ia_m0,is0,ik,ic) = aerosol_procs(ia_m0,is0,ik,ic) + m0t          
  1075          aerosol_procs(ia_m2,is0,ik,ic) = aerosol_procs(ia_m2,is0,ik,ic) + m2t          
  1076          aerosol_procs(ia_m3,is0,ik,ic) = aerosol_procs(ia_m3,is0,ik,ic) + m3t          
  1077          aerosol_procs(ia_ms,is0,ik,ic) = aerosol_procs(ia_ms,is0,ik,ic) + mst*1.e-9_rp 
  1078        elseif (dgt < d_lw(1,ic)) 
then  1079          aerosol_procs(ia_m0,1 ,ik,ic) = aerosol_procs(ia_m0,1 ,ik,ic) + m0t          
  1080          aerosol_procs(ia_m2,1 ,ik,ic) = aerosol_procs(ia_m2,1 ,ik,ic) + m2t          
  1081          aerosol_procs(ia_m3,1 ,ik,ic) = aerosol_procs(ia_m3,1 ,ik,ic) + m3t          
  1082          aerosol_procs(ia_ms,1 ,ik,ic) = aerosol_procs(ia_ms,1 ,ik,ic) + mst*1.e-9_rp 
  1083        elseif (dgt >= d_up(nsiz(ic),ic)) 
then  1084          aerosol_procs(ia_m0,nsiz(ic),ik,ic) = aerosol_procs(ia_m0,nsiz(ic),ik,ic) + m0t          
  1085          aerosol_procs(ia_m2,nsiz(ic),ik,ic) = aerosol_procs(ia_m2,nsiz(ic),ik,ic) + m2t          
  1086          aerosol_procs(ia_m3,nsiz(ic),ik,ic) = aerosol_procs(ia_m3,nsiz(ic),ik,ic) + m3t          
  1087          aerosol_procs(ia_ms,nsiz(ic),ik,ic) = aerosol_procs(ia_ms,nsiz(ic),ik,ic) + mst*1.e-9_rp 
  1093     conc_gas(:) = 0.0_rp
  1100        do is0 = 1, nsiz(ic) 
  1102           qtrc(k,i,j,it) = aerosol_procs(ia0,is0,ik,ic) / dens(k,i,j) 
  1110           qtrc(k,i,j,it) = conc_gas(ic) / dens(k,i,j) * 1.e-9_rp 
  1122        call saturation_pres2qsat_liq( temp(k,i,j), pres(k,i,j), qdry(k,i,j), & 
  1124        ssliq = qv(k,i,j) / qsat_tmp - 1.0_rp
  1126        call aerosol_activation( c_kappa, ssliq, temp(k,i,j), ia_m0, ia_m2, ia_m3,       &
  1127                                 n_atr, n_siz_max, n_kap_max, n_ctg, nsiz, n_kap, &
  1128                                 d_ct, aerosol_procs, aerosol_activ               )
  1130        do is0 = 1, nsiz(ic_mix)
  1131           ccn(k,i,j) = ccn(k,i,j) + aerosol_activ(ia_m0,is0,ik_out,ic_mix)
  1141        KA, KS, KE, IA, IS, IE, JA, JS, JE, QA_AE, &
  1144     integer, 
intent(in) :: IA, IS, IE
  1145     integer, 
intent(in) :: JA, JS, JE
  1146     integer, 
intent(in) :: QA_AE
  1148     real(RP), 
intent(inout) :: QTRC(
ka,ia,ja,qa_ae)
  1151     integer :: ic, ik, is0
  1158       do ik = 1, n_kap(ic)   
  1159       do is0 = 1, n_siz(ic)   
  1160         if (qtrc(k,i,j,it_procs2trans(ia_m0,is0,ik,ic)) < 0.0_rp .or. &
  1161             qtrc(k,i,j,it_procs2trans(ia_m2,is0,ik,ic)) < 0.0_rp .or. &
  1162             qtrc(k,i,j,it_procs2trans(ia_m3,is0,ik,ic)) < 0.0_rp .or. &
  1163             qtrc(k,i,j,it_procs2trans(ia_ms,is0,ik,ic)) < 0.0_rp .or. &
  1164             qtrc(k,i,j,it_procs2trans(ia_kp,is0,ik,ic)) < 0.0_rp ) 
then  1165           qtrc(k,i,j,it_procs2trans(1:n_atr,is0,ik,ic)) = 0.0_rp
  1182   subroutine aerosol_zerochem (        &
  1184         temp_k, pres_pa, super,        & !--- in
  1185 !        h2so4dt, c_ratio, c_kappa,     & !--- in
  1186         flag_npf, flag_cond, flag_coag,& !--- in
  1187         aerosol_procs,                 & !--- inout
  1188         conc_gas,                      & !--- inout
  1189         emis_procs,                    & !--- out
  1190         emis_gas,                      & !--- out
  1194     real(DP),
intent(in) :: deltt      
  1195     real(RP),
intent(in) :: temp_k     
  1196     real(RP),
intent(in) :: pres_pa    
  1197     real(RP),
intent(in) :: super      
  1201     logical, 
intent(in) :: flag_npf   
  1202     logical, 
intent(in) :: flag_cond  
  1203     logical, 
intent(in) :: flag_coag  
  1204     real(RP),
intent(inout) :: aerosol_procs(n_atr,n_siz_max,n_kap_max,n_ctg)
  1205     real(RP),
intent(inout) :: conc_gas(gas_ctg)
  1206     real(RP),
intent(out) :: aerosol_activ(n_atr,n_siz_max,n_kap_max,n_ctg)
  1207     real(RP),
intent(in) :: emis_procs(n_atr,n_siz_max,n_kap_max,n_ctg)
  1208     real(RP),
intent(in) :: emis_gas(gas_ctg)
  1215     real(RP)            :: chem_gas(gas_ctg)
  1222     do ik = 1, n_kap(ic)   
  1223     do is0 = 1, n_siz(ic)  
  1224       aerosol_procs(ia_ms,is0,ik,ic) = aerosol_procs(ia_ms,is0,ik,ic) * 1.e+9_rp 
  1231     do ik = 1, n_kap(ic)   
  1232     do is0 = 1, n_siz(ic)   
  1234       aerosol_procs(ia0,is0,ik,ic) = aerosol_procs(ia0,is0,ik,ic) &
  1235                                   + emis_procs(ia0,is0,ik,ic) * deltt
  1243     do ik = 1, n_kap(ic)   
  1244     do is0 = 1, n_siz(ic)   
  1245       if (aerosol_procs(ia_m0,is0,ik,ic) < cleannumber) 
then  1247           aerosol_procs(ia0,is0,ik,ic) = 0._rp 
  1256     conc_gas(ig_h2so4) = conc_gas(ig_h2so4) + chem_gas(ig_h2so4) * deltt  
  1257     conc_gas(ig_cgas)  = conc_gas(ig_cgas)  + chem_gas(ig_cgas) * deltt  
  1259     conc_gas(ig_h2so4) = conc_gas(ig_h2so4) + emis_gas(ig_h2so4) * deltt  
  1260     conc_gas(ig_cgas)  = conc_gas(ig_cgas)  + emis_gas(ig_cgas) * deltt  
  1262     if( chem_gas(ig_h2so4)+emis_gas(ig_h2so4) /= 0.0_rp ) 
then  1263       c_ratio = ( chem_gas(ig_h2so4)+emis_gas(ig_h2so4)+chem_gas(ig_cgas)+emis_gas(ig_cgas) ) &
  1264               / ( chem_gas(ig_h2so4)+emis_gas(ig_h2so4) )
  1275       call aerosol_nucleation(conc_gas(ig_h2so4), j_1nm)     
  1284       call aerosol_condensation(j_1nm, temp_k, pres_pa, deltt,     &
  1285              ic_nuc, ik_nuc, is_nuc, n_ctg, n_kap, n_siz, n_atr,   &
  1286              n_siz_max, n_kap_max, ia_m0, ia_m2, ia_m3, ia_ms,     &
  1289              conc_gas(ig_h2so4), c_ratio, aerosol_procs            )
  1296       call aerosol_coagulation(deltt, temp_k, pres_pa,                &
  1297              mcomb,is_i,is_j,is_k,ik_i,ik_j,ik_k,ic_i,ic_j,ic_k,      &
  1298              n_atr,n_siz_max,n_kap_max,n_ctg,ia_m0,ia_m2,ia_m3,ia_ms, &
  1299              n_siz,n_kap,d_lw,d_ct,d_up,aerosol_procs)
  1304                             super, temp_k, ia_m0, ia_m2, ia_m3, &
  1305                             n_atr,n_siz_max,n_kap_max,n_ctg,n_siz,n_kap, &
  1306                             d_ct,aerosol_procs, aerosol_activ)
  1312     do ik = 1, n_kap(ic)   
  1313     do is0 = 1, n_siz(ic)   
  1314       aerosol_procs(ia_ms,is0,ik,ic) = aerosol_procs(ia_ms,is0,ik,ic) * 1.e-9_rp 
  1320   end subroutine aerosol_zerochem
  1324   subroutine aerosol_nucleation(conc_h2so4, J_1nm)        
  1327     real(RP), 
intent(in)    :: conc_h2so4      
  1328     real(RP), 
intent(inout) :: J_1nm           
  1330     real(RP)                :: conc_num_h2so4  
  1334     conc_num_h2so4 = conc_h2so4 * conv_m2n     
  1340   end subroutine aerosol_nucleation
  1344   subroutine aerosol_condensation(J_1nm, temp_k, pres_pa, deltt,     &
  1345                ic_nuc, ik_nuc, is_nuc, n_ctg, n_kap, n_siz, n_atr,   &
  1346                n_siz_max, n_kap_max, ia_m0, ia_m2, ia_m3, ia_ms,     &
  1347   !              d_lw, d_ct, d_up, k_lw, k_ct, k_up,                   &
  1349                conc_h2so4, c_ratio, aerosol_procs )
  1353     real(RP),
intent(in)    :: J_1nm      
  1354     real(RP),
intent(in)    :: temp_k     
  1355     real(RP),
intent(in)    :: pres_pa    
  1356     real(DP),
intent(in)    :: deltt      
  1357     integer, 
intent(in)    :: ic_nuc     
  1358     integer, 
intent(in)    :: ik_nuc     
  1359     integer, 
intent(in)    :: is_nuc     
  1360     integer, 
intent(in)    :: n_ctg      
  1361     integer, 
intent(in)    :: n_kap(n_ctg)      
  1362     integer, 
intent(in)    :: n_siz(n_ctg)      
  1363     integer, 
intent(in)    :: n_atr      
  1364     integer, 
intent(in)    :: n_siz_max  
  1365     integer, 
intent(in)    :: n_kap_max  
  1366     integer, 
intent(in)    :: ia_m0      
  1367     integer, 
intent(in)    :: ia_m2      
  1368     integer, 
intent(in)    :: ia_m3      
  1369     integer, 
intent(in)    :: ia_ms      
  1370     real(RP),
intent(in)    :: c_ratio    
  1371     real(RP),
intent(in)   :: d_lw(n_siz_max,n_ctg), d_ct(n_siz_max,n_ctg), d_up(n_siz_max,n_ctg)
  1373     real(RP),
intent(inout) :: conc_h2so4 
  1374     real(RP),
intent(inout) :: aerosol_procs(n_atr,n_siz_max,n_kap_max,n_ctg)
  1376     real(RP), 
parameter  :: alpha = 0.1_rp       
  1377     real(RP), 
parameter  :: cour  = 0.5_rp       
  1380     real(RP) :: sq_cbar_h2so4                   
  1381     real(RP) :: cbar_h2so4                      
  1382     real(RP) :: dv_h2so4                        
  1384     real(RP) :: dm0dt_npf, dm2dt_npf, dm3dt_npf, dmsdt_npf     
  1385     real(RP) :: dm0dt_cnd(n_siz_max,n_kap_max,n_ctg) 
  1386     real(RP) :: dm2dt_cnd(n_siz_max,n_kap_max,n_ctg) 
  1387     real(RP) :: dm3dt_cnd(n_siz_max,n_kap_max,n_ctg) 
  1388     real(RP) :: dmsdt_cnd(n_siz_max,n_kap_max,n_ctg) 
  1389     real(RP) :: rm0_pls(n_siz_max), rm0_mns(n_siz_max) 
  1390     real(RP) :: rm2_pls(n_siz_max), rm2_mns(n_siz_max) 
  1391     real(RP) :: rm3_pls(n_siz_max), rm3_mns(n_siz_max) 
  1392     real(RP) :: rms_pls(n_siz_max), rms_mns(n_siz_max) 
  1393     real(RP) :: m0t,m1t,m2t,m3t,dgt,sgt,dm2
  1394     real(RP) :: gnc2,gnc3,gfm2,gfm3,harm2,harm3
  1395     real(RP) :: lossrate,tmps6
  1396     integer  :: ic, ik, is0, i, is1, is2
  1399     if (conc_h2so4 <= 0.0_rp) 
return  1402     drive         = conc_h2so4 * mwrat_s6 * conv_ms_vl  
  1403     sq_cbar_h2so4 = 8.0_rp*rgas*temp_k/(pi*mwh2so4*1.e-3_rp)
  1404     cbar_h2so4    = sqrt( sq_cbar_h2so4 )
  1405     dv_h2so4      = diffsulf*(stdatmpa/pres_pa)*(temp_k/stdtemp)**1.75_rp
  1417     dm0dt_npf = j_1nm * 1.e6_rp                  
  1418     dm2dt_npf = dm0dt_npf * 1.e-18_rp            
  1419     dm3dt_npf = dm0dt_npf * 1.e-27_rp            
  1420     dmsdt_npf = dm3dt_npf * pi6 * conv_vl_ms     
  1424     do ik = 1, n_kap(ic)   
  1425     do is0 = 1, n_siz(ic)   
  1427       dm2dt_cnd(is0,ik,ic) = 0._rp
  1428       dm3dt_cnd(is0,ik,ic) = 0._rp
  1430       if (aerosol_procs(ia_m0,is0,ik,ic) &
  1431          *aerosol_procs(ia_m2,is0,ik,ic) &
  1432          *aerosol_procs(ia_m3,is0,ik,ic) > 0._rp) 
then  1433         m0t = aerosol_procs(ia_m0,is0,ik,ic)
  1434         m2t = aerosol_procs(ia_m2,is0,ik,ic)
  1435         m3t = aerosol_procs(ia_m3,is0,ik,ic)
  1436         call diag_ds(m0t,m2t,m3t,dgt,sgt,dm2)
  1437         if (dgt <= 0._rp) dgt = d_ct(is0,ic)
  1438         if (sgt <= 0._rp) sgt = 1.3_rp
  1439         m1t = m0t*dgt*exp(0.5_rp*(log(sgt)**2._rp)) 
  1441         gnc2 =  8._rp * dv_h2so4 * m0t           
  1442         gnc3 = 12._rp * dv_h2so4 * m1t           
  1444         gfm2 =          alpha * cbar_h2so4 * m1t 
  1445         gfm3 = 1.5_rp * alpha * cbar_h2so4 * m2t 
  1447         harm2 = gnc2 * gfm2 / ( gnc2 + gfm2 )    
  1448         harm3 = gnc3 * gfm3 / ( gnc3 + gfm3 )    
  1450         dm2dt_cnd(is0,ik,ic) = harm2*drive                         
  1451         dm3dt_cnd(is0,ik,ic) = harm3*drive                         
  1452         dmsdt_cnd(is0,ik,ic) = dm3dt_cnd(is0,ik,ic)*pi6*conv_vl_ms 
  1462     lossrate = dmsdt_npf
  1465     do ik = 1, n_kap(ic)   
  1466     do is0 = 1, n_siz(ic)   
  1467       lossrate = lossrate + dmsdt_cnd(is0,ik,ic)
  1472     if (lossrate<=0._rp) 
return   1474     isplt  = int( (lossrate*deltt) / (cour*conc_h2so4) )+1
  1475     dtsplt = deltt/
float(isplt)
  1478     loop_split: 
do i = 1, isplt
  1480       conc_h2so4 = conc_h2so4 - (lossrate * dtsplt) * mwrat_s6_i 
  1482       if (conc_h2so4 < 0._rp) 
then  1483         dtsplt     = tmps6 * mwrat_s6 / lossrate
  1488       aerosol_procs(ia_m0,is_nuc,ik_nuc,ic_nuc) = & 
  1489           aerosol_procs(ia_m0,is_nuc,ik_nuc,ic_nuc) + dm0dt_npf * dtsplt
  1490       aerosol_procs(ia_m2,is_nuc,ik_nuc,ic_nuc) = & 
  1491           aerosol_procs(ia_m2,is_nuc,ik_nuc,ic_nuc) + dm2dt_npf * dtsplt
  1492       aerosol_procs(ia_m3,is_nuc,ik_nuc,ic_nuc) = & 
  1493           aerosol_procs(ia_m3,is_nuc,ik_nuc,ic_nuc) + dm3dt_npf * dtsplt
  1494       aerosol_procs(ia_ms,is_nuc,ik_nuc,ic_nuc) = & 
  1495           aerosol_procs(ia_ms,is_nuc,ik_nuc,ic_nuc) + dmsdt_npf * dtsplt
  1499       do ik = 1, n_kap(ic)   
  1500       do is0 = 1, n_siz(ic)   
  1502         aerosol_procs(ia_m2,is0,ik,ic) = & 
  1503             aerosol_procs(ia_m2,is0,ik,ic) + dm2dt_cnd(is0,ik,ic) * dtsplt & 
  1505         aerosol_procs(ia_m3,is0,ik,ic) = & 
  1506             aerosol_procs(ia_m3,is0,ik,ic) + dm3dt_cnd(is0,ik,ic) * dtsplt & 
  1508         aerosol_procs(ia_ms,is0,ik,ic) = & 
  1509             aerosol_procs(ia_ms,is0,ik,ic) + dmsdt_cnd(is0,ik,ic) * dtsplt & 
  1517       if (conc_h2so4 <= 0.0_rp) 
exit loop_split
  1523     conc_h2so4 = max(conc_h2so4, 0._rp)
  1527     do ik = 1, n_kap(ic)   
  1538       do is1 = 1, n_siz(ic)   
  1540         if (aerosol_procs(ia_m0,is1,ik,ic) &
  1541            *aerosol_procs(ia_m2,is1,ik,ic) &
  1542            *aerosol_procs(ia_m3,is1,ik,ic) > 0._rp) 
then  1543           m0t = aerosol_procs(ia_m0,is1,ik,ic)
  1544           m2t = aerosol_procs(ia_m2,is1,ik,ic)
  1545           m3t = aerosol_procs(ia_m3,is1,ik,ic)
  1546           call diag_ds(m0t,m2t,m3t,dgt,sgt,dm2)
  1547           if (dgt <= 0._rp) dgt = d_ct(is1,ic)
  1548           if (sgt <= 0._rp) sgt = 1.3_rp
  1550           do is2 = 1, n_siz(ic)
  1551             if (dgt >= d_lw(is2,ic) .AND. dgt < d_up(is2,ic)) 
then   1552               rm0_pls(is2) = rm0_pls(is2) + aerosol_procs(ia_m0,is1,ik,ic) 
  1553               rm0_mns(is1) =                aerosol_procs(ia_m0,is1,ik,ic) 
  1554               rm2_pls(is2) = rm2_pls(is2) + aerosol_procs(ia_m2,is1,ik,ic) 
  1555               rm2_mns(is1) =                aerosol_procs(ia_m2,is1,ik,ic) 
  1556               rm3_pls(is2) = rm3_pls(is2) + aerosol_procs(ia_m3,is1,ik,ic) 
  1557               rm3_mns(is1) =                aerosol_procs(ia_m3,is1,ik,ic) 
  1558               rms_pls(is2) = rms_pls(is2) + aerosol_procs(ia_ms,is1,ik,ic) 
  1559               rms_mns(is1) =                aerosol_procs(ia_ms,is1,ik,ic) 
  1568       do is0 = 1, n_siz(ic)
  1569         aerosol_procs(ia_m0,is0,ik,ic) = aerosol_procs(ia_m0,is0,ik,ic) &
  1570                                       + rm0_pls(is0) - rm0_mns(is0)
  1571         aerosol_procs(ia_m2,is0,ik,ic) = aerosol_procs(ia_m2,is0,ik,ic) &
  1572                                       + rm2_pls(is0) - rm2_mns(is0)
  1573         aerosol_procs(ia_m3,is0,ik,ic) = aerosol_procs(ia_m3,is0,ik,ic) &
  1574                                       + rm3_pls(is0) - rm3_mns(is0)
  1575         aerosol_procs(ia_ms,is0,ik,ic) = aerosol_procs(ia_ms,is0,ik,ic) &
  1576                                       + rms_pls(is0) - rms_mns(is0)
  1583   end subroutine aerosol_condensation
  1587   subroutine aerosol_coagulation(deltt, temp_k, pres_pa,                &
  1588               mcomb,is_i,is_j,is_k,ik_i,ik_j,ik_k,ic_i,ic_j,ic_k,      &
  1589               n_atr,n_siz_max,n_kap_max,n_ctg,ia_m0,ia_m2,ia_m3,ia_ms, &
  1590               n_siz,n_kap,d_lw,d_ct,d_up,aerosol_procs)
  1593       real(DP),
intent(in) :: deltt      
  1594       real(RP),
intent(in) :: temp_k     
  1595       real(RP),
intent(in) :: pres_pa    
  1596       integer, 
intent(in) :: mcomb 
  1597       integer, 
intent(in) :: is_i(mcomb), is_j(mcomb), is_k(mcomb)
  1598       integer, 
intent(in) :: ik_i(mcomb), ik_j(mcomb), ik_k(mcomb)
  1599       integer, 
intent(in) :: ic_i(mcomb), ic_j(mcomb), ic_k(mcomb)
  1600       integer, 
intent(in) :: n_ctg      
  1601       integer, 
intent(in) :: n_atr      
  1602       integer, 
intent(in) :: n_siz_max  
  1603       integer, 
intent(in) :: n_kap_max  
  1604       integer, 
intent(in) :: ia_m0      
  1605       integer, 
intent(in) :: ia_m2      
  1606       integer, 
intent(in) :: ia_m3      
  1607       integer, 
intent(in) :: ia_ms      
  1608       integer, 
intent(in) :: n_siz(n_ctg)      
  1609       integer, 
intent(in) :: n_kap(n_ctg)      
  1610       real(RP),
intent(in):: d_lw(n_siz_max,n_ctg), d_ct(n_siz_max,n_ctg), d_up(n_siz_max,n_ctg)
  1611       real(RP),
intent(inout) :: aerosol_procs(n_atr,n_siz_max,n_kap_max,n_ctg)
  1613       real(RP) :: dt_m0(n_siz_max,n_kap_max,n_ctg)
  1614       real(RP) :: dt_m3(n_siz_max,n_kap_max,n_ctg)
  1615       real(RP) :: dt_m6(n_siz_max,n_kap_max,n_ctg)
  1616       real(RP) :: dt_ms(n_siz_max,n_kap_max,n_ctg)
  1617       real(RP) :: sixth(n_siz_max,n_kap_max,n_ctg)
  1618       real(RP) :: rm0_pls(n_siz_max), rm0_mns(n_siz_max) 
  1619       real(RP) :: rm2_pls(n_siz_max), rm2_mns(n_siz_max) 
  1620       real(RP) :: rm3_pls(n_siz_max), rm3_mns(n_siz_max) 
  1621       real(RP) :: rms_pls(n_siz_max), rms_mns(n_siz_max) 
  1624       real(RP) :: rkfm,rknc
  1625       real(RP) :: m0i,m0j,m2i,m2j,m3i,m3j,msi,msj,dgi,dgj,sgi,sgj,dm2,m6i,m6j
  1626       real(RP) :: dm0i,dm0j,dm0k,dm3i,dm3j,dm3k,dm6i,dm6j,dm6k,dmsi,dmsj,dmsk
  1627       real(RP) :: m0t,m3t,m6t,dgt,sgt,m2t
  1628       integer  :: mc, ic, ik, is0, is1, is2
  1630       dt_m0(:,:,:) = 0._rp
  1631       dt_m3(:,:,:) = 0._rp
  1632       dt_m6(:,:,:) = 0._rp
  1633       dt_ms(:,:,:) = 0._rp
  1634       sixth(:,:,:) = 0._rp
  1636       visair=c1*temp_k**c2/(1._rp+c3/temp_k)                
  1637       rkfm=(3._rp*boltz*temp_k /rho_kg)**0.5_rp             
  1638       rknc= 2._rp*boltz*temp_k /(3._rp*visair)              
  1639       c4 = pi*boltz*temp_k/(8._rp*mwair/avo*1.e-3_rp)       
  1640       lambda = visair/(0.499_rp*pres_pa)*c4**0.5_rp*1.e2_rp 
  1645           m0i = aerosol_procs(ia_m0,is_i(mc),ik_i(mc),ic_i(mc))
  1646           m0j = aerosol_procs(ia_m0,is_j(mc),ik_j(mc),ic_j(mc))
  1647           m2i = aerosol_procs(ia_m2,is_i(mc),ik_i(mc),ic_i(mc))
  1648           m2j = aerosol_procs(ia_m2,is_j(mc),ik_j(mc),ic_j(mc))
  1649           m3i = aerosol_procs(ia_m3,is_i(mc),ik_i(mc),ic_i(mc))
  1650           m3j = aerosol_procs(ia_m3,is_j(mc),ik_j(mc),ic_j(mc))
  1651           msi = aerosol_procs(ia_ms,is_i(mc),ik_i(mc),ic_i(mc))
  1652           msj = aerosol_procs(ia_ms,is_j(mc),ik_j(mc),ic_j(mc))
  1654           if (m0i*m2i*m3i <= 0._rp .OR. m0j*m2j*m3j <= 0._rp) cycle
  1656           call diag_ds(m0i,m2i,m3i,dgi,sgi,dm2)
  1657           if (dgi <= 0._rp) dgi = d_ct(is_i(mc),ic_i(mc))
  1658           if (sgi <= 0._rp) sgi = 1.3_rp
  1659           call diag_ds(m0j,m2j,m3j,dgj,sgj,dm2)
  1660           if (dgj <= 0._rp) dgj = d_ct(is_j(mc),ic_j(mc))
  1661           if (sgj <= 0._rp) sgj = 1.3_rp
  1663           m6i = m0i*dgi**6._rp*exp(18._rp*(log(sgi)**2._rp))
  1664           m6j = m0j*dgj**6._rp*exp(18._rp*(log(sgj)**2._rp))
  1665           sixth(is_i(mc),ik_i(mc),ic_i(mc)) = m6i
  1666           sixth(is_j(mc),ik_j(mc),ic_j(mc)) = m6j
  1669           if ( is_i(mc) == is_j(mc) .AND. &
  1670                ik_i(mc) == ik_j(mc) .AND. &
  1671                ic_i(mc) == ic_j(mc) ) 
then  1672             call aero_intra(m0i,  m3i,              & 
  1673                             dm0i, dm3i, dm6i, dmsi, & 
  1686             call aero_inter(m0i ,m3i ,m6i ,      & 
  1688                             dm0i,dm3i,dm6i,dmsi, & 
  1689                             dm0j,dm3j,dm6j,dmsj, & 
  1690                             dm0k,dm3k,dm6k,dmsk, & 
  1691                             dgi ,sgi, rhod_ae,   & 
  1692                             dgj ,sgj, rhod_ae,   & 
  1693                             rknc, rkfm, lambda,  & 
  1697           dt_m0(is_i(mc),ik_i(mc),ic_i(mc)) = dt_m0(is_i(mc),ik_i(mc),ic_i(mc)) + dm0i
  1698           dt_m0(is_j(mc),ik_j(mc),ic_j(mc)) = dt_m0(is_j(mc),ik_j(mc),ic_j(mc)) + dm0j
  1699           dt_m0(is_k(mc),ik_k(mc),ic_k(mc)) = dt_m0(is_k(mc),ik_k(mc),ic_k(mc)) + dm0k
  1700           dt_m3(is_i(mc),ik_i(mc),ic_i(mc)) = dt_m3(is_i(mc),ik_i(mc),ic_i(mc)) + dm3i
  1701           dt_m3(is_j(mc),ik_j(mc),ic_j(mc)) = dt_m3(is_j(mc),ik_j(mc),ic_j(mc)) + dm3j
  1702           dt_m3(is_k(mc),ik_k(mc),ic_k(mc)) = dt_m3(is_k(mc),ik_k(mc),ic_k(mc)) + dm3k
  1703           dt_m6(is_i(mc),ik_i(mc),ic_i(mc)) = dt_m6(is_i(mc),ik_i(mc),ic_i(mc)) + dm6i
  1704           dt_m6(is_j(mc),ik_j(mc),ic_j(mc)) = dt_m6(is_j(mc),ik_j(mc),ic_j(mc)) + dm6j
  1705           dt_m6(is_k(mc),ik_k(mc),ic_k(mc)) = dt_m6(is_k(mc),ik_k(mc),ic_k(mc)) + dm6k
  1706           dt_ms(is_i(mc),ik_i(mc),ic_i(mc)) = dt_ms(is_i(mc),ik_i(mc),ic_i(mc)) + dmsi
  1707           dt_ms(is_j(mc),ik_j(mc),ic_j(mc)) = dt_ms(is_j(mc),ik_j(mc),ic_j(mc)) + dmsj
  1708           dt_ms(is_k(mc),ik_k(mc),ic_k(mc)) = dt_ms(is_k(mc),ik_k(mc),ic_k(mc)) + dmsk
  1715         do ik = 1, n_kap(ic)   
  1716         do is0 = 1, n_siz(ic)   
  1717           if (aerosol_procs(ia_m0,is0,ik,ic) > 0._rp .AND. &                 
  1718               dt_m0(is0,ik,ic)/aerosol_procs(ia_m0,is0,ik,ic) < 1._rp ) 
then   1719             aerosol_procs(ia_m0,is0,ik,ic)=aerosol_procs(ia_m0,is0,ik,ic) &
  1720                                           *exp( dt_m0(is0,ik,ic)/aerosol_procs(ia_m0,is0,ik,ic) )
  1722             aerosol_procs(ia_m0,is0,ik,ic)=aerosol_procs(ia_m0,is0,ik,ic) + dt_m0(is0,ik,ic)
  1724           if (aerosol_procs(ia_m3,is0,ik,ic) > 0._rp .AND. &                 
  1725               dt_m3(is0,ik,ic)/aerosol_procs(ia_m3,is0,ik,ic) < 1._rp ) 
then   1726             aerosol_procs(ia_m3,is0,ik,ic)=aerosol_procs(ia_m3,is0,ik,ic) &
  1727                                           *exp( dt_m3(is0,ik,ic)/aerosol_procs(ia_m3,is0,ik,ic) )
  1729             aerosol_procs(ia_m3,is0,ik,ic)=aerosol_procs(ia_m3,is0,ik,ic) + dt_m3(is0,ik,ic)
  1731           if (sixth(is0,ik,ic)               > 0._rp .AND. &                 
  1732               dt_m6(is0,ik,ic)/sixth(is0,ik,ic) < 1._rp               ) 
then   1733             sixth(is0,ik,ic)              =sixth(is0,ik,ic) &
  1734                                           *exp( dt_m6(is0,ik,ic)/sixth(is0,ik,ic)               )
  1736             sixth(is0,ik,ic)              =sixth(is0,ik,ic)               + dt_m6(is0,ik,ic)
  1738           if (aerosol_procs(ia_ms,is0,ik,ic) > 0._rp .AND. &                 
  1739               dt_ms(is0,ik,ic)/aerosol_procs(ia_ms,is0,ik,ic) < 1._rp ) 
then   1740             aerosol_procs(ia_ms,is0,ik,ic)=aerosol_procs(ia_ms,is0,ik,ic) &
  1741                                           *exp( dt_ms(is0,ik,ic)/aerosol_procs(ia_ms,is0,ik,ic) )
  1743             aerosol_procs(ia_ms,is0,ik,ic)=aerosol_procs(ia_ms,is0,ik,ic) + dt_ms(is0,ik,ic)
  1745           m0t = aerosol_procs(ia_m0,is0,ik,ic)
  1746           m3t = aerosol_procs(ia_m3,is0,ik,ic)
  1747           m6t = sixth(is0,ik,ic)
  1748           call diag_ds6(m0t,m3t,m6t,   &
  1750           if (dgt <= 0._rp) dgt = d_ct(is0,ic)
  1751           if (sgt <= 0._rp) sgt = 1.3_rp
  1752           if (m3t == 0._rp) aerosol_procs(ia_ms,is0,ik,ic)=0._rp
  1753           aerosol_procs(ia_m2,is0,ik,ic)=m2t
  1755           aerosol_procs(ia_m0,is0,ik,ic)=m0t
  1756           aerosol_procs(ia_m3,is0,ik,ic)=m3t
  1763         do ik = 1, n_kap(ic)   
  1774           do is1 = 1, n_siz(ic)   
  1776             if (aerosol_procs(ia_m0,is1,ik,ic) &
  1777                *aerosol_procs(ia_m2,is1,ik,ic) &
  1778                *aerosol_procs(ia_m3,is1,ik,ic) > 0._rp) 
then  1779               m0t = aerosol_procs(ia_m0,is1,ik,ic)
  1780               m2t = aerosol_procs(ia_m2,is1,ik,ic)
  1781               m3t = aerosol_procs(ia_m3,is1,ik,ic)
  1782               call diag_ds(m0t,m2t,m3t,dgt,sgt,dm2)
  1783               if (dgt <= 0._rp) dgt = d_ct(is1,ic)
  1784               if (sgt <= 0._rp) sgt = 1.3_rp
  1786               do is2 = 1, n_siz(ic)
  1787                 if (dgt >= d_lw(is2,ic) .AND. dgt < d_up(is2,ic)) 
then   1788                   rm0_pls(is2) = rm0_pls(is2) + aerosol_procs(ia_m0,is1,ik,ic) 
  1789                   rm0_mns(is1) =                aerosol_procs(ia_m0,is1,ik,ic) 
  1790                   rm2_pls(is2) = rm2_pls(is2) + aerosol_procs(ia_m2,is1,ik,ic) 
  1791                   rm2_mns(is1) =                aerosol_procs(ia_m2,is1,ik,ic) 
  1792                   rm3_pls(is2) = rm3_pls(is2) + aerosol_procs(ia_m3,is1,ik,ic) 
  1793                   rm3_mns(is1) =                aerosol_procs(ia_m3,is1,ik,ic) 
  1794                   rms_pls(is2) = rms_pls(is2) + aerosol_procs(ia_ms,is1,ik,ic) 
  1795                   rms_mns(is1) =                aerosol_procs(ia_ms,is1,ik,ic) 
  1804           do is0 = 1, n_siz(ic)
  1805             aerosol_procs(ia_m0,is0,ik,ic) = aerosol_procs(ia_m0,is0,ik,ic) &
  1806                                           + rm0_pls(is0) - rm0_mns(is0)
  1807             aerosol_procs(ia_m2,is0,ik,ic) = aerosol_procs(ia_m2,is0,ik,ic) &
  1808                                           + rm2_pls(is0) - rm2_mns(is0)
  1809             aerosol_procs(ia_m3,is0,ik,ic) = aerosol_procs(ia_m3,is0,ik,ic) &
  1810                                           + rm3_pls(is0) - rm3_mns(is0)
  1811             aerosol_procs(ia_ms,is0,ik,ic) = aerosol_procs(ia_ms,is0,ik,ic) &
  1812                                           + rms_pls(is0) - rms_mns(is0)
  1820   end subroutine aerosol_coagulation
  1826   subroutine aerosol_activation(c_kappa, super, temp_k, ia_m0, ia_m2, ia_m3, &
  1827                                 n_atr,n_siz_max,n_kap_max,n_ctg,n_siz,n_kap, &
  1828                                 d_ct, aerosol_procs, aerosol_activ)
  1832       real(RP),
intent(in) :: super      
  1833       real(RP),
intent(in) :: c_kappa    
  1834       real(RP),
intent(in) :: temp_k     
  1835       integer, 
intent(in) :: ia_m0, ia_m2, ia_m3
  1836       integer, 
intent(in) :: n_atr
  1837       integer, 
intent(in) :: n_siz_max
  1838       integer, 
intent(in) :: n_kap_max
  1839       integer, 
intent(in) :: n_ctg
  1840       integer, 
intent(in) :: n_siz(n_ctg), n_kap(n_ctg)
  1841       real(RP),
intent(in) :: d_ct(n_siz_max,n_ctg)
  1842       real(RP),
intent(in) :: aerosol_procs(n_atr,n_siz_max,n_kap_max,n_ctg)
  1844       real(RP), 
intent(out) :: aerosol_activ(n_atr,n_siz_max,n_kap_max,n_ctg)
  1847       real(RP),
parameter :: two3 = 2._rp/3._rp
  1848       real(RP),
parameter :: rt2  = sqrt(2._rp)
  1849       real(RP),
parameter :: twort2  = rt2       
  1850       real(RP),
parameter :: thrrt2  = 3._rp/rt2 
  1851       real(RP) :: smax_inv                
  1852       real(RP) :: am,scrit_am,aa,tc,st,bb,ac
  1853       real(RP) :: m0t,m2t,m3t,dgt,sgt,dm2
  1855       real(RP) :: tmp1, tmp2, tmp3        
  1856       real(RP) :: ccn_frc,cca_frc,ccv_frc 
  1857       integer  :: is0, ik, ic
  1859       aerosol_activ(:,:,:,:) = 0._rp
  1861       if (super<=0._rp) 
return  1863       smax_inv = 1._rp / super
  1866       tc = temp_k - stdtemp
  1867       if (tc >= 0._rp ) 
then  1868         st  = 75.94_rp-0.1365_rp*tc-0.3827e-3_rp*tc**2._rp 
  1870         st  = 75.93_rp                +0.115_rp*tc        &  
  1871             + 6.818e-2_rp*tc**2._rp+6.511e-3_rp*tc**3._rp &
  1872             + 2.933e-4_rp*tc**4._rp+6.283e-6_rp*tc**5._rp &
  1873             + 5.285e-8_rp*tc**6._rp
  1879       aa  = 2._rp * st * mwwat * 1.e-3_rp / (dnwat * rgas * temp_k ) 
  1882       do ik = 1, n_kap(ic)
  1883       do is0 = 1, n_siz(ic)
  1884         m0t = aerosol_procs(ia_m0,is0,ik,ic)
  1885         m2t = aerosol_procs(ia_m2,is0,ik,ic)
  1886         m3t = aerosol_procs(ia_m3,is0,ik,ic)
  1887         call diag_ds(m0t,m2t,m3t,dgt,sgt,dm2)
  1888         if (dgt <= 0._rp) dgt = d_ct(is0,ic)
  1889         if (sgt <= 0._rp) sgt = 1.3_rp
  1892         if (bb > 0._rp .AND. am > 0._rp ) 
then  1893           scrit_am = 2._rp/sqrt(bb)*(aa/(3._rp*am))**1.5_rp 
  1897         ac     = am * (scrit_am * smax_inv)**two3      
  1899         tmp1   = log(d_crit) - log(dgt)
  1900         tmp2   = 1._rp/(rt2*log(sgt))
  1902         ccn_frc= 0.5_rp*(1._rp-erf(tmp1*tmp2))
  1903         cca_frc= 0.5_rp*(1._rp-erf(tmp1*tmp2-twort2*tmp3))
  1904         ccv_frc= 0.5_rp*(1._rp-erf(tmp1*tmp2-thrrt2*tmp3))
  1905         ccn_frc=min(max(ccn_frc,0.0_rp),1.0_rp)
  1906         cca_frc=min(max(cca_frc,0.0_rp),1.0_rp)
  1907         ccv_frc=min(max(ccv_frc,0.0_rp),1.0_rp)
  1908         aerosol_activ(ia_m0,is0,ik,ic) = ccn_frc * aerosol_procs(ia_m0,is0,ik,ic)
  1909         aerosol_activ(ia_m2,is0,ik,ic) = cca_frc * aerosol_procs(ia_m2,is0,ik,ic)
  1910         aerosol_activ(ia_m3,is0,ik,ic) = ccv_frc * aerosol_procs(ia_m3,is0,ik,ic)
  1911         aerosol_activ(ia_ms,is0,ik,ic) = ccv_frc * aerosol_procs(ia_ms,is0,ik,ic)
  1917   end subroutine aerosol_activation
  1921   subroutine diag_ds(m0,m2,m3,  & !i
  1924     real(RP)            :: m0,m2,m3,dg,sg,m3_bar,m2_bar
  1925     real(RP)            :: m2_new,m2_old,dm2
  1926     real(RP), 
parameter :: sgmax=2.5_rp
  1927     real(RP), 
parameter :: rk1=2._rp
  1928     real(RP), 
parameter :: rk2=3._rp
  1929     real(RP), 
parameter :: ratio  =rk1/rk2
  1930     real(RP), 
parameter :: rk1_hat=1._rp/(ratio*(rk2-rk1))
  1931     real(RP), 
parameter :: rk2_hat=ratio/(rk1-rk2)
  1932     real(DP), 
parameter :: tiny=1.e-50_dp
  1936     if (m0 <= tiny .OR. m2 <= tiny .OR. m3 <= tiny) 
then  1948     dg     = m2_bar**rk1_hat*m3_bar**rk2_hat
  1950     if (m2_bar/m3_bar**ratio < 1._rp) 
then   1951       sg     = exp(sqrt(2._rp/(rk1*(rk1-rk2))  &
  1952              * log(m2_bar/m3_bar**ratio) ))
  1955     if (sg > sgmax) 
then  1958       call diag_d2(m0,m3,sg, & 
  1963     if (m2_bar/m3_bar**ratio >= 1._rp) 
then   1965       call diag_d2(m0,m3,sg, & 
  1971     dm2    = m2_old - m2_new 
  1974   end subroutine diag_ds
  1978   subroutine diag_d2(m0,m3,sg, & !i
  1981     real(RP)            :: dg,sg,m0,m2,m3,aaa
  1982     real(RP), 
parameter :: one3=1._rp/3._rp
  1983     aaa = m0               * exp( 4.5_rp * (log(sg)**2._rp) )
  1985     m2  = m0 * dg ** 2._rp * exp( 2.0_rp * (log(sg)**2._rp) )
  1988   end subroutine diag_d2
  1992   subroutine aero_intra(m0, m3,            & !input
  1993                           dm0,dm3,dm6,dmi, & !output
  1994                           dg, sg, lambda,  & !input
  2001     real(RP), 
intent(in)   :: m0         
  2002     real(RP), 
intent(in)   :: m3         
  2004     real(RP), 
intent(out)  :: dm0        
  2005     real(RP), 
intent(out)  :: dm3        
  2006     real(RP), 
intent(out)  :: dm6        
  2007     real(RP), 
intent(out)  :: dmi        
  2008     real(RP), 
intent(in)   :: dg         
  2009     real(RP), 
intent(in)   :: sg         
  2010     real(RP), 
intent(in)   :: lambda     
  2011     real(RP), 
intent(in)   :: rknc, rkfm 
  2012     real(DP), 
intent(in)   :: dt         
  2015     real(RP) :: mm2,mm1,m1,m4,m2,mm1p5,mm0p5,m0p5,m3p5,m1p5,m5,m2p5
  2016     real(RP) :: dm0dt,dm6dt,dm0dt_nc,dm6dt_nc,dm0dt_fm,dm6dt_fm,bbb
  2019     real(RP), 
parameter :: rk1=3._rp
  2020     real(RP), 
parameter :: rk2=6._rp
  2021     real(RP), 
parameter :: ratio=rk1/rk2
  2022     real(RP), 
parameter :: rk1_hat=1._rp/(ratio*(rk2-rk1))
  2023     real(RP), 
parameter :: rk2_hat=ratio/(rk1-rk2)
  2034     mm2  =m0*dg**(-2._rp) *exp(2.0_rp *(log(sg)**2._rp)) 
  2035     mm1  =m0*dg**(-1._rp) *exp(0.5_rp *(log(sg)**2._rp)) 
  2036     m1   =m0*dg**  1._rp  *exp(0.5_rp *(log(sg)**2._rp)) 
  2037     m2   =m0*dg**  2._rp  *exp(2.0_rp *(log(sg)**2._rp)) 
  2038     m4   =m0*dg**  4._rp  *exp(8.0_rp *(log(sg)**2._rp)) 
  2041     dm0dt_nc=-1._rp*rknc*(m0*m0+m1*mm1+2.492e-2_rp*lambda*(m0*mm1+m1*mm2))
  2042     dm6dt_nc= 2._rp*rknc*(m3*m3+m4*m2 +2.492e-2_rp*lambda*(m3*m2 +m4*m1 ))
  2046     mm1p5=m0*dg**(-1.5_rp)*exp(1.125_rp*(log(sg)**2._rp))
  2047     mm0p5=m0*dg**(-0.5_rp)*exp(0.125_rp*(log(sg)**2._rp))
  2048     m0p5 =m0*dg**  0.5_rp *exp(0.125_rp*(log(sg)**2._rp))
  2049     m3p5 =m0*dg**( 3.5_rp)*exp(6.125_rp*(log(sg)**2._rp))
  2050     m1p5 =m0*dg**( 1.5_rp)*exp(1.125_rp*(log(sg)**2._rp))
  2051     m5   =m0*dg**  5._rp  *exp(12.50_rp*(log(sg)**2._rp))
  2052     m2p5 =m0*dg**  2.5_rp *exp(3.125_rp*(log(sg)**2._rp))
  2054     bbb     =1._rp+1.2_rp *exp(-2._rp*sg)-0.646_rp*exp(-0.35_rp*sg**2._rp) 
  2056     dm0dt_fm= -1._rp*bbb*rkfm*(m0*m0p5+m2*mm1p5+2._rp*m1*mm0p5)
  2057     dm6dt_fm=  2._rp*bbb*rkfm*(m3p5*m3+m1p5*m5 +2._rp*m2p5*m4 )
  2060     if (dm0dt_fm+dm0dt_nc /= 0._rp) dm0dt = dm0dt_fm*dm0dt_nc/(dm0dt_fm+dm0dt_nc)
  2061     if (dm6dt_fm+dm6dt_nc /= 0._rp) dm6dt = dm6dt_fm*dm6dt_nc/(dm6dt_fm+dm6dt_nc)
  2070   end subroutine aero_intra
  2074   subroutine aero_inter(m0i ,m3i ,m6i ,      & !input unchanged
  2075                         m0j ,m3j ,m6j ,      & !input unchanged
  2076                         dm0i,dm3i,dm6i,dmsi, & !output
  2077                         dm0j,dm3j,dm6j,dmsj, & !output
  2078                         dm0k,dm3k,dm6k,dmsk, & !output
  2079                         dgi ,sgi, rhoi,      & !input unchanged
  2080                         dgj ,sgj, rhoj,      & !input unchanged
  2081                         rknc,  rkfm, lambda, & !input unchanged
  2088     real(RP), 
intent(in)  :: m0i,m0j                 
  2089     real(RP), 
intent(in)  :: m3i,m3j                 
  2090     real(RP), 
intent(in)  :: m6i,m6j                 
  2091     real(RP), 
intent(out) :: dm0i,dm0j,dm0k          
  2092     real(RP), 
intent(out) :: dm3i,dm3j,dm3k          
  2093     real(RP), 
intent(out) :: dm6i,dm6j,dm6k          
  2094     real(RP), 
intent(out) :: dmsi,dmsj,dmsk          
  2095     real(RP), 
intent(in)  :: dgi,sgi,rhoi            
  2096     real(RP), 
intent(in)  :: dgj,sgj,rhoj            
  2097     real(RP), 
intent(in)  :: rknc, rkfm              
  2098     real(RP), 
intent(in)  :: lambda     
  2102     real(RP) :: dm0dt_i,dm3dt_i,dm6dt_i
  2103     real(RP) :: dm0dt_j,dm3dt_j,dm6dt_j
  2104     real(RP) :: dm0dt_k,dm3dt_k,dm6dt_k
  2105     real(RP) :: mm2i,mm1i,m1i,m2i,m4i,m5i,m7i,m8i
  2106     real(RP) :: mm2j,mm1j,m1j,m2j,m4j,m5j,m7j,m8j
  2107     real(RP) :: mm1p5i,mm0p5i,m0p5i,m1p5i,m2p5i,m3p5i,m4p5i,m5p5i,m6p5i
  2108     real(RP) :: mm1p5j,mm0p5j,m0p5j,m1p5j,m2p5j,m3p5j,m4p5j,m5p5j,m6p5j
  2109     real(RP) :: dm6dt_nc_i,dm3dt_nc_i,dm0dt_nc_i
  2110     real(RP) :: dm6dt_nc_j,dm3dt_nc_j,dm0dt_nc_j
  2111     real(RP) :: dm6dt_nc_k,dm3dt_nc_k,dm0dt_nc_k
  2112     real(RP) :: dm6dt_fm_i,dm3dt_fm_i,dm0dt_fm_i
  2113     real(RP) :: dm6dt_fm_j,dm3dt_fm_j,dm0dt_fm_j
  2114     real(RP) :: dm6dt_fm_k,dm3dt_fm_k,dm0dt_fm_k
  2115     real(RP) :: bbb,gamma1,gamma2,alpha,beta
  2142     mm2i=m0i*dgi**(-2._rp)*exp( 2.0_rp*(log(sgi)**2._rp))
  2143     mm1i=m0i*dgi**(-1._rp)*exp( 0.5_rp*(log(sgi)**2._rp))
  2144     m1i =m0i*dgi**  1._rp *exp( 0.5_rp*(log(sgi)**2._rp))
  2145     m2i =m0i*dgi**  2._rp *exp( 2.0_rp*(log(sgi)**2._rp))
  2147     m4i =m0i*dgi**  4._rp *exp( 8.0_rp*(log(sgi)**2._rp))
  2148     m5i =m0i*dgi**  5._rp *exp(12.5_rp*(log(sgi)**2._rp))
  2150     m7i =m0i*dgi**  7._rp *exp(24.5_rp*(log(sgi)**2._rp))
  2152     mm2j=m0j*dgj**(-2._rp)*exp( 2.0_rp*(log(sgj)**2._rp))
  2153     mm1j=m0j*dgj**(-1._rp)*exp( 0.5_rp*(log(sgj)**2._rp))
  2154     m1j =m0j*dgj**  1._rp *exp( 0.5_rp*(log(sgj)**2._rp))
  2155     m2j =m0j*dgj**  2._rp *exp( 2.0_rp*(log(sgj)**2._rp))
  2157     m4j =m0j*dgj**  4._rp *exp( 8.0_rp*(log(sgj)**2._rp))
  2158     m5j =m0j*dgj**  5._rp *exp(12.5_rp*(log(sgj)**2._rp))
  2160     m7j =m0j*dgj**  7._rp *exp(24.5_rp*(log(sgj)**2._rp))
  2162     dm0dt_nc_i = -rknc*(2._rp*m0i*m0j+ mm1i*m1j             &
  2164                    +2.492e-2_rp*lambda*(m0i *mm1j+m1i*mm2j  &
  2165                                       + m0j *mm1i+m1j*mm2i) )
  2166     dm3dt_nc_i = -rknc*(2._rp*m3i*m0j+ m2i *m1j             &
  2168                    +2.492e-2_rp*lambda*(m3i *mm1j+m4i*mm2j  &
  2169                                       + m0j *m2i +m1j*m1i ) )
  2170     dm6dt_nc_i = -rknc*(2._rp*m6i*m0j+ m5i *m1j             &
  2172                    +2.492e-2_rp*lambda*(m6i *mm1j+m7i*mm2j  &
  2173                                       + m0j *m5i +m1j*m4i ) )
  2174     dm0dt_nc_j = -rknc*(2._rp*m0i*m0j+ mm1i*m1j             & 
  2176                    +2.492e-2_rp*lambda*(m0i *mm1j+m1i*mm2j  &
  2177                                       + m0j *mm1i+m1j*mm2i) )
  2178     dm3dt_nc_j = -rknc*(2._rp*m0i*m3j+ mm1i*m4j             &
  2180                    +2.492e-2_rp*lambda*(m0i *m2j +m1i*m1j   &
  2181                                       + m3j *mm1i+m4j*mm2i) )
  2182     dm6dt_nc_j = -rknc*(2._rp*m0i*m6j+ mm1i*m7j             &
  2184                    +2.492e-2_rp*lambda*(m0i *m5j +m1i*m4j   &
  2185                                       + m6j *mm1i+m7j*mm2i) )
  2186     dm0dt_nc_k = -dm0dt_nc_i
  2187     dm3dt_nc_k = -dm3dt_nc_i-dm3dt_nc_j
  2188     dm6dt_nc_k = -dm6dt_nc_i-dm6dt_nc_j                     &
  2189             +2._rp*rknc*(2._rp*m3i*m3j+ m2i *m4j            &
  2191                    +2.492e-2_rp*lambda*(m3i *m2j +m4i*m1j   &
  2192                                       + m3j *m2i +m4j*m1i)  )
  2196     mm1p5i=m0i*dgi**(-1.5_rp)*exp( 1.125_rp*(log(sgi)**2._rp) ) 
  2197     mm0p5i=m0i*dgi**(-0.5_rp)*exp( 0.125_rp*(log(sgi)**2._rp) ) 
  2198     m0p5i =m0i*dgi**  0.5_rp *exp( 0.125_rp*(log(sgi)**2._rp) ) 
  2199     m1p5i =m0i*dgi**  1.5_rp *exp( 1.125_rp*(log(sgi)**2._rp) ) 
  2200     m2p5i =m0i*dgi**  2.5_rp *exp( 3.125_rp*(log(sgi)**2._rp) ) 
  2201     m3p5i =m0i*dgi**  3.5_rp *exp( 6.125_rp*(log(sgi)**2._rp) ) 
  2202     m4p5i =m0i*dgi**  4.5_rp *exp(10.125_rp*(log(sgi)**2._rp) ) 
  2203     m5p5i =m0i*dgi**  5.5_rp *exp(15.125_rp*(log(sgi)**2._rp) ) 
  2204     m6p5i =m0i*dgi**  6.5_rp *exp(21.125_rp*(log(sgi)**2._rp) ) 
  2205     m8i   =m0i*dgi**   8._rp *exp(32.000_rp*(log(sgi)**2._rp) ) 
  2208     mm1p5j=m0j*dgj**(-1.5_rp)*exp( 1.125_rp*(log(sgj)**2._rp) ) 
  2209     mm0p5j=m0j*dgj**(-0.5_rp)*exp( 0.125_rp*(log(sgj)**2._rp) ) 
  2210     m0p5j =m0j*dgj**  0.5_rp *exp( 0.125_rp*(log(sgj)**2._rp) ) 
  2211     m1p5j =m0j*dgj**  1.5_rp *exp( 1.125_rp*(log(sgj)**2._rp) ) 
  2212     m2p5j =m0j*dgj**  2.5_rp *exp( 3.125_rp*(log(sgj)**2._rp) ) 
  2213     m3p5j =m0j*dgj**  3.5_rp *exp( 6.125_rp*(log(sgj)**2._rp) ) 
  2214     m4p5j =m0j*dgj**  4.5_rp *exp(10.125_rp*(log(sgj)**2._rp) ) 
  2215     m5p5j =m0j*dgj**  5.5_rp *exp(15.125_rp*(log(sgj)**2._rp) ) 
  2216     m6p5j =m0j*dgj**  6.5_rp *exp(21.125_rp*(log(sgj)**2._rp) ) 
  2217     m8j   =m0j*dgj**   8._rp *exp(32.000_rp*(log(sgj)**2._rp) ) 
  2221     beta    =(1._rp-(sqrt(1._rp+alpha**3._rp)/(1._rp+sqrt(alpha**3._rp))))/(1._rp-1._rp/sqrt(2._rp))
  2222     gamma1  =(sgi        +alpha*sgj       )/(1._rp+alpha)
  2223     gamma2  =(sgi**2._rp +alpha*sgj**2._rp)/(1._rp+alpha)
  2224     bbb     = 1._rp+1.2_rp*beta*dexp(
real(-2._rp*gamma1,kind=
dp))-0.646_rp*beta*dexp(
real(-0.35_rp*gamma2,kind=
dp))
  2226     dm0dt_fm_i=-bbb*rkfm*(                                 &
  2227                 m0i *m0p5j +m2i *mm1p5j +2._rp*m1i *mm0p5j &
  2228               + m0j *m0p5i +m2j *mm1p5i +2._rp*m1j *mm0p5i )
  2229     dm3dt_fm_i=-bbb*rkfm*(                                 &
  2230                 m3i *m0p5j +m5i *mm1p5j +2._rp*m4i *mm0p5j &
  2231               + m0j *m3p5i +m2j *m1p5i  +2._rp*m1j *m2p5i  )
  2232     dm6dt_fm_i=-bbb*rkfm*(                                 &
  2233                 m6i *m0p5j +m8i *mm1p5j +2._rp*m7i *mm0p5j &
  2234               + m0j *m6p5i +m2j *m4p5i  +2._rp*m1j *m5p5i  )
  2235     dm0dt_fm_j= dm0dt_fm_i
  2236     dm3dt_fm_j=-bbb*rkfm*(                                 &
  2237                 m0i *m3p5j +m2i *m1p5j  +2._rp*m1i *m2p5j  &
  2238               + m3j *m0p5i +m5j *mm1p5i +2._rp*m4j *mm0p5i )
  2239     dm6dt_fm_j=-bbb*rkfm*(                                 &
  2240                 m0i *m6p5j +m2i *m4p5j  +2._rp*m1i *m5p5j  &
  2241               + m6j *m0p5i +m8j *mm1p5i +2._rp*m7j *mm0p5i )
  2242     dm0dt_fm_k=-dm0dt_fm_i
  2243     dm3dt_fm_k=-dm3dt_fm_i-dm3dt_fm_j
  2244     dm6dt_fm_k=-dm6dt_fm_i-dm6dt_fm_j                      &
  2245         +2._rp * bbb*rkfm*(                                &
  2246                 m3i *m3p5j +m5i *m1p5j  +2._rp*m4i *m2p5j  &
  2247               + m3j *m3p5i +m5j *m1p5i  +2._rp*m4j *m2p5i  )
  2250     if (dm0dt_fm_i+dm0dt_nc_i/=0._rp) &
  2251     dm0dt_i  = dm0dt_fm_i*dm0dt_nc_i/(dm0dt_fm_i+dm0dt_nc_i)
  2252     if (dm3dt_fm_i+dm3dt_nc_i/=0._rp) &
  2253     dm3dt_i  = dm3dt_fm_i*dm3dt_nc_i/(dm3dt_fm_i+dm3dt_nc_i)
  2254     if (dm6dt_fm_i+dm6dt_nc_i/=0._rp) &
  2255     dm6dt_i  = dm6dt_fm_i*dm6dt_nc_i/(dm6dt_fm_i+dm6dt_nc_i)
  2256     if (dm0dt_fm_j+dm0dt_nc_j/=0._rp) &
  2257     dm0dt_j  = dm0dt_fm_j*dm0dt_nc_j/(dm0dt_fm_j+dm0dt_nc_j)
  2258     if (dm3dt_fm_j+dm3dt_nc_j/=0._rp) &
  2259     dm3dt_j  = dm3dt_fm_j*dm3dt_nc_j/(dm3dt_fm_j+dm3dt_nc_j)
  2260     if (dm6dt_fm_j+dm6dt_nc_j/=0._rp) &
  2261     dm6dt_j  = dm6dt_fm_j*dm6dt_nc_j/(dm6dt_fm_j+dm6dt_nc_j)
  2262     if (dm0dt_fm_k+dm0dt_nc_k/=0._rp) &
  2263     dm0dt_k  = dm0dt_fm_k*dm0dt_nc_k/(dm0dt_fm_k+dm0dt_nc_k)
  2266     if (dm6dt_fm_k+dm6dt_nc_k/=0._rp) &
  2267     dm6dt_k  = dm6dt_fm_k*dm6dt_nc_k/(dm6dt_fm_k+dm6dt_nc_k)
  2269     dm0i     = dm0dt_i * dtrest
  2270     dm3i     = dm3dt_i * dtrest
  2271     dm6i     = dm6dt_i * dtrest
  2272     dmsi     = dm3i    * rhoi *pi6*1.e12_rp 
  2273     dm0j     = dm0dt_j * dtrest
  2274     dm3j     = dm3dt_j * dtrest
  2275     dm6j     = dm6dt_j * dtrest
  2276     dmsj     = dm3j    * rhoj *pi6*1.e12_rp 
  2277     dm0k     = dm0dt_k * dtrest
  2279     dm6k     = dm6dt_k * dtrest
  2284   end subroutine aero_inter
  2288   subroutine diag_ds6(m0,m3,m6,m2,dg,sg,dm2)
  2290     real(RP)            :: m0,m2,m3,m6,dg,sg,m3_bar,m6_bar
  2291     real(RP)            :: m2_new,m2_old,dm2
  2292     real(RP), 
parameter :: sgmax=2.5_rp
  2293     real(RP), 
parameter :: rk1=3._rp
  2294     real(RP), 
parameter :: rk2=6._rp
  2295     real(RP), 
parameter :: ratio  =rk1/rk2
  2296     real(RP), 
parameter :: rk1_hat=1._rp/(ratio*(rk2-rk1))
  2297     real(RP), 
parameter :: rk2_hat=ratio/(rk1-rk2)
  2298     real(DP), 
parameter :: tiny=1.e-50_dp
  2301     if (m0 <= tiny .OR. m3 <= tiny .OR. m6 <= tiny) 
then  2313     dg     = m3_bar**rk1_hat*m6_bar**rk2_hat
  2315     if (m3_bar/m6_bar**ratio < 1._rp) 
then   2316       sg     = exp(dsqrt(
real(2._RP/(rk1*(rk1-rk2))  &
  2317                   * log(m3_bar/m6_bar**ratio), kind=
dp )))
  2318       m2     = m0*dg**2._rp*exp(2._rp*(log(sg)**2._rp))
  2322     if (sg > sgmax) 
then  2325       call diag_d2(m0,m3,sg, & 
  2330     if (m3_bar/m6_bar**ratio >= 1._rp) 
then   2332       call diag_d2(m0,m3,sg, & 
  2338     dm2    = m2_old - m2_new 
  2342   end subroutine diag_ds6
  2344   subroutine trans_ccn(aerosol_procs, aerosol_activ, t_ccn, t_cn,  &
  2345                     n_ctg, n_kap_max, n_siz_max, n_atr,         &
  2346                     ic_out, ia_m0, ia_m2, ia_m3, ik_out, n_siz, &
  2347 !                    rnum_out, nbins_out, dlog10d_out)
  2348                     rnum_out, nbins_out)
  2352     integer, 
intent(in) :: n_ctg, n_kap_max, n_siz_max, n_atr
  2353     integer, 
intent(in) :: ic_out, ik_out
  2354     integer, 
intent(in) :: ia_m0, ia_m2, ia_m3
  2355     integer, 
intent(in) :: n_siz(n_ctg)
  2356     real(RP),
intent(in) :: aerosol_procs(n_atr,n_siz_max,n_kap_max,n_ctg)
  2357     real(RP),
intent(in) :: aerosol_activ(n_atr,n_siz_max,n_kap_max,n_ctg)
  2358     integer, 
intent(in) :: nbins_out
  2359     real(RP),
intent(out):: rnum_out(nbins_out)
  2361     real(RP),
intent(out):: t_ccn 
  2362     real(RP),
intent(out):: t_cn  
  2384     do is0 = 1, n_siz(ic_out)
  2406       t_ccn = t_ccn + aerosol_activ(ia_m0,is0,ik_out,ic_out)
  2407       t_cn  = t_cn  + aerosol_procs(ia_m0,is0,ik_out,ic_out)
  2412   end subroutine trans_ccn
 real(rp), public atmos_phy_ae_kajino13_h2so4dt
 
module atmosphere / saturation 
 
logical, public atmos_phy_ae_kajino13_flag_coag
 
integer, parameter, public n_ae
 
character(len=h_short), dimension(:), allocatable, public atmos_phy_ae_kajino13_name
 
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. 
 
subroutine, public atmos_phy_ae_kajino13_tracer_setup(QA_AE)
Tracer setup. 
 
character(len=h_mid), dimension(:), allocatable, public atmos_phy_ae_kajino13_desc
 
subroutine, public atmos_phy_ae_kajino13_effective_radius(KA, IA, JA, QA_AE, QTRC, RH, Re)
Calculate Effective Radius. 
 
integer, public io_fid_conf
Config file ID. 
 
module atmosphere / aerosol 
 
character(len=h_short), dimension(:), allocatable, public atmos_phy_ae_kajino13_unit
 
logical, public atmos_phy_ae_kajino13_flag_npf
 
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)
 
logical, public atmos_phy_ae_kajino13_flag_ccn_interactive
 
integer, public ke
end point of inner domain: z, local 
 
real(rp), public atmos_phy_ae_kajino13_logk_aenucl
 
subroutine, public atmos_phy_ae_kajino13_setup
Setup. 
 
real(rp), public atmos_phy_ae_kajino13_ocgasdt
 
integer, public ks
start point of inner domain: z, local 
 
real(rp), public atmos_phy_ae_kajino13_c_kappa
 
subroutine, public prc_abort
Abort Process. 
 
subroutine, public atmos_phy_ae_kajino13_negative_fixer(KA, KS, KE, IA, IS, IE, JA, JS, JE, QA_AE, QTRC)
 
logical, public atmos_phy_ae_kajino13_flag_regeneration
 
subroutine, public prof_rapstart(rapname_base, level)
Start raptime. 
 
real(rp), public atmos_phy_ae_kajino13_sg_reg
 
module atmosphere / physics / aerosol / Kajino13 
 
integer, public ka
of whole cells: z, local, with HALO
 
real(rp), public const_pi
pi 
 
logical, public atmos_phy_ae_kajino13_flag_cond
 
real(rp), public const_mdry
mass weight (dry air) [g/mol] 
 
subroutine, public prof_rapend(rapname_base, level)
Save raptime. 
 
real(rp), public atmos_phy_ae_kajino13_dg_reg
 
integer, parameter, public rp
 
integer, parameter, public dp