source: palm/trunk/SOURCE/salsa_mod.f90 @ 3582

Last change on this file since 3582 was 3582, checked in by suehring, 3 years ago

Merge branch salsa with trunk

  • Property svn:keywords set to Id
File size: 459.4 KB
Line 
1!> @file salsa_mod.f90
2!--------------------------------------------------------------------------------!
3! This file is part of PALM-4U.
4!
5! PALM-4U is free software: you can redistribute it and/or modify it under the
6! terms of the GNU General Public License as published by the Free Software
7! Foundation, either version 3 of the License, or (at your option) any later
8! version.
9!
10! PALM-4U is distributed in the hope that it will be useful, but WITHOUT ANY
11! WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR
12! A PARTICULAR PURPOSE.  See the GNU General Public License for more details.
13!
14! You should have received a copy of the GNU General Public License along with
15! PALM. If not, see <http://www.gnu.org/licenses/>.
16!
17! Copyright 2018-2018 University of Helsinki
18! Copyright 1997-2018 Leibniz Universitaet Hannover
19!--------------------------------------------------------------------------------!
20!
21! Current revisions:
22! -----------------
23! - Moved the control parameter "salsa" from salsa_mod.f90 to control_parameters
24! - Updated salsa_rrd_local and salsa_wrd_local
25! - Add target attribute
26! - Revise initialization in case of restarts
27! - Revise masked data output
28!
29! Former revisions:
30! -----------------
31! $Id: salsa_mod.f90 3582 2018-11-29 19:16:36Z suehring $
32! missing comma separator inserted
33!
34! 3483 2018-11-02 14:19:26Z raasch
35! bugfix: directives added to allow compilation without netCDF
36!
37! 3481 2018-11-02 09:14:13Z raasch
38! temporary variable cc introduced to circumvent a possible Intel18 compiler bug
39! related to contiguous/non-contguous pointer/target attributes
40!
41! 3473 2018-10-30 20:50:15Z suehring
42! NetCDF input routine renamed
43!
44! 3467 2018-10-30 19:05:21Z suehring
45! Initial revision
46!
47! 3412 2018-10-24 07:25:57Z monakurppa
48!
49! Authors:
50! --------
51! @author Mona Kurppa (University of Helsinki)
52!
53!
54! Description:
55! ------------
56!> Sectional aerosol module for large scale applications SALSA
57!> (Kokkola et al., 2008, ACP 8, 2469-2483). Solves the aerosol number and mass
58!> concentration as well as chemical composition. Includes aerosol dynamic
59!> processes: nucleation, condensation/evaporation of vapours, coagulation and
60!> deposition on tree leaves, ground and roofs.
61!> Implementation is based on formulations implemented in UCLALES-SALSA except
62!> for deposition which is based on parametrisations by Zhang et al. (2001,
63!> Atmos. Environ. 35, 549-560) or Petroff&Zhang (2010, Geosci. Model Dev. 3,
64!> 753-769)
65!>
66!> @todo Implement turbulent inflow of aerosols in inflow_turbulence.
67!> @todo Deposition on subgrid scale vegetation
68!> @todo Deposition on vegetation calculated by default for deciduous broadleaf
69!>       trees
70!> @todo Revise masked data output. There is a potential bug in case of
71!>       terrain-following masked output, according to data_output_mask.
72!> @todo There are now improved interfaces for NetCDF data input which can be
73!>       used instead of get variable etc.
74!------------------------------------------------------------------------------!
75 MODULE salsa_mod
76
77    USE basic_constants_and_equations_mod,                                     &
78        ONLY:  c_p, g, p_0, pi, r_d
79 
80    USE chemistry_model_mod,                                                   &
81        ONLY:  chem_species, nspec, nvar, spc_names
82
83    USE chem_modules,                                                          &
84        ONLY:  call_chem_at_all_substeps, chem_gasphase_on
85
86    USE control_parameters
87       
88    USE indices,                                                               &
89        ONLY:  nbgp, nx, nxl, nxlg, nxr, nxrg, ny, nyn, nyng, nys, nysg, nzb,  &
90               nzb_s_inner, nz, nzt, wall_flags_0
91     
92    USE kinds
93   
94    USE pegrid
95   
96    USE salsa_util_mod
97
98    IMPLICIT NONE
99!
100!-- SALSA constants:
101!
102!-- Local constants:
103    INTEGER(iwp), PARAMETER ::  ngast   = 5 !< total number of gaseous tracers:
104                                            !< 1 = H2SO4, 2 = HNO3, 3 = NH3,
105                                            !< 4 = OCNV (non-volatile OC),
106                                            !< 5 = OCSV (semi-volatile) 
107    INTEGER(iwp), PARAMETER ::  nmod    = 7 !< number of modes for initialising
108                                            !< the aerosol size distribution                                             
109    INTEGER(iwp), PARAMETER ::  nreg    = 2 !< Number of main size subranges
110    INTEGER(iwp), PARAMETER ::  maxspec = 7 !< Max. number of aerosol species
111!   
112!-- Universal constants
113    REAL(wp), PARAMETER ::  abo    = 1.380662E-23_wp  !< Boltzmann constant (J/K)
114    REAL(wp), PARAMETER ::  alv    = 2.260E+6_wp      !< latent heat for H2O
115                                                      !< vaporisation (J/kg)
116    REAL(wp), PARAMETER ::  alv_d_rv  = 4896.96865_wp !< alv / rv
117    REAL(wp), PARAMETER ::  am_airmol = 4.8096E-26_wp !< Average mass of one air
118                                                      !< molecule (Jacobson,
119                                                      !< 2005, Eq. 2.3)                                                   
120    REAL(wp), PARAMETER ::  api6   = 0.5235988_wp     !< pi / 6   
121    REAL(wp), PARAMETER ::  argas  = 8.314409_wp      !< Gas constant (J/(mol K))
122    REAL(wp), PARAMETER ::  argas_d_cpd = 8.281283865E-3_wp !< argas per cpd
123    REAL(wp), PARAMETER ::  avo    = 6.02214E+23_wp   !< Avogadro constant (1/mol)
124    REAL(wp), PARAMETER ::  d_sa   = 5.539376964394570E-10_wp !< diameter of
125                                                      !< condensing sulphuric
126                                                      !< acid molecule (m) 
127    REAL(wp), PARAMETER ::  for_ppm_to_nconc =  7.243016311E+16_wp !<
128                                                 !< ppm * avo / R (K/(Pa*m3))
129    REAL(wp), PARAMETER ::  epsoc  = 0.15_wp          !< water uptake of organic
130                                                      !< material     
131    REAL(wp), PARAMETER ::  mclim  = 1.0E-23_wp    !< mass concentration min
132                                                   !< limit for aerosols (kg/m3)                                                   
133    REAL(wp), PARAMETER ::  n3     = 158.79_wp !< Number of H2SO4 molecules in
134                                               !< 3 nm cluster if d_sa=5.54e-10m
135    REAL(wp), PARAMETER ::  nclim  = 1.0_wp    !< number concentration min limit
136                                               !< for aerosols and gases (#/m3)
137    REAL(wp), PARAMETER ::  surfw0 = 0.073_wp  !< surface tension of pure water
138                                               !< at ~ 293 K (J/m2)   
139    REAL(wp), PARAMETER ::  vclim  = 1.0E-24_wp    !< volume concentration min
140                                                   !< limit for aerosols (m3/m3)                                           
141!-- Molar masses in kg/mol
142    REAL(wp), PARAMETER ::  ambc   = 12.0E-3_wp     !< black carbon (BC)
143    REAL(wp), PARAMETER ::  amdair = 28.970E-3_wp   !< dry air
144    REAL(wp), PARAMETER ::  amdu   = 100.E-3_wp     !< mineral dust
145    REAL(wp), PARAMETER ::  amh2o  = 18.0154E-3_wp  !< H2O
146    REAL(wp), PARAMETER ::  amh2so4  = 98.06E-3_wp  !< H2SO4
147    REAL(wp), PARAMETER ::  amhno3 = 63.01E-3_wp    !< HNO3
148    REAL(wp), PARAMETER ::  amn2o  = 44.013E-3_wp   !< N2O
149    REAL(wp), PARAMETER ::  amnh3  = 17.031E-3_wp   !< NH3
150    REAL(wp), PARAMETER ::  amo2   = 31.9988E-3_wp  !< O2
151    REAL(wp), PARAMETER ::  amo3   = 47.998E-3_wp   !< O3
152    REAL(wp), PARAMETER ::  amoc   = 150.E-3_wp     !< organic carbon (OC)
153    REAL(wp), PARAMETER ::  amss   = 58.44E-3_wp    !< sea salt (NaCl)
154!-- Densities in kg/m3
155    REAL(wp), PARAMETER ::  arhobc     = 2000.0_wp !< black carbon
156    REAL(wp), PARAMETER ::  arhodu     = 2650.0_wp !< mineral dust
157    REAL(wp), PARAMETER ::  arhoh2o    = 1000.0_wp !< H2O
158    REAL(wp), PARAMETER ::  arhoh2so4  = 1830.0_wp !< SO4
159    REAL(wp), PARAMETER ::  arhohno3   = 1479.0_wp !< HNO3
160    REAL(wp), PARAMETER ::  arhonh3    = 1530.0_wp !< NH3
161    REAL(wp), PARAMETER ::  arhooc     = 2000.0_wp !< organic carbon
162    REAL(wp), PARAMETER ::  arhoss     = 2165.0_wp !< sea salt (NaCl)
163!-- Volume of molecule in m3/#
164    REAL(wp), PARAMETER ::  amvh2o   = amh2o /avo / arhoh2o      !< H2O
165    REAL(wp), PARAMETER ::  amvh2so4 = amh2so4 / avo / arhoh2so4 !< SO4
166    REAL(wp), PARAMETER ::  amvhno3  = amhno3 / avo / arhohno3   !< HNO3
167    REAL(wp), PARAMETER ::  amvnh3   = amnh3 / avo / arhonh3     !< NH3 
168    REAL(wp), PARAMETER ::  amvoc    = amoc / avo / arhooc       !< OC
169    REAL(wp), PARAMETER ::  amvss    = amss / avo / arhoss       !< sea salt
170   
171!
172!-- SALSA switches:
173    INTEGER(iwp) ::  nj3 = 1 !< J3 parametrization (nucleation)
174                             !< 1 = condensational sink (Kerminen&Kulmala, 2002)
175                             !< 2 = coagulational sink (Lehtinen et al. 2007)
176                             !< 3 = coagS+self-coagulation (Anttila et al. 2010)                                       
177    INTEGER(iwp) ::  nsnucl = 0 !< Choice of the nucleation scheme:
178                                !< 0 = off   
179                                !< 1 = binary nucleation
180                                !< 2 = activation type nucleation
181                                !< 3 = kinetic nucleation
182                                !< 4 = ternary nucleation
183                                !< 5 = nucleation with ORGANICs
184                                !< 6 = activation type of nucleation with
185                                !<     H2SO4+ORG
186                                !< 7 = heteromolecular nucleation with H2SO4*ORG
187                                !< 8 = homomolecular nucleation of  H2SO4 +
188                                !<     heteromolecular nucleation with H2SO4*ORG
189                                !< 9 = homomolecular nucleation of  H2SO4 and ORG
190                                !<     +heteromolecular nucleation with H2SO4*ORG
191    LOGICAL ::  advect_particle_water = .TRUE.  !< advect water concentration of
192                                                !< particles                               
193    LOGICAL ::  decycle_lr            = .FALSE. !< Undo cyclic boundary
194                                                !< conditions: left and right
195    LOGICAL ::  decycle_ns            = .FALSE. !< north and south boundaries
196    LOGICAL ::  feedback_to_palm      = .FALSE. !< allow feedback due to
197                                                !< hydration and/or condensation
198                                                !< of H20
199    LOGICAL ::  no_insoluble          = .FALSE. !< Switch to exclude insoluble 
200                                                !< chemical components
201    LOGICAL ::  read_restart_data_salsa = .FALSE. !< read restart data for salsa
202    LOGICAL ::  salsa_gases_from_chem = .FALSE.   !< Transfer the gaseous
203                                                  !< components to SALSA from 
204                                                  !< from chemistry model
205    LOGICAL ::  van_der_waals_coagc   = .FALSE.   !< Enhancement of coagulation
206                                                  !< kernel by van der Waals and
207                                                  !< viscous forces
208    LOGICAL ::  write_binary_salsa    = .FALSE.   !< read binary for salsa
209!-- Process switches: nl* is read from the NAMELIST and is NOT changed.
210!--                   ls* is the switch used and will get the value of nl*
211!--                       except for special circumstances (spinup period etc.)
212    LOGICAL ::  nlcoag       = .FALSE. !< Coagulation master switch
213    LOGICAL ::  lscoag       = .FALSE. !<
214    LOGICAL ::  nlcnd        = .FALSE. !< Condensation master switch
215    LOGICAL ::  lscnd        = .FALSE. !<
216    LOGICAL ::  nlcndgas     = .FALSE. !< Condensation of precursor gases
217    LOGICAL ::  lscndgas     = .FALSE. !<
218    LOGICAL ::  nlcndh2oae   = .FALSE. !< Condensation of H2O on aerosol
219    LOGICAL ::  lscndh2oae   = .FALSE. !< particles (FALSE -> equilibrium calc.)
220    LOGICAL ::  nldepo       = .FALSE. !< Deposition master switch
221    LOGICAL ::  lsdepo       = .FALSE. !<
222    LOGICAL ::  nldepo_topo  = .FALSE. !< Deposition on vegetation master switch
223    LOGICAL ::  lsdepo_topo  = .FALSE. !<
224    LOGICAL ::  nldepo_vege  = .FALSE. !< Deposition on walls master switch
225    LOGICAL ::  lsdepo_vege  = .FALSE. !<
226    LOGICAL ::  nldistupdate = .TRUE.  !< Size distribution update master switch                                     
227    LOGICAL ::  lsdistupdate = .FALSE. !<                                     
228!
229!-- SALSA variables:
230    CHARACTER (LEN=20) ::  bc_salsa_b = 'neumann'   !< bottom boundary condition                                     
231    CHARACTER (LEN=20) ::  bc_salsa_t = 'neumann'   !< top boundary condition
232    CHARACTER (LEN=20) ::  depo_vege_type = 'zhang2001' !< or 'petroff2010'
233    CHARACTER (LEN=20) ::  depo_topo_type = 'zhang2001' !< or 'petroff2010'
234    CHARACTER (LEN=20), DIMENSION(4) ::  decycle_method = & 
235                             (/'dirichlet','dirichlet','dirichlet','dirichlet'/)
236                                 !< Decycling method at horizontal boundaries,
237                                 !< 1=left, 2=right, 3=south, 4=north
238                                 !< dirichlet = initial size distribution and
239                                 !< chemical composition set for the ghost and
240                                 !< first three layers
241                                 !< neumann = zero gradient
242    CHARACTER (LEN=3), DIMENSION(maxspec) ::  listspec = &  !< Active aerosols
243                                   (/'SO4','   ','   ','   ','   ','   ','   '/)
244    CHARACTER (LEN=20) ::  salsa_source_mode = 'no_source' 
245                                                    !< 'read_from_file',
246                                                    !< 'constant' or 'no_source'                                   
247    INTEGER(iwp) ::  dots_salsa = 0  !< starting index for salsa-timeseries
248    INTEGER(iwp) ::  fn1a = 1    !< last index for bin subranges:  subrange 1a
249    INTEGER(iwp) ::  fn2a = 1    !<                              subrange 2a
250    INTEGER(iwp) ::  fn2b = 1    !<                              subrange 2b
251    INTEGER(iwp), DIMENSION(ngast) ::  gas_index_chem = (/ 1, 1, 1, 1, 1/) !<
252                                 !< Index of gaseous compounds in the chemistry
253                                 !< model. In SALSA, 1 = H2SO4, 2 = HNO3,
254                                 !< 3 = NH3, 4 = OCNV, 5 = OCSV
255    INTEGER(iwp) ::  ibc_salsa_b !<
256    INTEGER(iwp) ::  ibc_salsa_t !<
257    INTEGER(iwp) ::  igctyp = 0  !< Initial gas concentration type
258                                 !< 0 = uniform (use H2SO4_init, HNO3_init,
259                                 !<     NH3_init, OCNV_init and OCSV_init)
260                                 !< 1 = read vertical profile from an input file 
261    INTEGER(iwp) ::  in1a = 1    !< start index for bin subranges: subrange 1a
262    INTEGER(iwp) ::  in2a = 1    !<                              subrange 2a
263    INTEGER(iwp) ::  in2b = 1    !<                              subrange 2b
264    INTEGER(iwp) ::  isdtyp = 0  !< Initial size distribution type
265                                 !< 0 = uniform
266                                 !< 1 = read vertical profile of the mode number
267                                 !<     concentration from an input file 
268    INTEGER(iwp) ::  ibc  = -1 !< Indice for: black carbon (BC)
269    INTEGER(iwp) ::  idu  = -1 !< dust
270    INTEGER(iwp) ::  inh  = -1 !< NH3
271    INTEGER(iwp) ::  ino  = -1 !< HNO3   
272    INTEGER(iwp) ::  ioc  = -1 !< organic carbon (OC)
273    INTEGER(iwp) ::  iso4 = -1 !< SO4 or H2SO4   
274    INTEGER(iwp) ::  iss  = -1 !< sea salt
275    INTEGER(iwp) ::  lod_aero = 0   !< level of detail for aerosol emissions
276    INTEGER(iwp) ::  lod_gases = 0  !< level of detail for gaseous emissions   
277    INTEGER(iwp), DIMENSION(nreg) ::  nbin = (/ 3, 7/)    !< Number of size bins
278                                               !< for each aerosol size subrange
279    INTEGER(iwp) ::  nbins = 1  !< total number of size bins
280    INTEGER(iwp) ::  ncc   = 1  !< number of chemical components used     
281    INTEGER(iwp) ::  ncc_tot = 1!< total number of chemical compounds (ncc+1
282                                !< if particle water is advected)
283    REAL(wp) ::  act_coeff = 1.0E-7_wp     !< Activation coefficient
284    REAL(wp) ::  aerosol_source = 0.0_wp   !< Constant aerosol flux (#/(m3*s))
285    REAL(wp), ALLOCATABLE, DIMENSION(:,:) ::  emission_mass_fracs  !< array for
286                                    !< aerosol composition per emission category
287                                    !< 1:SO4 2:OC 3:BC 4:DU 5:SS 6:NO 7:NH 
288    REAL(wp) ::  dt_salsa  = 0.00001_wp    !< Time step of SALSA
289    REAL(wp) ::  H2SO4_init = nclim        !< Init value for sulphuric acid gas
290    REAL(wp) ::  HNO3_init  = nclim        !< Init value for nitric acid gas
291    REAL(wp) ::  last_salsa_time = 0.0_wp  !< time of the previous salsa
292                                           !< timestep
293    REAL(wp) ::  nf2a = 1.0_wp             !< Number fraction allocated to a-
294                                           !< bins in subrange 2
295                                           !< (b-bins will get 1-nf2a)   
296    REAL(wp) ::  NH3_init  = nclim         !< Init value for ammonia gas
297    REAL(wp) ::  OCNV_init = nclim         !< Init value for non-volatile
298                                           !< organic gases
299    REAL(wp) ::  OCSV_init = nclim         !< Init value for semi-volatile
300                                           !< organic gases
301    REAL(wp), DIMENSION(nreg+1) ::  reglim = & !< Min&max diameters of size subranges
302                                 (/ 3.0E-9_wp, 5.0E-8_wp, 1.0E-5_wp/)
303    REAL(wp) ::  rhlim = 1.20_wp    !< RH limit in %/100. Prevents
304                                    !< unrealistically high RH in condensation                           
305    REAL(wp) ::  skip_time_do_salsa = 0.0_wp !< Starting time of SALSA (s)
306!-- Initial log-normal size distribution: mode diameter (dpg, micrometres),
307!-- standard deviation (sigmag) and concentration (n_lognorm, #/cm3)
308    REAL(wp), DIMENSION(nmod) ::  dpg   = (/0.013_wp, 0.054_wp, 0.86_wp,       &
309                                            0.2_wp, 0.2_wp, 0.2_wp, 0.2_wp/) 
310    REAL(wp), DIMENSION(nmod) ::  sigmag  = (/1.8_wp, 2.16_wp, 2.21_wp,        &
311                                              2.0_wp, 2.0_wp, 2.0_wp, 2.0_wp/) 
312    REAL(wp), DIMENSION(nmod) ::  n_lognorm = (/1.04e+5_wp, 3.23E+4_wp, 5.4_wp,&
313                                                0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp/)
314!-- Initial mass fractions / chemical composition of the size distribution   
315    REAL(wp), DIMENSION(maxspec) ::  mass_fracs_a = & !< mass fractions between
316             (/1.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0/) !< aerosol species for A bins
317    REAL(wp), DIMENSION(maxspec) ::  mass_fracs_b = & !< mass fractions between
318             (/0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0/) !< aerosol species for B bins
319             
320    REAL(wp), ALLOCATABLE, DIMENSION(:) ::  bin_low_limits  !< to deliver
321                                                            !< information about
322                                                            !< the lower
323                                                            !< diameters per bin                                       
324    REAL(wp), ALLOCATABLE, DIMENSION(:) ::  nsect     !< Background number
325                                                      !< concentration per bin
326    REAL(wp), ALLOCATABLE, DIMENSION(:) ::  massacc   !< Mass accomodation
327                                                      !< coefficients per bin                                             
328!
329!-- SALSA derived datatypes:
330!
331!-- Prognostic variable: Aerosol size bin information (number (#/m3) and
332!-- mass (kg/m3) concentration) and the concentration of gaseous tracers (#/m3).
333!-- Gas tracers are contained sequentially in dimension 4 as:
334!-- 1. H2SO4, 2. HNO3, 3. NH3, 4. OCNV (non-volatile organics),
335!-- 5. OCSV (semi-volatile)
336    TYPE salsa_variable
337       REAL(wp), POINTER, DIMENSION(:,:,:), CONTIGUOUS     ::  conc
338       REAL(wp), POINTER, DIMENSION(:,:,:), CONTIGUOUS     ::  conc_p
339       REAL(wp), POINTER, DIMENSION(:,:,:), CONTIGUOUS     ::  tconc_m
340       REAL(wp), ALLOCATABLE, DIMENSION(:,:)   ::  flux_s, diss_s
341       REAL(wp), ALLOCATABLE, DIMENSION(:,:,:) ::  flux_l, diss_l
342       REAL(wp), ALLOCATABLE, DIMENSION(:)     ::  init
343       REAL(wp), ALLOCATABLE, DIMENSION(:,:,:) ::  source
344       REAL(wp), ALLOCATABLE, DIMENSION(:,:)   ::  sums_ws_l
345    END TYPE salsa_variable
346   
347!-- Map bin indices between parallel size distributions   
348    TYPE t_parallelbin
349       INTEGER(iwp) ::  cur  ! Index for current distribution
350       INTEGER(iwp) ::  par  ! Index for corresponding parallel distribution
351    END TYPE t_parallelbin
352   
353!-- Datatype used to store information about the binned size distributions of
354!-- aerosols
355    TYPE t_section
356       REAL(wp) ::  vhilim   !< bin volume at the high limit
357       REAL(wp) ::  vlolim   !< bin volume at the low limit
358       REAL(wp) ::  vratiohi !< volume ratio between the center and high limit
359       REAL(wp) ::  vratiolo !< volume ratio between the center and low limit
360       REAL(wp) ::  dmid     !< bin middle diameter (m)
361       !******************************************************
362       ! ^ Do NOT change the stuff above after initialization !
363       !******************************************************
364       REAL(wp) ::  dwet    !< Wet diameter or mean droplet diameter (m)
365       REAL(wp), DIMENSION(maxspec+1) ::  volc !< Volume concentrations
366                            !< (m^3/m^3) of aerosols + water. Since most of
367                            !< the stuff in SALSA is hard coded, these *have to
368                            !< be* in the order
369                            !< 1:SO4, 2:OC, 3:BC, 4:DU, 5:SS, 6:NO, 7:NH, 8:H2O
370       REAL(wp) ::  veqh2o  !< Equilibrium H2O concentration for each particle
371       REAL(wp) ::  numc    !< Number concentration of particles/droplets (#/m3)
372       REAL(wp) ::  core    !< Volume of dry particle
373    END TYPE t_section 
374!
375!-- Local aerosol properties in SALSA
376    TYPE(t_section), ALLOCATABLE ::  aero(:)
377!
378!-- SALSA tracers:
379!-- Tracers as x = x(k,j,i,bin). The 4th dimension contains all the size bins
380!-- sequentially for each aerosol species  + water.
381!
382!-- Prognostic tracers:
383!
384!-- Number concentration (#/m3)
385    TYPE(salsa_variable), ALLOCATABLE, DIMENSION(:), TARGET ::  aerosol_number
386    REAL(wp), ALLOCATABLE, DIMENSION(:,:,:,:), TARGET ::  nconc_1
387    REAL(wp), ALLOCATABLE, DIMENSION(:,:,:,:), TARGET ::  nconc_2
388    REAL(wp), ALLOCATABLE, DIMENSION(:,:,:,:), TARGET ::  nconc_3
389!
390!-- Mass concentration (kg/m3)
391    TYPE(salsa_variable), ALLOCATABLE, DIMENSION(:), TARGET ::  aerosol_mass
392    REAL(wp), ALLOCATABLE, DIMENSION(:,:,:,:), TARGET ::  mconc_1
393    REAL(wp), ALLOCATABLE, DIMENSION(:,:,:,:), TARGET ::  mconc_2
394    REAL(wp), ALLOCATABLE, DIMENSION(:,:,:,:), TARGET ::  mconc_3
395!
396!-- Gaseous tracers (#/m3)
397    TYPE(salsa_variable), ALLOCATABLE, DIMENSION(:), TARGET ::  salsa_gas
398    REAL(wp), ALLOCATABLE, DIMENSION(:,:,:,:), TARGET ::  gconc_1
399    REAL(wp), ALLOCATABLE, DIMENSION(:,:,:,:), TARGET ::  gconc_2
400    REAL(wp), ALLOCATABLE, DIMENSION(:,:,:,:), TARGET ::  gconc_3
401!
402!-- Diagnostic tracers
403    REAL(wp), ALLOCATABLE, DIMENSION(:,:,:,:) ::  sedim_vd !< sedimentation
404                                                           !< velocity per size
405                                                           !< bin (m/s)
406    REAL(wp), ALLOCATABLE, DIMENSION(:,:,:,:) ::  Ra_dry !< dry radius (m)
407   
408!-- Particle component index tables
409    TYPE(component_index) :: prtcl !< Contains "getIndex" which gives the index
410                                   !< for a given aerosol component name, i.e.
411                                   !< 1:SO4, 2:OC, 3:BC, 4:DU,
412                                   !< 5:SS, 6:NO, 7:NH, 8:H2O 
413!                                   
414!-- Data output arrays:
415!-- Gases:
416    REAL(wp), ALLOCATABLE, DIMENSION(:,:,:), TARGET ::  g_H2SO4_av  !< H2SO4
417    REAL(wp), ALLOCATABLE, DIMENSION(:,:,:), TARGET ::  g_HNO3_av   !< HNO3
418    REAL(wp), ALLOCATABLE, DIMENSION(:,:,:), TARGET ::  g_NH3_av    !< NH3
419    REAL(wp), ALLOCATABLE, DIMENSION(:,:,:), TARGET ::  g_OCNV_av   !< non-vola-
420                                                                    !< tile OC
421    REAL(wp), ALLOCATABLE, DIMENSION(:,:,:), TARGET ::  g_OCSV_av   !< semi-vol.
422                                                                    !< OC
423!-- Integrated:                                                                   
424    REAL(wp), ALLOCATABLE, DIMENSION(:,:,:), TARGET ::  LDSA_av  !< lung-
425                                                                 !< deposited
426                                                                 !< surface area                                                   
427    REAL(wp), ALLOCATABLE, DIMENSION(:,:,:), TARGET ::  Ntot_av  !< total number
428                                                                 !< conc.
429    REAL(wp), ALLOCATABLE, DIMENSION(:,:,:), TARGET ::  PM25_av  !< PM2.5
430    REAL(wp), ALLOCATABLE, DIMENSION(:,:,:), TARGET ::  PM10_av  !< PM10
431!-- In the particle phase:   
432    REAL(wp), ALLOCATABLE, DIMENSION(:,:,:), TARGET ::  s_BC_av  !< black carbon
433    REAL(wp), ALLOCATABLE, DIMENSION(:,:,:), TARGET ::  s_DU_av  !< dust
434    REAL(wp), ALLOCATABLE, DIMENSION(:,:,:), TARGET ::  s_H2O_av !< liquid water
435    REAL(wp), ALLOCATABLE, DIMENSION(:,:,:), TARGET ::  s_NH_av  !< ammonia
436    REAL(wp), ALLOCATABLE, DIMENSION(:,:,:), TARGET ::  s_NO_av  !< nitrates
437    REAL(wp), ALLOCATABLE, DIMENSION(:,:,:), TARGET ::  s_OC_av  !< org. carbon
438    REAL(wp), ALLOCATABLE, DIMENSION(:,:,:), TARGET ::  s_SO4_av !< sulphates
439    REAL(wp), ALLOCATABLE, DIMENSION(:,:,:), TARGET ::  s_SS_av  !< sea salt
440!-- Bins:   
441    REAL(wp), ALLOCATABLE, DIMENSION(:,:,:,:), TARGET ::  mbins_av !< bin mass 
442    REAL(wp), ALLOCATABLE, DIMENSION(:,:,:,:), TARGET ::  Nbins_av !< bin number
443
444   
445!
446!-- PALM interfaces:
447!
448!-- Boundary conditions:
449    INTERFACE salsa_boundary_conds
450       MODULE PROCEDURE salsa_boundary_conds
451       MODULE PROCEDURE salsa_boundary_conds_decycle
452    END INTERFACE salsa_boundary_conds
453!   
454!-- Data output checks for 2D/3D data to be done in check_parameters
455    INTERFACE salsa_check_data_output
456       MODULE PROCEDURE salsa_check_data_output
457    END INTERFACE salsa_check_data_output
458   
459!
460!-- Input parameter checks to be done in check_parameters
461    INTERFACE salsa_check_parameters
462       MODULE PROCEDURE salsa_check_parameters
463    END INTERFACE salsa_check_parameters
464
465!
466!-- Averaging of 3D data for output
467    INTERFACE salsa_3d_data_averaging
468       MODULE PROCEDURE salsa_3d_data_averaging
469    END INTERFACE salsa_3d_data_averaging
470
471!
472!-- Data output of 2D quantities
473    INTERFACE salsa_data_output_2d
474       MODULE PROCEDURE salsa_data_output_2d
475    END INTERFACE salsa_data_output_2d
476
477!
478!-- Data output of 3D data
479    INTERFACE salsa_data_output_3d
480       MODULE PROCEDURE salsa_data_output_3d
481    END INTERFACE salsa_data_output_3d
482   
483!
484!-- Data output of 3D data
485    INTERFACE salsa_data_output_mask
486       MODULE PROCEDURE salsa_data_output_mask
487    END INTERFACE salsa_data_output_mask
488
489!
490!-- Definition of data output quantities
491    INTERFACE salsa_define_netcdf_grid
492       MODULE PROCEDURE salsa_define_netcdf_grid
493    END INTERFACE salsa_define_netcdf_grid
494   
495!
496!-- Output of information to the header file
497    INTERFACE salsa_header
498       MODULE PROCEDURE salsa_header
499    END INTERFACE salsa_header
500 
501!
502!-- Initialization actions 
503    INTERFACE salsa_init
504       MODULE PROCEDURE salsa_init
505    END INTERFACE salsa_init
506 
507!
508!-- Initialization of arrays
509    INTERFACE salsa_init_arrays
510       MODULE PROCEDURE salsa_init_arrays
511    END INTERFACE salsa_init_arrays
512
513!
514!-- Writing of binary output for restart runs  !!! renaming?!
515    INTERFACE salsa_wrd_local
516       MODULE PROCEDURE salsa_wrd_local
517    END INTERFACE salsa_wrd_local
518   
519!
520!-- Reading of NAMELIST parameters
521    INTERFACE salsa_parin
522       MODULE PROCEDURE salsa_parin
523    END INTERFACE salsa_parin
524
525!
526!-- Reading of parameters for restart runs
527    INTERFACE salsa_rrd_local
528       MODULE PROCEDURE salsa_rrd_local
529    END INTERFACE salsa_rrd_local
530   
531!
532!-- Swapping of time levels (required for prognostic variables)
533    INTERFACE salsa_swap_timelevel
534       MODULE PROCEDURE salsa_swap_timelevel
535    END INTERFACE salsa_swap_timelevel
536
537    INTERFACE salsa_driver
538       MODULE PROCEDURE salsa_driver
539    END INTERFACE salsa_driver
540
541    INTERFACE salsa_tendency
542       MODULE PROCEDURE salsa_tendency
543       MODULE PROCEDURE salsa_tendency_ij
544    END INTERFACE salsa_tendency
545   
546   
547   
548    SAVE
549
550    PRIVATE
551!
552!-- Public functions:
553    PUBLIC salsa_boundary_conds, salsa_check_data_output,                      &
554           salsa_check_parameters, salsa_3d_data_averaging,                    &
555           salsa_data_output_2d, salsa_data_output_3d, salsa_data_output_mask, &
556           salsa_define_netcdf_grid, salsa_diagnostics, salsa_driver,          &
557           salsa_header, salsa_init, salsa_init_arrays, salsa_parin,           &
558           salsa_rrd_local, salsa_swap_timelevel, salsa_tendency,              &
559           salsa_wrd_local
560!
561!-- Public parameters, constants and initial values
562    PUBLIC dots_salsa, dt_salsa, last_salsa_time, lsdepo, salsa,               &
563           salsa_gases_from_chem, skip_time_do_salsa
564!
565!-- Public prognostic variables
566    PUBLIC aerosol_mass, aerosol_number, fn2a, fn2b, gconc_2, in1a, in2b,      &
567           mconc_2, nbins, ncc, ncc_tot, nclim, nconc_2, ngast, prtcl, Ra_dry, &
568           salsa_gas, sedim_vd
569           
570
571 CONTAINS
572
573!------------------------------------------------------------------------------!
574! Description:
575! ------------
576!> Parin for &salsa_par for new modules
577!------------------------------------------------------------------------------!
578 SUBROUTINE salsa_parin
579
580    IMPLICIT NONE
581
582    CHARACTER (LEN=80) ::  line   !< dummy string that contains the current line
583                                  !< of the parameter file
584                                 
585    NAMELIST /salsa_parameters/             &
586                          advect_particle_water, & ! Switch for advecting
587                                                ! particle water. If .FALSE.,
588                                                ! equilibration is called at
589                                                ! each time step.       
590                          bc_salsa_b,       &   ! bottom boundary condition
591                          bc_salsa_t,       &   ! top boundary condition
592                          decycle_lr,       &   ! decycle SALSA components
593                          decycle_method,   &   ! decycle method applied:
594                                                ! 1=left 2=right 3=south 4=north
595                          decycle_ns,       &   ! decycle SALSA components
596                          depo_vege_type,   &   ! Parametrisation type
597                          depo_topo_type,   &   ! Parametrisation type
598                          dpg,              &   ! Mean diameter for the initial
599                                                ! log-normal modes
600                          dt_salsa,         &   ! SALSA timestep in seconds
601                          feedback_to_palm, &   ! allow feedback due to
602                                                ! hydration / condensation
603                          H2SO4_init,       &   ! Init value for sulphuric acid
604                          HNO3_init,        &   ! Init value for nitric acid
605                          igctyp,           &   ! Initial gas concentration type
606                          isdtyp,           &   ! Initial size distribution type                                               
607                          listspec,         &   ! List of actived aerosols
608                                                ! (string list)
609                          mass_fracs_a,     &   ! Initial relative contribution 
610                                                ! of each species to particle 
611                                                ! volume in a-bins, 0 for unused
612                          mass_fracs_b,     &   ! Initial relative contribution 
613                                                ! of each species to particle
614                                                ! volume in b-bins, 0 for unused
615                          n_lognorm,        &   ! Number concentration for the
616                                                ! log-normal modes                                               
617                          nbin,             &   ! Number of size bins for
618                                                ! aerosol size subranges 1 & 2
619                          nf2a,             &   ! Number fraction of particles
620                                                ! allocated to a-bins in
621                                                ! subrange 2 b-bins will get
622                                                ! 1-nf2a                         
623                          NH3_init,         &   ! Init value for ammonia
624                          nj3,              &   ! J3 parametrization
625                                                ! 1 = condensational sink
626                                                !     (Kerminen&Kulmala, 2002)
627                                                ! 2 = coagulational sink
628                                                !     (Lehtinen et al. 2007)
629                                                ! 3 = coagS+self-coagulation
630                                                !     (Anttila et al. 2010)                                                   
631                          nlcnd,            &   ! Condensation master switch
632                          nlcndgas,         &   ! Condensation of gases
633                          nlcndh2oae,       &   ! Condensation of H2O                           
634                          nlcoag,           &   ! Coagulation master switch
635                          nldepo,           &   ! Deposition master switch
636                          nldepo_vege,      &   ! Deposition on vegetation
637                                                ! master switch
638                          nldepo_topo,      &   ! Deposition on topo master
639                                                ! switch                         
640                          nldistupdate,     &   ! Size distribution update
641                                                ! master switch
642                          nsnucl,           &   ! Nucleation scheme:
643                                                ! 0 = off,
644                                                ! 1 = binary nucleation
645                                                ! 2 = activation type nucleation
646                                                ! 3 = kinetic nucleation
647                                                ! 4 = ternary nucleation
648                                                ! 5 = nucleation with organics
649                                                ! 6 = activation type of
650                                                !     nucleation with H2SO4+ORG
651                                                ! 7 = heteromolecular nucleation
652                                                !     with H2SO4*ORG
653                                                ! 8 = homomolecular nucleation 
654                                                !     of H2SO4 + heteromolecular
655                                                !     nucleation with H2SO4*ORG
656                                                ! 9 = homomolecular nucleation
657                                                !     of H2SO4 and ORG + hetero-
658                                                !     molecular nucleation with
659                                                !     H2SO4*ORG
660                          OCNV_init,        &   ! Init value for non-volatile
661                                                ! organic gases
662                          OCSV_init,        &   ! Init value for semi-volatile
663                                                ! organic gases
664                          read_restart_data_salsa, & ! read restart data for
665                                                     ! salsa
666                          reglim,           &   ! Min&max diameter limits of
667                                                ! size subranges
668                          salsa,            &   ! Master switch for SALSA
669                          salsa_source_mode,&   ! 'read_from_file' or 'constant'
670                                                ! or 'no_source'
671                          sigmag,           &   ! stdev for the initial log-
672                                                ! normal modes                                               
673                          skip_time_do_salsa, & ! Starting time of SALSA (s)
674                          van_der_waals_coagc,& ! include van der Waals forces
675                          write_binary_salsa    ! Write binary for salsa
676                           
677       
678    line = ' '
679       
680!
681!-- Try to find salsa package
682    REWIND ( 11 )
683    line = ' '
684    DO WHILE ( INDEX( line, '&salsa_parameters' ) == 0 )
685       READ ( 11, '(A)', END=10 )  line
686    ENDDO
687    BACKSPACE ( 11 )
688
689!
690!-- Read user-defined namelist
691    READ ( 11, salsa_parameters )
692
693!
694!-- Enable salsa (salsa switch in modules.f90)
695    salsa = .TRUE.
696
697 10 CONTINUE
698       
699 END SUBROUTINE salsa_parin
700
701 
702!------------------------------------------------------------------------------!
703! Description:
704! ------------
705!> Check parameters routine for salsa.
706!------------------------------------------------------------------------------!
707 SUBROUTINE salsa_check_parameters
708
709    USE control_parameters,                                                    &
710        ONLY:  message_string
711       
712    IMPLICIT NONE
713   
714!
715!-- Checks go here (cf. check_parameters.f90).
716    IF ( salsa  .AND.  .NOT.  humidity )  THEN
717       WRITE( message_string, * ) 'salsa = ', salsa, ' is ',                   &
718              'not allowed with humidity = ', humidity
719       CALL message( 'check_parameters', 'SA0009', 1, 2, 0, 6, 0 )
720    ENDIF
721   
722    IF ( bc_salsa_b == 'dirichlet' )  THEN
723       ibc_salsa_b = 0
724    ELSEIF ( bc_salsa_b == 'neumann' )  THEN
725       ibc_salsa_b = 1
726    ELSE
727       message_string = 'unknown boundary condition: bc_salsa_b = "'           &
728                         // TRIM( bc_salsa_t ) // '"'
729       CALL message( 'check_parameters', 'SA0011', 1, 2, 0, 6, 0 )                 
730    ENDIF
731   
732    IF ( bc_salsa_t == 'dirichlet' )  THEN
733       ibc_salsa_t = 0
734    ELSEIF ( bc_salsa_t == 'neumann' )  THEN
735       ibc_salsa_t = 1
736    ELSE
737       message_string = 'unknown boundary condition: bc_salsa_t = "'           &
738                         // TRIM( bc_salsa_t ) // '"'
739       CALL message( 'check_parameters', 'SA0012', 1, 2, 0, 6, 0 )                 
740    ENDIF
741   
742    IF ( nj3 < 1  .OR.  nj3 > 3 )  THEN
743       message_string = 'unknown nj3 (must be 1-3)'
744       CALL message( 'check_parameters', 'SA0044', 1, 2, 0, 6, 0 )
745    ENDIF
746           
747 END SUBROUTINE salsa_check_parameters
748
749!------------------------------------------------------------------------------!
750!
751! Description:
752! ------------
753!> Subroutine defining appropriate grid for netcdf variables.
754!> It is called out from subroutine netcdf.
755!> Same grid as for other scalars (see netcdf_interface_mod.f90)
756!------------------------------------------------------------------------------!
757 SUBROUTINE salsa_define_netcdf_grid( var, found, grid_x, grid_y, grid_z )
758   
759    IMPLICIT NONE
760
761    CHARACTER (LEN=*), INTENT(OUT) ::  grid_x   !<
762    CHARACTER (LEN=*), INTENT(OUT) ::  grid_y   !<
763    CHARACTER (LEN=*), INTENT(OUT) ::  grid_z   !<
764    CHARACTER (LEN=*), INTENT(IN)  ::  var      !<
765   
766    LOGICAL, INTENT(OUT) ::  found   !<
767   
768    found  = .TRUE.
769!
770!-- Check for the grid
771
772    IF ( var(1:2) == 'g_' )  THEN
773       grid_x = 'x' 
774       grid_y = 'y' 
775       grid_z = 'zu'   
776    ELSEIF ( var(1:4) == 'LDSA' )  THEN
777       grid_x = 'x' 
778       grid_y = 'y' 
779       grid_z = 'zu'
780    ELSEIF ( var(1:5) == 'm_bin' )  THEN
781       grid_x = 'x' 
782       grid_y = 'y' 
783       grid_z = 'zu'
784    ELSEIF ( var(1:5) == 'N_bin' )  THEN
785       grid_x = 'x' 
786       grid_y = 'y' 
787       grid_z = 'zu'
788    ELSEIF ( var(1:4) == 'Ntot' ) THEN
789       grid_x = 'x' 
790       grid_y = 'y' 
791       grid_z = 'zu'
792    ELSEIF ( var(1:2) == 'PM' )  THEN
793       grid_x = 'x' 
794       grid_y = 'y' 
795       grid_z = 'zu'
796    ELSEIF ( var(1:2) == 's_' )  THEN
797       grid_x = 'x' 
798       grid_y = 'y' 
799       grid_z = 'zu'
800    ELSE
801       found  = .FALSE.
802       grid_x = 'none'
803       grid_y = 'none'
804       grid_z = 'none'
805    ENDIF
806
807 END SUBROUTINE salsa_define_netcdf_grid
808
809 
810!------------------------------------------------------------------------------!
811! Description:
812! ------------
813!> Header output for new module
814!------------------------------------------------------------------------------!
815 SUBROUTINE salsa_header( io )
816
817    IMPLICIT NONE
818 
819    INTEGER(iwp), INTENT(IN) ::  io   !< Unit of the output file
820!
821!-- Write SALSA header
822    WRITE( io, 1 )
823    WRITE( io, 2 ) skip_time_do_salsa
824    WRITE( io, 3 ) dt_salsa
825    WRITE( io, 12 )  SHAPE( aerosol_number(1)%conc ), nbins
826    IF ( advect_particle_water )  THEN
827       WRITE( io, 16 )  SHAPE( aerosol_mass(1)%conc ), ncc_tot*nbins,          &
828                        advect_particle_water
829    ELSE
830       WRITE( io, 16 )  SHAPE( aerosol_mass(1)%conc ), ncc*nbins,              &
831                        advect_particle_water
832    ENDIF
833    IF ( .NOT. salsa_gases_from_chem )  THEN
834       WRITE( io, 17 )  SHAPE( aerosol_mass(1)%conc ), ngast,                  &
835                        salsa_gases_from_chem
836    ENDIF
837    WRITE( io, 4 ) 
838    IF ( nsnucl > 0 )  THEN
839       WRITE( io, 5 ) nsnucl, nj3
840    ENDIF
841    IF ( nlcoag )  THEN
842       WRITE( io, 6 ) 
843    ENDIF
844    IF ( nlcnd )  THEN
845       WRITE( io, 7 ) nlcndgas, nlcndh2oae
846    ENDIF
847    IF ( nldepo )  THEN
848       WRITE( io, 14 ) nldepo_vege, nldepo_topo
849    ENDIF
850    WRITE( io, 8 )  reglim, nbin, bin_low_limits
851    WRITE( io, 15 ) nsect
852    WRITE( io, 13 ) ncc, listspec, mass_fracs_a, mass_fracs_b
853    IF ( .NOT. salsa_gases_from_chem )  THEN
854       WRITE( io, 18 ) ngast, H2SO4_init, HNO3_init, NH3_init, OCNV_init,      &
855                       OCSV_init
856    ENDIF
857    WRITE( io, 9 )  isdtyp, igctyp
858    IF ( isdtyp == 0 )  THEN
859       WRITE( io, 10 )  dpg, sigmag, n_lognorm
860    ELSE
861       WRITE( io, 11 )
862    ENDIF
863   
864
8651   FORMAT (//' SALSA information:'/                                           &
866              ' ------------------------------'/)
8672   FORMAT   ('    Starts at: skip_time_do_salsa = ', F10.2, '  s')
8683   FORMAT  (/'    Timestep: dt_salsa = ', F6.2, '  s')
86912  FORMAT  (/'    Array shape (z,y,x,bins):'/                                 &
870              '       aerosol_number:  ', 4(I3)) 
87116  FORMAT  (/'       aerosol_mass:    ', 4(I3),/                              &
872              '       (advect_particle_water = ', L1, ')')
87317  FORMAT   ('       salsa_gas: ', 4(I3),/                                    &
874              '       (salsa_gases_from_chem = ', L1, ')')
8754   FORMAT  (/'    Aerosol dynamic processes included: ')
8765   FORMAT  (/'       nucleation (scheme = ', I1, ' and J3 parametrization = ',&
877               I1, ')')
8786   FORMAT  (/'       coagulation')
8797   FORMAT  (/'       condensation (of precursor gases = ', L1,                &
880              '          and water vapour = ', L1, ')' )
88114  FORMAT  (/'       dry deposition (on vegetation = ', L1,                   &
882              '          and on topography = ', L1, ')')             
8838   FORMAT  (/'    Aerosol bin subrange limits (in metres): ',  3(ES10.2E3), / &
884              '    Number of size bins for each aerosol subrange: ', 2I3,/     &
885              '    Aerosol bin limits (in metres): ', *(ES10.2E3))
88615  FORMAT   ('    Initial number concentration in bins at the lowest level',  &
887              ' (#/m**3):', *(ES10.2E3))       
88813  FORMAT  (/'    Number of chemical components used: ', I1,/                 &
889              '       Species: ',7(A6),/                                       &
890              '    Initial relative contribution of each species to particle', & 
891              ' volume in:',/                                                  &
892              '       a-bins: ', 7(F6.3),/                                     &
893              '       b-bins: ', 7(F6.3))
89418  FORMAT  (/'    Number of gaseous tracers used: ', I1,/                     &
895              '    Initial gas concentrations:',/                              &
896              '       H2SO4: ',ES12.4E3, ' #/m**3',/                           &
897              '       HNO3:  ',ES12.4E3, ' #/m**3',/                           &
898              '       NH3:   ',ES12.4E3, ' #/m**3',/                           &
899              '       OCNV:  ',ES12.4E3, ' #/m**3',/                           &
900              '       OCSV:  ',ES12.4E3, ' #/m**3')
9019    FORMAT (/'   Initialising concentrations: ', /                            &
902              '      Aerosol size distribution: isdtyp = ', I1,/               &
903              '      Gas concentrations: igctyp = ', I1 )
90410   FORMAT ( '      Mode diametres: dpg(nmod) = ', 7(F7.3),/                  &
905              '      Standard deviation: sigmag(nmod) = ', 7(F7.2),/           &
906              '      Number concentration: n_lognorm(nmod) = ', 7(ES12.4E3) )
90711   FORMAT (/'      Size distribution read from a file.')
908
909 END SUBROUTINE salsa_header
910
911!------------------------------------------------------------------------------!
912! Description:
913! ------------
914!> Allocate SALSA arrays and define pointers if required
915!------------------------------------------------------------------------------!
916 SUBROUTINE salsa_init_arrays
917 
918    USE surface_mod,                                                           &
919        ONLY:  surf_def_h, surf_def_v, surf_lsm_h, surf_lsm_v, surf_usm_h,     &
920               surf_usm_v
921
922    IMPLICIT NONE
923   
924    INTEGER(iwp) ::  gases_available !< Number of available gas components in
925                                     !< the chemistry model
926    INTEGER(iwp) ::  i   !< loop index for allocating
927    INTEGER(iwp) ::  l   !< loop index for allocating: surfaces
928    INTEGER(iwp) ::  lsp !< loop index for chem species in the chemistry model
929   
930    gases_available = 0
931
932!
933!-- Allocate prognostic variables (see salsa_swap_timelevel)
934#if defined( __nopointer )
935    message_string = 'SALSA runs only with POINTER Version'
936    CALL message( 'salsa_mod: salsa_init_arrays', 'SA0023', 1, 2, 0, 6, 0 )
937#else         
938!
939!-- Set derived indices:
940!-- (This does the same as the subroutine salsa_initialize in SALSA/
941!-- UCLALES-SALSA)       
942    in1a = 1                ! 1st index of subrange 1a
943    in2a = in1a + nbin(1)   ! 1st index of subrange 2a
944    fn1a = in2a - 1         ! last index of subrange 1a
945    fn2a = fn1a + nbin(2)   ! last index of subrange 2a
946   
947!   
948!-- If the fraction of insoluble aerosols in subrange 2 is zero: do not allocate
949!-- arrays for them
950    IF ( nf2a > 0.999999_wp  .AND.  SUM( mass_fracs_b ) < 0.00001_wp )  THEN
951       no_insoluble = .TRUE.
952       in2b = fn2a+1    ! 1st index of subrange 2b
953       fn2b = fn2a      ! last index of subrange 2b
954    ELSE
955       in2b = in2a + nbin(2)   ! 1st index of subrange 2b
956       fn2b = fn2a + nbin(2)   ! last index of subrange 2b
957    ENDIF
958   
959   
960    nbins = fn2b   ! total number of aerosol size bins
961!   
962!-- Create index tables for different aerosol components
963    CALL component_index_constructor( prtcl, ncc, maxspec, listspec )
964   
965    ncc_tot = ncc
966    IF ( advect_particle_water )  ncc_tot = ncc + 1  ! Add water
967   
968!
969!-- Allocate:
970    ALLOCATE( aero(nbins), bin_low_limits(nbins), nsect(nbins), massacc(nbins) )
971    IF ( nldepo ) ALLOCATE( sedim_vd(nzb:nzt+1,nysg:nyng,nxlg:nxrg,nbins) )         
972    ALLOCATE( Ra_dry(nzb:nzt+1,nysg:nyng,nxlg:nxrg,nbins) )
973   
974!   
975!-- Aerosol number concentration
976    ALLOCATE( aerosol_number(nbins) )
977    ALLOCATE( nconc_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg,nbins),                    &
978              nconc_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg,nbins),                    &
979              nconc_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg,nbins) )
980    nconc_1 = 0.0_wp
981    nconc_2 = 0.0_wp
982    nconc_3 = 0.0_wp
983   
984    DO i = 1, nbins
985       aerosol_number(i)%conc(nzb:nzt+1,nysg:nyng,nxlg:nxrg)    => nconc_1(:,:,:,i)
986       aerosol_number(i)%conc_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg)  => nconc_2(:,:,:,i)
987       aerosol_number(i)%tconc_m(nzb:nzt+1,nysg:nyng,nxlg:nxrg) => nconc_3(:,:,:,i)
988       ALLOCATE( aerosol_number(i)%flux_s(nzb+1:nzt,0:threads_per_task-1),     &
989                 aerosol_number(i)%diss_s(nzb+1:nzt,0:threads_per_task-1),     &
990                 aerosol_number(i)%flux_l(nzb+1:nzt,nys:nyn,0:threads_per_task-1),&
991                 aerosol_number(i)%diss_l(nzb+1:nzt,nys:nyn,0:threads_per_task-1),&
992                 aerosol_number(i)%init(nzb:nzt+1),                            &
993                 aerosol_number(i)%sums_ws_l(nzb:nzt+1,0:threads_per_task-1) )
994    ENDDO     
995   
996!   
997!-- Aerosol mass concentration   
998    ALLOCATE( aerosol_mass(ncc_tot*nbins) ) 
999    ALLOCATE( mconc_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg,ncc_tot*nbins),            &
1000              mconc_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg,ncc_tot*nbins),            &
1001              mconc_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg,ncc_tot*nbins) )
1002    mconc_1 = 0.0_wp
1003    mconc_2 = 0.0_wp
1004    mconc_3 = 0.0_wp
1005   
1006    DO i = 1, ncc_tot*nbins
1007       aerosol_mass(i)%conc(nzb:nzt+1,nysg:nyng,nxlg:nxrg)    => mconc_1(:,:,:,i)
1008       aerosol_mass(i)%conc_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg)  => mconc_2(:,:,:,i)
1009       aerosol_mass(i)%tconc_m(nzb:nzt+1,nysg:nyng,nxlg:nxrg) => mconc_3(:,:,:,i)       
1010       ALLOCATE( aerosol_mass(i)%flux_s(nzb+1:nzt,0:threads_per_task-1),       &
1011                 aerosol_mass(i)%diss_s(nzb+1:nzt,0:threads_per_task-1),       &
1012                 aerosol_mass(i)%flux_l(nzb+1:nzt,nys:nyn,0:threads_per_task-1),&
1013                 aerosol_mass(i)%diss_l(nzb+1:nzt,nys:nyn,0:threads_per_task-1),&
1014                 aerosol_mass(i)%init(nzb:nzt+1),                              &
1015                 aerosol_mass(i)%sums_ws_l(nzb:nzt+1,0:threads_per_task-1)  )
1016    ENDDO
1017   
1018!
1019!-- Surface fluxes: answs = aerosol number, amsws = aerosol mass
1020!
1021!-- Horizontal surfaces: default type
1022    DO  l = 0, 2   ! upward (l=0), downward (l=1) and model top (l=2)
1023       ALLOCATE( surf_def_h(l)%answs( 1:surf_def_h(l)%ns, nbins ) )
1024       ALLOCATE( surf_def_h(l)%amsws( 1:surf_def_h(l)%ns, nbins*ncc_tot ) )
1025       surf_def_h(l)%answs = 0.0_wp
1026       surf_def_h(l)%amsws = 0.0_wp
1027    ENDDO
1028!-- Horizontal surfaces: natural type   
1029    IF ( land_surface )  THEN
1030       ALLOCATE( surf_lsm_h%answs( 1:surf_lsm_h%ns, nbins ) )
1031       ALLOCATE( surf_lsm_h%amsws( 1:surf_lsm_h%ns, nbins*ncc_tot ) )
1032       surf_lsm_h%answs = 0.0_wp
1033       surf_lsm_h%amsws = 0.0_wp
1034    ENDIF
1035!-- Horizontal surfaces: urban type
1036    IF ( urban_surface )  THEN
1037       ALLOCATE( surf_usm_h%answs( 1:surf_usm_h%ns, nbins ) )
1038       ALLOCATE( surf_usm_h%amsws( 1:surf_usm_h%ns, nbins*ncc_tot ) )
1039       surf_usm_h%answs = 0.0_wp
1040       surf_usm_h%amsws = 0.0_wp
1041    ENDIF
1042!
1043!-- Vertical surfaces: northward (l=0), southward (l=1), eastward (l=2) and
1044!-- westward (l=3) facing
1045    DO  l = 0, 3   
1046       ALLOCATE( surf_def_v(l)%answs( 1:surf_def_v(l)%ns, nbins ) )
1047       surf_def_v(l)%answs = 0.0_wp
1048       ALLOCATE( surf_def_v(l)%amsws( 1:surf_def_v(l)%ns, nbins*ncc_tot ) )
1049       surf_def_v(l)%amsws = 0.0_wp
1050       
1051       IF ( land_surface)  THEN
1052          ALLOCATE( surf_lsm_v(l)%answs( 1:surf_lsm_v(l)%ns, nbins ) )
1053          surf_lsm_v(l)%answs = 0.0_wp
1054          ALLOCATE( surf_lsm_v(l)%amsws( 1:surf_lsm_v(l)%ns, nbins*ncc_tot ) )
1055          surf_lsm_v(l)%amsws = 0.0_wp
1056       ENDIF
1057       
1058       IF ( urban_surface )  THEN
1059          ALLOCATE( surf_usm_v(l)%answs( 1:surf_usm_v(l)%ns, nbins ) )
1060          surf_usm_v(l)%answs = 0.0_wp
1061          ALLOCATE( surf_usm_v(l)%amsws( 1:surf_usm_v(l)%ns, nbins*ncc_tot ) )
1062          surf_usm_v(l)%amsws = 0.0_wp
1063       ENDIF
1064    ENDDO   
1065   
1066!
1067!-- Concentration of gaseous tracers (1. SO4, 2. HNO3, 3. NH3, 4. OCNV, 5. OCSV)
1068!-- (number concentration (#/m3) )
1069!
1070!-- If chemistry is on, read gas phase concentrations from there. Otherwise,
1071!-- allocate salsa_gas array.
1072
1073    IF ( air_chemistry )  THEN   
1074       DO  lsp = 1, nvar
1075          IF ( TRIM( chem_species(lsp)%name ) == 'H2SO4' )  THEN
1076             gases_available = gases_available + 1
1077             gas_index_chem(1) = lsp
1078          ELSEIF ( TRIM( chem_species(lsp)%name ) == 'HNO3' )  THEN
1079             gases_available = gases_available + 1 
1080             gas_index_chem(2) = lsp
1081          ELSEIF ( TRIM( chem_species(lsp)%name ) == 'NH3' )  THEN
1082             gases_available = gases_available + 1
1083             gas_index_chem(3) = lsp
1084          ELSEIF ( TRIM( chem_species(lsp)%name ) == 'OCNV' )  THEN
1085             gases_available = gases_available + 1
1086             gas_index_chem(4) = lsp
1087          ELSEIF ( TRIM( chem_species(lsp)%name ) == 'OCSV' )  THEN
1088             gases_available = gases_available + 1
1089             gas_index_chem(5) = lsp
1090          ENDIF
1091       ENDDO
1092
1093       IF ( gases_available == ngast )  THEN
1094          salsa_gases_from_chem = .TRUE.
1095       ELSE
1096          WRITE( message_string, * ) 'SALSA is run together with chemistry '// &
1097                                     'but not all gaseous components are '//   &
1098                                     'provided by kpp (H2SO4, HNO3, NH3, '//   &
1099                                     'OCNV, OCSC)'
1100       CALL message( 'check_parameters', 'SA0024', 1, 2, 0, 6, 0 )
1101       ENDIF
1102
1103    ELSE
1104
1105       ALLOCATE( salsa_gas(ngast) ) 
1106       ALLOCATE( gconc_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg,ngast),                 &
1107                 gconc_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg,ngast),                 &
1108                 gconc_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg,ngast) )
1109       gconc_1 = 0.0_wp
1110       gconc_2 = 0.0_wp
1111       gconc_3 = 0.0_wp
1112       
1113       DO i = 1, ngast
1114          salsa_gas(i)%conc(nzb:nzt+1,nysg:nyng,nxlg:nxrg)    => gconc_1(:,:,:,i)
1115          salsa_gas(i)%conc_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg)  => gconc_2(:,:,:,i)
1116          salsa_gas(i)%tconc_m(nzb:nzt+1,nysg:nyng,nxlg:nxrg) => gconc_3(:,:,:,i)
1117          ALLOCATE( salsa_gas(i)%flux_s(nzb+1:nzt,0:threads_per_task-1),       &
1118                    salsa_gas(i)%diss_s(nzb+1:nzt,0:threads_per_task-1),       &
1119                    salsa_gas(i)%flux_l(nzb+1:nzt,nys:nyn,0:threads_per_task-1),&
1120                    salsa_gas(i)%diss_l(nzb+1:nzt,nys:nyn,0:threads_per_task-1),&
1121                    salsa_gas(i)%init(nzb:nzt+1),                              &
1122                    salsa_gas(i)%sums_ws_l(nzb:nzt+1,0:threads_per_task-1) )
1123       ENDDO       
1124!
1125!--    Surface fluxes: gtsws = gaseous tracer flux
1126!
1127!--    Horizontal surfaces: default type
1128       DO  l = 0, 2   ! upward (l=0), downward (l=1) and model top (l=2)
1129          ALLOCATE( surf_def_h(l)%gtsws( 1:surf_def_h(l)%ns, ngast ) )
1130          surf_def_h(l)%gtsws = 0.0_wp
1131       ENDDO
1132!--    Horizontal surfaces: natural type   
1133       IF ( land_surface )  THEN
1134          ALLOCATE( surf_lsm_h%gtsws( 1:surf_lsm_h%ns, ngast ) )
1135          surf_lsm_h%gtsws = 0.0_wp
1136       ENDIF
1137!--    Horizontal surfaces: urban type         
1138       IF ( urban_surface )  THEN
1139          ALLOCATE( surf_usm_h%gtsws( 1:surf_usm_h%ns, ngast ) )
1140          surf_usm_h%gtsws = 0.0_wp
1141       ENDIF
1142!
1143!--    Vertical surfaces: northward (l=0), southward (l=1), eastward (l=2) and
1144!--    westward (l=3) facing
1145       DO  l = 0, 3     
1146          ALLOCATE( surf_def_v(l)%gtsws( 1:surf_def_v(l)%ns, ngast ) )
1147          surf_def_v(l)%gtsws = 0.0_wp
1148          IF ( land_surface )  THEN
1149             ALLOCATE( surf_lsm_v(l)%gtsws( 1:surf_lsm_v(l)%ns, ngast ) )
1150             surf_lsm_v(l)%gtsws = 0.0_wp
1151          ENDIF
1152          IF ( urban_surface )  THEN
1153             ALLOCATE( surf_usm_v(l)%gtsws( 1:surf_usm_v(l)%ns, ngast ) )
1154             surf_usm_v(l)%gtsws = 0.0_wp
1155          ENDIF
1156       ENDDO
1157    ENDIF
1158   
1159#endif
1160
1161 END SUBROUTINE salsa_init_arrays
1162
1163!------------------------------------------------------------------------------!
1164! Description:
1165! ------------
1166!> Initialization of SALSA. Based on salsa_initialize in UCLALES-SALSA.
1167!> Subroutines salsa_initialize, SALSAinit and DiagInitAero in UCLALES-SALSA are
1168!> also merged here.
1169!------------------------------------------------------------------------------!
1170 SUBROUTINE salsa_init
1171
1172    IMPLICIT NONE
1173   
1174    INTEGER(iwp) :: b
1175    INTEGER(iwp) :: c
1176    INTEGER(iwp) :: g
1177    INTEGER(iwp) :: i
1178    INTEGER(iwp) :: j
1179   
1180    bin_low_limits = 0.0_wp
1181    nsect          = 0.0_wp
1182    massacc        = 1.0_wp 
1183   
1184!
1185!-- Indices for chemical components used (-1 = not used)
1186    i = 0
1187    IF ( is_used( prtcl, 'SO4' ) )  THEN
1188       iso4 = get_index( prtcl,'SO4' )
1189       i = i + 1
1190    ENDIF
1191    IF ( is_used( prtcl,'OC' ) )  THEN
1192       ioc = get_index(prtcl, 'OC')
1193       i = i + 1
1194    ENDIF
1195    IF ( is_used( prtcl, 'BC' ) )  THEN
1196       ibc = get_index( prtcl, 'BC' )
1197       i = i + 1
1198    ENDIF
1199    IF ( is_used( prtcl, 'DU' ) )  THEN
1200       idu = get_index( prtcl, 'DU' )
1201       i = i + 1
1202    ENDIF
1203    IF ( is_used( prtcl, 'SS' ) )  THEN
1204       iss = get_index( prtcl, 'SS' )
1205       i = i + 1
1206    ENDIF
1207    IF ( is_used( prtcl, 'NO' ) )  THEN
1208       ino = get_index( prtcl, 'NO' )
1209       i = i + 1
1210    ENDIF
1211    IF ( is_used( prtcl, 'NH' ) )  THEN
1212       inh = get_index( prtcl, 'NH' )
1213       i = i + 1
1214    ENDIF
1215!   
1216!-- All species must be known
1217    IF ( i /= ncc )  THEN
1218       message_string = 'Unknown aerosol species/component(s) given in the' // &
1219                        ' initialization'
1220       CALL message( 'salsa_mod: salsa_init', 'SA0020', 1, 2, 0, 6, 0 )
1221    ENDIF
1222   
1223!
1224!-- Initialise
1225!
1226!-- Aerosol size distribution (TYPE t_section)
1227    aero(:)%dwet     = 1.0E-10_wp
1228    aero(:)%veqh2o   = 1.0E-10_wp
1229    aero(:)%numc     = nclim
1230    aero(:)%core     = 1.0E-10_wp
1231    DO c = 1, maxspec+1    ! 1:SO4, 2:OC, 3:BC, 4:DU, 5:SS, 6:NO, 7:NH, 8:H2O
1232       aero(:)%volc(c) = 0.0_wp
1233    ENDDO
1234   
1235    IF ( nldepo )  sedim_vd = 0.0_wp 
1236   
1237    DO  b = 1, nbins
1238       IF ( .NOT. read_restart_data_salsa )  aerosol_number(b)%conc = nclim
1239       aerosol_number(b)%conc_p    = 0.0_wp
1240       aerosol_number(b)%tconc_m   = 0.0_wp
1241       aerosol_number(b)%flux_s    = 0.0_wp
1242       aerosol_number(b)%diss_s    = 0.0_wp
1243       aerosol_number(b)%flux_l    = 0.0_wp
1244       aerosol_number(b)%diss_l    = 0.0_wp
1245       aerosol_number(b)%init      = nclim
1246       aerosol_number(b)%sums_ws_l = 0.0_wp
1247    ENDDO
1248    DO  c = 1, ncc_tot*nbins
1249       IF ( .NOT. read_restart_data_salsa )  aerosol_mass(c)%conc = mclim
1250       aerosol_mass(c)%conc_p    = 0.0_wp
1251       aerosol_mass(c)%tconc_m   = 0.0_wp
1252       aerosol_mass(c)%flux_s    = 0.0_wp
1253       aerosol_mass(c)%diss_s    = 0.0_wp
1254       aerosol_mass(c)%flux_l    = 0.0_wp
1255       aerosol_mass(c)%diss_l    = 0.0_wp
1256       aerosol_mass(c)%init      = mclim
1257       aerosol_mass(c)%sums_ws_l = 0.0_wp
1258    ENDDO
1259   
1260    IF ( .NOT. salsa_gases_from_chem )  THEN
1261       DO  g = 1, ngast
1262          salsa_gas(g)%conc_p    = 0.0_wp
1263          salsa_gas(g)%tconc_m   = 0.0_wp
1264          salsa_gas(g)%flux_s    = 0.0_wp
1265          salsa_gas(g)%diss_s    = 0.0_wp
1266          salsa_gas(g)%flux_l    = 0.0_wp
1267          salsa_gas(g)%diss_l    = 0.0_wp
1268          salsa_gas(g)%sums_ws_l = 0.0_wp
1269       ENDDO
1270       IF ( .NOT. read_restart_data_salsa )  THEN
1271          salsa_gas(1)%conc = H2SO4_init
1272          salsa_gas(2)%conc = HNO3_init
1273          salsa_gas(3)%conc = NH3_init
1274          salsa_gas(4)%conc = OCNV_init
1275          salsa_gas(5)%conc = OCSV_init
1276       ENDIF
1277!
1278!--    Set initial value for gas compound tracers and initial values
1279       salsa_gas(1)%init = H2SO4_init
1280       salsa_gas(2)%init = HNO3_init
1281       salsa_gas(3)%init = NH3_init
1282       salsa_gas(4)%init = OCNV_init
1283       salsa_gas(5)%init = OCSV_init     
1284    ENDIF
1285!
1286!-- Aerosol radius in each bin: dry and wet (m)
1287    Ra_dry = 1.0E-10_wp
1288!   
1289!-- Initialise aerosol tracers   
1290    aero(:)%vhilim   = 0.0_wp
1291    aero(:)%vlolim   = 0.0_wp
1292    aero(:)%vratiohi = 0.0_wp
1293    aero(:)%vratiolo = 0.0_wp
1294    aero(:)%dmid     = 0.0_wp
1295!
1296!-- Initialise the sectional particle size distribution
1297    CALL set_sizebins()
1298!
1299!-- Initialise location-dependent aerosol size distributions and
1300!-- chemical compositions:
1301    CALL aerosol_init
1302!
1303!-- Initalisation run of SALSA
1304    DO  i = nxl, nxr
1305       DO  j = nys, nyn
1306          CALL salsa_driver( i, j, 1 )
1307          CALL salsa_diagnostics( i, j )
1308       ENDDO
1309    ENDDO
1310!
1311!-- Set the aerosol and gas sources
1312    IF ( salsa_source_mode == 'read_from_file' )  THEN
1313       CALL salsa_set_source
1314    ENDIF
1315   
1316 END SUBROUTINE salsa_init
1317
1318!------------------------------------------------------------------------------!
1319! Description:
1320! ------------
1321!> Initializes particle size distribution grid by calculating size bin limits
1322!> and mid-size for *dry* particles in each bin. Called from salsa_initialize
1323!> (only at the beginning of simulation).
1324!> Size distribution described using:
1325!>   1) moving center method (subranges 1 and 2)
1326!>      (Jacobson, Atmos. Env., 31, 131-144, 1997)
1327!>   2) fixed sectional method (subrange 3)
1328!> Size bins in each subrange are spaced logarithmically
1329!> based on given subrange size limits and bin number.
1330!
1331!> Mona changed 06/2017: Use geometric mean diameter to describe the mean
1332!> particle diameter in a size bin, not the arithmeric mean which clearly
1333!> overestimates the total particle volume concentration.
1334!
1335!> Coded by:
1336!> Hannele Korhonen (FMI) 2005
1337!> Harri Kokkola (FMI) 2006
1338!
1339!> Bug fixes for box model + updated for the new aerosol datatype:
1340!> Juha Tonttila (FMI) 2014
1341!------------------------------------------------------------------------------!
1342 SUBROUTINE set_sizebins
1343               
1344    IMPLICIT NONE
1345!   
1346!-- Local variables
1347    INTEGER(iwp) ::  cc
1348    INTEGER(iwp) ::  dd
1349    REAL(wp) ::  ratio_d !< ratio of the upper and lower diameter of subranges
1350!
1351!-- vlolim&vhilim: min & max *dry* volumes [fxm]
1352!-- dmid: bin mid *dry* diameter (m)
1353!-- vratiolo&vratiohi: volume ratio between the center and low/high limit
1354!
1355!-- 1) Size subrange 1:
1356    ratio_d = reglim(2) / reglim(1)   ! section spacing (m)
1357    DO  cc = in1a,fn1a
1358       aero(cc)%vlolim = api6 * ( reglim(1) * ratio_d **                       &
1359                                ( REAL( cc-1 ) / nbin(1) ) ) ** 3.0_wp
1360       aero(cc)%vhilim = api6 * ( reglim(1) * ratio_d **                       &
1361                                ( REAL( cc ) / nbin(1) ) ) ** 3.0_wp
1362       aero(cc)%dmid = SQRT( ( aero(cc)%vhilim / api6 ) ** ( 1.0_wp / 3.0_wp ) &
1363                           * ( aero(cc)%vlolim / api6 ) ** ( 1.0_wp / 3.0_wp ) )
1364       aero(cc)%vratiohi = aero(cc)%vhilim / ( api6 * aero(cc)%dmid ** 3.0_wp )
1365       aero(cc)%vratiolo = aero(cc)%vlolim / ( api6 * aero(cc)%dmid ** 3.0_wp )
1366    ENDDO
1367!
1368!-- 2) Size subrange 2:
1369!-- 2.1) Sub-subrange 2a: high hygroscopicity
1370    ratio_d = reglim(3) / reglim(2)   ! section spacing
1371    DO  dd = in2a, fn2a
1372       cc = dd - in2a
1373       aero(dd)%vlolim = api6 * ( reglim(2) * ratio_d **                       &
1374                                  ( REAL( cc ) / nbin(2) ) ) ** 3.0_wp
1375       aero(dd)%vhilim = api6 * ( reglim(2) * ratio_d **                       &
1376                                  ( REAL( cc+1 ) / nbin(2) ) ) ** 3.0_wp
1377       aero(dd)%dmid = SQRT( ( aero(dd)%vhilim / api6 ) ** ( 1.0_wp / 3.0_wp ) &
1378                           * ( aero(dd)%vlolim / api6 ) ** ( 1.0_wp / 3.0_wp ) )
1379       aero(dd)%vratiohi = aero(dd)%vhilim / ( api6 * aero(dd)%dmid ** 3.0_wp )
1380       aero(dd)%vratiolo = aero(dd)%vlolim / ( api6 * aero(dd)%dmid ** 3.0_wp )
1381    ENDDO
1382!         
1383!-- 2.2) Sub-subrange 2b: low hygroscopicity
1384    IF ( .NOT. no_insoluble )  THEN
1385       aero(in2b:fn2b)%vlolim   = aero(in2a:fn2a)%vlolim
1386       aero(in2b:fn2b)%vhilim   = aero(in2a:fn2a)%vhilim
1387       aero(in2b:fn2b)%dmid     = aero(in2a:fn2a)%dmid
1388       aero(in2b:fn2b)%vratiohi = aero(in2a:fn2a)%vratiohi
1389       aero(in2b:fn2b)%vratiolo = aero(in2a:fn2a)%vratiolo
1390    ENDIF
1391!         
1392!-- Initialize the wet diameter with the bin dry diameter to avoid numerical
1393!-- problems later
1394    aero(:)%dwet = aero(:)%dmid
1395!
1396!-- Save bin limits (lower diameter) to be delivered to the host model if needed
1397    DO cc = 1, nbins
1398       bin_low_limits(cc) = ( aero(cc)%vlolim / api6 )**( 1.0_wp / 3.0_wp )
1399    ENDDO   
1400   
1401 END SUBROUTINE set_sizebins
1402 
1403!------------------------------------------------------------------------------!
1404! Description:
1405! ------------
1406!> Initilize altitude-dependent aerosol size distributions and compositions.
1407!>
1408!> Mona added 06/2017: Correct the number and mass concentrations by normalizing
1409!< by the given total number and mass concentration.
1410!>
1411!> Tomi Raatikainen, FMI, 29.2.2016
1412!------------------------------------------------------------------------------!
1413 SUBROUTINE aerosol_init
1414 
1415    USE arrays_3d,                                                             &
1416        ONLY:  zu
1417 
1418!    USE NETCDF
1419   
1420    USE netcdf_data_input_mod,                                                 &
1421        ONLY:  get_attribute, get_variable,                                    &
1422               netcdf_data_input_get_dimension_length, open_read_file
1423   
1424    IMPLICIT NONE
1425   
1426    INTEGER(iwp) ::  b          !< loop index: size bins
1427    INTEGER(iwp) ::  c          !< loop index: chemical components
1428    INTEGER(iwp) ::  ee         !< index: end
1429    INTEGER(iwp) ::  g          !< loop index: gases
1430    INTEGER(iwp) ::  i          !< loop index: x-direction
1431    INTEGER(iwp) ::  id_faero   !< NetCDF id of PIDS_SALSA
1432    INTEGER(iwp) ::  id_fchem   !< NetCDF id of PIDS_CHEM
1433    INTEGER(iwp) ::  j          !< loop index: y-direction
1434    INTEGER(iwp) ::  k          !< loop index: z-direction
1435    INTEGER(iwp) ::  kk         !< loop index: z-direction
1436    INTEGER(iwp) ::  nz_file    !< Number of grid-points in file (heights)                           
1437    INTEGER(iwp) ::  prunmode
1438    INTEGER(iwp) ::  ss !< index: start
1439    LOGICAL  ::  netcdf_extend = .FALSE. !< Flag indicating wether netcdf
1440                                         !< topography input file or not
1441    REAL(wp), DIMENSION(nbins) ::  core  !< size of the bin mid aerosol particle,
1442    REAL(wp) ::  flag           !< flag to mask topography grid points
1443    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  pr_gas !< gas profiles
1444    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  pr_mass_fracs_a !< mass fraction
1445                                                              !< profiles: a
1446    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  pr_mass_fracs_b !< and b
1447    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  pr_nsect !< sectional size
1448                                                       !< distribution profile
1449    REAL(wp), DIMENSION(nbins)            ::  nsect  !< size distribution (#/m3)
1450    REAL(wp), DIMENSION(0:nz+1,nbins)     ::  pndist !< size dist as a function
1451                                                     !< of height (#/m3)
1452    REAL(wp), DIMENSION(0:nz+1)           ::  pnf2a  !< number fraction: bins 2a
1453    REAL(wp), DIMENSION(0:nz+1,maxspec)   ::  pvf2a  !< mass distributions of 
1454                                                     !< aerosol species for a 
1455    REAL(wp), DIMENSION(0:nz+1,maxspec)   ::  pvf2b  !< and b-bins     
1456    REAL(wp), DIMENSION(0:nz+1)           ::  pvfOC1a !< mass fraction between
1457                                                     !< SO4 and OC in 1a
1458    REAL(wp), DIMENSION(:), ALLOCATABLE   ::  pr_z
1459
1460    prunmode = 1
1461!
1462!-- Bin mean aerosol particle volume (m3)
1463    core(:) = 0.0_wp
1464    core(1:nbins) = api6 * aero(1:nbins)%dmid ** 3.0_wp
1465!   
1466!-- Set concentrations to zero
1467    nsect(:)     = 0.0_wp
1468    pndist(:,:)  = 0.0_wp
1469    pnf2a(:)     = nf2a   
1470    pvf2a(:,:)   = 0.0_wp
1471    pvf2b(:,:)   = 0.0_wp
1472    pvfOC1a(:)   = 0.0_wp
1473
1474    IF ( isdtyp == 1 )  THEN
1475!
1476!--    Read input profiles from PIDS_SALSA   
1477#if defined( __netcdf )
1478!   
1479!--    Location-dependent size distributions and compositions.     
1480       INQUIRE( FILE='PIDS_SALSA'// TRIM( coupling_char ), EXIST=netcdf_extend )
1481       IF ( netcdf_extend )  THEN
1482!
1483!--       Open file in read-only mode 
1484          CALL open_read_file( 'PIDS_SALSA' // TRIM( coupling_char ), id_faero )
1485!
1486!--       Input heights   
1487          CALL netcdf_data_input_get_dimension_length( id_faero, nz_file,      &
1488                                                       "profile_z" ) 
1489         
1490          ALLOCATE( pr_z(nz_file), pr_mass_fracs_a(maxspec,nz_file),           &
1491                    pr_mass_fracs_b(maxspec,nz_file), pr_nsect(nbins,nz_file) ) 
1492          CALL get_variable( id_faero, 'profile_z', pr_z ) 
1493!       
1494!--       Mass fracs profile: 1: H2SO4 (sulphuric acid), 2: OC (organic carbon),
1495!--                           3: BC (black carbon),      4: DU (dust), 
1496!--                           5: SS (sea salt),          6: HNO3 (nitric acid),
1497!--                           7: NH3 (ammonia)         
1498          CALL get_variable( id_faero, "profile_mass_fracs_a", pr_mass_fracs_a,&
1499                             0, nz_file-1, 0, maxspec-1 )
1500          CALL get_variable( id_faero, "profile_mass_fracs_b", pr_mass_fracs_b,&
1501                             0, nz_file-1, 0, maxspec-1 )
1502          CALL get_variable( id_faero, "profile_nsect", pr_nsect, 0, nz_file-1,&
1503                             0, nbins-1 )                   
1504         
1505          kk = 1
1506          DO  k = nzb, nz+1
1507             IF ( kk < nz_file )  THEN
1508                DO  WHILE ( pr_z(kk+1) <= zu(k) )
1509                   kk = kk + 1
1510                   IF ( kk == nz_file )  EXIT
1511                ENDDO
1512             ENDIF
1513             IF ( kk < nz_file )  THEN
1514!             
1515!--             Set initial value for gas compound tracers and initial values
1516                pvf2a(k,:) = pr_mass_fracs_a(:,kk) + ( zu(k) - pr_z(kk) ) / (  &
1517                            pr_z(kk+1) - pr_z(kk) ) * ( pr_mass_fracs_a(:,kk+1)&
1518                            - pr_mass_fracs_a(:,kk) )   
1519                pvf2b(k,:) = pr_mass_fracs_b(:,kk) + ( zu(k) - pr_z(kk) ) / (  &
1520                            pr_z(kk+1) - pr_z(kk) ) * ( pr_mass_fracs_b(:,kk+1)&
1521                            - pr_mass_fracs_b(:,kk) )             
1522                pndist(k,:) = pr_nsect(:,kk) + ( zu(k) - pr_z(kk) ) / (        &
1523                              pr_z(kk+1) - pr_z(kk) ) * ( pr_nsect(:,kk+1) -   &
1524                              pr_nsect(:,kk) )
1525             ELSE
1526                pvf2a(k,:) = pr_mass_fracs_a(:,kk)       
1527                pvf2b(k,:) = pr_mass_fracs_b(:,kk)
1528                pndist(k,:) = pr_nsect(:,kk)
1529             ENDIF
1530             IF ( iso4 < 0 )  THEN
1531                pvf2a(k,1) = 0.0_wp
1532                pvf2b(k,1) = 0.0_wp
1533             ENDIF
1534             IF ( ioc < 0 )  THEN
1535                pvf2a(k,2) = 0.0_wp
1536                pvf2b(k,2) = 0.0_wp
1537             ENDIF
1538             IF ( ibc < 0 )  THEN
1539                pvf2a(k,3) = 0.0_wp
1540                pvf2b(k,3) = 0.0_wp
1541             ENDIF
1542             IF ( idu < 0 )  THEN
1543                pvf2a(k,4) = 0.0_wp
1544                pvf2b(k,4) = 0.0_wp
1545             ENDIF
1546             IF ( iss < 0 )  THEN
1547                pvf2a(k,5) = 0.0_wp
1548                pvf2b(k,5) = 0.0_wp
1549             ENDIF
1550             IF ( ino < 0 )  THEN
1551                pvf2a(k,6) = 0.0_wp
1552                pvf2b(k,6) = 0.0_wp
1553             ENDIF
1554             IF ( inh < 0 )  THEN
1555                pvf2a(k,7) = 0.0_wp
1556                pvf2b(k,7) = 0.0_wp
1557             ENDIF
1558!
1559!--          Then normalise the mass fraction so that SUM = 1
1560             pvf2a(k,:) = pvf2a(k,:) / SUM( pvf2a(k,:) )
1561             IF ( SUM( pvf2b(k,:) ) > 0.0_wp ) pvf2b(k,:) = pvf2b(k,:) /       &
1562                                                            SUM( pvf2b(k,:) )
1563          ENDDO         
1564          DEALLOCATE( pr_z, pr_mass_fracs_a, pr_mass_fracs_b, pr_nsect )
1565       ELSE
1566          message_string = 'Input file '// TRIM( 'PIDS_SALSA' ) //             &
1567                           TRIM( coupling_char ) // ' for SALSA missing!'
1568          CALL message( 'salsa_mod: aerosol_init', 'SA0032', 1, 2, 0, 6, 0 )               
1569       ENDIF   ! netcdf_extend   
1570#endif
1571 
1572    ELSEIF ( isdtyp == 0 )  THEN
1573!
1574!--    Mass fractions for species in a and b-bins
1575       IF ( iso4 > 0 )  THEN
1576          pvf2a(:,1) = mass_fracs_a(iso4) 
1577          pvf2b(:,1) = mass_fracs_b(iso4)
1578       ENDIF
1579       IF ( ioc > 0 )  THEN
1580          pvf2a(:,2) = mass_fracs_a(ioc)
1581          pvf2b(:,2) = mass_fracs_b(ioc) 
1582       ENDIF
1583       IF ( ibc > 0 )  THEN
1584          pvf2a(:,3) = mass_fracs_a(ibc) 
1585          pvf2b(:,3) = mass_fracs_b(ibc)
1586       ENDIF
1587       IF ( idu > 0 )  THEN
1588          pvf2a(:,4) = mass_fracs_a(idu)
1589          pvf2b(:,4) = mass_fracs_b(idu) 
1590       ENDIF
1591       IF ( iss > 0 )  THEN
1592          pvf2a(:,5) = mass_fracs_a(iss)
1593          pvf2b(:,5) = mass_fracs_b(iss) 
1594       ENDIF
1595       IF ( ino > 0 )  THEN
1596          pvf2a(:,6) = mass_fracs_a(ino)
1597          pvf2b(:,6) = mass_fracs_b(ino)
1598       ENDIF
1599       IF ( inh > 0 )  THEN
1600          pvf2a(:,7) = mass_fracs_a(inh)
1601          pvf2b(:,7) = mass_fracs_b(inh)
1602       ENDIF
1603       DO  k = nzb, nz+1
1604          pvf2a(k,:) = pvf2a(k,:) / SUM( pvf2a(k,:) )
1605          IF ( SUM( pvf2b(k,:) ) > 0.0_wp ) pvf2b(k,:) = pvf2b(k,:) /          &
1606                                                         SUM( pvf2b(k,:) )
1607       ENDDO
1608       
1609       CALL size_distribution( n_lognorm, dpg, sigmag, nsect )
1610!
1611!--    Normalize by the given total number concentration
1612       nsect = nsect * SUM( n_lognorm ) * 1.0E+6_wp / SUM( nsect )     
1613       DO  b = in1a, fn2b
1614          pndist(:,b) = nsect(b)
1615       ENDDO
1616    ENDIF
1617   
1618    IF ( igctyp == 1 )  THEN
1619!
1620!--    Read input profiles from PIDS_CHEM   
1621#if defined( __netcdf )
1622!   
1623!--    Location-dependent size distributions and compositions.     
1624       INQUIRE( FILE='PIDS_CHEM' // TRIM( coupling_char ), EXIST=netcdf_extend )
1625       IF ( netcdf_extend  .AND.  .NOT. salsa_gases_from_chem )  THEN
1626!
1627!--       Open file in read-only mode     
1628          CALL open_read_file( 'PIDS_CHEM' // TRIM( coupling_char ), id_fchem )
1629!
1630!--       Input heights   
1631          CALL netcdf_data_input_get_dimension_length( id_fchem, nz_file,      &
1632                                                       "profile_z" ) 
1633          ALLOCATE( pr_z(nz_file), pr_gas(ngast,nz_file) ) 
1634          CALL get_variable( id_fchem, 'profile_z', pr_z ) 
1635!       
1636!--       Gases:
1637          CALL get_variable( id_fchem, "profile_H2SO4", pr_gas(1,:) )
1638          CALL get_variable( id_fchem, "profile_HNO3", pr_gas(2,:) )
1639          CALL get_variable( id_fchem, "profile_NH3", pr_gas(3,:) )
1640          CALL get_variable( id_fchem, "profile_OCNV", pr_gas(4,:) )
1641          CALL get_variable( id_fchem, "profile_OCSV", pr_gas(5,:) )
1642         
1643          kk = 1
1644          DO  k = nzb, nz+1
1645             IF ( kk < nz_file )  THEN
1646                DO  WHILE ( pr_z(kk+1) <= zu(k) )
1647                   kk = kk + 1
1648                   IF ( kk == nz_file )  EXIT
1649                ENDDO
1650             ENDIF
1651             IF ( kk < nz_file )  THEN
1652!             
1653!--             Set initial value for gas compound tracers and initial values
1654                DO  g = 1, ngast
1655                   salsa_gas(g)%init(k) =  pr_gas(g,kk) + ( zu(k) - pr_z(kk) ) &
1656                                           / ( pr_z(kk+1) - pr_z(kk) ) *       &
1657                                           ( pr_gas(g,kk+1) - pr_gas(g,kk) )
1658                   salsa_gas(g)%conc(k,:,:) = salsa_gas(g)%init(k)
1659                ENDDO
1660             ELSE
1661                DO  g = 1, ngast
1662                   salsa_gas(g)%init(k) =  pr_gas(g,kk) 
1663                   salsa_gas(g)%conc(k,:,:) = salsa_gas(g)%init(k)
1664                ENDDO
1665             ENDIF
1666          ENDDO
1667         
1668          DEALLOCATE( pr_z, pr_gas )
1669       ELSEIF ( .NOT. netcdf_extend  .AND.  .NOT.  salsa_gases_from_chem )  THEN
1670          message_string = 'Input file '// TRIM( 'PIDS_CHEM' ) //              &
1671                           TRIM( coupling_char ) // ' for SALSA missing!'
1672          CALL message( 'salsa_mod: aerosol_init', 'SA0033', 1, 2, 0, 6, 0 )               
1673       ENDIF   ! netcdf_extend     
1674#endif
1675
1676    ENDIF
1677
1678    IF ( ioc > 0  .AND.  iso4 > 0 )  THEN     
1679!--    Both are there, so use the given "massDistrA"
1680       pvfOC1a(:) = pvf2a(:,2) / ( pvf2a(:,2) + pvf2a(:,1) )  ! Normalize
1681    ELSEIF ( ioc > 0 )  THEN
1682!--    Pure organic carbon
1683       pvfOC1a(:) = 1.0_wp
1684    ELSEIF ( iso4 > 0 )  THEN
1685!--    Pure SO4
1686       pvfOC1a(:) = 0.0_wp   
1687    ELSE
1688       message_string = 'Either OC or SO4 must be active for aerosol region 1a!'
1689       CALL message( 'salsa_mod: aerosol_init', 'SA0021', 1, 2, 0, 6, 0 )
1690    ENDIF   
1691   
1692!
1693!-- Initialize concentrations
1694    DO  i = nxlg, nxrg
1695       DO  j = nysg, nyng
1696          DO  k = nzb, nzt+1
1697!
1698!--          Predetermine flag to mask topography         
1699             flag = MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_0(k,j,i), 0 ) )
1700!         
1701!--          a) Number concentrations
1702!--           Region 1:
1703             DO  b = in1a, fn1a
1704                aerosol_number(b)%conc(k,j,i) = pndist(k,b) * flag
1705                IF ( prunmode == 1 )  THEN
1706                   aerosol_number(b)%init = pndist(:,b)
1707                ENDIF
1708             ENDDO
1709!             
1710!--           Region 2:
1711             IF ( nreg > 1 )  THEN
1712                DO  b = in2a, fn2a
1713                   aerosol_number(b)%conc(k,j,i) = MAX( 0.0_wp, pnf2a(k) ) *   &
1714                                                    pndist(k,b) * flag
1715                   IF ( prunmode == 1 )  THEN
1716                      aerosol_number(b)%init = MAX( 0.0_wp, nf2a ) * pndist(:,b)
1717                   ENDIF
1718                ENDDO
1719                IF ( .NOT. no_insoluble )  THEN
1720                   DO  b = in2b, fn2b
1721                      IF ( pnf2a(k) < 1.0_wp )  THEN             
1722                         aerosol_number(b)%conc(k,j,i) = MAX( 0.0_wp, 1.0_wp   &
1723                                               - pnf2a(k) ) * pndist(k,b) * flag
1724                         IF ( prunmode == 1 )  THEN
1725                            aerosol_number(b)%init = MAX( 0.0_wp, 1.0_wp -     &
1726                                                          nf2a ) * pndist(:,b)
1727                         ENDIF
1728                      ENDIF
1729                   ENDDO
1730                ENDIF
1731             ENDIF
1732!
1733!--          b) Aerosol mass concentrations
1734!--             bin subrange 1: done here separately due to the SO4/OC convention
1735!--          SO4:
1736             IF ( iso4 > 0 )  THEN
1737                ss = ( iso4 - 1 ) * nbins + in1a !< start
1738                ee = ( iso4 - 1 ) * nbins + fn1a !< end
1739                b = in1a
1740                DO  c = ss, ee
1741                   aerosol_mass(c)%conc(k,j,i) = MAX( 0.0_wp, 1.0_wp -         &
1742                                                  pvfOC1a(k) ) * pndist(k,b) * &
1743                                                  core(b) * arhoh2so4 * flag
1744                   IF ( prunmode == 1 )  THEN
1745                      aerosol_mass(c)%init = MAX( 0.0_wp, 1.0_wp - MAXVAL(     &
1746                                             pvfOC1a ) ) * pndist(:,b) *       &
1747                                             core(b) * arhoh2so4
1748                   ENDIF
1749                   b = b+1
1750                ENDDO
1751             ENDIF
1752!--          OC:
1753             IF ( ioc > 0 ) THEN
1754                ss = ( ioc - 1 ) * nbins + in1a !< start
1755                ee = ( ioc - 1 ) * nbins + fn1a !< end
1756                b = in1a
1757                DO  c = ss, ee
1758                   aerosol_mass(c)%conc(k,j,i) = MAX( 0.0_wp, pvfOC1a(k) ) *   &
1759                                           pndist(k,b) * core(b) * arhooc * flag
1760                   IF ( prunmode == 1 )  THEN
1761                      aerosol_mass(c)%init = MAX( 0.0_wp, MAXVAL( pvfOC1a ) )  &
1762                                             * pndist(:,b) *  core(b) * arhooc
1763                   ENDIF
1764                   b = b+1
1765                ENDDO
1766             ENDIF
1767             
1768             prunmode = 3  ! Init only once
1769 
1770          ENDDO !< k
1771       ENDDO !< j
1772    ENDDO !< i
1773   
1774!
1775!-- c) Aerosol mass concentrations
1776!--    bin subrange 2:
1777    IF ( nreg > 1 ) THEN
1778   
1779       IF ( iso4 > 0 ) THEN
1780          CALL set_aero_mass( iso4, pvf2a(:,1), pvf2b(:,1), pnf2a, pndist,     &
1781                              core, arhoh2so4 )
1782       ENDIF
1783       IF ( ioc > 0 ) THEN
1784          CALL set_aero_mass( ioc, pvf2a(:,2), pvf2b(:,2), pnf2a, pndist, core,&
1785                              arhooc )
1786       ENDIF
1787       IF ( ibc > 0 ) THEN
1788          CALL set_aero_mass( ibc, pvf2a(:,3), pvf2b(:,3), pnf2a, pndist, core,&
1789                              arhobc )
1790       ENDIF
1791       IF ( idu > 0 ) THEN
1792          CALL set_aero_mass( idu, pvf2a(:,4), pvf2b(:,4), pnf2a, pndist, core,&
1793                              arhodu )
1794       ENDIF
1795       IF ( iss > 0 ) THEN
1796          CALL set_aero_mass( iss, pvf2a(:,5), pvf2b(:,5), pnf2a, pndist, core,&
1797                              arhoss )
1798       ENDIF
1799       IF ( ino > 0 ) THEN
1800          CALL set_aero_mass( ino, pvf2a(:,6), pvf2b(:,6), pnf2a, pndist, core,&
1801                              arhohno3 )
1802       ENDIF
1803       IF ( inh > 0 ) THEN
1804          CALL set_aero_mass( inh, pvf2a(:,7), pvf2b(:,7), pnf2a, pndist, core,&
1805                              arhonh3 )
1806       ENDIF
1807
1808    ENDIF
1809   
1810 END SUBROUTINE aerosol_init
1811 
1812!------------------------------------------------------------------------------!
1813! Description:
1814! ------------
1815!> Create a lognormal size distribution and discretise to a sectional
1816!> representation.
1817!------------------------------------------------------------------------------!
1818 SUBROUTINE size_distribution( in_ntot, in_dpg, in_sigma, psd_sect )
1819   
1820    IMPLICIT NONE
1821   
1822!-- Log-normal size distribution: modes   
1823    REAL(wp), DIMENSION(:), INTENT(in) ::  in_dpg    !< geometric mean diameter
1824                                                     !< (micrometres)
1825    REAL(wp), DIMENSION(:), INTENT(in) ::  in_ntot   !< number conc. (#/cm3)
1826    REAL(wp), DIMENSION(:), INTENT(in) ::  in_sigma  !< standard deviation
1827    REAL(wp), DIMENSION(:), INTENT(inout) ::  psd_sect !< sectional size
1828                                                       !< distribution
1829    INTEGER(iwp) ::  b          !< running index: bin
1830    INTEGER(iwp) ::  ib         !< running index: iteration
1831    REAL(wp) ::  d1             !< particle diameter (m, dummy)
1832    REAL(wp) ::  d2             !< particle diameter (m, dummy)
1833    REAL(wp) ::  delta_d        !< (d2-d1)/10                                                     
1834    REAL(wp) ::  deltadp        !< bin width
1835    REAL(wp) ::  dmidi          !< ( d1 + d2 ) / 2
1836   
1837    DO  b = in1a, fn2b !< aerosol size bins
1838       psd_sect(b) = 0.0_wp
1839!--    Particle diameter at the low limit (largest in the bin) (m)
1840       d1 = ( aero(b)%vlolim / api6 ) ** ( 1.0_wp / 3.0_wp )
1841!--    Particle diameter at the high limit (smallest in the bin) (m)
1842       d2 = ( aero(b)%vhilim / api6 ) ** ( 1.0_wp / 3.0_wp )
1843!--    Span of particle diameter in a bin (m)
1844       delta_d = ( d2 - d1 ) / 10.0_wp
1845!--    Iterate:             
1846       DO  ib = 1, 10
1847          d1 = ( aero(b)%vlolim / api6 ) ** ( 1.0_wp / 3.0_wp ) + ( ib - 1)    &
1848               * delta_d
1849          d2 = d1 + delta_d
1850          dmidi = ( d1 + d2 ) / 2.0_wp
1851          deltadp = LOG10( d2 / d1 )
1852         
1853!--       Size distribution
1854!--       in_ntot = total number, total area, or total volume concentration
1855!--       in_dpg = geometric-mean number, area, or volume diameter
1856!--       n(k) = number, area, or volume concentration in a bin
1857!--       n_lognorm and dpg converted to units of #/m3 and m
1858          psd_sect(b) = psd_sect(b) + SUM( in_ntot * 1.0E+6_wp * deltadp /     &
1859                     ( SQRT( 2.0_wp * pi ) * LOG10( in_sigma ) ) *             &
1860                     EXP( -LOG10( dmidi / ( 1.0E-6_wp * in_dpg ) )**2.0_wp /   &
1861                     ( 2.0_wp * LOG10( in_sigma ) ** 2.0_wp ) ) )
1862 
1863       ENDDO
1864    ENDDO
1865   
1866 END SUBROUTINE size_distribution
1867
1868!------------------------------------------------------------------------------!
1869! Description:
1870! ------------
1871!> Sets the mass concentrations to aerosol arrays in 2a and 2b.
1872!>
1873!> Tomi Raatikainen, FMI, 29.2.2016
1874!------------------------------------------------------------------------------!
1875 SUBROUTINE set_aero_mass( ispec, ppvf2a, ppvf2b, ppnf2a, ppndist, pcore, prho )
1876   
1877    IMPLICIT NONE
1878
1879    INTEGER(iwp), INTENT(in) :: ispec  !< Aerosol species index
1880    REAL(wp), INTENT(in) ::  pcore(nbins) !< Aerosol bin mid core volume   
1881    REAL(wp), INTENT(in) ::  ppndist(0:nz+1,nbins) !< Aerosol size distribution
1882    REAL(wp), INTENT(in) ::  ppnf2a(0:nz+1) !< Number fraction for 2a   
1883    REAL(wp), INTENT(in) ::  ppvf2a(0:nz+1) !< Mass distributions for a
1884    REAL(wp), INTENT(in) ::  ppvf2b(0:nz+1) !< and b bins   
1885    REAL(wp), INTENT(in) ::  prho !< Aerosol density
1886    INTEGER(iwp) ::  b  !< loop index
1887    INTEGER(iwp) ::  c  !< loop index       
1888    INTEGER(iwp) ::  ee !< index: end
1889    INTEGER(iwp) ::  i  !< loop index
1890    INTEGER(iwp) ::  j  !< loop index
1891    INTEGER(iwp) ::  k  !< loop index
1892    INTEGER(iwp) ::  prunmode  !< 1 = initialise
1893    INTEGER(iwp) ::  ss !< index: start
1894    REAL(wp) ::  flag   !< flag to mask topography grid points
1895   
1896    prunmode = 1
1897   
1898    DO i = nxlg, nxrg
1899       DO j = nysg, nyng
1900          DO k = nzb, nzt+1 
1901!
1902!--          Predetermine flag to mask topography
1903             flag = MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_0(k,j,i), 0 ) ) 
1904!             
1905!--          Regime 2a:
1906             ss = ( ispec - 1 ) * nbins + in2a
1907             ee = ( ispec - 1 ) * nbins + fn2a
1908             b = in2a
1909             DO c = ss, ee
1910                aerosol_mass(c)%conc(k,j,i) = MAX( 0.0_wp, ppvf2a(k) ) *       &
1911                               ppnf2a(k) * ppndist(k,b) * pcore(b) * prho * flag
1912                IF ( prunmode == 1 )  THEN
1913                   aerosol_mass(c)%init = MAX( 0.0_wp, MAXVAL( ppvf2a(:) ) ) * &
1914                                          MAXVAL( ppnf2a ) * pcore(b) * prho * &
1915                                          MAXVAL( ppndist(:,b) ) 
1916                ENDIF
1917                b = b+1
1918             ENDDO
1919!--          Regime 2b:
1920             IF ( .NOT. no_insoluble )  THEN
1921                ss = ( ispec - 1 ) * nbins + in2b
1922                ee = ( ispec - 1 ) * nbins + fn2b
1923                b = in2a
1924                DO c = ss, ee
1925                   aerosol_mass(c)%conc(k,j,i) = MAX( 0.0_wp, ppvf2b(k) ) * (  &
1926                                         1.0_wp - ppnf2a(k) ) * ppndist(k,b) * &
1927                                         pcore(b) * prho * flag
1928                   IF ( prunmode == 1 )  THEN
1929                      aerosol_mass(c)%init = MAX( 0.0_wp, MAXVAL( ppvf2b(:) ) )&
1930                                        * ( 1.0_wp - MAXVAL( ppnf2a ) ) *      &
1931                                        MAXVAL( ppndist(:,b) ) * pcore(b) * prho
1932                   ENDIF
1933                   b = b+1
1934                ENDDO
1935             ENDIF
1936             prunmode = 3  ! Init only once
1937          ENDDO
1938       ENDDO
1939    ENDDO
1940 END SUBROUTINE set_aero_mass
1941
1942!------------------------------------------------------------------------------!
1943! Description:
1944! ------------
1945!> Swapping of timelevels
1946!------------------------------------------------------------------------------!
1947 SUBROUTINE salsa_swap_timelevel( mod_count )
1948
1949    IMPLICIT NONE
1950
1951    INTEGER(iwp), INTENT(IN) ::  mod_count  !<
1952    INTEGER(iwp) ::  b  !<   
1953    INTEGER(iwp) ::  c  !<   
1954    INTEGER(iwp) ::  cc !<
1955    INTEGER(iwp) ::  g  !<
1956
1957!
1958!-- Example for prognostic variable "prog_var"
1959#if defined( __nopointer )
1960    IF ( myid == 0 )  THEN
1961       message_string =  ' SALSA runs only with POINTER Version'
1962       CALL message( 'salsa_swap_timelevel', 'SA0022', 1, 2, 0, 6, 0 )
1963    ENDIF
1964#else
1965   
1966    SELECT CASE ( mod_count )
1967
1968       CASE ( 0 )
1969
1970          DO  b = 1, nbins
1971             aerosol_number(b)%conc(nzb:nzt+1,nysg:nyng,nxlg:nxrg)   =>        &
1972                nconc_1(:,:,:,b)
1973             aerosol_number(b)%conc_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg) =>        &
1974                nconc_2(:,:,:,b)
1975             DO  c = 1, ncc_tot
1976                cc = ( c-1 ) * nbins + b  ! required due to possible Intel18 bug
1977                aerosol_mass(cc)%conc(nzb:nzt+1,nysg:nyng,nxlg:nxrg)   =>      &
1978                   mconc_1(:,:,:,cc)
1979                aerosol_mass(cc)%conc_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg) =>      &
1980                   mconc_2(:,:,:,cc)
1981             ENDDO
1982          ENDDO
1983         
1984          IF ( .NOT. salsa_gases_from_chem )  THEN
1985             DO  g = 1, ngast
1986                salsa_gas(g)%conc(nzb:nzt+1,nysg:nyng,nxlg:nxrg)   =>          &
1987                   gconc_1(:,:,:,g)
1988                salsa_gas(g)%conc_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg) =>          &
1989                   gconc_2(:,:,:,g)
1990             ENDDO
1991          ENDIF
1992
1993       CASE ( 1 )
1994
1995          DO  b = 1, nbins
1996             aerosol_number(b)%conc(nzb:nzt+1,nysg:nyng,nxlg:nxrg)   =>        &
1997                nconc_2(:,:,:,b)
1998             aerosol_number(b)%conc_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg) =>        &
1999                nconc_1(:,:,:,b)
2000             DO  c = 1, ncc_tot
2001                cc = ( c-1 ) * nbins + b  ! required due to possible Intel18 bug
2002                aerosol_mass(cc)%conc(nzb:nzt+1,nysg:nyng,nxlg:nxrg)   =>      &
2003                   mconc_2(:,:,:,cc)
2004                aerosol_mass(cc)%conc_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg) =>      &
2005                   mconc_1(:,:,:,cc)
2006             ENDDO
2007          ENDDO
2008         
2009          IF ( .NOT. salsa_gases_from_chem )  THEN
2010             DO  g = 1, ngast
2011                salsa_gas(g)%conc(nzb:nzt+1,nysg:nyng,nxlg:nxrg)   =>          &
2012                   gconc_2(:,:,:,g)
2013                salsa_gas(g)%conc_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg) =>          &
2014                   gconc_1(:,:,:,g)
2015             ENDDO
2016          ENDIF
2017
2018    END SELECT
2019#endif
2020
2021 END SUBROUTINE salsa_swap_timelevel
2022
2023
2024!------------------------------------------------------------------------------!
2025! Description:
2026! ------------
2027!> This routine reads the respective restart data.
2028!------------------------------------------------------------------------------!
2029 SUBROUTINE salsa_rrd_local( i, k, nxlf, nxlc, nxl_on_file, nxrf, nxrc,        &
2030                             nxr_on_file, nynf, nync, nyn_on_file, nysf,       &
2031                             nysc, nys_on_file, tmp_3d, found )
2032
2033   
2034    IMPLICIT NONE
2035   
2036    CHARACTER (LEN=20) :: field_char   !<
2037    INTEGER(iwp) ::  b  !<   
2038    INTEGER(iwp) ::  c  !<
2039    INTEGER(iwp) ::  g  !<
2040    INTEGER(iwp) ::  i  !<
2041    INTEGER(iwp) ::  k  !<
2042    INTEGER(iwp) ::  nxlc            !<
2043    INTEGER(iwp) ::  nxlf            !<
2044    INTEGER(iwp) ::  nxl_on_file     !<
2045    INTEGER(iwp) ::  nxrc            !<
2046    INTEGER(iwp) ::  nxrf            !<
2047    INTEGER(iwp) ::  nxr_on_file     !<
2048    INTEGER(iwp) ::  nync            !<
2049    INTEGER(iwp) ::  nynf            !<
2050    INTEGER(iwp) ::  nyn_on_file     !<
2051    INTEGER(iwp) ::  nysc            !<
2052    INTEGER(iwp) ::  nysf            !<
2053    INTEGER(iwp) ::  nys_on_file     !<
2054
2055    LOGICAL, INTENT(OUT)  ::  found
2056
2057    REAL(wp), &
2058       DIMENSION(nzb:nzt+1,nys_on_file-nbgp:nyn_on_file+nbgp,nxl_on_file-nbgp:nxr_on_file+nbgp) :: tmp_3d   !<
2059       
2060    found = .FALSE.
2061   
2062    IF ( read_restart_data_salsa )  THEN
2063   
2064       SELECT CASE ( restart_string(1:length) )
2065       
2066          CASE ( 'aerosol_number' )
2067             DO  b = 1, nbins
2068                IF ( k == 1 )  READ ( 13 ) tmp_3d
2069                aerosol_number(b)%conc(:,nysc-nbgp:nync+nbgp,nxlc-nbgp:nxrc+nbgp) = & 
2070                               tmp_3d(:,nysf-nbgp:nynf+nbgp,nxlf-nbgp:nxrf+nbgp)
2071                found = .TRUE.
2072             ENDDO
2073       
2074          CASE ( 'aerosol_mass' )
2075             DO  c = 1, ncc_tot * nbins
2076                IF ( k == 1 )  READ ( 13 ) tmp_3d
2077                aerosol_mass(c)%conc(:,nysc-nbgp:nync+nbgp,nxlc-nbgp:nxrc+nbgp) = & 
2078                               tmp_3d(:,nysf-nbgp:nynf+nbgp,nxlf-nbgp:nxrf+nbgp)
2079                found = .TRUE.
2080             ENDDO
2081         
2082          CASE ( 'salsa_gas' )
2083             DO  g = 1, ngast
2084                IF ( k == 1 )  READ ( 13 ) tmp_3d
2085                salsa_gas(g)%conc(:,nysc-nbgp:nync+nbgp,nxlc-nbgp:nxrc+nbgp) =  & 
2086                               tmp_3d(:,nysf-nbgp:nynf+nbgp,nxlf-nbgp:nxrf+nbgp)
2087                found = .TRUE.
2088             ENDDO
2089             
2090          CASE DEFAULT
2091             found = .FALSE.
2092             
2093       END SELECT
2094    ENDIF
2095
2096 END SUBROUTINE salsa_rrd_local
2097   
2098
2099!------------------------------------------------------------------------------!
2100! Description:
2101! ------------
2102!> This routine writes the respective restart data.
2103!> Note that the following input variables in PARIN have to be equal between
2104!> restart runs:
2105!>    listspec, nbin, nbin2, nf2a, ncc, mass_fracs_a, mass_fracs_b
2106!------------------------------------------------------------------------------!
2107 SUBROUTINE salsa_wrd_local
2108
2109    IMPLICIT NONE
2110   
2111    INTEGER(iwp) ::  b  !<   
2112    INTEGER(iwp) ::  c  !<
2113    INTEGER(iwp) ::  g  !<
2114   
2115    IF ( write_binary  .AND.  write_binary_salsa )  THEN
2116   
2117       CALL wrd_write_string( 'aerosol_number' )
2118       DO  b = 1, nbins
2119          WRITE ( 14 )  aerosol_number(b)%conc
2120       ENDDO
2121       
2122       CALL wrd_write_string( 'aerosol_mass' )
2123       DO  c = 1, nbins*ncc_tot
2124          WRITE ( 14 )  aerosol_mass(c)%conc
2125       ENDDO
2126       
2127       CALL wrd_write_string( 'salsa_gas' )
2128       DO  g = 1, ngast
2129          WRITE ( 14 )  salsa_gas(g)%conc
2130       ENDDO
2131         
2132    ENDIF
2133       
2134 END SUBROUTINE salsa_wrd_local   
2135
2136
2137!------------------------------------------------------------------------------!
2138! Description:
2139! ------------
2140!> Performs necessary unit and dimension conversion between the host model and
2141!> SALSA module, and calls the main SALSA routine.
2142!> Partially adobted form the original SALSA boxmodel version.
2143!> Now takes masses in as kg/kg from LES!! Converted to m3/m3 for SALSA
2144!> 05/2016 Juha: This routine is still pretty much in its original shape.
2145!>               It's dumb as a mule and twice as ugly, so implementation of
2146!>               an improved solution is necessary sooner or later.
2147!> Juha Tonttila, FMI, 2014
2148!> Jaakko Ahola, FMI, 2016
2149!> Only aerosol processes included, Mona Kurppa, UHel, 2017
2150!------------------------------------------------------------------------------!
2151 SUBROUTINE salsa_driver( i, j, prunmode )
2152
2153    USE arrays_3d,                                                             &
2154        ONLY: pt_p, q_p, rho_air_zw, u, v, w
2155       
2156    USE plant_canopy_model_mod,                                                &
2157        ONLY: lad_s
2158       
2159    USE surface_mod,                                                           &
2160        ONLY:  surf_def_h, surf_def_v, surf_lsm_h, surf_lsm_v, surf_usm_h,     &
2161               surf_usm_v
2162 
2163    IMPLICIT NONE
2164   
2165    INTEGER(iwp), INTENT(in) ::  i   !< loop index
2166    INTEGER(iwp), INTENT(in) ::  j   !< loop index
2167    INTEGER(iwp), INTENT(in) ::  prunmode !< 1: Initialization call
2168                                          !< 2: Spinup period call
2169                                          !< 3: Regular runtime call
2170!-- Local variables
2171    TYPE(t_section), DIMENSION(fn2b) ::  aero_old !< helper array
2172    INTEGER(iwp) ::  bb     !< loop index
2173    INTEGER(iwp) ::  cc     !< loop index
2174    INTEGER(iwp) ::  endi   !< end index
2175    INTEGER(iwp) ::  k_wall !< vertical index of topography top
2176    INTEGER(iwp) ::  k      !< loop index
2177    INTEGER(iwp) ::  l      !< loop index
2178    INTEGER(iwp) ::  nc_h2o !< index of H2O in the prtcl index table
2179    INTEGER(iwp) ::  ss     !< loop index
2180    INTEGER(iwp) ::  str    !< start index
2181    INTEGER(iwp) ::  vc     !< default index in prtcl
2182    REAL(wp) ::  cw_old     !< previous H2O mixing ratio
2183    REAL(wp) ::  flag       !< flag to mask topography grid points
2184    REAL(wp), DIMENSION(nzb:nzt+1) ::  in_adn !< air density (kg/m3)   
2185    REAL(wp), DIMENSION(nzb:nzt+1) ::  in_cs  !< H2O sat. vapour conc.
2186    REAL(wp), DIMENSION(nzb:nzt+1) ::  in_cw  !< H2O vapour concentration
2187    REAL(wp) ::  in_lad                       !< leaf area density (m2/m3)
2188    REAL(wp), DIMENSION(nzb:nzt+1) ::  in_p   !< pressure (Pa)     
2189    REAL(wp) ::  in_rh                        !< relative humidity                     
2190    REAL(wp), DIMENSION(nzb:nzt+1) ::  in_t   !< temperature (K)
2191    REAL(wp), DIMENSION(nzb:nzt+1) ::  in_u   !< wind magnitude (m/s)
2192    REAL(wp), DIMENSION(nzb:nzt+1) ::  kvis   !< kinematic viscosity of air(m2/s)                                           
2193    REAL(wp), DIMENSION(nzb:nzt+1,fn2b) ::  Sc      !< particle Schmidt number   
2194    REAL(wp), DIMENSION(nzb:nzt+1,fn2b) ::  vd      !< particle fall seed (m/s,
2195                                                    !< sedimentation velocity)
2196    REAL(wp), DIMENSION(nzb:nzt+1) ::  ppm_to_nconc !< Conversion factor
2197                                                    !< from ppm to #/m3                                                     
2198    REAL(wp) ::  zgso4  !< SO4
2199    REAL(wp) ::  zghno3 !< HNO3
2200    REAL(wp) ::  zgnh3  !< NH3
2201    REAL(wp) ::  zgocnv !< non-volatile OC
2202    REAL(wp) ::  zgocsv !< semi-volatile OC
2203   
2204    aero_old(:)%numc = 0.0_wp
2205    in_adn           = 0.0_wp   
2206    in_cs            = 0.0_wp
2207    in_cw            = 0.0_wp 
2208    in_lad           = 0.0_wp
2209    in_rh            = 0.0_wp
2210    in_p             = 0.0_wp 
2211    in_t             = 0.0_wp 
2212    in_u             = 0.0_wp
2213    kvis             = 0.0_wp
2214    Sc               = 0.0_wp
2215    vd               = 0.0_wp
2216    ppm_to_nconc     = 1.0_wp
2217    zgso4            = nclim
2218    zghno3           = nclim
2219    zgnh3            = nclim
2220    zgocnv           = nclim
2221    zgocsv           = nclim
2222   
2223!       
2224!-- Aerosol number is always set, but mass can be uninitialized
2225    DO cc = 1, nbins
2226       aero(cc)%volc     = 0.0_wp
2227       aero_old(cc)%volc = 0.0_wp
2228    ENDDO
2229!   
2230!-- Set the salsa runtime config (How to make this more efficient?)
2231    CALL set_salsa_runtime( prunmode )
2232!             
2233!-- Calculate thermodynamic quantities needed in SALSA
2234    CALL salsa_thrm_ij( i, j, p_ij=in_p, temp_ij=in_t, cw_ij=in_cw,            &
2235                        cs_ij=in_cs, adn_ij=in_adn )
2236!
2237!-- Magnitude of wind: needed for deposition
2238    IF ( lsdepo )  THEN
2239       in_u(nzb+1:nzt) = SQRT(                                                 &
2240                   ( 0.5_wp * ( u(nzb+1:nzt,j,i) + u(nzb+1:nzt,j,i+1) ) )**2 + & 
2241                   ( 0.5_wp * ( v(nzb+1:nzt,j,i) + v(nzb+1:nzt,j+1,i) ) )**2 + &
2242                   ( 0.5_wp * ( w(nzb:nzt-1,j,i) + w(nzb+1:nzt,j,  i) ) )**2 )
2243    ENDIF
2244!
2245!-- Calculate conversion factors for gas concentrations
2246    ppm_to_nconc = for_ppm_to_nconc * in_p / in_t
2247!
2248!-- Determine topography-top index on scalar grid
2249    k_wall = MAXLOC( MERGE( 1, 0, BTEST( wall_flags_0(:,j,i), 12 ) ),          &
2250                     DIM = 1 ) - 1     
2251               
2252    DO k = nzb+1, nzt
2253!
2254!--    Predetermine flag to mask topography
2255       flag = MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_0(k,j,i), 0 ) )
2256!       
2257!--    Do not run inside buildings       
2258       IF ( flag == 0.0_wp )  CYCLE   
2259!
2260!--    Wind velocity for dry depositon on vegetation   
2261       IF ( lsdepo_vege  .AND.  plant_canopy  )  THEN
2262          in_lad = lad_s(k-k_wall,j,i)
2263       ENDIF       
2264!
2265!--    For initialization and spinup, limit the RH with the parameter rhlim
2266       IF ( prunmode < 3 ) THEN
2267          in_cw(k) = MIN( in_cw(k), in_cs(k) * rhlim )
2268       ELSE
2269          in_cw(k) = in_cw(k)
2270       ENDIF
2271       cw_old = in_cw(k) !* in_adn(k)
2272!               
2273!--    Set volume concentrations:
2274!--    Sulphate (SO4) or sulphuric acid H2SO4
2275       IF ( iso4 > 0 )  THEN
2276          vc = 1
2277          str = ( iso4-1 ) * nbins + 1    ! start index
2278          endi = iso4 * nbins             ! end index
2279          cc = 1
2280          DO ss = str, endi
2281             aero(cc)%volc(vc) = aerosol_mass(ss)%conc(k,j,i) / arhoh2so4
2282             cc = cc+1
2283          ENDDO
2284          aero_old(1:nbins)%volc(vc) = aero(1:nbins)%volc(vc)
2285       ENDIF
2286       
2287!--    Organic carbon (OC) compounds
2288       IF ( ioc > 0 )  THEN
2289          vc = 2
2290          str = ( ioc-1 ) * nbins + 1
2291          endi = ioc * nbins
2292          cc = 1
2293          DO ss = str, endi
2294             aero(cc)%volc(vc) = aerosol_mass(ss)%conc(k,j,i) / arhooc
2295             cc = cc+1
2296          ENDDO
2297          aero_old(1:nbins)%volc(vc) = aero(1:nbins)%volc(vc)
2298       ENDIF
2299       
2300!--    Black carbon (BC)
2301       IF ( ibc > 0 )  THEN
2302          vc = 3
2303          str = ( ibc-1 ) * nbins + 1 + fn1a
2304          endi = ibc * nbins
2305          cc = 1 + fn1a
2306          DO ss = str, endi
2307             aero(cc)%volc(vc) = aerosol_mass(ss)%conc(k,j,i) / arhobc
2308             cc = cc+1
2309          ENDDO                   
2310          aero_old(1:nbins)%volc(vc) = aero(1:nbins)%volc(vc)
2311       ENDIF
2312
2313!--    Dust (DU)
2314       IF ( idu > 0 )  THEN
2315          vc = 4
2316          str = ( idu-1 ) * nbins + 1 + fn1a
2317          endi = idu * nbins
2318          cc = 1 + fn1a
2319          DO ss = str, endi
2320             aero(cc)%volc(vc) = aerosol_mass(ss)%conc(k,j,i) / arhodu
2321             cc = cc+1
2322          ENDDO
2323          aero_old(1:nbins)%volc(vc) = aero(1:nbins)%volc(vc)
2324       ENDIF
2325
2326!--    Sea salt (SS)
2327       IF ( iss > 0 )  THEN
2328          vc = 5
2329          str = ( iss-1 ) * nbins + 1 + fn1a
2330          endi = iss * nbins
2331          cc = 1 + fn1a
2332          DO ss = str, endi
2333             aero(cc)%volc(vc) = aerosol_mass(ss)%conc(k,j,i) / arhoss 
2334             cc = cc+1
2335          ENDDO
2336          aero_old(1:nbins)%volc(vc) = aero(1:nbins)%volc(vc)
2337       ENDIF
2338
2339!--    Nitrate (NO(3-)) or nitric acid HNO3
2340       IF ( ino > 0 )  THEN
2341          vc = 6
2342          str = ( ino-1 ) * nbins + 1 
2343          endi = ino * nbins
2344          cc = 1
2345          DO ss = str, endi
2346             aero(cc)%volc(vc) = aerosol_mass(ss)%conc(k,j,i) / arhohno3
2347             cc = cc+1
2348          ENDDO
2349          aero_old(1:nbins)%volc(vc) = aero(1:nbins)%volc(vc)
2350       ENDIF
2351
2352!--    Ammonium (NH(4+)) or ammonia NH3
2353       IF ( inh > 0 )  THEN
2354          vc = 7
2355          str = ( inh-1 ) * nbins + 1
2356          endi = inh * nbins
2357          cc = 1
2358          DO ss = str, endi
2359             aero(cc)%volc(vc) = aerosol_mass(ss)%conc(k,j,i) / arhonh3 
2360             cc = cc+1
2361          ENDDO
2362          aero_old(1:nbins)%volc(vc) = aero(1:nbins)%volc(vc)
2363       ENDIF
2364
2365!--    Water (always used)
2366       nc_h2o = get_index( prtcl,'H2O' )
2367       vc = 8
2368       str = ( nc_h2o-1 ) * nbins + 1
2369       endi = nc_h2o * nbins
2370       cc = 1
2371       IF ( advect_particle_water )  THEN
2372          DO ss = str, endi
2373             aero(cc)%volc(vc) = aerosol_mass(ss)%conc(k,j,i) / arhoh2o
2374             cc = cc+1
2375          ENDDO
2376       ELSE
2377         aero(1:nbins)%volc(vc) = mclim
2378       ENDIF
2379       aero_old(1:nbins)%volc(vc) = aero(1:nbins)%volc(vc)
2380!
2381!--    Number concentrations (numc) and particle sizes
2382!--    (dwet = wet diameter, core = dry volume)
2383       DO  bb = 1, nbins
2384          aero(bb)%numc = aerosol_number(bb)%conc(k,j,i) 
2385          aero_old(bb)%numc = aero(bb)%numc
2386          IF ( aero(bb)%numc > nclim )  THEN
2387             aero(bb)%dwet = ( SUM( aero(bb)%volc(:) ) / aero(bb)%numc / api6 )&
2388                                **( 1.0_wp / 3.0_wp )
2389             aero(bb)%core = SUM( aero(bb)%volc(1:7) ) / aero(bb)%numc
2390          ELSE
2391             aero(bb)%dwet = aero(bb)%dmid
2392             aero(bb)%core = api6 * ( aero(bb)%dwet ) ** 3.0_wp
2393          ENDIF
2394       ENDDO
2395!       
2396!--    On EACH call of salsa_driver, calculate the ambient sizes of
2397!--    particles by equilibrating soluble fraction of particles with water
2398!--    using the ZSR method.
2399       in_rh = in_cw(k) / in_cs(k)
2400       IF ( prunmode==1  .OR.  .NOT. advect_particle_water )  THEN
2401          CALL equilibration( in_rh, in_t(k), aero, .TRUE. )
2402       ENDIF
2403!
2404!--    Gaseous tracer concentrations in #/m3
2405       IF ( salsa_gases_from_chem )  THEN       
2406!       
2407!--       Convert concentrations in ppm to #/m3
2408          zgso4  = chem_species(gas_index_chem(1))%conc(k,j,i) * ppm_to_nconc(k)
2409          zghno3 = chem_species(gas_index_chem(2))%conc(k,j,i) * ppm_to_nconc(k)
2410          zgnh3  = chem_species(gas_index_chem(3))%conc(k,j,i) * ppm_to_nconc(k)
2411          zgocnv = chem_species(gas_index_chem(4))%conc(k,j,i) * ppm_to_nconc(k)     
2412          zgocsv = chem_species(gas_index_chem(5))%conc(k,j,i) * ppm_to_nconc(k)                 
2413       ELSE
2414          zgso4  = salsa_gas(1)%conc(k,j,i) 
2415          zghno3 = salsa_gas(2)%conc(k,j,i) 
2416          zgnh3  = salsa_gas(3)%conc(k,j,i) 
2417          zgocnv = salsa_gas(4)%conc(k,j,i) 
2418          zgocsv = salsa_gas(5)%conc(k,j,i)
2419       ENDIF   
2420!
2421!--    ***************************************!
2422!--                   Run SALSA               !
2423!--    ***************************************!
2424       CALL run_salsa( in_p(k), in_cw(k), in_cs(k), in_t(k), in_u(k),          &
2425                       in_adn(k), in_lad, zgso4, zgocnv, zgocsv, zghno3, zgnh3,&
2426                       aero, prtcl, kvis(k), Sc(k,:), vd(k,:), dt_salsa )
2427!--    ***************************************!
2428       IF ( lsdepo ) sedim_vd(k,j,i,:) = vd(k,:)
2429!                           
2430!--    Calculate changes in concentrations
2431       DO bb = 1, nbins
2432          aerosol_number(bb)%conc(k,j,i) = aerosol_number(bb)%conc(k,j,i)      &
2433                                 +  ( aero(bb)%numc - aero_old(bb)%numc ) * flag
2434       ENDDO
2435       
2436       IF ( iso4 > 0 )  THEN
2437          vc = 1
2438          str = ( iso4-1 ) * nbins + 1
2439          endi = iso4 * nbins
2440          cc = 1
2441          DO ss = str, endi
2442             aerosol_mass(ss)%conc(k,j,i) = aerosol_mass(ss)%conc(k,j,i)       &
2443                               + ( aero(cc)%volc(vc) - aero_old(cc)%volc(vc) ) &
2444                               * arhoh2so4 * flag
2445             cc = cc+1
2446          ENDDO
2447       ENDIF
2448       
2449       IF ( ioc > 0 )  THEN
2450          vc = 2
2451          str = ( ioc-1 ) * nbins + 1
2452          endi = ioc * nbins
2453          cc = 1
2454          DO ss = str, endi
2455             aerosol_mass(ss)%conc(k,j,i) = aerosol_mass(ss)%conc(k,j,i)       &
2456                               + ( aero(cc)%volc(vc) - aero_old(cc)%volc(vc) ) &
2457                               * arhooc * flag
2458             cc = cc+1
2459          ENDDO
2460       ENDIF
2461       
2462       IF ( ibc > 0 )  THEN
2463          vc = 3
2464          str = ( ibc-1 ) * nbins + 1 + fn1a
2465          endi = ibc * nbins
2466          cc = 1 + fn1a
2467          DO ss = str, endi
2468             aerosol_mass(ss)%conc(k,j,i) = aerosol_mass(ss)%conc(k,j,i)       &
2469                               + ( aero(cc)%volc(vc) - aero_old(cc)%volc(vc) ) &
2470                               * arhobc * flag
2471             cc = cc+1
2472          ENDDO
2473       ENDIF
2474       
2475       IF ( idu > 0 )  THEN
2476          vc = 4
2477          str = ( idu-1 ) * nbins + 1 + fn1a
2478          endi = idu * nbins
2479          cc = 1 + fn1a
2480          DO ss = str, endi
2481             aerosol_mass(ss)%conc(k,j,i) = aerosol_mass(ss)%conc(k,j,i)       &
2482                               + ( aero(cc)%volc(vc) - aero_old(cc)%volc(vc) ) &
2483                               * arhodu * flag
2484             cc = cc+1
2485          ENDDO
2486       ENDIF
2487       
2488       IF ( iss > 0 )  THEN
2489          vc = 5
2490          str = ( iss-1 ) * nbins + 1 + fn1a
2491          endi = iss * nbins
2492          cc = 1 + fn1a
2493          DO ss = str, endi
2494             aerosol_mass(ss)%conc(k,j,i) = aerosol_mass(ss)%conc(k,j,i)       &
2495                               + ( aero(cc)%volc(vc) - aero_old(cc)%volc(vc) ) &
2496                               * arhoss * flag
2497             cc = cc+1
2498          ENDDO
2499       ENDIF
2500       
2501       IF ( ino > 0 )  THEN
2502          vc = 6
2503          str = ( ino-1 ) * nbins + 1
2504          endi = ino * nbins
2505          cc = 1
2506          DO ss = str, endi
2507             aerosol_mass(ss)%conc(k,j,i) = aerosol_mass(ss)%conc(k,j,i)       &
2508                               + ( aero(cc)%volc(vc) - aero_old(cc)%volc(vc) ) &
2509                               * arhohno3 * flag
2510             cc = cc+1
2511          ENDDO
2512       ENDIF
2513       
2514       IF ( inh > 0 )  THEN
2515          vc = 7
2516          str = ( ino-1 ) * nbins + 1
2517          endi = ino * nbins
2518          cc = 1
2519          DO ss = str, endi
2520             aerosol_mass(ss)%conc(k,j,i) = aerosol_mass(ss)%conc(k,j,i)       &
2521                               + ( aero(cc)%volc(vc) - aero_old(cc)%volc(vc) ) &
2522                               * arhonh3 * flag
2523             cc = cc+1
2524          ENDDO
2525       ENDIF
2526       
2527       IF ( advect_particle_water )  THEN
2528          nc_h2o = get_index( prtcl,'H2O' )
2529          vc = 8
2530          str = ( nc_h2o-1 ) * nbins + 1
2531          endi = nc_h2o * nbins
2532          cc = 1
2533          DO ss = str, endi
2534             aerosol_mass(ss)%conc(k,j,i) = aerosol_mass(ss)%conc(k,j,i)       &
2535                               + ( aero(cc)%volc(vc) - aero_old(cc)%volc(vc) ) &
2536                               * arhoh2o * flag
2537             IF ( prunmode == 1 )  THEN
2538                aerosol_mass(ss)%init(k) = MAX( aerosol_mass(ss)%init(k),      &
2539                                               aerosol_mass(ss)%conc(k,j,i) )
2540             ENDIF
2541             cc = cc+1                             
2542          ENDDO
2543       ENDIF
2544
2545!--    Condensation of precursor gases
2546       IF ( lscndgas )  THEN
2547          IF ( salsa_gases_from_chem )  THEN         
2548!         
2549!--          SO4 (or H2SO4)
2550             chem_species( gas_index_chem(1) )%conc(k,j,i) =                &
2551                            chem_species( gas_index_chem(1) )%conc(k,j,i) + &
2552                                                  ( zgso4 / ppm_to_nconc(k) - &
2553                       chem_species( gas_index_chem(1) )%conc(k,j,i) ) * flag
2554!                           
2555!--          HNO3
2556             chem_species( gas_index_chem(2) )%conc(k,j,i) =                &
2557                            chem_species( gas_index_chem(2) )%conc(k,j,i) + &
2558                                                 ( zghno3 / ppm_to_nconc(k) - &
2559                       chem_species( gas_index_chem(2) )%conc(k,j,i) ) * flag
2560!                           
2561!--          NH3
2562             chem_species( gas_index_chem(3) )%conc(k,j,i) =                &
2563                            chem_species( gas_index_chem(3) )%conc(k,j,i) + &
2564                                                  ( zgnh3 / ppm_to_nconc(k) - &
2565                       chem_species( gas_index_chem(3) )%conc(k,j,i) ) * flag
2566!                           
2567!--          non-volatile OC
2568             chem_species( gas_index_chem(4) )%conc(k,j,i) =                &
2569                            chem_species( gas_index_chem(4) )%conc(k,j,i) + &
2570                                                 ( zgocnv / ppm_to_nconc(k) - &
2571                       chem_species( gas_index_chem(4) )%conc(k,j,i) ) * flag
2572!                           
2573!--          semi-volatile OC
2574             chem_species( gas_index_chem(5) )%conc(k,j,i) =                &
2575                            chem_species( gas_index_chem(5) )%conc(k,j,i) + &
2576                                                 ( zgocsv / ppm_to_nconc(k) - &
2577                       chem_species( gas_index_chem(5) )%conc(k,j,i) ) * flag                 
2578         
2579          ELSE
2580!         
2581!--          SO4 (or H2SO4)
2582             salsa_gas(1)%conc(k,j,i) = salsa_gas(1)%conc(k,j,i) + ( zgso4 -   &
2583                                          salsa_gas(1)%conc(k,j,i) ) * flag
2584!                           
2585!--          HNO3
2586             salsa_gas(2)%conc(k,j,i) = salsa_gas(2)%conc(k,j,i) + ( zghno3 -  &
2587                                          salsa_gas(2)%conc(k,j,i) ) * flag
2588!                           
2589!--          NH3
2590             salsa_gas(3)%conc(k,j,i) = salsa_gas(3)%conc(k,j,i) + ( zgnh3 -   &
2591                                          salsa_gas(3)%conc(k,j,i) ) * flag
2592!                           
2593!--          non-volatile OC
2594             salsa_gas(4)%conc(k,j,i) = salsa_gas(4)%conc(k,j,i) + ( zgocnv -  &
2595                                          salsa_gas(4)%conc(k,j,i) ) * flag
2596!                           
2597!--          semi-volatile OC
2598             salsa_gas(5)%conc(k,j,i) = salsa_gas(5)%conc(k,j,i) + ( zgocsv -  &
2599                                          salsa_gas(5)%conc(k,j,i) ) * flag
2600          ENDIF
2601       ENDIF
2602!               
2603!--    Tendency of water vapour mixing ratio is obtained from the
2604!--    change in RH during SALSA run. This releases heat and changes pt.
2605!--    Assumes no temperature change during SALSA run.
2606!--    q = r / (1+r), Euler method for integration
2607!
2608       IF ( feedback_to_palm )  THEN
2609          q_p(k,j,i) = q_p(k,j,i) + 1.0_wp / ( in_cw(k) * in_adn(k) + 1.0_wp ) &
2610                       ** 2.0_wp * ( in_cw(k) - cw_old ) * in_adn(k) 
2611          pt_p(k,j,i) = pt_p(k,j,i) + alv / c_p * ( in_cw(k) - cw_old ) *      &
2612                        in_adn(k) / ( in_cw(k) / in_adn(k) + 1.0_wp ) ** 2.0_wp&
2613                        * pt_p(k,j,i) / in_t(k)
2614       ENDIF
2615                         
2616    ENDDO   ! k
2617!   
2618!-- Set surfaces and wall fluxes due to deposition 
2619    IF ( lsdepo_topo  .AND.  prunmode == 3 )  THEN
2620       IF ( .NOT. land_surface  .AND.  .NOT. urban_surface )  THEN
2621          CALL depo_topo( i, j, surf_def_h(0), vd, Sc, kvis, in_u, rho_air_zw )
2622          DO  l = 0, 3
2623             CALL depo_topo( i, j, surf_def_v(l), vd, Sc, kvis, in_u,          &
2624                             rho_air_zw**0.0_wp )
2625          ENDDO
2626       ELSE
2627          CALL depo_topo( i, j, surf_usm_h, vd, Sc, kvis, in_u, rho_air_zw )
2628          DO  l = 0, 3
2629             CALL depo_topo( i, j, surf_usm_v(l), vd, Sc, kvis, in_u,          &
2630                             rho_air_zw**0.0_wp )
2631          ENDDO
2632          CALL depo_topo( i, j, surf_lsm_h, vd, Sc, kvis, in_u, rho_air_zw )
2633          DO  l = 0, 3
2634             CALL depo_topo( i, j, surf_lsm_v(l), vd, Sc, kvis, in_u,          &
2635                             rho_air_zw**0.0_wp )
2636          ENDDO
2637       ENDIF
2638    ENDIF
2639   
2640 END SUBROUTINE salsa_driver
2641
2642!------------------------------------------------------------------------------!
2643! Description:
2644! ------------
2645!> The SALSA subroutine
2646!> Modified for the new aerosol datatype,
2647!> Juha Tonttila, FMI, 2014.
2648!> Only aerosol processes included, Mona Kurppa, UHel, 2017
2649!------------------------------------------------------------------------------!   
2650 SUBROUTINE run_salsa( ppres, pcw, pcs, ptemp, mag_u, adn, lad, pc_h2so4,      &
2651                       pc_ocnv, pc_ocsv, pc_hno3, pc_nh3, paero, prtcl, kvis,  &
2652                       Sc, vc, ptstep )
2653
2654    IMPLICIT NONE
2655!
2656!-- Input parameters and variables
2657    REAL(wp), INTENT(in) ::  adn    !< air density (kg/m3)
2658    REAL(wp), INTENT(in) ::  lad    !< leaf area density (m2/m3)
2659    REAL(wp), INTENT(in) ::  mag_u  !< magnitude of wind (m/s)
2660    REAL(wp), INTENT(in) ::  ppres  !< atmospheric pressure at each grid
2661                                    !< point (Pa)
2662    REAL(wp), INTENT(in) ::  ptemp  !< temperature at each grid point (K)
2663    REAL(wp), INTENT(in) ::  ptstep !< time step of salsa processes (s)
2664    TYPE(component_index), INTENT(in) :: prtcl  !< part. component index table
2665!       
2666!-- Input variables that are changed within:
2667    REAL(wp), INTENT(inout) ::  kvis     !< kinematic viscosity of air (m2/s)
2668    REAL(wp), INTENT(inout) ::  Sc(:)    !< particle Schmidt number
2669    REAL(wp), INTENT(inout) ::  vc(:)    !< particle fall speed (m/s,
2670                                         !< sedimentation velocity)
2671!-- Gas phase concentrations at each grid point (#/m3)
2672    REAL(wp), INTENT(inout) ::  pc_h2so4 !< sulphuric acid
2673    REAL(wp), INTENT(inout) ::  pc_hno3  !< nitric acid
2674    REAL(wp), INTENT(inout) ::  pc_nh3   !< ammonia
2675    REAL(wp), INTENT(inout) ::  pc_ocnv  !< nonvolatile OC
2676    REAL(wp), INTENT(inout) ::  pc_ocsv  !< semivolatile OC
2677    REAL(wp), INTENT(inout) ::  pcs      !< Saturation concentration of water
2678                                         !< vapour (kg/m3)
2679    REAL(wp), INTENT(inout) ::  pcw      !< Water vapour concentration (kg/m3)                                                   
2680    TYPE(t_section), INTENT(inout) ::  paero(fn2b) 
2681!
2682!-- Coagulation
2683    IF ( lscoag )   THEN
2684       CALL coagulation( paero, ptstep, ptemp, ppres )
2685    ENDIF
2686!
2687!-- Condensation
2688    IF ( lscnd )   THEN
2689       CALL condensation( paero, pc_h2so4, pc_ocnv, pc_ocsv,  pc_hno3, pc_nh3, &
2690                          pcw, pcs, ptemp, ppres, ptstep, prtcl )
2691    ENDIF   
2692!   
2693!-- Deposition
2694    IF ( lsdepo )  THEN
2695       CALL deposition( paero, ptemp, adn, mag_u, lad, kvis, Sc, vc ) 
2696    ENDIF       
2697!
2698!-- Size distribution bin update
2699!-- Mona: why done 3 times in SALSA-standalone?
2700    IF ( lsdistupdate )   THEN
2701       CALL distr_update( paero )
2702    ENDIF
2703   
2704  END SUBROUTINE run_salsa
2705 
2706!------------------------------------------------------------------------------!
2707! Description:
2708! ------------
2709!> Set logical switches according to the host model state and user-specified
2710!> NAMELIST options.
2711!> Juha Tonttila, FMI, 2014
2712!> Only aerosol processes included, Mona Kurppa, UHel, 2017
2713!------------------------------------------------------------------------------!
2714 SUBROUTINE set_salsa_runtime( prunmode )
2715 
2716    IMPLICIT NONE
2717   
2718    INTEGER(iwp), INTENT(in) ::  prunmode
2719   
2720    SELECT CASE(prunmode)
2721
2722       CASE(1) !< Initialization
2723          lscoag       = .FALSE.
2724          lscnd        = .FALSE.
2725          lscndgas     = .FALSE.
2726          lscndh2oae   = .FALSE.
2727          lsdepo       = .FALSE.
2728          lsdepo_vege  = .FALSE.
2729          lsdepo_topo  = .FALSE.
2730          lsdistupdate = .TRUE.
2731
2732       CASE(2)  !< Spinup period
2733          lscoag      = ( .FALSE. .AND. nlcoag   )
2734          lscnd       = ( .TRUE.  .AND. nlcnd    )
2735          lscndgas    = ( .TRUE.  .AND. nlcndgas )
2736          lscndh2oae  = ( .TRUE.  .AND. nlcndh2oae )
2737
2738       CASE(3)  !< Run
2739          lscoag       = nlcoag
2740          lscnd        = nlcnd
2741          lscndgas     = nlcndgas
2742          lscndh2oae   = nlcndh2oae
2743          lsdepo       = nldepo
2744          lsdepo_vege  = nldepo_vege
2745          lsdepo_topo  = nldepo_topo
2746          lsdistupdate = nldistupdate
2747
2748    END SELECT
2749
2750
2751 END SUBROUTINE set_salsa_runtime
2752 
2753!------------------------------------------------------------------------------!
2754! Description:
2755! ------------
2756!> Calculates the absolute temperature (using hydrostatic pressure), saturation
2757!> vapour pressure and mixing ratio over water, relative humidity and air
2758!> density needed in the SALSA model.
2759!> NOTE, no saturation adjustment takes place -> the resulting water vapour
2760!> mixing ratio can be supersaturated, allowing the microphysical calculations
2761!> in SALSA.
2762!
2763!> Juha Tonttila, FMI, 2014 (original SALSAthrm)
2764!> Mona Kurppa, UHel, 2017 (adjustment for PALM and only aerosol processes)
2765!------------------------------------------------------------------------------!
2766 SUBROUTINE salsa_thrm_ij( i, j, p_ij, temp_ij, cw_ij, cs_ij, adn_ij )
2767 
2768    USE arrays_3d,                                                             &
2769        ONLY: p, pt, q, zu
2770       
2771    USE basic_constants_and_equations_mod,                                     &
2772        ONLY:  barometric_formula, exner_function, ideal_gas_law_rho, magnus
2773       
2774    USE control_parameters,                                                    &
2775        ONLY: pt_surface, surface_pressure
2776       
2777    IMPLICIT NONE
2778   
2779    INTEGER(iwp), INTENT(in) ::  i
2780    INTEGER(iwp), INTENT(in) ::  j
2781    REAL(wp), DIMENSION(:), INTENT(inout) ::  adn_ij
2782    REAL(wp), DIMENSION(:), INTENT(inout) ::  p_ij       
2783    REAL(wp), DIMENSION(:), INTENT(inout) ::  temp_ij
2784    REAL(wp), DIMENSION(:), INTENT(inout), OPTIONAL ::  cw_ij
2785    REAL(wp), DIMENSION(:), INTENT(inout), OPTIONAL ::  cs_ij
2786    REAL(wp), DIMENSION(nzb:nzt+1) ::  e_s !< saturation vapour pressure
2787                                           !< over water (Pa)
2788    REAL(wp) ::  t_surface !< absolute surface temperature (K)
2789!
2790!-- Pressure p_ijk (Pa) = hydrostatic pressure + perturbation pressure (p)
2791    t_surface = pt_surface * exner_function( surface_pressure * 100.0_wp )
2792    p_ij(:) = barometric_formula( zu, t_surface, surface_pressure * 100.0_wp ) &
2793              + p(:,j,i)
2794!             
2795!-- Absolute ambient temperature (K)
2796    temp_ij(:) = pt(:,j,i) * exner_function( p_ij(:) )       
2797!
2798!-- Air density
2799    adn_ij(:) = ideal_gas_law_rho( p_ij(:), temp_ij(:) )
2800!
2801!-- Water vapour concentration r_v (kg/m3)
2802    IF ( PRESENT( cw_ij ) )  THEN
2803       cw_ij(:) = ( q(:,j,i) / ( 1.0_wp - q(:,j,i) ) ) * adn_ij(:) 
2804    ENDIF
2805!
2806!-- Saturation mixing ratio r_s (kg/kg) from vapour pressure at temp (Pa)
2807    IF ( PRESENT( cs_ij ) )  THEN
2808       e_s(:) = 611.0_wp * EXP( alv_d_rv * ( 3.6609E-3_wp - 1.0_wp /           &
2809                temp_ij(:) ) )! magnus( temp_ij(:) )
2810       cs_ij(:) = ( 0.622_wp * e_s / ( p_ij(:) - e_s(:) ) ) * adn_ij(:) 
2811    ENDIF
2812
2813 END SUBROUTINE salsa_thrm_ij
2814
2815!------------------------------------------------------------------------------!
2816! Description:
2817! ------------
2818!> Calculates ambient sizes of particles by equilibrating soluble fraction of
2819!> particles with water using the ZSR method (Stokes and Robinson, 1966).
2820!> Method:
2821!> Following chemical components are assumed water-soluble
2822!> - (ammonium) sulphate (100%)
2823!> - sea salt (100 %)
2824!> - organic carbon (epsoc * 100%)
2825!> Exact thermodynamic considerations neglected.
2826!> - If particles contain no sea salt, calculation according to sulphate
2827!>   properties
2828!> - If contain sea salt but no sulphate, calculation according to sea salt
2829!>   properties
2830!> - If contain both sulphate and sea salt -> the molar fraction of these
2831!>   compounds determines which one of them is used as the basis of calculation.
2832!> If sulphate and sea salt coexist in a particle, it is assumed that the Cl is
2833!> replaced by sulphate; thus only either sulphate + organics or sea salt +
2834!> organics is included in the calculation of soluble fraction.
2835!> Molality parameterizations taken from Table 1 of Tang: Thermodynamic and
2836!> optical properties of mixed-salt aerosols of atmospheric importance,
2837!> J. Geophys. Res., 102 (D2), 1883-1893 (1997)
2838!
2839!> Coded by:
2840!> Hannele Korhonen (FMI) 2005
2841!> Harri Kokkola (FMI) 2006
2842!> Matti Niskanen(FMI) 2012
2843!> Anton Laakso  (FMI) 2013
2844!> Modified for the new aerosol datatype, Juha Tonttila (FMI) 2014
2845!
2846!> fxm: should sea salt form a solid particle when prh is very low (even though
2847!> it could be mixed with e.g. sulphate)?
2848!> fxm: crashes if no sulphate or sea salt
2849!> fxm: do we really need to consider Kelvin effect for subrange 2
2850!------------------------------------------------------------------------------!     
2851 SUBROUTINE equilibration( prh, ptemp, paero, init )
2852     
2853    IMPLICIT NONE
2854!
2855!-- Input variables
2856    LOGICAL, INTENT(in) ::  init   !< TRUE: Initialization call
2857                                   !< FALSE: Normal runtime: update water
2858                                   !<        content only for 1a
2859    REAL(wp), INTENT(in) ::  prh   !< relative humidity [0-1]
2860    REAL(wp), INTENT(in) ::  ptemp !< temperature (K)
2861!
2862!-- Output variables
2863    TYPE(t_section), INTENT(inout) ::  paero(fn2b)     
2864!
2865!-- Local
2866    INTEGER(iwp) :: b      !< loop index
2867    INTEGER(iwp) :: counti  !< loop index
2868    REAL(wp) ::  zaw        !< water activity [0-1]       
2869    REAL(wp) ::  zbinmol(7) !< binary molality of each components (mol/kg)
2870    REAL(wp) ::  zcore      !< Volume of dry particle   
2871    REAL(wp) ::  zdold      !< Old diameter
2872    REAL(wp) ::  zdwet      !< Wet diameter or mean droplet diameter
2873    REAL(wp) ::  zke        !< Kelvin term in the Köhler equation
2874    REAL(wp) ::  zlwc       !< liquid water content [kg/m3-air]
2875    REAL(wp) ::  zrh        !< Relative humidity
2876    REAL(wp) ::  zvpart(7)  !< volume of chem. compounds in one particle
2877   
2878    zaw       = 0.0_wp
2879    zbinmol   = 0.0_wp
2880    zcore     = 0.0_wp
2881    zdold     = 0.0_wp
2882    zdwet     = 0.0_wp
2883    zlwc      = 0.0_wp
2884    zrh       = 0.0_wp
2885   
2886!               
2887!-- Relative humidity:
2888    zrh = prh
2889    zrh = MAX( zrh, 0.05_wp )
2890    zrh = MIN( zrh, 0.98_wp)   
2891!
2892!-- 1) Regime 1: sulphate and partly water-soluble OC. Done for every CALL
2893    DO  b = in1a, fn1a   ! size bin
2894         
2895       zbinmol = 0.0_wp
2896       zdold   = 1.0_wp 
2897       zke     = 1.02_wp
2898       
2899       IF ( paero(b)%numc > nclim )  THEN
2900!
2901!--       Volume in one particle
2902          zvpart = 0.0_wp
2903          zvpart(1:2) = paero(b)%volc(1:2) / paero(b)%numc
2904          zvpart(6:7) = paero(b)%volc(6:7) / paero(b)%numc
2905!               
2906!--       Total volume and wet diameter of one dry particle
2907          zcore = SUM( zvpart(1:2) )
2908          zdwet = paero(b)%dwet
2909         
2910          counti = 0
2911          DO  WHILE ( ABS( zdwet / zdold - 1.0_wp ) > 1.0E-2_wp ) 
2912         
2913             zdold = MAX( zdwet, 1.0E-20_wp )
2914             zaw = MAX( 1.0E-3_wp, zrh / zke ) ! To avoid underflow
2915!                   
2916!--          Binary molalities (mol/kg):
2917!--          Sulphate
2918             zbinmol(1) = 1.1065495E+2_wp - 3.6759197E+2_wp * zaw              &
2919                                          + 5.0462934E+2_wp * zaw**2.0_wp      &
2920                                          - 3.1543839E+2_wp * zaw**3.0_wp      &
2921                                          + 6.770824E+1_wp  * zaw**4.0_wp 
2922!--          Organic carbon                     
2923             zbinmol(2) = 1.0_wp / ( zaw * amh2o ) - 1.0_wp / amh2o
2924!--          Nitric acid                             
2925             zbinmol(6) = 2.306844303E+1_wp - 3.563608869E+1_wp * zaw          &
2926                                            - 6.210577919E+1_wp * zaw**2.0_wp  &
2927                                            + 5.510176187E+2_wp * zaw**3.0_wp  &
2928                                            - 1.460055286E+3_wp * zaw**4.0_wp  &
2929                                            + 1.894467542E+3_wp * zaw**5.0_wp  &
2930                                            - 1.220611402E+3_wp * zaw**6.0_wp  &
2931                                            + 3.098597737E+2_wp * zaw**7.0_wp 
2932!
2933!--          Calculate the liquid water content (kg/m3-air) using ZSR (see e.g.
2934!--          Eq. 10.98 in Seinfeld and Pandis (2006))
2935             zlwc = ( paero(b)%volc(1) * ( arhoh2so4 / amh2so4 ) ) /           &
2936                    zbinmol(1) + epsoc * paero(b)%volc(2) * ( arhooc / amoc )  &
2937                    / zbinmol(2) + ( paero(b)%volc(6) * ( arhohno3/amhno3 ) )  &
2938                    / zbinmol(6)
2939!                           
2940!--          Particle wet diameter (m)
2941             zdwet = ( zlwc / paero(b)%numc / arhoh2o / api6 +                 &
2942                     ( SUM( zvpart(6:7) ) / api6 ) +      &
2943                       zcore / api6 )**( 1.0_wp / 3.0_wp )
2944!                             
2945!--          Kelvin effect (Eq. 10.85 in in Seinfeld and Pandis (2006)). Avoid
2946!--          overflow.
2947             zke = EXP( MIN( 50.0_wp,                                          &
2948                       4.0_wp * surfw0 * amvh2so4 / ( abo * ptemp *  zdwet ) ) )
2949             
2950             counti = counti + 1
2951             IF ( counti > 1000 )  THEN
2952                message_string = 'Subrange 1: no convergence!'
2953                CALL message( 'salsa_mod: equilibration', 'SA0042',            &
2954                              1, 2, 0, 6, 0 )
2955             ENDIF
2956          ENDDO
2957!               
2958!--       Instead of lwc, use the volume concentration of water from now on
2959!--       (easy to convert...)
2960          paero(b)%volc(8) = zlwc / arhoh2o
2961!               
2962!--       If this is initialization, update the core and wet diameter
2963          IF ( init )  THEN
2964             paero(b)%dwet = zdwet
2965             paero(b)%core = zcore
2966          ENDIF
2967         
2968       ELSE
2969!--       If initialization
2970!--       1.2) empty bins given bin average values 
2971          IF ( init )  THEN
2972             paero(b)%dwet = paero(b)%dmid
2973             paero(b)%core = api6 * paero(b)%dmid ** 3.0_wp
2974          ENDIF
2975         
2976       ENDIF
2977             
2978    ENDDO !< b
2979!
2980!-- 2) Regime 2a: sulphate, OC, BC and sea salt
2981!--    This is done only for initialization call, otherwise the water contents
2982!--    are computed via condensation
2983    IF ( init )  THEN
2984       DO  b = in2a, fn2b
2985             
2986!--       Initialize
2987          zke     = 1.02_wp
2988          zbinmol = 0.0_wp
2989          zdold   = 1.0_wp
2990!               
2991!--       1) Particle properties calculated for non-empty bins
2992          IF ( paero(b)%numc > nclim )  THEN
2993!               
2994!--          Volume in one particle [fxm]
2995             zvpart = 0.0_wp
2996             zvpart(1:7) = paero(b)%volc(1:7) / paero(b)%numc
2997!
2998!--          Total volume and wet diameter of one dry particle [fxm]
2999             zcore = SUM( zvpart(1:5) )
3000             zdwet = paero(b)%dwet
3001
3002             counti = 0
3003             DO  WHILE ( ABS( zdwet / zdold - 1.0_wp ) > 1.0E-12_wp )
3004             
3005                zdold = MAX( zdwet, 1.0E-20_wp )
3006                zaw = zrh / zke
3007!                     
3008!--             Binary molalities (mol/kg):
3009!--             Sulphate
3010                zbinmol(1) = 1.1065495E+2_wp - 3.6759197E+2_wp * zaw           & 
3011                        + 5.0462934E+2_wp * zaw**2 - 3.1543839E+2_wp * zaw**3  &
3012                        + 6.770824E+1_wp  * zaw**4 
3013!--             Organic carbon                       
3014                zbinmol(2) = 1.0_wp / ( zaw * amh2o ) - 1.0_wp / amh2o
3015!--             Nitric acid
3016                zbinmol(6) = 2.306844303E+1_wp - 3.563608869E+1_wp * zaw       &
3017                     - 6.210577919E+1_wp * zaw**2 + 5.510176187E+2_wp * zaw**3 &
3018                     - 1.460055286E+3_wp * zaw**4 + 1.894467542E+3_wp * zaw**5 &
3019                     - 1.220611402E+3_wp * zaw**6 + 3.098597737E+2_wp * zaw**7 
3020!--             Sea salt (natrium chloride)                                 
3021                zbinmol(5) = 5.875248E+1_wp - 1.8781997E+2_wp * zaw            &
3022                         + 2.7211377E+2_wp * zaw**2 - 1.8458287E+2_wp * zaw**3 &
3023                         + 4.153689E+1_wp  * zaw**4 
3024!                                 
3025!--             Calculate the liquid water content (kg/m3-air)
3026                zlwc = ( paero(b)%volc(1) * ( arhoh2so4 / amh2so4 ) ) /        &
3027                       zbinmol(1) + epsoc * ( paero(b)%volc(2) * ( arhooc /    &
3028                       amoc ) ) / zbinmol(2) + ( paero(b)%volc(6) * ( arhohno3 &
3029                       / amhno3 ) ) / zbinmol(6) + ( paero(b)%volc(5) *        &
3030                       ( arhoss / amss ) ) / zbinmol(5)
3031                       
3032!--             Particle wet radius (m)
3033                zdwet = ( zlwc / paero(b)%numc / arhoh2o / api6 +              &
3034                          ( SUM( zvpart(6:7) ) / api6 )  + &
3035                           zcore / api6 ) ** ( 1.0_wp / 3.0_wp )
3036!                               
3037!--             Kelvin effect (Eq. 10.85 in Seinfeld and Pandis (2006))
3038                zke = EXP( MIN( 50.0_wp,                                       &
3039                        4.0_wp * surfw0 * amvh2so4 / ( abo * zdwet * ptemp ) ) )
3040                         
3041                counti = counti + 1
3042                IF ( counti > 1000 )  THEN
3043                   message_string = 'Subrange 2: no convergence!'
3044                CALL message( 'salsa_mod: equilibration', 'SA0043',            &
3045                              1, 2, 0, 6, 0 )
3046                ENDIF
3047             ENDDO
3048!                   
3049!--          Liquid water content; instead of LWC use the volume concentration
3050             paero(b)%volc(8) = zlwc / arhoh2o
3051             paero(b)%dwet    = zdwet
3052             paero(b)%core    = zcore
3053             
3054          ELSE
3055!--          2.2) empty bins given bin average values
3056             paero(b)%dwet = paero(b)%dmid
3057             paero(b)%core = api6 * paero(b)%dmid ** 3.0_wp
3058          ENDIF
3059               
3060       ENDDO   ! b
3061    ENDIF
3062
3063 END SUBROUTINE equilibration
3064 
3065!------------------------------------------------------------------------------!
3066!> Description:
3067!> ------------
3068!> Calculation of the settling velocity vc (m/s) per aerosol size bin and
3069!> deposition on plant canopy (lsdepo_vege).
3070!
3071!> Deposition is based on either the scheme presented in:
3072!> Zhang et al. (2001), Atmos. Environ. 35, 549-560 (includes collection due to
3073!> Brownian diffusion, impaction, interception and sedimentation)
3074!> OR
3075!> Petroff & Zhang (2010), Geosci. Model Dev. 3, 753-769 (includes also
3076!> collection due to turbulent impaction)
3077!
3078!> Equation numbers refer to equation in Jacobson (2005): Fundamentals of
3079!> Atmospheric Modeling, 2nd Edition.
3080!
3081!> Subroutine follows closely sedim_SALSA in UCLALES-SALSA written by Juha
3082!> Tonttila (KIT/FMI) and Zubair Maalick (UEF).
3083!> Rewritten to PALM by Mona Kurppa (UH), 2017.
3084!
3085!> Call for grid point i,j,k
3086!------------------------------------------------------------------------------!
3087
3088 SUBROUTINE deposition( paero, tk, adn, mag_u, lad, kvis, Sc, vc )
3089 
3090    USE plant_canopy_model_mod,                                                &
3091        ONLY: cdc
3092 
3093    IMPLICIT NONE
3094   
3095    REAL(wp), INTENT(in)    ::  adn    !< air density (kg/m3) 
3096    REAL(wp), INTENT(out)   ::  kvis   !< kinematic viscosity of air (m2/s)
3097    REAL(wp), INTENT(in) ::     lad    !< leaf area density (m2/m3)
3098    REAL(wp), INTENT(in)    ::  mag_u  !< wind velocity (m/s)
3099    REAL(wp), INTENT(out)   ::  Sc(:)  !< particle Schmidt number 
3100    REAL(wp), INTENT(in)    ::  tk     !< abs.temperature (K)   
3101    REAL(wp), INTENT(out)   ::  vc(:)  !< critical fall speed i.e. settling
3102                                       !< velocity of an aerosol particle (m/s)
3103    TYPE(t_section), INTENT(inout) ::  paero(fn2b)       
3104   
3105    INTEGER(iwp) ::  b      !< loop index
3106    INTEGER(iwp) ::  c      !< loop index
3107    REAL(wp) ::  avis       !< molecular viscocity of air (kg/(m*s))
3108    REAL(wp), PARAMETER ::  c_A = 1.249_wp !< Constants A, B and C for
3109    REAL(wp), PARAMETER ::  c_B = 0.42_wp  !< calculating  the Cunningham 
3110    REAL(wp), PARAMETER ::  c_C = 0.87_wp  !< slip-flow correction (Cc) 
3111                                           !< according to Jacobson (2005),
3112                                           !< Eq. 15.30
3113    REAL(wp) ::  Cc         !< Cunningham slip-flow correction factor     
3114    REAL(wp) ::  Kn         !< Knudsen number   
3115    REAL(wp) ::  lambda     !< molecular mean free path (m)
3116    REAL(wp) ::  mdiff      !< particle diffusivity coefficient   
3117    REAL(wp) ::  pdn        !< particle density (kg/m3)     
3118    REAL(wp) ::  ustar      !< friction velocity (m/s)   
3119    REAL(wp) ::  va         !< thermal speed of an air molecule (m/s)
3120    REAL(wp) ::  zdwet      !< wet diameter (m)                             
3121!
3122!-- Initialise
3123    Cc            = 0.0_wp
3124    Kn            = 0.0_wp
3125    mdiff         = 0.0_wp
3126    pdn           = 1500.0_wp    ! default value
3127    ustar         = 0.0_wp 
3128!
3129!-- Molecular viscosity of air (Eq. 4.54)
3130    avis = 1.8325E-5_wp * ( 416.16_wp / ( tk + 120.0_wp ) ) * ( tk /           &
3131           296.16_wp )**1.5_wp
3132!             
3133!-- Kinematic viscosity (Eq. 4.55)
3134    kvis =  avis / adn
3135!       
3136!-- Thermal velocity of an air molecule (Eq. 15.32)
3137    va = SQRT( 8.0_wp * abo * tk / ( pi * am_airmol ) ) 
3138!
3139!-- Mean free path (m) (Eq. 15.24)
3140    lambda = 2.0_wp * avis / ( adn * va )
3141   
3142    DO  b = 1, nbins
3143   
3144       IF ( paero(b)%numc < nclim )  CYCLE
3145       zdwet = paero(b)%dwet
3146!
3147!--    Knudsen number (Eq. 15.23)
3148       Kn = MAX( 1.0E-2_wp, lambda / ( zdwet * 0.5_wp ) ) ! To avoid underflow
3149!
3150!--    Cunningham slip-flow correction (Eq. 15.30)
3151       Cc = 1.0_wp + Kn * ( c_A + c_B * EXP( -c_C / Kn ) )
3152
3153!--    Particle diffusivity coefficient (Eq. 15.29)
3154       mdiff = ( abo * tk * Cc ) / ( 3.0_wp * pi * avis * zdwet )
3155!       
3156!--    Particle Schmidt number (Eq. 15.36)
3157       Sc(b) = kvis / mdiff       
3158!       
3159!--    Critical fall speed i.e. settling velocity  (Eq. 20.4)                 
3160       vc(b) = MIN( 1.0_wp, terminal_vel( 0.5_wp * zdwet, pdn, adn, avis, Cc) )
3161       
3162       IF ( lsdepo_vege  .AND.  plant_canopy  .AND.  lad > 0.0_wp )  THEN
3163!       
3164!--       Friction velocity calculated following Prandtl (1925):
3165          ustar = SQRT( cdc ) * mag_u
3166          CALL depo_vege( paero, b, vc(b), mag_u, ustar, kvis, Sc(b), lad )
3167       ENDIF
3168    ENDDO
3169 
3170 END SUBROUTINE deposition
3171 
3172!------------------------------------------------------------------------------!
3173! Description:
3174! ------------
3175!> Calculate change in number and volume concentrations due to deposition on
3176!> plant canopy.
3177!------------------------------------------------------------------------------!
3178 SUBROUTINE depo_vege( paero, b, vc, mag_u, ustar, kvis_a, Sc, lad )
3179 
3180    IMPLICIT NONE
3181   
3182    INTEGER(iwp), INTENT(in) ::  b  !< loop index
3183    REAL(wp), INTENT(in) ::  kvis_a !< kinematic viscosity of air (m2/s)
3184    REAL(wp), INTENT(in) ::  lad    !< leaf area density (m2/m3)
3185    REAL(wp), INTENT(in) ::  mag_u  !< wind velocity (m/s)   
3186    REAL(wp), INTENT(in) ::  Sc     !< particle Schmidt number
3187    REAL(wp), INTENT(in) ::  ustar  !< friction velocity (m/s)                                   
3188    REAL(wp), INTENT(in) ::  vc     !< terminal velocity (m/s) 
3189    TYPE(t_section), INTENT(inout) ::  paero(fn2b) 
3190   
3191    INTEGER(iwp) ::  c      !< loop index
3192    REAL(wp), PARAMETER ::  c_A = 1.249_wp !< Constants A, B and C for
3193    REAL(wp), PARAMETER ::  c_B = 0.42_wp  !< calculating  the Cunningham 
3194    REAL(wp), PARAMETER ::  c_C = 0.87_wp  !< slip-flow correction (Cc) 
3195                                           !< according to Jacobson (2005),
3196                                           !< Eq. 15.30
3197    REAL(wp) ::  alpha       !< parameter, Table 3 in Zhang et al. (2001) 
3198    REAL(wp) ::  depo        !< deposition efficiency
3199    REAL(wp) ::  C_Br        !< coefficient for Brownian diffusion
3200    REAL(wp) ::  C_IM        !< coefficient for inertial impaction
3201    REAL(wp) ::  C_IN        !< coefficient for interception
3202    REAL(wp) ::  C_IT        !< coefficient for turbulent impaction   
3203    REAL(wp) ::  gamma       !< parameter, Table 3 in Zhang et al. (2001)   
3204    REAL(wp) ::  par_A       !< parameter A for the characteristic radius of
3205                             !< collectors, Table 3 in Zhang et al. (2001)   
3206    REAL(wp) ::  rt          !< the overall quasi-laminar resistance for
3207                             !< particles
3208    REAL(wp) ::  St          !< Stokes number for smooth surfaces or bluff
3209                             !< surface elements                                 
3210    REAL(wp) ::  tau_plus    !< dimensionless particle relaxation time   
3211    REAL(wp) ::  v_bd        !< deposition velocity due to Brownian diffusion
3212    REAL(wp) ::  v_im        !< deposition velocity due to impaction
3213    REAL(wp) ::  v_in        !< deposition velocity due to interception
3214    REAL(wp) ::  v_it        !< deposition velocity due to turbulent impaction                               
3215!
3216!-- Initialise
3217    depo     = 0.0_wp 
3218    rt       = 0.0_wp
3219    St       = 0.0_wp
3220    tau_plus = 0.0_wp
3221    v_bd     = 0.0_wp     
3222    v_im     = 0.0_wp       
3223    v_in     = 0.0_wp       
3224    v_it     = 0.0_wp         
3225       
3226    IF ( depo_vege_type == 'zhang2001' )  THEN
3227!       
3228!--    Parameters for the land use category 'deciduous broadleaf trees'(Table 3)     
3229       par_A = 5.0E-3_wp
3230       alpha = 0.8_wp
3231       gamma = 0.56_wp 
3232!       
3233!--    Stokes number for vegetated surfaces (Seinfeld & Pandis (2006): Eq.19.24) 
3234       St = vc * ustar / ( g * par_A )         
3235!         
3236!--    The overall quasi-laminar resistance for particles (Zhang et al., Eq. 5)       
3237       rt = MAX( EPSILON( 1.0_wp ), ( 3.0_wp * ustar * EXP( -St**0.5_wp ) *    &
3238                         ( Sc**( -gamma ) + ( St / ( alpha + St ) )**2.0_wp +  &
3239                           0.5_wp * ( paero(b)%dwet / par_A )**2.0_wp ) ) )
3240       depo = ( rt + vc ) * lad
3241       paero(b)%numc = paero(b)%numc - depo * paero(b)%numc * dt_salsa
3242       DO  c = 1, maxspec+1
3243          paero(b)%volc(c) = paero(b)%volc(c) - depo * paero(b)%volc(c) *      &
3244                             dt_salsa
3245       ENDDO
3246       
3247    ELSEIF ( depo_vege_type == 'petroff2010' )  THEN
3248!
3249!--    vd = v_BD + v_IN + v_IM + v_IT + vc
3250!--    Deposition efficiencies from Table 1. Constants from Table 2.
3251       C_Br  = 1.262_wp
3252       C_IM  = 0.130_wp
3253       C_IN  = 0.216_wp
3254       C_IT  = 0.056_wp
3255       par_A = 0.03_wp   ! Here: leaf width (m)     
3256!       
3257!--    Stokes number for vegetated surfaces (Seinfeld & Pandis (2006): Eq.19.24) 
3258       St = vc * ustar / ( g * par_A )         
3259!
3260!--    Non-dimensional relexation time of the particle on top of canopy
3261       tau_plus = vc * ustar**2.0_wp / ( kvis_a * g ) 
3262!
3263!--    Brownian diffusion
3264       v_bd = mag_u * C_Br * Sc**( -2.0_wp / 3.0_wp ) *                        &
3265              ( mag_u * par_A / kvis_a )**( -0.5_wp )
3266!
3267!--    Interception
3268       v_in = mag_u * C_IN * paero(b)%dwet / par_A * ( 2.0_wp + LOG( 2.0_wp *  &
3269              par_A / paero(b)%dwet ) )                     
3270!
3271!--    Impaction: Petroff (2009) Eq. 18
3272       v_im = mag_u * C_IM * ( St / ( St + 0.47_wp ) )**2.0_wp
3273       
3274       IF ( tau_plus < 20.0_wp )  THEN
3275          v_it = 2.5E-3_wp * C_IT * tau_plus**2.0_wp
3276       ELSE
3277          v_it = C_IT
3278       ENDIF
3279       depo = ( v_bd + v_in + v_im + v_it + vc ) * lad     
3280       paero(b)%numc = paero(b)%numc - depo * paero(b)%numc * dt_salsa     
3281       DO  c = 1, maxspec+1
3282          paero(b)%volc(c) = paero(b)%volc(c) - depo * paero(b)%volc(c) *      &
3283                             dt_salsa
3284       ENDDO
3285    ENDIF 
3286 
3287 END SUBROUTINE depo_vege
3288 
3289!------------------------------------------------------------------------------!
3290! Description:
3291! ------------ 
3292!> Calculate deposition on horizontal and vertical surfaces. Implement as
3293!> surface flux.
3294!------------------------------------------------------------------------------!
3295
3296 SUBROUTINE depo_topo( i, j, surf, vc, Sc, kvis, mag_u, norm )
3297 
3298    USE surface_mod,                                                           &
3299        ONLY:  surf_type
3300 
3301    IMPLICIT NONE
3302   
3303    INTEGER(iwp), INTENT(in) ::  i     !< loop index
3304    INTEGER(iwp), INTENT(in) ::  j     !< loop index
3305    REAL(wp), INTENT(in) ::  kvis(:)   !< kinematic viscosity of air (m2/s)
3306    REAL(wp), INTENT(in) ::  mag_u(:)  !< wind velocity (m/s)                                                 
3307    REAL(wp), INTENT(in) ::  norm(:)   !< normalisation (usually air density)
3308    REAL(wp), INTENT(in) ::  Sc(:,:)  !< particle Schmidt number
3309    REAL(wp), INTENT(in) ::  vc(:,:)  !< terminal velocity (m/s)   
3310    TYPE(surf_type), INTENT(inout) :: surf  !< respective surface type
3311    INTEGER(iwp) ::  b      !< loop index
3312    INTEGER(iwp) ::  c      !< loop index
3313    INTEGER(iwp) ::  k      !< loop index
3314    INTEGER(iwp) ::  m      !< loop index
3315    INTEGER(iwp) ::  surf_e !< End index of surface elements at (j,i)-gridpoint
3316    INTEGER(iwp) ::  surf_s !< Start index of surface elements at (j,i)-gridpoint
3317    REAL(wp) ::  alpha      !< parameter, Table 3 in Zhang et al. (2001)
3318    REAL(wp) ::  C_Br       !< coefficient for Brownian diffusion
3319    REAL(wp) ::  C_IM       !< coefficient for inertial impaction
3320    REAL(wp) ::  C_IN       !< coefficient for interception
3321    REAL(wp) ::  C_IT       !< coefficient for turbulent impaction
3322    REAL(wp) ::  depo       !< deposition efficiency
3323    REAL(wp) ::  gamma      !< parameter, Table 3 in Zhang et al. (2001)
3324    REAL(wp) ::  par_A      !< parameter A for the characteristic radius of
3325                            !< collectors, Table 3 in Zhang et al. (2001)
3326    REAL(wp) ::  rt         !< the overall quasi-laminar resistance for
3327                            !< particles
3328    REAL(wp) ::  St         !< Stokes number for bluff surface elements 
3329    REAL(wp) ::  tau_plus   !< dimensionless particle relaxation time   
3330    REAL(wp) ::  v_bd       !< deposition velocity due to Brownian diffusion
3331    REAL(wp) ::  v_im       !< deposition velocity due to impaction
3332    REAL(wp) ::  v_in       !< deposition velocity due to interception
3333    REAL(wp) ::  v_it       !< deposition velocity due to turbulent impaction 
3334!
3335!-- Initialise
3336    rt       = 0.0_wp
3337    St       = 0.0_wp
3338    tau_plus = 0.0_wp
3339    v_bd     = 0.0_wp     
3340    v_im     = 0.0_wp       
3341    v_in     = 0.0_wp       
3342    v_it     = 0.0_wp                                 
3343    surf_s   = surf%start_index(j,i)
3344    surf_e   = surf%end_index(j,i) 
3345   
3346    DO  m = surf_s, surf_e
3347       k = surf%k(m)       
3348       DO  b = 1, nbins
3349          IF ( aerosol_number(b)%conc(k,j,i) <= nclim  .OR.                    &
3350               Sc(k+1,b) < 1.0_wp )  CYCLE   
3351                   
3352          IF ( depo_topo_type == 'zhang2001' )  THEN
3353!       
3354!--          Parameters for the land use category 'urban' in Table 3
3355             alpha = 1.5_wp
3356             gamma = 0.56_wp 
3357             par_A = 10.0E-3_wp
3358!       
3359!--          Stokes number for smooth surfaces or surfaces with bluff roughness
3360!--          elements (Seinfeld and Pandis, 2nd edition (2006): Eq. 19.23)       
3361             St = MAX( 0.01_wp, vc(k+1,b) * surf%us(m) ** 2.0_wp /             &
3362                       ( g * kvis(k+1)  ) ) 
3363!         
3364!--          The overall quasi-laminar resistance for particles (Eq. 5)       
3365             rt = MAX( EPSILON( 1.0_wp ), ( 3.0_wp * surf%us(m) * (            &
3366                       Sc(k+1,b)**( -gamma ) + ( St / ( alpha + St ) )**2.0_wp &
3367                        + 0.5_wp * ( Ra_dry(k,j,i,b) / par_A )**2.0_wp ) *     &
3368                       EXP( -St**0.5_wp ) ) ) 
3369             depo = vc(k+1,b) + rt
3370             
3371          ELSEIF ( depo_topo_type == 'petroff2010' )  THEN 
3372!
3373!--          vd = v_BD + v_IN + v_IM + v_IT + vc
3374!--          Deposition efficiencies from Table 1. Constants from Table 2.
3375             C_Br  = 1.262_wp
3376             C_IM  = 0.130_wp
3377             C_IN  = 0.216_wp
3378             C_IT  = 0.056_wp
3379             par_A = 0.03_wp   ! Here: leaf width (m) 
3380!       
3381!--          Stokes number for smooth surfaces or surfaces with bluff roughness
3382!--          elements (Seinfeld and Pandis, 2nd edition (2006): Eq. 19.23)       
3383             St = MAX( 0.01_wp, vc(k+1,b) * surf%us(m) ** 2.0_wp /             &
3384                       ( g *  kvis(k+1) ) )             
3385!
3386!--          Non-dimensional relexation time of the particle on top of canopy
3387             tau_plus = vc(k+1,b) * surf%us(m)**2.0_wp / ( kvis(k+1) * g ) 
3388!
3389!--          Brownian diffusion
3390             v_bd = mag_u(k+1) * C_Br * Sc(k+1,b)**( -2.0_wp / 3.0_wp ) *      &
3391                    ( mag_u(k+1) * par_A / kvis(k+1) )**( -0.5_wp )
3392!
3393!--          Interception
3394             v_in = mag_u(k+1) * C_IN * Ra_dry(k,j,i,b)/ par_A * ( 2.0_wp +    &
3395                    LOG( 2.0_wp * par_A / Ra_dry(k,j,i,b) ) )                     
3396!
3397!--          Impaction: Petroff (2009) Eq. 18
3398             v_im = mag_u(k+1) * C_IM * ( St / ( St + 0.47_wp ) )**2.0_wp
3399             
3400             IF ( tau_plus < 20.0_wp )  THEN
3401                v_it = 2.5E-3_wp * C_IT * tau_plus**2.0_wp
3402             ELSE
3403                v_it = C_IT
3404             ENDIF
3405             depo =  v_bd + v_in + v_im + v_it + vc(k+1,b)       
3406         
3407          ENDIF
3408          IF ( lod_aero == 3  .OR.  salsa_source_mode ==  'no_source' )  THEN
3409             surf%answs(m,b) = -depo * norm(k) * aerosol_number(b)%conc(k,j,i) 
3410             DO  c = 1, ncc_tot   
3411                surf%amsws(m,(c-1)*nbins+b) = -depo *  norm(k) *               &
3412                                         aerosol_mass((c-1)*nbins+b)%conc(k,j,i)
3413             ENDDO    ! c
3414          ELSE
3415             surf%answs(m,b) = SUM( aerosol_number(b)%source(:,j,i) ) -        &
3416                               MAX( 0.0_wp, depo * norm(k) *                   &
3417                               aerosol_number(b)%conc(k,j,i) )
3418             DO  c = 1, ncc_tot   
3419                surf%amsws(m,(c-1)*nbins+b) = SUM(                             &
3420                               aerosol_mass((c-1)*nbins+b)%source(:,j,i) ) -   &
3421                               MAX(  0.0_wp, depo *  norm(k) *                 &
3422                               aerosol_mass((c-1)*nbins+b)%conc(k,j,i) )
3423             ENDDO
3424          ENDIF
3425       ENDDO    ! b
3426    ENDDO    ! m     
3427     
3428 END SUBROUTINE depo_topo
3429 
3430!------------------------------------------------------------------------------!
3431! Description:
3432! ------------
3433! Function for calculating terminal velocities for different particles sizes.
3434!------------------------------------------------------------------------------!
3435 REAL(wp) FUNCTION terminal_vel( radius, rhop, rhoa, visc, beta )
3436 
3437    IMPLICIT NONE
3438   
3439    REAL(wp), INTENT(in) ::  beta    !< Cunningham correction factor
3440    REAL(wp), INTENT(in) ::  radius  !< particle radius (m)
3441    REAL(wp), INTENT(in) ::  rhop    !< particle density (kg/m3)
3442    REAL(wp), INTENT(in) ::  rhoa    !< air density (kg/m3)
3443    REAL(wp), INTENT(in) ::  visc    !< molecular viscosity of air (kg/(m*s))
3444   
3445    REAL(wp), PARAMETER ::  rhoa_ref = 1.225_wp ! reference air density (kg/m3)
3446!
3447!-- Stokes law with Cunningham slip correction factor
3448    terminal_vel = ( 4.0_wp * radius**2.0_wp ) * ( rhop - rhoa ) * g * beta /  &
3449                   ( 18.0_wp * visc ) ! (m/s)
3450       
3451 END FUNCTION terminal_vel
3452 
3453!------------------------------------------------------------------------------!
3454! Description:
3455! ------------
3456!> Calculates particle loss and change in size distribution due to (Brownian)
3457!> coagulation. Only for particles with dwet < 30 micrometres.
3458!
3459!> Method:
3460!> Semi-implicit, non-iterative method: (Jacobson, 1994)
3461!> Volume concentrations of the smaller colliding particles added to the bin of
3462!> the larger colliding particles. Start from first bin and use the updated
3463!> number and volume for calculation of following bins. NB! Our bin numbering
3464!> does not follow particle size in subrange 2.
3465!
3466!> Schematic for bin numbers in different subranges:
3467!>             1                            2
3468!>    +-------------------------------------------+
3469!>  a | 1 | 2 | 3 || 4 | 5 | 6 | 7 |  8 |  9 | 10||
3470!>  b |           ||11 |12 |13 |14 | 15 | 16 | 17||
3471!>    +-------------------------------------------+
3472!
3473!> Exact coagulation coefficients for each pressure level are scaled according
3474!> to current particle wet size (linear scaling).
3475!> Bins are organized in terms of the dry size of the condensation nucleus,
3476!> while coagulation kernell is calculated with the actual hydrometeor
3477!> size.
3478!
3479!> Called from salsa_driver
3480!> fxm: Process selection should be made smarter - now just lots of IFs inside
3481!>      loops
3482!
3483!> Coded by:
3484!> Hannele Korhonen (FMI) 2005
3485!> Harri Kokkola (FMI) 2006
3486!> Tommi Bergman (FMI) 2012
3487!> Matti Niskanen(FMI) 2012
3488!> Anton Laakso  (FMI) 2013
3489!> Juha Tonttila (FMI) 2014
3490!------------------------------------------------------------------------------!
3491 SUBROUTINE coagulation( paero, ptstep, ptemp, ppres )
3492               
3493    IMPLICIT NONE
3494   
3495!-- Input and output variables
3496    TYPE(t_section), INTENT(inout) ::  paero(fn2b) !< Aerosol properties
3497    REAL(wp), INTENT(in) ::  ppres  !< ambient pressure (Pa)
3498    REAL(wp), INTENT(in) ::  ptemp  !< ambient temperature (K)
3499    REAL(wp), INTENT(in) ::  ptstep !< time step (s)
3500!-- Local variables
3501    INTEGER(iwp) ::  index_2a !< corresponding bin in subrange 2a
3502    INTEGER(iwp) ::  index_2b !< corresponding bin in subrange 2b
3503    INTEGER(iwp) ::  b !< loop index
3504    INTEGER(iwp) ::  ll !< loop index
3505    INTEGER(iwp) ::  mm !< loop index
3506    INTEGER(iwp) ::  nn !< loop index
3507    REAL(wp) ::  pressi !< pressure
3508    REAL(wp) ::  temppi !< temperature
3509    REAL(wp) ::  zcc(fn2b,fn2b)   !< updated coagulation coefficients (m3/s) 
3510    REAL(wp) ::  zdpart_mm        !< diameter of particle (m)
3511    REAL(wp) ::  zdpart_nn        !< diameter of particle (m)   
3512    REAL(wp) ::  zminusterm       !< coagulation loss in a bin (1/s)
3513    REAL(wp) ::  zplusterm(8)     !< coagulation gain in a bin (fxm/s)
3514                                  !< (for each chemical compound)
3515    REAL(wp) ::  zmpart(fn2b)     !< approximate mass of particles (kg)
3516   
3517    zcc       = 0.0_wp
3518    zmpart    = 0.0_wp
3519    zdpart_mm = 0.0_wp
3520    zdpart_nn = 0.0_wp
3521!
3522!-- 1) Coagulation to coarse mode calculated in a simplified way:
3523!--    CoagSink ~ Dp in continuum subrange, thus we calculate 'effective'
3524!--    number concentration of coarse particles
3525
3526!-- 2) Updating coagulation coefficients
3527!   
3528!-- Aerosol mass (kg). Density of 1500 kg/m3 assumed
3529    zmpart(1:fn2b) = api6 * ( MIN( paero(1:fn2b)%dwet, 30.0E-6_wp )**3.0_wp  ) &
3530                     * 1500.0_wp 
3531    temppi = ptemp
3532    pressi = ppres
3533    zcc    = 0.0_wp
3534!
3535!-- Aero-aero coagulation
3536    DO  mm = 1, fn2b   ! smaller colliding particle
3537       IF ( paero(mm)%numc < nclim )  CYCLE
3538       DO  nn = mm, fn2b   ! larger colliding particle
3539          IF ( paero(nn)%numc < nclim )  CYCLE
3540         
3541          zdpart_mm = MIN( paero(mm)%dwet, 30.0E-6_wp )     ! Limit to 30 um
3542          zdpart_nn = MIN( paero(nn)%dwet, 30.0E-6_wp )     ! Limit to 30 um
3543!             
3544!--       Coagulation coefficient of particles (m3/s)
3545          zcc(mm,nn) = coagc( zdpart_mm, zdpart_nn, zmpart(mm), zmpart(nn),    &
3546                              temppi, pressi )
3547          zcc(nn,mm) = zcc(mm,nn)
3548       ENDDO
3549    ENDDO
3550       
3551!   
3552!-- 3) New particle and volume concentrations after coagulation:
3553!--    Calculated according to Jacobson (2005) eq. 15.9
3554!
3555!-- Aerosols in subrange 1a:
3556    DO  b = in1a, fn1a
3557       IF ( paero(b)%numc < nclim )  CYCLE
3558       zminusterm   = 0.0_wp
3559       zplusterm(:) = 0.0_wp
3560!       
3561!--    Particles lost by coagulation with larger aerosols
3562       DO  ll = b+1, fn2b
3563          zminusterm = zminusterm + zcc(b,ll) * paero(ll)%numc
3564       ENDDO
3565!       
3566!--    Coagulation gain in a bin: change in volume conc. (cm3/cm3):
3567       DO ll = in1a, b-1
3568          zplusterm(1:2) = zplusterm(1:2) + zcc(ll,b) * paero(ll)%volc(1:2)
3569          zplusterm(6:7) = zplusterm(6:7) + zcc(ll,b) * paero(ll)%volc(6:7)
3570          zplusterm(8)   = zplusterm(8)   + zcc(ll,b) * paero(ll)%volc(8)
3571       ENDDO
3572!       
3573!--    Volume and number concentrations after coagulation update [fxm]
3574       paero(b)%volc(1:2) = ( paero(b)%volc(1:2) + ptstep * zplusterm(1:2) * &
3575                             paero(b)%numc ) / ( 1.0_wp + ptstep * zminusterm )
3576       paero(b)%volc(6:7) = ( paero(b)%volc(6:7) + ptstep * zplusterm(6:7) * &
3577                             paero(b)%numc ) / ( 1.0_wp + ptstep * zminusterm )
3578       paero(b)%volc(8)   = ( paero(b)%volc(8)   + ptstep * zplusterm(8) *   &
3579                             paero(b)%numc ) / ( 1.0_wp + ptstep * zminusterm )
3580       paero(b)%numc = paero(b)%numc / ( 1.0_wp + ptstep * zminusterm  +     &
3581                        0.5_wp * ptstep * zcc(b,b) * paero(b)%numc )               
3582    ENDDO
3583!             
3584!-- Aerosols in subrange 2a:
3585    DO  b = in2a, fn2a
3586       IF ( paero(b)%numc < nclim )  CYCLE
3587       zminusterm   = 0.0_wp
3588       zplusterm(:) = 0.0_wp
3589!       
3590!--    Find corresponding size bin in subrange 2b
3591       index_2b = b - in2a + in2b
3592!       
3593!--    Particles lost by larger particles in 2a
3594       DO  ll = b+1, fn2a
3595          zminusterm = zminusterm + zcc(b,ll) * paero(ll)%numc
3596       ENDDO
3597!       
3598!--    Particles lost by larger particles in 2b
3599       IF ( .NOT. no_insoluble )  THEN
3600          DO  ll = index_2b+1, fn2b
3601             zminusterm = zminusterm + zcc(b,ll) * paero(ll)%numc
3602          ENDDO
3603       ENDIF
3604!       
3605!--    Particle volume gained from smaller particles in subranges 1, 2a and 2b
3606       DO  ll = in1a, b-1
3607          zplusterm(1:2) = zplusterm(1:2) + zcc(ll,b) * paero(ll)%volc(1:2)
3608          zplusterm(6:7) = zplusterm(6:7) + zcc(ll,b) * paero(ll)%volc(6:7)
3609          zplusterm(8)   = zplusterm(8)   + zcc(ll,b) * paero(ll)%volc(8)
3610       ENDDO 
3611!       
3612!--    Particle volume gained from smaller particles in 2a
3613!--    (Note, for components not included in the previous loop!)
3614       DO  ll = in2a, b-1
3615          zplusterm(3:5) = zplusterm(3:5) + zcc(ll,b)*paero(ll)%volc(3:5)             
3616       ENDDO
3617       
3618!       
3619!--    Particle volume gained from smaller (and equal) particles in 2b
3620       IF ( .NOT. no_insoluble )  THEN
3621          DO  ll = in2b, index_2b
3622             zplusterm(1:8) = zplusterm(1:8) + zcc(ll,b) * paero(ll)%volc(1:8)
3623          ENDDO
3624       ENDIF
3625!       
3626!--    Volume and number concentrations after coagulation update [fxm]
3627       paero(b)%volc(1:8) = ( paero(b)%volc(1:8) + ptstep * zplusterm(1:8) * &
3628                             paero(b)%numc ) / ( 1.0_wp + ptstep * zminusterm )
3629       paero(b)%numc = paero(b)%numc / ( 1.0_wp + ptstep * zminusterm +      &
3630                        0.5_wp * ptstep * zcc(b,b) * paero(b)%numc )
3631    ENDDO
3632!             
3633!-- Aerosols in subrange 2b:
3634    IF ( .NOT. no_insoluble )  THEN
3635       DO  b = in2b, fn2b
3636          IF ( paero(b)%numc < nclim )  CYCLE
3637          zminusterm   = 0.0_wp
3638          zplusterm(:) = 0.0_wp
3639!       
3640!--       Find corresponding size bin in subsubrange 2a
3641          index_2a = b - in2b + in2a
3642!       
3643!--       Particles lost to larger particles in subranges 2b
3644          DO  ll = b+1, fn2b
3645             zminusterm = zminusterm + zcc(b,ll) * paero(ll)%numc
3646          ENDDO
3647!       
3648!--       Particles lost to larger and equal particles in 2a
3649          DO  ll = index_2a, fn2a
3650             zminusterm = zminusterm + zcc(b,ll) * paero(ll)%numc
3651          ENDDO
3652!       
3653!--       Particle volume gained from smaller particles in subranges 1 & 2a
3654          DO  ll = in1a, index_2a-1
3655             zplusterm(1:8) = zplusterm(1:8) + zcc(ll,b) * paero(ll)%volc(1:8)
3656          ENDDO
3657!       
3658!--       Particle volume gained from smaller particles in 2b
3659          DO  ll = in2b, b-1
3660             zplusterm(1:8) = zplusterm(1:8) + zcc(ll,b) * paero(ll)%volc(1:8)
3661          ENDDO
3662!       
3663!--       Volume and number concentrations after coagulation update [fxm]
3664          paero(b)%volc(1:8) = ( paero(b)%volc(1:8) + ptstep * zplusterm(1:8)&
3665                           * paero(b)%numc ) / ( 1.0_wp + ptstep * zminusterm )
3666          paero(b)%numc = paero(b)%numc / ( 1.0_wp + ptstep * zminusterm  +  &
3667                           0.5_wp * ptstep * zcc(b,b) * paero(b)%numc )
3668       ENDDO
3669    ENDIF
3670
3671 END SUBROUTINE coagulation
3672
3673!------------------------------------------------------------------------------!
3674! Description:
3675! ------------
3676!> Calculation of coagulation coefficients. Extended version of the function
3677!> originally found in mo_salsa_init.
3678!
3679!> J. Tonttila, FMI, 05/2014
3680!------------------------------------------------------------------------------!
3681 REAL(wp) FUNCTION coagc( diam1, diam2, mass1, mass2, temp, pres )
3682 
3683    IMPLICIT NONE
3684!       
3685!-- Input and output variables
3686    REAL(wp), INTENT(in) ::  diam1 !< diameter of colliding particle 1 (m)
3687    REAL(wp), INTENT(in) ::  diam2 !< diameter of colliding particle 2 (m)
3688    REAL(wp), INTENT(in) ::  mass1 !< mass of colliding particle 1 (kg)
3689    REAL(wp), INTENT(in) ::  mass2 !< mass of colliding particle 2 (kg)
3690    REAL(wp), INTENT(in) ::  pres  !< ambient pressure (Pa?) [fxm]
3691    REAL(wp), INTENT(in) ::  temp  !< ambient temperature (K)       
3692!
3693!-- Local variables
3694    REAL(wp) ::  fmdist !< distance of flux matching (m)   
3695    REAL(wp) ::  knud_p !< particle Knudsen number
3696    REAL(wp) ::  mdiam  !< mean diameter of colliding particles (m) 
3697    REAL(wp) ::  mfp    !< mean free path of air molecules (m)   
3698    REAL(wp) ::  visc   !< viscosity of air (kg/(m s))                   
3699    REAL(wp), DIMENSION (2) ::  beta   !< Cunningham correction factor
3700    REAL(wp), DIMENSION (2) ::  dfpart !< particle diffusion coefficient
3701                                       !< (m2/s)       
3702    REAL(wp), DIMENSION (2) ::  diam   !< diameters of particles (m)
3703    REAL(wp), DIMENSION (2) ::  flux   !< flux in continuum and free molec.
3704                                       !< regime (m/s)       
3705    REAL(wp), DIMENSION (2) ::  knud   !< particle Knudsen number       
3706    REAL(wp), DIMENSION (2) ::  mpart  !< masses of particles (kg)
3707    REAL(wp), DIMENSION (2) ::  mtvel  !< particle mean thermal velocity (m/s)
3708    REAL(wp), DIMENSION (2) ::  omega  !< particle mean free path             
3709    REAL(wp), DIMENSION (2) ::  tva    !< temporary variable (m)       
3710!
3711!-- Initialisation
3712    coagc   = 0.0_wp
3713!
3714!-- 1) Initializing particle and ambient air variables
3715    diam  = (/ diam1, diam2 /) !< particle diameters (m)
3716    mpart = (/ mass1, mass2 /) !< particle masses (kg)
3717!-- Viscosity of air (kg/(m s))       
3718    visc = ( 7.44523E-3_wp * temp ** 1.5_wp ) /                                &
3719           ( 5093.0_wp * ( temp + 110.4_wp ) ) 
3720!-- Mean free path of air (m)           
3721    mfp = ( 1.656E-10_wp * temp + 1.828E-8_wp ) * ( p_0 + 1325.0_wp ) / pres
3722!
3723!-- 2) Slip correction factor for small particles
3724    knud = 2.0_wp * EXP( LOG(mfp) - LOG(diam) )! Knudsen number for air (15.23)
3725!-- Cunningham correction factor (Allen and Raabe, Aerosol Sci. Tech. 4, 269)       
3726    beta = 1.0_wp + knud * ( 1.142_wp + 0.558_wp * EXP( -0.999_wp / knud ) )
3727!
3728!-- 3) Particle properties
3729!-- Diffusion coefficient (m2/s) (Jacobson (2005) eq. 15.29)
3730    dfpart = beta * abo * temp / ( 3.0_wp * pi * visc * diam ) 
3731!-- Mean thermal velocity (m/s) (Jacobson (2005) eq. 15.32)
3732    mtvel = SQRT( ( 8.0_wp * abo * temp ) / ( pi * mpart ) )
3733!-- Particle mean free path (m) (Jacobson (2005) eq. 15.34 )
3734    omega = 8.0_wp * dfpart / ( pi * mtvel ) 
3735!-- Mean diameter (m)
3736    mdiam = 0.5_wp * ( diam(1) + diam(2) )
3737!
3738!-- 4) Calculation of fluxes (Brownian collision kernels) and flux matching
3739!-- following Jacobson (2005):
3740!-- Flux in continuum regime (m3/s) (eq. 15.28)
3741    flux(1) = 4.0_wp * pi * mdiam * ( dfpart(1) + dfpart(2) )
3742!-- Flux in free molec. regime (m3/s) (eq. 15.31)
3743    flux(2) = pi * SQRT( ( mtvel(1)**2.0_wp ) + ( mtvel(2)**2.0_wp ) ) *      &
3744              ( mdiam**2.0_wp )
3745!-- temporary variables (m) to calculate flux matching distance (m)
3746    tva(1) = ( ( mdiam + omega(1) )**3.0_wp - ( mdiam**2.0_wp +                &
3747               omega(1)**2.0_wp ) * SQRT( ( mdiam**2.0_wp + omega(1)**2.0_wp ) &
3748               ) ) / ( 3.0_wp * mdiam * omega(1) ) - mdiam
3749    tva(2) = ( ( mdiam + omega(2) )**3.0_wp - ( mdiam**2.0_wp +                &
3750               omega(2)**2.0_wp ) * SQRT( ( mdiam**2 + omega(2)**2 ) ) ) /     &
3751             ( 3.0_wp * mdiam * omega(2) ) - mdiam
3752!-- Flux matching distance (m) i.e. the mean distance from the centre of a
3753!-- sphere reached by particles leaving sphere's surface and travelling a
3754!-- distance of particle mean free path mfp (eq. 15 34)                 
3755    fmdist = SQRT( tva(1)**2 + tva(2)**2.0_wp) 
3756!
3757!-- 5) Coagulation coefficient (m3/s) (eq. 15.33). Here assumed
3758!-- coalescence efficiency 1!!
3759    coagc = flux(1) / ( mdiam / ( mdiam + fmdist) + flux(1) / flux(2) ) 
3760!-- coagulation coefficient = coalescence efficiency * collision kernel
3761!
3762!-- Corrected collision kernel following Karl et al., 2016 (ACP):
3763!-- Inclusion of van der Waals and viscous forces
3764    IF ( van_der_waals_coagc )  THEN
3765       knud_p = SQRT( omega(1)**2 + omega(2)**2 ) / mdiam   
3766       IF ( knud_p >= 0.1_wp  .AND.  knud_p <= 10.0_wp )  THEN
3767          coagc = coagc * ( 2.0_wp + 0.4_wp * LOG( knud_p ) )
3768       ELSE
3769          coagc = coagc * 3.0_wp
3770       ENDIF
3771    ENDIF
3772   
3773 END FUNCTION coagc
3774 
3775!------------------------------------------------------------------------------!   
3776! Description:
3777! ------------
3778!> Calculates the change in particle volume and gas phase
3779!> concentrations due to nucleation, condensation and dissolutional growth.
3780!
3781!> Sulphuric acid and organic vapour: only condensation and no evaporation.
3782!
3783!> New gas and aerosol phase concentrations calculated according to Jacobson
3784!> (1997): Numerical techniques to solve condensational and dissolutional growth
3785!> equations when growth is coupled to reversible reactions, Aerosol Sci. Tech.,
3786!> 27, pp 491-498.
3787!
3788!> Following parameterization has been used:
3789!> Molecular diffusion coefficient of condensing vapour (m2/s)
3790!> (Reid et al. (1987): Properties of gases and liquids, McGraw-Hill, New York.)
3791!> D = {1.d-7*sqrt(1/M_air + 1/M_gas)*T^1.75} / &
3792!      {p_atm/p_stand * (d_air^(1/3) + d_gas^(1/3))^2 }
3793! M_air = 28.965 : molar mass of air (g/mol)
3794! d_air = 19.70  : diffusion volume of air
3795! M_h2so4 = 98.08 : molar mass of h2so4 (g/mol)
3796! d_h2so4 = 51.96  : diffusion volume of h2so4
3797!
3798!> Called from main aerosol model
3799!
3800!> fxm: calculated for empty bins too
3801!> fxm: same diffusion coefficients and mean free paths used for sulphuric acid
3802!>      and organic vapours (average values? 'real' values for each?)
3803!> fxm: one should really couple with vapour production and loss terms as well
3804!>      should nucleation be coupled here as well????
3805!
3806! Coded by:
3807! Hannele Korhonen (FMI) 2005
3808! Harri Kokkola (FMI) 2006
3809! Juha Tonttila (FMI) 2014
3810! Rewritten to PALM by Mona Kurppa (UHel) 2017
3811!------------------------------------------------------------------------------!
3812 SUBROUTINE condensation( paero, pcsa, pcocnv, pcocsv, pchno3, pcnh3, pcw, pcs,&
3813                          ptemp, ppres, ptstep, prtcl )
3814       
3815    IMPLICIT NONE
3816   
3817!-- Input and output variables
3818    REAL(wp), INTENT(IN) ::  ppres !< ambient pressure (Pa)
3819    REAL(wp), INTENT(IN) ::  pcs   !< Water vapour saturation concentration
3820                                   !< (kg/m3)     
3821    REAL(wp), INTENT(IN) ::  ptemp !< ambient temperature (K)
3822    REAL(wp), INTENT(IN) ::  ptstep            !< timestep (s) 
3823    TYPE(component_index), INTENT(in) :: prtcl !< Keeps track which substances
3824                                               !< are used                                               
3825    REAL(wp), INTENT(INOUT) ::  pchno3 !< Gas concentrations (#/m3):
3826                                       !< nitric acid HNO3
3827    REAL(wp), INTENT(INOUT) ::  pcnh3  !< ammonia NH3
3828    REAL(wp), INTENT(INOUT) ::  pcocnv !< non-volatile organics
3829    REAL(wp), INTENT(INOUT) ::  pcocsv !< semi-volatile organics
3830    REAL(wp), INTENT(INOUT) ::  pcsa   !< sulphuric acid H2SO4
3831    REAL(wp), INTENT(INOUT) ::  pcw    !< Water vapor concentration (kg/m3)
3832    TYPE(t_section), INTENT(inout) ::  paero(fn2b) !< Aerosol properties                                     
3833!-- Local variables
3834    REAL(wp) ::  zbeta(fn2b) !< transitional correction factor for aerosols
3835    REAL(wp) ::  zcolrate(fn2b) !< collision rate of molecules to particles
3836                                !< (1/s)
3837    REAL(wp) ::  zcolrate_ocnv(fn2b) !< collision rate of organic molecules
3838                                     !< to particles (1/s)
3839    REAL(wp) ::  zcs_ocnv !< condensation sink of nonvolatile organics (1/s)       
3840    REAL(wp) ::  zcs_ocsv !< condensation sink of semivolatile organics (1/s)
3841    REAL(wp) ::  zcs_su !< condensation sink of sulfate (1/s)
3842    REAL(wp) ::  zcs_tot!< total condensation sink (1/s) (gases)
3843!-- vapour concentration after time step (#/m3)
3844    REAL(wp) ::  zcvap_new1 !< sulphuric acid
3845    REAL(wp) ::  zcvap_new2 !< nonvolatile organics
3846    REAL(wp) ::  zcvap_new3 !< semivolatile organics
3847    REAL(wp) ::  zdfpart(in1a+1) !< particle diffusion coefficient (m2/s)     
3848    REAL(wp) ::  zdfvap !< air diffusion coefficient (m2/s)
3849!-- change in vapour concentration (#/m3)
3850    REAL(wp) ::  zdvap1 !< sulphuric acid
3851    REAL(wp) ::  zdvap2 !< nonvolatile organics
3852    REAL(wp) ::  zdvap3 !< semivolatile organics
3853    REAL(wp) ::  zdvoloc(fn2b) !< change of organics volume in each bin [fxm]   
3854    REAL(wp) ::  zdvolsa(fn2b) !< change of sulphate volume in each bin [fxm]
3855    REAL(wp) ::  zj3n3(2)      !< Formation massrate of molecules in
3856                               !< nucleation, (molec/m3s). 1: H2SO4
3857                               !< and 2: organic vapor       
3858    REAL(wp) ::  zknud(fn2b) !< particle Knudsen number       
3859    REAL(wp) ::  zmfp    !< mean free path of condensing vapour (m)
3860    REAL(wp) ::  zrh     !< Relative humidity [0-1]         
3861    REAL(wp) ::  zvisc   !< viscosity of air (kg/(m s))     
3862    REAL(wp) ::  zn_vs_c !< ratio of nucleation of all mass transfer in the
3863                         !< smallest bin
3864    REAL(wp) ::  zxocnv  !< ratio of organic vapour in 3nm particles
3865    REAL(wp) ::  zxsa    !< Ratio in 3nm particles: sulphuric acid
3866   
3867    zj3n3  = 0.0_wp
3868    zrh    = pcw / pcs   
3869    zxocnv = 0.0_wp
3870    zxsa   = 0.0_wp
3871!
3872!-- Nucleation
3873    IF ( nsnucl > 0 )  THEN
3874       CALL nucleation( paero, ptemp, zrh, ppres, pcsa, pcocnv, pcnh3, ptstep, &
3875                        zj3n3, zxsa, zxocnv )
3876    ENDIF
3877!
3878!-- Condensation on pre-existing particles
3879    IF ( lscndgas )  THEN
3880!
3881!--    Initialise:
3882       zdvolsa = 0.0_wp 
3883       zdvoloc = 0.0_wp
3884       zcolrate = 0.0_wp
3885!             
3886!--    1) Properties of air and condensing gases:
3887!--    Viscosity of air (kg/(m s)) (Eq. 4.54 in Jabonson (2005))
3888       zvisc = ( 7.44523E-3_wp * ptemp ** 1.5_wp ) / ( 5093.0_wp *             &
3889                 ( ptemp + 110.4_wp ) )
3890!--    Diffusion coefficient of air (m2/s)
3891       zdfvap = 5.1111E-10_wp * ptemp ** 1.75_wp * ( p_0 + 1325.0_wp ) / ppres
3892!--    Mean free path (m): same for H2SO4 and organic compounds
3893       zmfp = 3.0_wp * zdfvap * SQRT( pi * amh2so4 / ( 8.0_wp * argas * ptemp ) )
3894!                   
3895!--    2) Transition regime correction factor zbeta for particles:
3896!--       Fuchs and Sutugin (1971), In: Hidy et al. (ed.) Topics in current
3897!--       aerosol research, Pergamon. Size of condensing molecule considered 
3898!--       only for nucleation mode (3 - 20 nm)
3899!
3900!--    Particle Knudsen number: condensation of gases on aerosols
3901       zknud(in1a:in1a+1) = 2.0_wp * zmfp / ( paero(in1a:in1a+1)%dwet + d_sa )
3902       zknud(in1a+2:fn2b) = 2.0_wp * zmfp / paero(in1a+2:fn2b)%dwet
3903!   
3904!--    Transitional correction factor: aerosol + gas (the semi-empirical Fuchs-
3905!--    Sutugin interpolation function (Fuchs and Sutugin, 1971))
3906       zbeta = ( zknud + 1.0_wp ) / ( 0.377_wp * zknud + 1.0_wp + 4.0_wp /     &
3907               ( 3.0_wp * massacc ) * ( zknud + zknud ** 2.0_wp ) )
3908!                   
3909!--    3) Collision rate of molecules to particles
3910!--       Particle diffusion coefficient considered only for nucleation mode
3911!--       (3 - 20 nm)
3912!
3913!--    Particle diffusion coefficient (m2/s) (e.g. Eq. 15.29 in Jacobson (2005))
3914       zdfpart = abo * ptemp * zbeta(in1a:in1a+1) / ( 3.0_wp * pi * zvisc *    &
3915                 paero(in1a:in1a+1)%dwet )
3916!             
3917!--    Collision rate (mass-transfer coefficient): gases on aerosols (1/s)
3918!--    (Eq. 16.64 in Jacobson (2005))
3919       zcolrate(in1a:in1a+1) = MERGE( 2.0_wp * pi *                            &
3920                                      ( paero(in1a:in1a+1)%dwet + d_sa ) *     &
3921                                      ( zdfvap + zdfpart ) * zbeta(in1a:in1a+1)& 
3922                                        * paero(in1a:in1a+1)%numc, 0.0_wp,     &
3923                                      paero(in1a:in1a+1)%numc > nclim )
3924       zcolrate(in1a+2:fn2b) = MERGE( 2.0_wp * pi * paero(in1a+2:fn2b)%dwet *  &
3925                                      zdfvap * zbeta(in1a+2:fn2b) *            &
3926                                      paero(in1a+2:fn2b)%numc, 0.0_wp,         &
3927                                      paero(in1a+2:fn2b)%numc > nclim )
3928!                 
3929!-- 4) Condensation sink (1/s)
3930       zcs_tot = SUM( zcolrate )   ! total sink
3931!
3932!--    5) Changes in gas-phase concentrations and particle volume
3933!
3934!--    5.1) Organic vapours
3935!
3936!--    5.1.1) Non-volatile organic compound: condenses onto all bins
3937       IF ( pcocnv > 1.0E+10_wp  .AND.  zcs_tot > 1.0E-30_wp  .AND.            &
3938            is_used( prtcl,'OC' ) )                                            &
3939       THEN
3940!--       Ratio of nucleation vs. condensation rates in the smallest bin   
3941          zn_vs_c = 0.0_wp 
3942          IF ( zj3n3(2) > 1.0_wp )  THEN
3943             zn_vs_c = ( zj3n3(2) ) / ( zj3n3(2) + pcocnv * zcolrate(in1a) )
3944          ENDIF
3945!       
3946!--       Collision rate in the smallest bin, including nucleation and
3947!--       condensation(see Jacobson, Fundamentals of Atmospheric Modeling, 2nd
3948!--       Edition (2005), equation (16.73) )
3949          zcolrate_ocnv = zcolrate
3950          zcolrate_ocnv(in1a) = zcolrate_ocnv(in1a) + zj3n3(2) / pcocnv
3951!       
3952!--       Total sink for organic vapor
3953          zcs_ocnv = zcs_tot + zj3n3(2) / pcocnv
3954!       
3955!--       New gas phase concentration (#/m3)
3956          zcvap_new2 = pcocnv / ( 1.0_wp + ptstep * zcs_ocnv )
3957!       
3958!--       Change in gas concentration (#/m3)
3959          zdvap2 = pcocnv - zcvap_new2
3960!
3961!--       Updated vapour concentration (#/m3)               
3962          pcocnv = zcvap_new2
3963!       
3964!--       Volume change of particles (m3(OC)/m3(air))
3965          zdvoloc = zcolrate_ocnv(in1a:fn2b) / zcs_ocnv * amvoc * zdvap2
3966!       
3967!--       Change of volume due to condensation in 1a-2b
3968          paero(in1a:fn2b)%volc(2) = paero(in1a:fn2b)%volc(2) + zdvoloc 
3969!       
3970!--       Change of number concentration in the smallest bin caused by
3971!--       nucleation (Jacobson (2005), equation (16.75)). If zxocnv = 0, then 
3972!--       the chosen nucleation mechanism doesn't take into account the non-
3973!--       volatile organic vapors and thus the paero doesn't have to be updated.
3974          IF ( zxocnv > 0.0_wp )  THEN
3975             paero(in1a)%numc = paero(in1a)%numc + zn_vs_c * zdvoloc(in1a) /   &
3976                                amvoc / ( n3 * zxocnv )
3977          ENDIF
3978       ENDIF
3979!   
3980!--    5.1.2) Semivolatile organic compound: all bins except subrange 1
3981       zcs_ocsv = SUM( zcolrate(in2a:fn2b) ) !< sink for semi-volatile organics
3982       IF ( pcocsv > 1.0E+10_wp  .AND.  zcs_ocsv > 1.0E-30  .AND.              &
3983            is_used( prtcl,'OC') )                                             &
3984       THEN
3985!
3986!--       New gas phase concentration (#/m3)
3987          zcvap_new3 = pcocsv / ( 1.0_wp + ptstep * zcs_ocsv )
3988!       
3989!--       Change in gas concentration (#/m3)
3990          zdvap3 = pcocsv - zcvap_new3
3991!       
3992!--       Updated gas concentration (#/m3)               
3993          pcocsv = zcvap_new3
3994!       
3995!--       Volume change of particles (m3(OC)/m3(air))
3996          zdvoloc(in2a:fn2b) = zdvoloc(in2a:fn2b) + zcolrate(in2a:fn2b) /      &
3997                               zcs_ocsv * amvoc * zdvap3
3998!                           
3999!--       Change of volume due to condensation in 1a-2b
4000          paero(in1a:fn2b)%volc(2) = paero(in1a:fn2b)%volc(2) + zdvoloc
4001       ENDIF
4002!
4003!-- 5.2) Sulphate: condensed on all bins
4004       IF ( pcsa > 1.0E+10_wp  .AND.  zcs_tot > 1.0E-30_wp  .AND.              &
4005            is_used( prtcl,'SO4' ) )                                           &
4006       THEN
4007!   
4008!--    Ratio of mass transfer between nucleation and condensation
4009          zn_vs_c = 0.0_wp
4010          IF ( zj3n3(1) > 1.0_wp )  THEN
4011             zn_vs_c = ( zj3n3(1) ) / ( zj3n3(1) + pcsa * zcolrate(in1a) )
4012          ENDIF
4013!       
4014!--       Collision rate in the smallest bin, including nucleation and
4015!--       condensation (see Jacobson, Fundamentals of Atmospheric Modeling, 2nd
4016!--       Edition (2005), equation (16.73))
4017          zcolrate(in1a) = zcolrate(in1a) + zj3n3(1) / pcsa     
4018!       
4019!--       Total sink for sulfate (1/s)
4020          zcs_su = zcs_tot + zj3n3(1) / pcsa
4021!       
4022!--       Sulphuric acid:
4023!--       New gas phase concentration (#/m3)
4024          zcvap_new1 = pcsa / ( 1.0_wp + ptstep * zcs_su )
4025!       
4026!--       Change in gas concentration (#/m3)
4027          zdvap1 = pcsa - zcvap_new1
4028!       
4029!--       Updating vapour concentration (#/m3)
4030          pcsa = zcvap_new1
4031!       
4032!--       Volume change of particles (m3(SO4)/m3(air)) by condensation
4033          zdvolsa = zcolrate(in1a:fn2b) / zcs_su * amvh2so4 * zdvap1
4034!--       For validation: zdvolsa = 5.5 mum3/cm3 per 12 h       
4035       !   zdvolsa = zdvolsa / SUM( zdvolsa ) * 5.5E-12_wp * dt_salsa / 43200.0_wp 
4036          !0.3E-12_wp, 0.6E-12_wp, 11.0E-12_wp, 4.6E-12_wp, 9.2E-12_wp   
4037!       
4038!--       Change of volume concentration of sulphate in aerosol [fxm]
4039          paero(in1a:fn2b)%volc(1) = paero(in1a:fn2b)%volc(1) + zdvolsa
4040!       
4041!--       Change of number concentration in the smallest bin caused by nucleation
4042!--       (Jacobson (2005), equation (16.75))
4043          IF ( zxsa > 0.0_wp )  THEN
4044             paero(in1a)%numc = paero(in1a)%numc + zn_vs_c * zdvolsa(in1a) /   &
4045                                amvh2so4 / ( n3 * zxsa )
4046          ENDIF
4047       ENDIF
4048    ENDIF
4049!
4050!
4051!-- Condensation of water vapour
4052    IF ( lscndh2oae )  THEN
4053       CALL gpparth2o( paero, ptemp, ppres, pcs, pcw, ptstep )
4054    ENDIF
4055!   
4056!
4057!-- Partitioning of H2O, HNO3, and NH3: Dissolutional growth
4058    IF ( lscndgas  .AND.  ino > 0  .AND.  inh > 0  .AND.                       &
4059         ( pchno3 > 1.0E+10_wp  .OR.  pcnh3 > 1.0E+10_wp ) )                   &
4060    THEN
4061       CALL gpparthno3( ppres, ptemp, paero, pchno3, pcnh3, pcw, pcs, zbeta,   &
4062                        ptstep )
4063    ENDIF
4064   
4065 END SUBROUTINE condensation
4066 
4067!------------------------------------------------------------------------------!
4068! Description:
4069! ------------
4070!> Calculates the particle number and volume increase, and gas-phase
4071!> concentration decrease due to nucleation subsequent growth to detectable size
4072!> of 3 nm.
4073!
4074!> Method:
4075!> When the formed clusters grow by condensation (possibly also by self-
4076!> coagulation), their number is reduced due to scavenging to pre-existing
4077!> particles. Thus, the apparent nucleation rate at 3 nm is significantly lower
4078!> than the real nucleation rate (at ~1 nm).
4079!
4080!> Calculation of the formation rate of detectable particles at 3 nm (i.e. J3):
4081!> nj3 = 1: Kerminen, V.-M. and Kulmala, M. (2002), J. Aerosol Sci.,33, 609-622.
4082!> nj3 = 2: Lehtinen et al. (2007), J. Aerosol Sci., 38(9), 988-994.
4083!> nj3 = 3: Anttila et al. (2010), J. Aerosol Sci., 41(7), 621-636.
4084!
4085!> Called from subroutine condensation (in module salsa_dynamics_mod.f90)
4086!
4087!> Calls one of the following subroutines:
4088!>  - binnucl
4089!>  - ternucl
4090!>  - kinnucl
4091!>  - actnucl
4092!
4093!> fxm: currently only sulphuric acid grows particles from 1 to 3 nm
4094!>  (if asked from Markku, this is terribly wrong!!!)
4095!
4096!> Coded by:
4097!> Hannele Korhonen (FMI) 2005
4098!> Harri Kokkola (FMI) 2006
4099!> Matti Niskanen(FMI) 2012
4100!> Anton Laakso  (FMI) 2013
4101!------------------------------------------------------------------------------!
4102
4103 SUBROUTINE nucleation( paero, ptemp, prh, ppres, pcsa, pcocnv, pcnh3, ptstep, &
4104                        pj3n3, pxsa, pxocnv )
4105    IMPLICIT NONE
4106!       
4107!-- Input and output variables
4108    REAL(wp), INTENT(in) ::  pcnh3    !< ammonia concentration (#/m3)
4109    REAL(wp), INTENT(in) ::  pcocnv   !< conc. of non-volatile OC (#/m3)     
4110    REAL(wp), INTENT(in) ::  pcsa     !< sulphuric acid conc. (#/m3)
4111    REAL(wp), INTENT(in) ::  ppres    !< ambient air pressure (Pa)
4112    REAL(wp), INTENT(in) ::  prh      !< ambient rel. humidity [0-1]       
4113    REAL(wp), INTENT(in) ::  ptemp    !< ambient temperature (K)
4114    REAL(wp), INTENT(in) ::  ptstep   !< time step (s) of SALSA
4115    TYPE(t_section), INTENT(inout) ::  paero(fn2b) !< aerosol properties                                                 
4116    REAL(wp), INTENT(inout) ::  pj3n3(2) !< formation mass rate of molecules
4117                                         !< (molec/m3s) for 1: H2SO4 and
4118                                         !< 2: organic vapour
4119    REAL(wp), INTENT(out) ::  pxocnv !< ratio of non-volatile organic vapours in
4120                                     !< 3nm aerosol particles
4121    REAL(wp), INTENT(out) ::  pxsa   !< ratio of H2SO4 in 3nm aerosol particles
4122!-- Local variables
4123    INTEGER(iwp) ::  iteration
4124    REAL(wp) ::  zbeta(fn2b)  !< transitional correction factor                                         
4125    REAL(wp) ::  zc_h2so4     !< H2SO4 conc. (#/cm3) !UNITS!
4126    REAL(wp) ::  zc_org       !< organic vapour conc. (#/cm3)
4127    REAL(wp) ::  zCoagStot    !< total losses due to coagulation, including
4128                              !< condensation and self-coagulation       
4129    REAL(wp) ::  zcocnv_local !< organic vapour conc. (#/m3)
4130    REAL(wp) ::  zcsink       !< condensational sink (#/m2)       
4131    REAL(wp) ::  zcsa_local   !< H2SO4 conc. (#/m3)       
4132    REAL(wp) ::  zdcrit       !< diameter of critical cluster (m)
4133    REAL(wp) ::  zdelta_vap   !< change of H2SO4 and organic vapour
4134                              !< concentration (#/m3)       
4135    REAL(wp) ::  zdfvap       !< air diffusion coefficient (m2/s)
4136    REAL(wp) ::  zdmean       !< mean diameter of existing particles (m)
4137    REAL(wp) ::  zeta         !< constant: proportional to ratio of CS/GR (m)
4138                              !< (condensation sink / growth rate)                                   
4139    REAL(wp) ::  zgamma       !< proportionality factor ((nm2*m2)/h)                                       
4140    REAL(wp) ::  zGRclust     !< growth rate of formed clusters (nm/h)
4141    REAL(wp) ::  zGRtot       !< total growth rate       
4142    REAL(wp) ::  zj3          !< number conc. of formed 3nm particles (#/m3)       
4143    REAL(wp) ::  zjnuc        !< nucleation rate at ~1nm (#/m3s)
4144    REAL(wp) ::  zKeff        !< effective cogulation coefficient between
4145                              !< freshly nucleated particles       
4146    REAL(wp) ::  zknud(fn2b)  !< particle Knudsen number       
4147    REAL(wp) ::  zkocnv       !< lever: zkocnv=1 --> organic compounds involved
4148                              !< in nucleation   
4149    REAL(wp) ::  zksa         !< lever: zksa=1 --> H2SO4 involved in nucleation
4150    REAL(wp) ::  zlambda      !< parameter for adjusting the growth rate due to
4151                              !< self-coagulation                                 
4152    REAL(wp) ::  zmfp         !< mean free path of condesing vapour(m)                                       
4153    REAL(wp) ::  zmixnh3      !< ammonia mixing ratio (ppt)
4154    REAL(wp) ::  zNnuc        !< number of clusters/particles at the size range
4155                              !< d1-dx (#/m3) 
4156    REAL(wp) ::  znoc         !< number of organic molecules in critical cluster
4157    REAL(wp) ::  znsa         !< number of H2SO4 molecules in critical cluster                                           
4158!
4159!-- Variable determined for the m-parameter
4160    REAL(wp) ::  zCc_2(fn2b) !<
4161    REAL(wp) ::  zCc_c !<
4162    REAL(wp) ::  zCc_x !<
4163    REAL(wp) ::  zCoagS_c !<
4164    REAL(wp) ::  zCoagS_x !<
4165    REAL(wp) ::  zcv_2(fn2b) !<
4166    REAL(wp) ::  zcv_c !<
4167    REAL(wp) ::  zcv_c2(fn2b) !<
4168    REAL(wp) ::  zcv_x !<
4169    REAL(wp) ::  zcv_x2(fn2b) !<
4170    REAL(wp) ::  zDc_2(fn2b) !<
4171    REAL(wp) ::  zDc_c(fn2b) !<
4172    REAL(wp) ::  zDc_c2(fn2b) !<
4173    REAL(wp) ::  zDc_x(fn2b) !<
4174    REAL(wp) ::  zDc_x2(fn2b) !<
4175    REAL(wp) ::  zgammaF_2(fn2b) !<
4176    REAL(wp) ::  zgammaF_c(fn2b) !<
4177    REAL(wp) ::  zgammaF_x(fn2b) !<
4178    REAL(wp) ::  zK_c2(fn2b) !<
4179    REAL(wp) ::  zK_x2(fn2b) !<
4180    REAL(wp) ::  zknud_2(fn2b) !<
4181    REAL(wp) ::  zknud_c !<
4182    REAL(wp) ::  zknud_x !<       
4183    REAL(wp) ::  zm_2(fn2b) !<
4184    REAL(wp) ::  zm_c !<
4185    REAL(wp) ::  zm_para !<
4186    REAL(wp) ::  zm_x !<
4187    REAL(wp) ::  zmyy !<
4188    REAL(wp) ::  zomega_2c(fn2b) !<
4189    REAL(wp) ::  zomega_2x(fn2b) !<
4190    REAL(wp) ::  zomega_c(fn2b) !<
4191    REAL(wp) ::  zomega_x(fn2b) !<
4192    REAL(wp) ::  zRc2(fn2b) !<
4193    REAL(wp) ::  zRx2(fn2b) !<
4194    REAL(wp) ::  zsigma_c2(fn2b) !<
4195    REAL(wp) ::  zsigma_x2(fn2b) !<
4196!
4197!-- 1) Nucleation rate (zjnuc) and diameter of critical cluster (zdcrit)
4198    zjnuc  = 0.0_wp
4199    znsa   = 0.0_wp
4200    znoc   = 0.0_wp
4201    zdcrit = 0.0_wp
4202    zksa   = 0.0_wp
4203    zkocnv = 0.0_wp
4204   
4205    SELECT CASE ( nsnucl )
4206   
4207    CASE(1)   ! Binary H2SO4-H2O nucleation
4208       
4209       zc_h2so4 = pcsa * 1.0E-6_wp   ! sulphuric acid conc. to #/cm3
4210       CALL binnucl( zc_h2so4, ptemp, prh, zjnuc, znsa, znoc, zdcrit,  zksa, &
4211                     zkocnv )     
4212   
4213    CASE(2)   ! Activation type nucleation
4214   
4215       zc_h2so4 = pcsa * 1.0E-6_wp   ! sulphuric acid conc. to #/cm3
4216       CALL binnucl( zc_h2so4, ptemp, prh, zjnuc, znsa,  znoc, zdcrit, zksa,  &
4217                     zkocnv )
4218       CALL actnucl( pcsa, zjnuc, zdcrit, znsa, znoc, zksa, zkocnv, act_coeff )
4219   
4220    CASE(3)   ! Kinetically limited nucleation of (NH4)HSO4 clusters
4221       
4222       zc_h2so4 = pcsa * 1.0E-6_wp   ! sulphuric acid conc. to #/cm3
4223       CALL binnucl( zc_h2so4, ptemp, prh, zjnuc, znsa, znoc, zdcrit, zksa,    &
4224                     zkocnv )
4225
4226       CALL kinnucl( zc_h2so4, zjnuc, zdcrit, znsa, znoc, zksa, zkocnv )
4227   
4228    CASE(4)   ! Ternary H2SO4-H2O-NH3 nucleation
4229   
4230       zmixnh3 = pcnh3 * ptemp * argas / ( ppres * avo )
4231       zc_h2so4 = pcsa * 1.0E-6_wp   ! sulphuric acid conc. to #/cm3
4232       CALL ternucl( zc_h2so4, zmixnh3, ptemp, prh, zjnuc, znsa, znoc, zdcrit, &
4233                     zksa, zkocnv ) 
4234   
4235    CASE(5)   ! Organic nucleation, J~[ORG] or J~[ORG]**2
4236   
4237       zc_org = pcocnv * 1.0E-6_wp   ! conc. of non-volatile OC to #/cm3
4238       zc_h2so4 = pcsa * 1.0E-6_wp   ! sulphuric acid conc. to #/cm3
4239       CALL binnucl( zc_h2so4, ptemp, prh, zjnuc, znsa, znoc, zdcrit, zksa,    &
4240                     zkocnv ) 
4241       CALL orgnucl( pcocnv, zjnuc, zdcrit, znsa, znoc, zksa, zkocnv )
4242   
4243    CASE(6)   ! Sum of H2SO4 and organic activation type nucleation,
4244              ! J~[H2SO4]+[ORG]
4245       
4246       zc_h2so4 = pcsa * 1.0E-6_wp   ! sulphuric acid conc. to #/cm3
4247       CALL binnucl( zc_h2so4, ptemp, prh, zjnuc, znsa, znoc, zdcrit, zksa,    &
4248                     zkocnv ) 
4249       CALL sumnucl( pcsa, pcocnv, zjnuc, zdcrit, znsa, znoc, zksa, zkocnv )
4250
4251           
4252    CASE(7)   ! Heteromolecular nucleation, J~[H2SO4]*[ORG]
4253       
4254       zc_h2so4 = pcsa * 1.0E-6_wp   ! sulphuric acid conc. to #/cm3
4255       zc_org = pcocnv * 1.0E-6_wp   ! conc. of non-volatile OC to #/cm3
4256       CALL binnucl( zc_h2so4, ptemp, prh, zjnuc, znsa, znoc, zdcrit, zksa,    &
4257                     zkocnv ) 
4258       CALL hetnucl( zc_h2so4, zc_org, zjnuc, zdcrit, znsa, znoc, zksa, zkocnv )
4259   
4260    CASE(8)   ! Homomolecular nucleation of H2SO4 and heteromolecular
4261              ! nucleation of H2SO4 and organic vapour,
4262              ! J~[H2SO4]**2 + [H2SO4]*[ORG] (EUCAARI project)
4263       zc_h2so4 = pcsa * 1.0E-6_wp   ! sulphuric acid conc. to #/cm3
4264       zc_org = pcocnv * 1.0E-6_wp   ! conc. of non-volatile OC to #/cm3
4265       CALL binnucl( zc_h2so4, ptemp, prh, zjnuc, znsa, znoc, zdcrit, zksa,    &
4266                     zkocnv ) 
4267       CALL SAnucl( zc_h2so4, zc_org, zjnuc, zdcrit, znsa, znoc, zksa, zkocnv )
4268   
4269    CASE(9)   ! Homomolecular nucleation of H2SO4 and organic vapour and
4270              ! heteromolecular nucleation of H2SO4 and organic vapour,
4271              ! J~[H2SO4]**2 + [H2SO4]*[ORG]+[ORG]**2 (EUCAARI project)
4272   
4273       zc_h2so4 = pcsa * 1.0E-6_wp   ! sulphuric acid conc. to #/cm3
4274       zc_org = pcocnv * 1.0E-6_wp   ! conc. of non-volatile OC to #/cm3
4275       CALL binnucl( zc_h2so4, ptemp, prh, zjnuc, znsa, znoc, zdcrit, zksa,    &
4276                     zkocnv ) 
4277
4278       CALL SAORGnucl( zc_h2so4, zc_org, zjnuc, zdcrit, znsa, znoc, zksa,      &
4279                       zkocnv )
4280    END SELECT
4281   
4282    zcsa_local = pcsa
4283    zcocnv_local = pcocnv
4284!
4285!-- 2) Change of particle and gas concentrations due to nucleation
4286!         
4287!-- 2.1) Check that there is enough H2SO4 and organic vapour to produce the
4288!--      nucleation 
4289    IF ( nsnucl <= 4 )  THEN 
4290!--    If the chosen nucleation scheme is 1-4, nucleation occurs only due to
4291!--    H2SO4. All of the total vapour concentration that is taking part to the
4292!--    nucleation is there for sulphuric acid (sa = H2SO4) and non-volatile
4293!--    organic vapour is zero.
4294       pxsa   = 1.0_wp   ! ratio of sulphuric acid in 3nm particles
4295       pxocnv = 0.0_wp   ! ratio of non-volatile origanic vapour
4296                                ! in 3nm particles
4297    ELSEIF ( nsnucl > 4 )  THEN
4298!--    If the chosen nucleation scheme is 5-9, nucleation occurs due to organic
4299!--    vapour or the combination of organic vapour and H2SO4. The number of
4300!--    needed molecules depends on the chosen nucleation type and it has an
4301!--    effect also on the minimum ratio of the molecules present.
4302       IF ( pcsa * znsa + pcocnv * znoc < 1.E-14_wp )  THEN
4303          pxsa   = 0.0_wp
4304          pxocnv = 0.0_wp             
4305       ELSE
4306          pxsa   = pcsa * znsa / ( pcsa * znsa + pcocnv * znoc ) 
4307          pxocnv = pcocnv * znoc / ( pcsa * znsa + pcocnv * znoc )
4308       ENDIF 
4309    ENDIF
4310!   
4311!-- The change in total vapour concentration is the sum of the concentrations
4312!-- of the vapours taking part to the nucleation (depends on the chosen
4313!-- nucleation scheme)
4314    zdelta_vap = MIN( zjnuc * ( znoc + znsa ), ( pcocnv * zkocnv + pcsa *      &
4315                      zksa ) / ptstep ) 
4316!                     
4317!-- Nucleation rate J at ~1nm (#/m3s)                           
4318    zjnuc = zdelta_vap / ( znoc + znsa )
4319!   
4320!-- H2SO4 concentration after nucleation in #/m3           
4321    zcsa_local = MAX( 1.0_wp, pcsa - zdelta_vap * pxsa ) 
4322!   
4323!-- Non-volative organic vapour concentration after nucleation (#/m3)
4324    zcocnv_local = MAX( 1.0_wp, pcocnv - zdelta_vap * pxocnv )
4325!
4326!-- 2.2) Formation rate of 3 nm particles (Kerminen & Kulmala, 2002)
4327!
4328!-- 2.2.1) Growth rate of clusters formed by H2SO4
4329!
4330!-- GR = 3.0e-15 / dens_clus * sum( molecspeed * molarmass * conc )
4331
4332!-- dens_clus  = density of the clusters (here 1830 kg/m3)
4333!-- molarmass  = molar mass of condensing species (here 98.08 g/mol)
4334!-- conc       = concentration of condensing species [#/m3]
4335!-- molecspeed = molecular speed of condensing species [m/s]
4336!--            = sqrt( 8.0 * R * ptemp / ( pi * molarmass ) )
4337!-- (Seinfeld & Pandis, 1998)
4338!
4339!-- Growth rate by H2SO4 and organic vapour in nm/h (Eq. 21)
4340    zGRclust = 2.3623E-15_wp * SQRT( ptemp ) * ( zcsa_local + zcocnv_local )
4341!   
4342!-- 2.2.2) Condensational sink of pre-existing particle population
4343!
4344!-- Diffusion coefficient (m2/s)
4345    zdfvap = 5.1111E-10_wp * ptemp ** 1.75_wp * ( p_0 + 1325.0_wp ) / ppres
4346!-- Mean free path of condensing vapour (m) (Jacobson (2005), Eq. 15.25 and
4347!-- 16.29)
4348    zmfp = 3.0_wp * zdfvap * SQRT( pi * amh2so4 / ( 8.0_wp * argas * ptemp ) )
4349!-- Knudsen number           
4350    zknud = 2.0_wp * zmfp / ( paero(:)%dwet + d_sa )                     
4351!-- Transitional regime correction factor (zbeta) according to Fuchs and
4352!-- Sutugin (1971), In: Hidy et al. (ed.), Topics in current  aerosol research,
4353!-- Pergamon. (Eq. 4 in Kerminen and Kulmala, 2002)
4354    zbeta = ( zknud + 1.0_wp) / ( 0.377_wp * zknud + 1.0_wp + 4.0_wp /         &
4355            ( 3.0_wp * massacc ) * ( zknud + zknud ** 2 ) ) 
4356!-- Condensational sink (#/m2) (Eq. 3)
4357    zcsink = SUM( paero(:)%dwet * zbeta * paero(:)%numc )
4358!
4359!-- Parameterised formation rate of detectable 3 nm particles (i.e. J3)
4360    IF ( nj3 == 1 )  THEN   ! Kerminen and Kulmala (2002)
4361!--    2.2.3) Parameterised formation rate of detectable 3 nm particles
4362!--    Constants needed for the parameterisation:
4363!--    dapp = 3 nm and dens_nuc = 1830 kg/m3
4364       IF ( zcsink < 1.0E-30_wp )  THEN
4365          zeta = 0._dp
4366       ELSE
4367!--       Mean diameter of backgroud population (nm)
4368          zdmean = 1.0_wp / SUM( paero(:)%numc ) * SUM( paero(:)%numc *        &
4369                   paero(:)%dwet ) * 1.0E+9_wp
4370!--       Proportionality factor (nm2*m2/h) (Eq. 22)
4371          zgamma = 0.23_wp * ( zdcrit * 1.0E+9_wp ) ** 0.2_wp * ( zdmean /     &
4372                 150.0_wp ) ** 0.048_wp * ( ptemp / 293.0_wp ) ** ( -0.75_wp ) &
4373                 * ( arhoh2so4 / 1000.0_wp ) ** ( -0.33_wp )
4374!--       Factor eta (nm) (Eq. 11)
4375          zeta = MIN( zgamma * zcsink / zGRclust, zdcrit * 1.0E11_wp ) 
4376       ENDIF
4377!       
4378!--    Number conc. of clusters surviving to 3 nm in a time step (#/m3) (Eq.14)
4379       zj3 = zjnuc * EXP( MIN( 0.0_wp, zeta / 3.0_wp - zeta /                  &
4380                               ( zdcrit * 1.0E9_wp ) ) )                   
4381
4382    ELSEIF ( nj3 > 1 )  THEN
4383!--    Defining the value for zm_para. The growth is investigated between
4384!--    [d1,reglim(1)] = [zdcrit,3nm]   
4385!--    m = LOG( CoagS_dx / CoagX_zdcrit ) / LOG( reglim / zdcrit )
4386!--    (Lehtinen et al. 2007, Eq. 5)
4387!--    The steps for the coagulation sink for reglim = 3nm and zdcrit ~= 1nm are
4388!--    explained in article of Kulmala et al. (2001). The particles of diameter
4389!--    zdcrit ~1.14 nm  and reglim = 3nm are both in turn the "number 1"
4390!--    variables (Kulmala et al. 2001).             
4391!--    c = critical (1nm), x = 3nm, 2 = wet or mean droplet
4392!--    Sum of the radii, R12 = R1 + zR2 (m) of two particles 1 and 2
4393       zRc2 = zdcrit / 2.0_wp + paero(:)%dwet / 2.0_wp
4394       zRx2 = reglim(1) / 2.0_wp + paero(:)%dwet / 2.0_wp
4395!       
4396!--    The mass of particle (kg) (comes only from H2SO4)
4397       zm_c = 4.0_wp / 3.0_wp * pi * ( zdcrit / 2.0_wp ) ** 3.0_wp * arhoh2so4                     
4398       zm_x = 4.0_wp / 3.0_wp * pi * ( reglim(1) / 2.0_wp ) ** 3.0_wp *        &
4399              arhoh2so4                 
4400       zm_2 = 4.0_wp / 3.0_wp * pi * ( paero(:)%dwet / 2.0_wp )** 3.0_wp *     &
4401              arhoh2so4
4402!             
4403!--    Mean relative thermal velocity between the particles (m/s)
4404       zcv_c = SQRT( 8.0_wp * abo * ptemp / ( pi * zm_c ) )
4405       zcv_x = SQRT( 8.0_wp * abo * ptemp / ( pi * zm_x ) )
4406       zcv_2 = SQRT( 8.0_wp * abo * ptemp / ( pi * zm_2 ) )
4407!       
4408!--    Average velocity after coagulation               
4409       zcv_c2 = SQRT( zcv_c ** 2.0_wp + zcv_2 ** 2.0_wp )
4410       zcv_x2 = SQRT( zcv_x ** 2.0_wp + zcv_2 ** 2.0_wp )
4411!       
4412!--    Knudsen number (zmfp = mean free path of condensing vapour)
4413       zknud_c = 2.0_wp * zmfp / zdcrit
4414       zknud_x = 2.0_wp * zmfp / reglim(1)
4415       zknud_2 = MAX( 0.0_wp, 2.0_wp * zmfp / paero(:)%dwet )
4416!
4417!--    Cunningham correction factor               
4418       zCc_c = 1.0_wp + zknud_c * ( 1.142_wp + 0.558_wp *                      &
4419               EXP( -0.999_wp / zknud_c ) ) 
4420       zCc_x = 1.0_wp + zknud_x * ( 1.142_wp + 0.558_wp *                      &
4421               EXP( -0.999_wp / zknud_x ) )
4422       zCc_2 = 1.0_wp + zknud_2 * ( 1.142_wp + 0.558_wp *                      &
4423               EXP( -0.999_wp / zknud_2 ) )
4424!                     
4425!--    Gas dynamic viscosity (N*s/m2).
4426!--    Viscocity(air @20C) = 1.81e-5_dp N/m2 *s (Hinds, p. 25)                     
4427       zmyy = 1.81E-5_wp * ( ptemp / 293.0_wp) ** ( 0.74_wp ) 
4428!       
4429!--    Particle diffusion coefficient (m2/s)               
4430       zDc_c = abo * ptemp * zCc_c / ( 3.0_wp * pi * zmyy * zdcrit ) 
4431       zDc_x = abo * ptemp * zCc_x / ( 3.0_wp * pi * zmyy * reglim(1) )
4432       zDc_2 = abo * ptemp * zCc_2 / ( 3.0_wp * pi * zmyy * paero(:)%dwet )
4433!       
4434!--    D12 = D1+D2 (Seinfield and Pandis, 2nd ed. Eq. 13.38)
4435       zDc_c2 = zDc_c + zDc_2   
4436       zDc_x2 = zDc_x + zDc_2
4437!       
4438!--    zgammaF = 8*D/pi/zcv (m) for calculating zomega
4439       zgammaF_c = 8.0_wp * zDc_c / pi / zcv_c
4440       zgammaF_x = 8.0_wp * zDc_x / pi / zcv_x
4441       zgammaF_2 = 8.0_wp * zDc_2 / pi / zcv_2
4442!       
4443!--    zomega (m) for calculating zsigma             
4444       zomega_c = ( ( zRc2 + zgammaF_c ) ** 3 - ( zRc2 ** 2 +                  &
4445                      zgammaF_c ) ** ( 3.0_wp / 2.0_wp ) ) / ( 3.0_wp *        &
4446                      zRc2 * zgammaF_c ) - zRc2 
4447       zomega_x = ( ( zRx2 + zgammaF_x ) ** 3.0_wp - ( zRx2 ** 2.0_wp +        &
4448                      zgammaF_x ) ** ( 3.0_wp / 2.0_wp ) ) / ( 3.0_wp *        &
4449                      zRx2 * zgammaF_x ) - zRx2
4450       zomega_2c = ( ( zRc2 + zgammaF_2 ) ** 3.0_wp - ( zRc2 ** 2.0_wp +       &
4451                       zgammaF_2 ) ** ( 3.0_wp / 2.0_wp ) ) / ( 3.0_wp *       &
4452                       zRc2 * zgammaF_2 ) - zRc2
4453       zomega_2x = ( ( zRx2 + zgammaF_2 ) ** 3.0_wp - ( zRx2 ** 2.0_wp +       &
4454                       zgammaF_2 ) ** ( 3.0_wp / 2.0_wp ) ) / ( 3.0_wp *       &
4455                       zRx2 * zgammaF_2 ) - zRx2
4456!                       
4457!--    The distance (m) at which the two fluxes are matched (condensation and
4458!--    coagulation sinks?)           
4459       zsigma_c2 = SQRT( zomega_c ** 2.0_wp + zomega_2c ** 2.0_wp ) 
4460       zsigma_x2 = SQRT( zomega_x ** 2.0_wp + zomega_2x ** 2.0_wp ) 
4461!       
4462!--    Coagulation coefficient in the continuum regime (m*m2/s)
4463       zK_c2 = 4.0_wp * pi * zRc2 * zDc_c2 / ( zRc2 / ( zRc2 + zsigma_c2 ) +   &
4464               4.0_wp * zDc_c2 / ( zcv_c2 * zRc2 ) ) 
4465       zK_x2 = 4.0_wp * pi * zRx2 * zDc_x2 / ( zRx2 / ( zRx2 + zsigma_x2 ) +   &
4466               4.0_wp * zDc_x2 / ( zcv_x2 * zRx2 ) )
4467!               
4468!--    Coagulation sink (1/s)
4469       zCoagS_c = MAX( 1.0E-20_wp, SUM( zK_c2 * paero(:)%numc ) )         
4470       zCoagS_x = MAX( 1.0E-20_wp, SUM( zK_x2 * paero(:)%numc ) ) 
4471!       
4472!--    Parameter m for calculating the coagulation sink onto background
4473!--    particles (Eq. 5&6 in Lehtinen et al. 2007)             
4474       zm_para = LOG( zCoagS_x / zCoagS_c ) / LOG( reglim(1) / zdcrit )
4475!       
4476!--    Parameter gamma for calculating the formation rate J of particles having
4477!--    a diameter zdcrit < d < reglim(1) (Anttila et al. 2010, eq. 5)
4478       zgamma = ( ( ( reglim(1) / zdcrit ) ** ( zm_para + 1.0_wp ) ) - 1.0_wp )&
4479                / ( zm_para + 1.0_wp )     
4480               
4481       IF ( nj3 == 2 )  THEN   ! Coagulation sink
4482!       
4483!--       Formation rate J before iteration (#/m3s)               
4484          zj3 = zjnuc * EXP( MIN( 0.0_wp, -zgamma * zdcrit * zCoagS_c /        &
4485                ( zGRclust * 1.0E-9_wp / ( 60.0_wp ** 2.0_wp ) ) ) )
4486               
4487       ELSEIF ( nj3 == 3 )  THEN  ! Coagulation sink and self-coag.
4488!--       IF polluted air... then the self-coagulation becomes important.
4489!--       Self-coagulation of small particles < 3 nm.
4490!
4491!--       "Effective" coagulation coefficient between freshly-nucleated
4492!--       particles:
4493          zKeff = 5.0E-16_wp   ! cm3/s
4494!         
4495!--       zlambda parameter for "adjusting" the growth rate due to the
4496!--       self-coagulation
4497          zlambda = 6.0_wp 
4498          IF ( reglim(1) >= 10.0E-9_wp )  THEN   ! for particles >10 nm:
4499             zKeff   = 5.0E-17_wp
4500             zlambda = 3.0_wp
4501          ENDIF
4502!         
4503!--       Initial values for coagulation sink and growth rate  (m/s)
4504          zCoagStot = zCoagS_c
4505          zGRtot = zGRclust * 1.0E-9_wp / 60.0_wp ** 2.0_wp 
4506!         
4507!--       Number of clusters/particles at the size range [d1,dx] (#/m3):
4508          zNnuc = zjnuc / zCoagStot !< Initial guess
4509!         
4510!--       Coagulation sink and growth rate due to self-coagulation:
4511          DO  iteration = 1, 5
4512             zCoagStot = zCoagS_c + zKeff * zNnuc * 1.0E-6_wp   ! (1/s) 
4513             zGRtot = zGRclust * 1.0E-9_wp / ( 3600.0_wp ) +  1.5708E-6_wp *   &
4514                      zlambda * zdcrit ** 3.0_wp * ( zNnuc * 1.0E-6_wp ) *     &
4515                      zcv_c * avo * 1.0E-9_wp / 3600.0_wp 
4516             zeta = - zCoagStot / ( ( zm_para + 1.0_wp ) * zGRtot * ( zdcrit **&
4517                      zm_para ) )   ! Eq. 7b (Anttila)
4518             zNnuc =  zNnuc_tayl( zdcrit, reglim(1), zm_para, zjnuc, zeta,     &
4519                      zGRtot )
4520          ENDDO
4521!         
4522!--       Calculate the final values with new zNnuc:   
4523          zCoagStot = zCoagS_c + zKeff * zNnuc * 1.0E-6_wp   ! (1/s)
4524          zGRtot = zGRclust * 1.0E-9_wp / 3600.0_wp + 1.5708E-6_wp *  zlambda  &
4525                   * zdcrit ** 3.0_wp * ( zNnuc * 1.0E-6_wp ) * zcv_c * avo *  &
4526                   1.0E-9_wp / 3600.0_wp !< (m/s)
4527          zj3 = zjnuc * EXP( MIN( 0.0_wp, -zgamma * zdcrit * zCoagStot /       &
4528                zGRtot ) )   ! (Eq. 5a) (#/m3s)
4529               
4530       ENDIF
4531       
4532    ENDIF
4533!-- If J3 very small (< 1 #/cm3), neglect particle formation. In real atmosphere
4534!-- this would mean that clusters form but coagulate to pre-existing particles
4535!-- who gain sulphate. Since CoagS ~ CS (4piD*CS'), we do *not* update H2SO4
4536!-- concentration here but let condensation take care of it.
4537!-- Formation mass rate of molecules (molec/m3s) for 1: H2SO4 and 2: organic
4538!-- vapour
4539    pj3n3(1) = zj3 * n3 * pxsa
4540    pj3n3(2) = zj3 * n3 * pxocnv
4541                                 
4542                         
4543 END SUBROUTINE nucleation
4544
4545!------------------------------------------------------------------------------!
4546! Description:
4547! ------------
4548!> Calculate the nucleation rate and the size of critical clusters assuming
4549!> binary nucleation.
4550!> Parametrisation according to Vehkamaki et al. (2002), J. Geophys. Res.,
4551!> 107(D22), 4622. Called from subroutine nucleation.
4552!------------------------------------------------------------------------------!
4553 SUBROUTINE binnucl( pc_sa, ptemp, prh, pnuc_rate, pn_crit_sa, pn_crit_ocnv,   &
4554                     pd_crit, pk_sa, pk_ocnv )
4555                   
4556    IMPLICIT NONE
4557!       
4558!-- Input and output variables       
4559    REAL(wp), INTENT(in) ::   pc_sa        !< H2SO4 conc. (#/cm3)
4560    REAL(wp), INTENT(in) ::   prh          !< relative humidity [0-1]       
4561    REAL(wp), INTENT(in) ::   ptemp        !< ambient temperature (K)
4562    REAL(wp), INTENT(out) ::  pnuc_rate    !< nucleation rate (#/(m3 s))
4563    REAL(wp), INTENT(out) ::  pn_crit_sa   !< number of H2SO4 molecules in
4564                                           !< cluster (#)
4565    REAL(wp), INTENT(out) ::  pn_crit_ocnv !< number of organic molecules in
4566                                           !< cluster (#)
4567    REAL(wp), INTENT(out) ::  pd_crit      !< diameter of critical cluster (m)
4568    REAL(wp), INTENT(out) ::  pk_sa        !< Lever: if pk_sa = 1, H2SO4 is
4569                                           !< involved in nucleation.
4570    REAL(wp), INTENT(out) ::  pk_ocnv      !< Lever: if pk_ocnv = 1, organic
4571                                           !< compounds are involved in
4572                                           !< nucleation.
4573!-- Local variables
4574    REAL(wp) ::  zx    !< mole fraction of sulphate in critical cluster
4575    REAL(wp) ::  zntot !< number of molecules in critical cluster
4576    REAL(wp) ::  zt    !< temperature
4577    REAL(wp) ::  zpcsa !< sulfuric acid concentration
4578    REAL(wp) ::  zrh   !< relative humidity
4579    REAL(wp) ::  zma   !<
4580    REAL(wp) ::  zmw   !<
4581    REAL(wp) ::  zxmass!<
4582    REAL(wp) ::  za    !<
4583    REAL(wp) ::  zb    !<
4584    REAL(wp) ::  zc    !<
4585    REAL(wp) ::  zroo  !<
4586    REAL(wp) ::  zm1   !<
4587    REAL(wp) ::  zm2   !<
4588    REAL(wp) ::  zv1   !<
4589    REAL(wp) ::  zv2   !<
4590    REAL(wp) ::  zcoll !<
4591   
4592    pnuc_rate = 0.0_wp
4593    pd_crit   = 1.0E-9_wp
4594
4595!             
4596!-- 1) Checking that we are in the validity range of the parameterization 
4597    zt    = MAX( ptemp, 190.15_wp )
4598    zt    = MIN( zt,    300.15_wp )
4599    zpcsa = MAX( pc_sa, 1.0E4_wp  )
4600    zpcsa = MIN( zpcsa, 1.0E11_wp ) 
4601    zrh   = MAX( prh,   0.0001_wp )
4602    zrh   = MIN( zrh,   1.0_wp    )
4603!               
4604!-- 2) Mole fraction of sulphate in a critical cluster (Eq. 11)
4605    zx = 0.7409967177282139_wp                                           &
4606         - 0.002663785665140117_wp * zt                                  &
4607         + 0.002010478847383187_wp * LOG( zrh )                          &
4608         - 0.0001832894131464668_wp* zt * LOG( zrh )                     &
4609         + 0.001574072538464286_wp * LOG( zrh ) ** 2                     &
4610         - 0.00001790589121766952_wp * zt * LOG( zrh ) ** 2              &
4611         + 0.0001844027436573778_wp * LOG( zrh ) ** 3                    &
4612         - 1.503452308794887E-6_wp * zt * LOG( zrh ) ** 3                &
4613         - 0.003499978417957668_wp * LOG( zpcsa )                        &
4614         + 0.0000504021689382576_wp * zt * LOG( zpcsa )
4615!                   
4616!-- 3) Nucleation rate (Eq. 12)
4617    pnuc_rate = 0.1430901615568665_wp                                    &
4618        + 2.219563673425199_wp * zt                                      &
4619        - 0.02739106114964264_wp * zt ** 2                               &
4620        + 0.00007228107239317088_wp * zt ** 3                            &
4621        + 5.91822263375044_wp / zx                                       &
4622        + 0.1174886643003278_wp * LOG( zrh )                             &
4623        + 0.4625315047693772_wp * zt * LOG( zrh )                        &
4624        - 0.01180591129059253_wp * zt ** 2 * LOG( zrh )                  &
4625        + 0.0000404196487152575_wp * zt ** 3 * LOG( zrh )                &
4626        + ( 15.79628615047088_wp * LOG( zrh ) ) / zx                     &
4627        - 0.215553951893509_wp * LOG( zrh ) ** 2                         &
4628        - 0.0810269192332194_wp * zt * LOG( zrh ) ** 2                   &
4629        + 0.001435808434184642_wp * zt ** 2 * LOG( zrh ) ** 2            &
4630        - 4.775796947178588E-6_wp * zt ** 3 * LOG( zrh ) ** 2            &
4631        - (2.912974063702185_wp * LOG( zrh ) ** 2 ) / zx                 &
4632        - 3.588557942822751_wp * LOG( zrh ) ** 3                         &
4633        + 0.04950795302831703_wp * zt * LOG( zrh ) ** 3                  &
4634        - 0.0002138195118737068_wp * zt ** 2 * LOG( zrh ) ** 3           &
4635        + 3.108005107949533E-7_wp * zt ** 3 * LOG( zrh ) ** 3            &
4636        - ( 0.02933332747098296_wp * LOG( zrh ) ** 3 ) / zx              &
4637        + 1.145983818561277_wp * LOG( zpcsa )                            &
4638        - 0.6007956227856778_wp * zt * LOG( zpcsa )                      &
4639        + 0.00864244733283759_wp * zt ** 2 * LOG( zpcsa )                &
4640        - 0.00002289467254710888_wp * zt ** 3 * LOG( zpcsa )             &
4641        - ( 8.44984513869014_wp * LOG( zpcsa ) ) / zx                    &
4642        + 2.158548369286559_wp * LOG( zrh ) * LOG( zpcsa )               &
4643        + 0.0808121412840917_wp * zt * LOG( zrh ) * LOG( zpcsa )         &
4644        - 0.0004073815255395214_wp * zt ** 2 * LOG( zrh ) * LOG( zpcsa ) &
4645        - 4.019572560156515E-7_wp * zt ** 3 * LOG( zrh ) * LOG( zpcsa )  & 
4646        + ( 0.7213255852557236_wp * LOG( zrh ) * LOG( zpcsa ) ) / zx     &
4647        + 1.62409850488771_wp * LOG( zrh ) ** 2 * LOG( zpcsa )           &
4648        - 0.01601062035325362_wp * zt * LOG( zrh ) ** 2 * LOG( zpcsa )   &
4649        + 0.00003771238979714162_wp*zt**2* LOG( zrh )**2 * LOG( zpcsa )  &
4650        + 3.217942606371182E-8_wp * zt**3 * LOG( zrh )**2 * LOG( zpcsa ) &
4651        - (0.01132550810022116_wp * LOG( zrh )**2 * LOG( zpcsa ) ) / zx  &
4652        + 9.71681713056504_wp * LOG( zpcsa ) ** 2                        &
4653        - 0.1150478558347306_wp * zt * LOG( zpcsa ) ** 2                 &
4654        + 0.0001570982486038294_wp * zt ** 2 * LOG( zpcsa ) ** 2         &
4655        + 4.009144680125015E-7_wp * zt ** 3 * LOG( zpcsa ) ** 2          &
4656        + ( 0.7118597859976135_wp * LOG( zpcsa ) ** 2 ) / zx             &
4657        - 1.056105824379897_wp * LOG( zrh ) * LOG( zpcsa ) ** 2          &
4658        + 0.00903377584628419_wp * zt * LOG( zrh ) * LOG( zpcsa )**2     &
4659        - 0.00001984167387090606_wp*zt**2*LOG( zrh )*LOG( zpcsa )**2     &
4660        + 2.460478196482179E-8_wp * zt**3 * LOG( zrh ) * LOG( zpcsa )**2 &
4661        - ( 0.05790872906645181_wp * LOG( zrh ) * LOG( zpcsa )**2 ) / zx &
4662        - 0.1487119673397459_wp * LOG( zpcsa ) ** 3                      &
4663        + 0.002835082097822667_wp * zt * LOG( zpcsa ) ** 3               &
4664        - 9.24618825471694E-6_wp * zt ** 2 * LOG( zpcsa ) ** 3           &
4665        + 5.004267665960894E-9_wp * zt ** 3 * LOG( zpcsa ) ** 3          &
4666        - ( 0.01270805101481648_wp * LOG( zpcsa ) ** 3 ) / zx
4667!           
4668!-- Nucleation rate in #/(cm3 s)
4669    pnuc_rate = EXP( pnuc_rate ) 
4670!       
4671!-- Check the validity of parameterization
4672    IF ( pnuc_rate < 1.0E-7_wp )  THEN
4673       pnuc_rate = 0.0_wp
4674       pd_crit   = 1.0E-9_wp
4675    ENDIF
4676!               
4677!-- 4) Total number of molecules in the critical cluster (Eq. 13)
4678    zntot = - 0.002954125078716302_wp                                    &
4679      - 0.0976834264241286_wp * zt                                       &
4680      + 0.001024847927067835_wp * zt ** 2                                &
4681      - 2.186459697726116E-6_wp * zt ** 3                                &
4682      - 0.1017165718716887_wp / zx                                       &
4683      - 0.002050640345231486_wp * LOG( zrh )                             &
4684      - 0.007585041382707174_wp * zt * LOG( zrh )                        &
4685      + 0.0001926539658089536_wp * zt ** 2 * LOG( zrh )                  &
4686      - 6.70429719683894E-7_wp * zt ** 3 * LOG( zrh )                    &
4687      - ( 0.2557744774673163_wp * LOG( zrh ) ) / zx                      &
4688      + 0.003223076552477191_wp * LOG( zrh ) ** 2                        &
4689      + 0.000852636632240633_wp * zt * LOG( zrh ) ** 2                   &
4690      - 0.00001547571354871789_wp * zt ** 2 * LOG( zrh ) ** 2            &
4691      + 5.666608424980593E-8_wp * zt ** 3 * LOG( zrh ) ** 2              &
4692      + ( 0.03384437400744206_wp * LOG( zrh ) ** 2 ) / zx                &
4693      + 0.04743226764572505_wp * LOG( zrh ) ** 3                         &
4694      - 0.0006251042204583412_wp * zt * LOG( zrh ) ** 3                  &
4695      + 2.650663328519478E-6_wp * zt ** 2 * LOG( zrh ) ** 3              &
4696      - 3.674710848763778E-9_wp * zt ** 3 * LOG( zrh ) ** 3              &
4697      - ( 0.0002672510825259393_wp * LOG( zrh ) ** 3 ) / zx              &
4698      - 0.01252108546759328_wp * LOG( zpcsa )                            &
4699      + 0.005806550506277202_wp * zt * LOG( zpcsa )                      &
4700      - 0.0001016735312443444_wp * zt ** 2 * LOG( zpcsa )                &
4701      + 2.881946187214505E-7_wp * zt ** 3 * LOG( zpcsa )                 &
4702      + ( 0.0942243379396279_wp * LOG( zpcsa ) ) / zx                    &
4703      - 0.0385459592773097_wp * LOG( zrh ) * LOG( zpcsa )                &
4704      - 0.0006723156277391984_wp * zt * LOG( zrh ) * LOG( zpcsa )        &
4705      + 2.602884877659698E-6_wp * zt ** 2 * LOG( zrh ) * LOG( zpcsa )    &
4706      + 1.194163699688297E-8_wp * zt ** 3 * LOG( zrh ) * LOG( zpcsa )    &
4707      - ( 0.00851515345806281_wp * LOG( zrh ) * LOG( zpcsa ) ) / zx      &
4708      - 0.01837488495738111_wp * LOG( zrh ) ** 2 * LOG( zpcsa )          &
4709      + 0.0001720723574407498_wp * zt * LOG( zrh ) ** 2 * LOG( zpcsa )   &
4710      - 3.717657974086814E-7_wp * zt**2 * LOG( zrh )**2 * LOG( zpcsa )   &
4711      - 5.148746022615196E-10_wp * zt**3 * LOG( zrh )**2 * LOG( zpcsa )  &
4712      + ( 0.0002686602132926594_wp * LOG(zrh)**2 * LOG(zpcsa) ) / zx     &
4713      - 0.06199739728812199_wp * LOG( zpcsa ) ** 2                       &
4714      + 0.000906958053583576_wp * zt * LOG( zpcsa ) ** 2                 &
4715      - 9.11727926129757E-7_wp * zt ** 2 * LOG( zpcsa ) ** 2             &
4716      - 5.367963396508457E-9_wp * zt ** 3 * LOG( zpcsa ) ** 2            &
4717      - ( 0.007742343393937707_wp * LOG( zpcsa ) ** 2 ) / zx             &
4718      + 0.0121827103101659_wp * LOG( zrh ) * LOG( zpcsa ) ** 2           &
4719      - 0.0001066499571188091_wp * zt * LOG( zrh ) * LOG( zpcsa ) ** 2   &
4720      + 2.534598655067518E-7_wp * zt**2 * LOG( zrh ) * LOG( zpcsa )**2   &
4721      - 3.635186504599571E-10_wp * zt**3 * LOG( zrh ) * LOG( zpcsa )**2  &
4722      + ( 0.0006100650851863252_wp * LOG( zrh ) * LOG( zpcsa ) **2 )/ zx &
4723      + 0.0003201836700403512_wp * LOG( zpcsa ) ** 3                     &
4724      - 0.0000174761713262546_wp * zt * LOG( zpcsa ) ** 3                &
4725      + 6.065037668052182E-8_wp * zt ** 2 * LOG( zpcsa ) ** 3            &
4726      - 1.421771723004557E-11_wp * zt ** 3 * LOG( zpcsa ) ** 3           &
4727      + ( 0.0001357509859501723_wp * LOG( zpcsa ) ** 3 ) / zx
4728    zntot = EXP( zntot )  ! in #
4729!
4730!-- 5) Size of the critical cluster pd_crit (m) (diameter) (Eq. 14)
4731    pn_crit_sa = zx * zntot
4732    pd_crit    = 2.0E-9_wp * EXP( -1.6524245_wp + 0.42316402_wp  * zx +        &
4733                 0.33466487_wp * LOG( zntot ) )
4734!
4735!-- 6) Organic compounds not involved when binary nucleation is assumed
4736    pn_crit_ocnv = 0.0_wp   ! number of organic molecules
4737    pk_sa        = 1.0_wp   ! if = 1, H2SO4 involved in nucleation
4738    pk_ocnv      = 0.0_wp   ! if = 1, organic compounds involved
4739!               
4740!-- Set nucleation rate to collision rate               
4741    IF ( pn_crit_sa < 4.0_wp ) THEN
4742!       
4743!--    Volumes of the colliding objects
4744       zma    = 96.0_wp   ! molar mass of SO4 in g/mol
4745       zmw    = 18.0_wp   ! molar mass of water in g/mol
4746       zxmass = 1.0_wp    ! mass fraction of H2SO4
4747       za = 0.7681724_wp + zxmass * ( 2.1847140_wp + zxmass * (     &
4748            7.1630022_wp + zxmass * ( -44.31447_wp + zxmass * (     &
4749            88.75606 + zxmass * ( -75.73729_wp + zxmass *           &
4750            23.43228_wp ) ) ) ) )<