!> @file salsa_mod.f90 !--------------------------------------------------------------------------------! ! This file is part of PALM-4U. ! ! PALM-4U is free software: you can redistribute it and/or modify it under the ! terms of the GNU General Public License as published by the Free Software ! Foundation, either version 3 of the License, or (at your option) any later ! version. ! ! PALM-4U is distributed in the hope that it will be useful, but WITHOUT ANY ! WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR ! A PARTICULAR PURPOSE. See the GNU General Public License for more details. ! ! You should have received a copy of the GNU General Public License along with ! PALM. If not, see . ! ! Copyright 2018-2018 University of Helsinki ! Copyright 1997-2018 Leibniz Universitaet Hannover !--------------------------------------------------------------------------------! ! ! Current revisions: ! ----------------- ! ! ! Former revisions: ! ----------------- ! $Id: salsa_mod.f90 3637 2018-12-20 01:51:36Z kanani $ ! Implementation of the PALM module interface ! ! 3636 2018-12-19 13:48:34Z raasch ! nopointer option removed ! ! 3630 2018-12-17 11:04:17Z knoop ! - Moved the control parameter "salsa" from salsa_mod.f90 to control_parameters ! - Updated salsa_rrd_local and salsa_wrd_local ! - Add target attribute ! - Revise initialization in case of restarts ! - Revise masked data output ! ! 3582 2018-11-29 19:16:36Z suehring ! missing comma separator inserted ! ! 3483 2018-11-02 14:19:26Z raasch ! bugfix: directives added to allow compilation without netCDF ! ! 3481 2018-11-02 09:14:13Z raasch ! temporary variable cc introduced to circumvent a possible Intel18 compiler bug ! related to contiguous/non-contguous pointer/target attributes ! ! 3473 2018-10-30 20:50:15Z suehring ! NetCDF input routine renamed ! ! 3467 2018-10-30 19:05:21Z suehring ! Initial revision ! ! 3412 2018-10-24 07:25:57Z monakurppa ! ! Authors: ! -------- ! @author Mona Kurppa (University of Helsinki) ! ! ! Description: ! ------------ !> Sectional aerosol module for large scale applications SALSA !> (Kokkola et al., 2008, ACP 8, 2469-2483). Solves the aerosol number and mass !> concentration as well as chemical composition. Includes aerosol dynamic !> processes: nucleation, condensation/evaporation of vapours, coagulation and !> deposition on tree leaves, ground and roofs. !> Implementation is based on formulations implemented in UCLALES-SALSA except !> for deposition which is based on parametrisations by Zhang et al. (2001, !> Atmos. Environ. 35, 549-560) or Petroff&Zhang (2010, Geosci. Model Dev. 3, !> 753-769) !> !> @todo Implement turbulent inflow of aerosols in inflow_turbulence. !> @todo Deposition on subgrid scale vegetation !> @todo Deposition on vegetation calculated by default for deciduous broadleaf !> trees !> @todo Revise masked data output. There is a potential bug in case of !> terrain-following masked output, according to data_output_mask. !> @todo There are now improved interfaces for NetCDF data input which can be !> used instead of get variable etc. !------------------------------------------------------------------------------! MODULE salsa_mod USE basic_constants_and_equations_mod, & ONLY: c_p, g, p_0, pi, r_d USE chemistry_model_mod, & ONLY: chem_species, nspec, nvar, spc_names USE chem_modules, & ONLY: call_chem_at_all_substeps, chem_gasphase_on USE control_parameters USE indices, & ONLY: nbgp, nx, nxl, nxlg, nxr, nxrg, ny, nyn, nyng, nys, nysg, nzb, & nzb_s_inner, nz, nzt, wall_flags_0 USE kinds USE pegrid USE salsa_util_mod IMPLICIT NONE ! !-- SALSA constants: ! !-- Local constants: INTEGER(iwp), PARAMETER :: ngast = 5 !< total number of gaseous tracers: !< 1 = H2SO4, 2 = HNO3, 3 = NH3, !< 4 = OCNV (non-volatile OC), !< 5 = OCSV (semi-volatile) INTEGER(iwp), PARAMETER :: nmod = 7 !< number of modes for initialising !< the aerosol size distribution INTEGER(iwp), PARAMETER :: nreg = 2 !< Number of main size subranges INTEGER(iwp), PARAMETER :: maxspec = 7 !< Max. number of aerosol species ! !-- Universal constants REAL(wp), PARAMETER :: abo = 1.380662E-23_wp !< Boltzmann constant (J/K) REAL(wp), PARAMETER :: alv = 2.260E+6_wp !< latent heat for H2O !< vaporisation (J/kg) REAL(wp), PARAMETER :: alv_d_rv = 4896.96865_wp !< alv / rv REAL(wp), PARAMETER :: am_airmol = 4.8096E-26_wp !< Average mass of one air !< molecule (Jacobson, !< 2005, Eq. 2.3) REAL(wp), PARAMETER :: api6 = 0.5235988_wp !< pi / 6 REAL(wp), PARAMETER :: argas = 8.314409_wp !< Gas constant (J/(mol K)) REAL(wp), PARAMETER :: argas_d_cpd = 8.281283865E-3_wp !< argas per cpd REAL(wp), PARAMETER :: avo = 6.02214E+23_wp !< Avogadro constant (1/mol) REAL(wp), PARAMETER :: d_sa = 5.539376964394570E-10_wp !< diameter of !< condensing sulphuric !< acid molecule (m) REAL(wp), PARAMETER :: for_ppm_to_nconc = 7.243016311E+16_wp !< !< ppm * avo / R (K/(Pa*m3)) REAL(wp), PARAMETER :: epsoc = 0.15_wp !< water uptake of organic !< material REAL(wp), PARAMETER :: mclim = 1.0E-23_wp !< mass concentration min !< limit for aerosols (kg/m3) REAL(wp), PARAMETER :: n3 = 158.79_wp !< Number of H2SO4 molecules in !< 3 nm cluster if d_sa=5.54e-10m REAL(wp), PARAMETER :: nclim = 1.0_wp !< number concentration min limit !< for aerosols and gases (#/m3) REAL(wp), PARAMETER :: surfw0 = 0.073_wp !< surface tension of pure water !< at ~ 293 K (J/m2) REAL(wp), PARAMETER :: vclim = 1.0E-24_wp !< volume concentration min !< limit for aerosols (m3/m3) !-- Molar masses in kg/mol REAL(wp), PARAMETER :: ambc = 12.0E-3_wp !< black carbon (BC) REAL(wp), PARAMETER :: amdair = 28.970E-3_wp !< dry air REAL(wp), PARAMETER :: amdu = 100.E-3_wp !< mineral dust REAL(wp), PARAMETER :: amh2o = 18.0154E-3_wp !< H2O REAL(wp), PARAMETER :: amh2so4 = 98.06E-3_wp !< H2SO4 REAL(wp), PARAMETER :: amhno3 = 63.01E-3_wp !< HNO3 REAL(wp), PARAMETER :: amn2o = 44.013E-3_wp !< N2O REAL(wp), PARAMETER :: amnh3 = 17.031E-3_wp !< NH3 REAL(wp), PARAMETER :: amo2 = 31.9988E-3_wp !< O2 REAL(wp), PARAMETER :: amo3 = 47.998E-3_wp !< O3 REAL(wp), PARAMETER :: amoc = 150.E-3_wp !< organic carbon (OC) REAL(wp), PARAMETER :: amss = 58.44E-3_wp !< sea salt (NaCl) !-- Densities in kg/m3 REAL(wp), PARAMETER :: arhobc = 2000.0_wp !< black carbon REAL(wp), PARAMETER :: arhodu = 2650.0_wp !< mineral dust REAL(wp), PARAMETER :: arhoh2o = 1000.0_wp !< H2O REAL(wp), PARAMETER :: arhoh2so4 = 1830.0_wp !< SO4 REAL(wp), PARAMETER :: arhohno3 = 1479.0_wp !< HNO3 REAL(wp), PARAMETER :: arhonh3 = 1530.0_wp !< NH3 REAL(wp), PARAMETER :: arhooc = 2000.0_wp !< organic carbon REAL(wp), PARAMETER :: arhoss = 2165.0_wp !< sea salt (NaCl) !-- Volume of molecule in m3/# REAL(wp), PARAMETER :: amvh2o = amh2o /avo / arhoh2o !< H2O REAL(wp), PARAMETER :: amvh2so4 = amh2so4 / avo / arhoh2so4 !< SO4 REAL(wp), PARAMETER :: amvhno3 = amhno3 / avo / arhohno3 !< HNO3 REAL(wp), PARAMETER :: amvnh3 = amnh3 / avo / arhonh3 !< NH3 REAL(wp), PARAMETER :: amvoc = amoc / avo / arhooc !< OC REAL(wp), PARAMETER :: amvss = amss / avo / arhoss !< sea salt ! !-- SALSA switches: INTEGER(iwp) :: nj3 = 1 !< J3 parametrization (nucleation) !< 1 = condensational sink (Kerminen&Kulmala, 2002) !< 2 = coagulational sink (Lehtinen et al. 2007) !< 3 = coagS+self-coagulation (Anttila et al. 2010) INTEGER(iwp) :: nsnucl = 0 !< Choice of the nucleation scheme: !< 0 = off !< 1 = binary nucleation !< 2 = activation type nucleation !< 3 = kinetic nucleation !< 4 = ternary nucleation !< 5 = nucleation with ORGANICs !< 6 = activation type of nucleation with !< H2SO4+ORG !< 7 = heteromolecular nucleation with H2SO4*ORG !< 8 = homomolecular nucleation of H2SO4 + !< heteromolecular nucleation with H2SO4*ORG !< 9 = homomolecular nucleation of H2SO4 and ORG !< +heteromolecular nucleation with H2SO4*ORG LOGICAL :: advect_particle_water = .TRUE. !< advect water concentration of !< particles LOGICAL :: decycle_lr = .FALSE. !< Undo cyclic boundary !< conditions: left and right LOGICAL :: decycle_ns = .FALSE. !< north and south boundaries LOGICAL :: feedback_to_palm = .FALSE. !< allow feedback due to !< hydration and/or condensation !< of H20 LOGICAL :: no_insoluble = .FALSE. !< Switch to exclude insoluble !< chemical components LOGICAL :: read_restart_data_salsa = .FALSE. !< read restart data for salsa LOGICAL :: salsa_gases_from_chem = .FALSE. !< Transfer the gaseous !< components to SALSA from !< from chemistry model LOGICAL :: van_der_waals_coagc = .FALSE. !< Enhancement of coagulation !< kernel by van der Waals and !< viscous forces LOGICAL :: write_binary_salsa = .FALSE. !< read binary for salsa !-- Process switches: nl* is read from the NAMELIST and is NOT changed. !-- ls* is the switch used and will get the value of nl* !-- except for special circumstances (spinup period etc.) LOGICAL :: nlcoag = .FALSE. !< Coagulation master switch LOGICAL :: lscoag = .FALSE. !< LOGICAL :: nlcnd = .FALSE. !< Condensation master switch LOGICAL :: lscnd = .FALSE. !< LOGICAL :: nlcndgas = .FALSE. !< Condensation of precursor gases LOGICAL :: lscndgas = .FALSE. !< LOGICAL :: nlcndh2oae = .FALSE. !< Condensation of H2O on aerosol LOGICAL :: lscndh2oae = .FALSE. !< particles (FALSE -> equilibrium calc.) LOGICAL :: nldepo = .FALSE. !< Deposition master switch LOGICAL :: lsdepo = .FALSE. !< LOGICAL :: nldepo_topo = .FALSE. !< Deposition on vegetation master switch LOGICAL :: lsdepo_topo = .FALSE. !< LOGICAL :: nldepo_vege = .FALSE. !< Deposition on walls master switch LOGICAL :: lsdepo_vege = .FALSE. !< LOGICAL :: nldistupdate = .TRUE. !< Size distribution update master switch LOGICAL :: lsdistupdate = .FALSE. !< ! !-- SALSA variables: CHARACTER (LEN=20) :: bc_salsa_b = 'neumann' !< bottom boundary condition CHARACTER (LEN=20) :: bc_salsa_t = 'neumann' !< top boundary condition CHARACTER (LEN=20) :: depo_vege_type = 'zhang2001' !< or 'petroff2010' CHARACTER (LEN=20) :: depo_topo_type = 'zhang2001' !< or 'petroff2010' CHARACTER (LEN=20), DIMENSION(4) :: decycle_method = & (/'dirichlet','dirichlet','dirichlet','dirichlet'/) !< Decycling method at horizontal boundaries, !< 1=left, 2=right, 3=south, 4=north !< dirichlet = initial size distribution and !< chemical composition set for the ghost and !< first three layers !< neumann = zero gradient CHARACTER (LEN=3), DIMENSION(maxspec) :: listspec = & !< Active aerosols (/'SO4',' ',' ',' ',' ',' ',' '/) CHARACTER (LEN=20) :: salsa_source_mode = 'no_source' !< 'read_from_file', !< 'constant' or 'no_source' INTEGER(iwp) :: dots_salsa = 0 !< starting index for salsa-timeseries INTEGER(iwp) :: fn1a = 1 !< last index for bin subranges: subrange 1a INTEGER(iwp) :: fn2a = 1 !< subrange 2a INTEGER(iwp) :: fn2b = 1 !< subrange 2b INTEGER(iwp), DIMENSION(ngast) :: gas_index_chem = (/ 1, 1, 1, 1, 1/) !< !< Index of gaseous compounds in the chemistry !< model. In SALSA, 1 = H2SO4, 2 = HNO3, !< 3 = NH3, 4 = OCNV, 5 = OCSV INTEGER(iwp) :: ibc_salsa_b !< INTEGER(iwp) :: ibc_salsa_t !< INTEGER(iwp) :: igctyp = 0 !< Initial gas concentration type !< 0 = uniform (use H2SO4_init, HNO3_init, !< NH3_init, OCNV_init and OCSV_init) !< 1 = read vertical profile from an input file INTEGER(iwp) :: in1a = 1 !< start index for bin subranges: subrange 1a INTEGER(iwp) :: in2a = 1 !< subrange 2a INTEGER(iwp) :: in2b = 1 !< subrange 2b INTEGER(iwp) :: isdtyp = 0 !< Initial size distribution type !< 0 = uniform !< 1 = read vertical profile of the mode number !< concentration from an input file INTEGER(iwp) :: ibc = -1 !< Indice for: black carbon (BC) INTEGER(iwp) :: idu = -1 !< dust INTEGER(iwp) :: inh = -1 !< NH3 INTEGER(iwp) :: ino = -1 !< HNO3 INTEGER(iwp) :: ioc = -1 !< organic carbon (OC) INTEGER(iwp) :: iso4 = -1 !< SO4 or H2SO4 INTEGER(iwp) :: iss = -1 !< sea salt INTEGER(iwp) :: lod_aero = 0 !< level of detail for aerosol emissions INTEGER(iwp) :: lod_gases = 0 !< level of detail for gaseous emissions INTEGER(iwp), DIMENSION(nreg) :: nbin = (/ 3, 7/) !< Number of size bins !< for each aerosol size subrange INTEGER(iwp) :: nbins = 1 !< total number of size bins INTEGER(iwp) :: ncc = 1 !< number of chemical components used INTEGER(iwp) :: ncc_tot = 1!< total number of chemical compounds (ncc+1 !< if particle water is advected) REAL(wp) :: act_coeff = 1.0E-7_wp !< Activation coefficient REAL(wp) :: aerosol_source = 0.0_wp !< Constant aerosol flux (#/(m3*s)) REAL(wp), ALLOCATABLE, DIMENSION(:,:) :: emission_mass_fracs !< array for !< aerosol composition per emission category !< 1:SO4 2:OC 3:BC 4:DU 5:SS 6:NO 7:NH REAL(wp) :: dt_salsa = 0.00001_wp !< Time step of SALSA REAL(wp) :: H2SO4_init = nclim !< Init value for sulphuric acid gas REAL(wp) :: HNO3_init = nclim !< Init value for nitric acid gas REAL(wp) :: last_salsa_time = 0.0_wp !< time of the previous salsa !< timestep REAL(wp) :: nf2a = 1.0_wp !< Number fraction allocated to a- !< bins in subrange 2 !< (b-bins will get 1-nf2a) REAL(wp) :: NH3_init = nclim !< Init value for ammonia gas REAL(wp) :: OCNV_init = nclim !< Init value for non-volatile !< organic gases REAL(wp) :: OCSV_init = nclim !< Init value for semi-volatile !< organic gases REAL(wp), DIMENSION(nreg+1) :: reglim = & !< Min&max diameters of size subranges (/ 3.0E-9_wp, 5.0E-8_wp, 1.0E-5_wp/) REAL(wp) :: rhlim = 1.20_wp !< RH limit in %/100. Prevents !< unrealistically high RH in condensation REAL(wp) :: skip_time_do_salsa = 0.0_wp !< Starting time of SALSA (s) !-- Initial log-normal size distribution: mode diameter (dpg, micrometres), !-- standard deviation (sigmag) and concentration (n_lognorm, #/cm3) REAL(wp), DIMENSION(nmod) :: dpg = (/0.013_wp, 0.054_wp, 0.86_wp, & 0.2_wp, 0.2_wp, 0.2_wp, 0.2_wp/) REAL(wp), DIMENSION(nmod) :: sigmag = (/1.8_wp, 2.16_wp, 2.21_wp, & 2.0_wp, 2.0_wp, 2.0_wp, 2.0_wp/) REAL(wp), DIMENSION(nmod) :: n_lognorm = (/1.04e+5_wp, 3.23E+4_wp, 5.4_wp,& 0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp/) !-- Initial mass fractions / chemical composition of the size distribution REAL(wp), DIMENSION(maxspec) :: mass_fracs_a = & !< mass fractions between (/1.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0/) !< aerosol species for A bins REAL(wp), DIMENSION(maxspec) :: mass_fracs_b = & !< mass fractions between (/0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0/) !< aerosol species for B bins REAL(wp), ALLOCATABLE, DIMENSION(:) :: bin_low_limits !< to deliver !< information about !< the lower !< diameters per bin REAL(wp), ALLOCATABLE, DIMENSION(:) :: nsect !< Background number !< concentration per bin REAL(wp), ALLOCATABLE, DIMENSION(:) :: massacc !< Mass accomodation !< coefficients per bin ! !-- SALSA derived datatypes: ! !-- Prognostic variable: Aerosol size bin information (number (#/m3) and !-- mass (kg/m3) concentration) and the concentration of gaseous tracers (#/m3). !-- Gas tracers are contained sequentially in dimension 4 as: !-- 1. H2SO4, 2. HNO3, 3. NH3, 4. OCNV (non-volatile organics), !-- 5. OCSV (semi-volatile) TYPE salsa_variable REAL(wp), POINTER, DIMENSION(:,:,:), CONTIGUOUS :: conc REAL(wp), POINTER, DIMENSION(:,:,:), CONTIGUOUS :: conc_p REAL(wp), POINTER, DIMENSION(:,:,:), CONTIGUOUS :: tconc_m REAL(wp), ALLOCATABLE, DIMENSION(:,:) :: flux_s, diss_s REAL(wp), ALLOCATABLE, DIMENSION(:,:,:) :: flux_l, diss_l REAL(wp), ALLOCATABLE, DIMENSION(:) :: init REAL(wp), ALLOCATABLE, DIMENSION(:,:,:) :: source REAL(wp), ALLOCATABLE, DIMENSION(:,:) :: sums_ws_l END TYPE salsa_variable !-- Map bin indices between parallel size distributions TYPE t_parallelbin INTEGER(iwp) :: cur ! Index for current distribution INTEGER(iwp) :: par ! Index for corresponding parallel distribution END TYPE t_parallelbin !-- Datatype used to store information about the binned size distributions of !-- aerosols TYPE t_section REAL(wp) :: vhilim !< bin volume at the high limit REAL(wp) :: vlolim !< bin volume at the low limit REAL(wp) :: vratiohi !< volume ratio between the center and high limit REAL(wp) :: vratiolo !< volume ratio between the center and low limit REAL(wp) :: dmid !< bin middle diameter (m) !****************************************************** ! ^ Do NOT change the stuff above after initialization ! !****************************************************** REAL(wp) :: dwet !< Wet diameter or mean droplet diameter (m) REAL(wp), DIMENSION(maxspec+1) :: volc !< Volume concentrations !< (m^3/m^3) of aerosols + water. Since most of !< the stuff in SALSA is hard coded, these *have to !< be* in the order !< 1:SO4, 2:OC, 3:BC, 4:DU, 5:SS, 6:NO, 7:NH, 8:H2O REAL(wp) :: veqh2o !< Equilibrium H2O concentration for each particle REAL(wp) :: numc !< Number concentration of particles/droplets (#/m3) REAL(wp) :: core !< Volume of dry particle END TYPE t_section ! !-- Local aerosol properties in SALSA TYPE(t_section), ALLOCATABLE :: aero(:) ! !-- SALSA tracers: !-- Tracers as x = x(k,j,i,bin). The 4th dimension contains all the size bins !-- sequentially for each aerosol species + water. ! !-- Prognostic tracers: ! !-- Number concentration (#/m3) TYPE(salsa_variable), ALLOCATABLE, DIMENSION(:), TARGET :: aerosol_number REAL(wp), ALLOCATABLE, DIMENSION(:,:,:,:), TARGET :: nconc_1 REAL(wp), ALLOCATABLE, DIMENSION(:,:,:,:), TARGET :: nconc_2 REAL(wp), ALLOCATABLE, DIMENSION(:,:,:,:), TARGET :: nconc_3 ! !-- Mass concentration (kg/m3) TYPE(salsa_variable), ALLOCATABLE, DIMENSION(:), TARGET :: aerosol_mass REAL(wp), ALLOCATABLE, DIMENSION(:,:,:,:), TARGET :: mconc_1 REAL(wp), ALLOCATABLE, DIMENSION(:,:,:,:), TARGET :: mconc_2 REAL(wp), ALLOCATABLE, DIMENSION(:,:,:,:), TARGET :: mconc_3 ! !-- Gaseous tracers (#/m3) TYPE(salsa_variable), ALLOCATABLE, DIMENSION(:), TARGET :: salsa_gas REAL(wp), ALLOCATABLE, DIMENSION(:,:,:,:), TARGET :: gconc_1 REAL(wp), ALLOCATABLE, DIMENSION(:,:,:,:), TARGET :: gconc_2 REAL(wp), ALLOCATABLE, DIMENSION(:,:,:,:), TARGET :: gconc_3 ! !-- Diagnostic tracers REAL(wp), ALLOCATABLE, DIMENSION(:,:,:,:) :: sedim_vd !< sedimentation !< velocity per size !< bin (m/s) REAL(wp), ALLOCATABLE, DIMENSION(:,:,:,:) :: Ra_dry !< dry radius (m) !-- Particle component index tables TYPE(component_index) :: prtcl !< Contains "getIndex" which gives the index !< for a given aerosol component name, i.e. !< 1:SO4, 2:OC, 3:BC, 4:DU, !< 5:SS, 6:NO, 7:NH, 8:H2O ! !-- Data output arrays: !-- Gases: REAL(wp), ALLOCATABLE, DIMENSION(:,:,:), TARGET :: g_H2SO4_av !< H2SO4 REAL(wp), ALLOCATABLE, DIMENSION(:,:,:), TARGET :: g_HNO3_av !< HNO3 REAL(wp), ALLOCATABLE, DIMENSION(:,:,:), TARGET :: g_NH3_av !< NH3 REAL(wp), ALLOCATABLE, DIMENSION(:,:,:), TARGET :: g_OCNV_av !< non-vola- !< tile OC REAL(wp), ALLOCATABLE, DIMENSION(:,:,:), TARGET :: g_OCSV_av !< semi-vol. !< OC !-- Integrated: REAL(wp), ALLOCATABLE, DIMENSION(:,:,:), TARGET :: LDSA_av !< lung- !< deposited !< surface area REAL(wp), ALLOCATABLE, DIMENSION(:,:,:), TARGET :: Ntot_av !< total number !< conc. REAL(wp), ALLOCATABLE, DIMENSION(:,:,:), TARGET :: PM25_av !< PM2.5 REAL(wp), ALLOCATABLE, DIMENSION(:,:,:), TARGET :: PM10_av !< PM10 !-- In the particle phase: REAL(wp), ALLOCATABLE, DIMENSION(:,:,:), TARGET :: s_BC_av !< black carbon REAL(wp), ALLOCATABLE, DIMENSION(:,:,:), TARGET :: s_DU_av !< dust REAL(wp), ALLOCATABLE, DIMENSION(:,:,:), TARGET :: s_H2O_av !< liquid water REAL(wp), ALLOCATABLE, DIMENSION(:,:,:), TARGET :: s_NH_av !< ammonia REAL(wp), ALLOCATABLE, DIMENSION(:,:,:), TARGET :: s_NO_av !< nitrates REAL(wp), ALLOCATABLE, DIMENSION(:,:,:), TARGET :: s_OC_av !< org. carbon REAL(wp), ALLOCATABLE, DIMENSION(:,:,:), TARGET :: s_SO4_av !< sulphates REAL(wp), ALLOCATABLE, DIMENSION(:,:,:), TARGET :: s_SS_av !< sea salt !-- Bins: REAL(wp), ALLOCATABLE, DIMENSION(:,:,:,:), TARGET :: mbins_av !< bin mass REAL(wp), ALLOCATABLE, DIMENSION(:,:,:,:), TARGET :: Nbins_av !< bin number ! !-- PALM interfaces: ! !-- Boundary conditions: INTERFACE salsa_boundary_conds MODULE PROCEDURE salsa_boundary_conds MODULE PROCEDURE salsa_boundary_conds_decycle END INTERFACE salsa_boundary_conds ! !-- Data output checks for 2D/3D data to be done in check_parameters INTERFACE salsa_check_data_output MODULE PROCEDURE salsa_check_data_output END INTERFACE salsa_check_data_output ! !-- Input parameter checks to be done in check_parameters INTERFACE salsa_check_parameters MODULE PROCEDURE salsa_check_parameters END INTERFACE salsa_check_parameters ! !-- Averaging of 3D data for output INTERFACE salsa_3d_data_averaging MODULE PROCEDURE salsa_3d_data_averaging END INTERFACE salsa_3d_data_averaging ! !-- Data output of 2D quantities INTERFACE salsa_data_output_2d MODULE PROCEDURE salsa_data_output_2d END INTERFACE salsa_data_output_2d ! !-- Data output of 3D data INTERFACE salsa_data_output_3d MODULE PROCEDURE salsa_data_output_3d END INTERFACE salsa_data_output_3d ! !-- Data output of 3D data INTERFACE salsa_data_output_mask MODULE PROCEDURE salsa_data_output_mask END INTERFACE salsa_data_output_mask ! !-- Definition of data output quantities INTERFACE salsa_define_netcdf_grid MODULE PROCEDURE salsa_define_netcdf_grid END INTERFACE salsa_define_netcdf_grid ! !-- Output of information to the header file INTERFACE salsa_header MODULE PROCEDURE salsa_header END INTERFACE salsa_header ! !-- Initialization actions INTERFACE salsa_init MODULE PROCEDURE salsa_init END INTERFACE salsa_init ! !-- Initialization of arrays INTERFACE salsa_init_arrays MODULE PROCEDURE salsa_init_arrays END INTERFACE salsa_init_arrays ! !-- Writing of binary output for restart runs !!! renaming?! INTERFACE salsa_wrd_local MODULE PROCEDURE salsa_wrd_local END INTERFACE salsa_wrd_local ! !-- Reading of NAMELIST parameters INTERFACE salsa_parin MODULE PROCEDURE salsa_parin END INTERFACE salsa_parin ! !-- Reading of parameters for restart runs INTERFACE salsa_rrd_local MODULE PROCEDURE salsa_rrd_local END INTERFACE salsa_rrd_local ! !-- Swapping of time levels (required for prognostic variables) INTERFACE salsa_swap_timelevel MODULE PROCEDURE salsa_swap_timelevel END INTERFACE salsa_swap_timelevel INTERFACE salsa_driver MODULE PROCEDURE salsa_driver END INTERFACE salsa_driver INTERFACE salsa_tendency MODULE PROCEDURE salsa_tendency MODULE PROCEDURE salsa_tendency_ij END INTERFACE salsa_tendency SAVE PRIVATE ! !-- Public functions: PUBLIC salsa_boundary_conds, salsa_check_data_output, & salsa_check_parameters, salsa_3d_data_averaging, & salsa_data_output_2d, salsa_data_output_3d, salsa_data_output_mask, & salsa_define_netcdf_grid, salsa_diagnostics, salsa_driver, & salsa_header, salsa_init, salsa_init_arrays, salsa_parin, & salsa_rrd_local, salsa_swap_timelevel, salsa_tendency, & salsa_wrd_local ! !-- Public parameters, constants and initial values PUBLIC dots_salsa, dt_salsa, last_salsa_time, lsdepo, salsa, & salsa_gases_from_chem, skip_time_do_salsa ! !-- Public prognostic variables PUBLIC aerosol_mass, aerosol_number, fn2a, fn2b, gconc_2, in1a, in2b, & mconc_2, nbins, ncc, ncc_tot, nclim, nconc_2, ngast, prtcl, Ra_dry, & salsa_gas, sedim_vd CONTAINS !------------------------------------------------------------------------------! ! Description: ! ------------ !> Parin for &salsa_par for new modules !------------------------------------------------------------------------------! SUBROUTINE salsa_parin IMPLICIT NONE CHARACTER (LEN=80) :: line !< dummy string that contains the current line !< of the parameter file NAMELIST /salsa_parameters/ & advect_particle_water, & ! Switch for advecting ! particle water. If .FALSE., ! equilibration is called at ! each time step. bc_salsa_b, & ! bottom boundary condition bc_salsa_t, & ! top boundary condition decycle_lr, & ! decycle SALSA components decycle_method, & ! decycle method applied: ! 1=left 2=right 3=south 4=north decycle_ns, & ! decycle SALSA components depo_vege_type, & ! Parametrisation type depo_topo_type, & ! Parametrisation type dpg, & ! Mean diameter for the initial ! log-normal modes dt_salsa, & ! SALSA timestep in seconds feedback_to_palm, & ! allow feedback due to ! hydration / condensation H2SO4_init, & ! Init value for sulphuric acid HNO3_init, & ! Init value for nitric acid igctyp, & ! Initial gas concentration type isdtyp, & ! Initial size distribution type listspec, & ! List of actived aerosols ! (string list) mass_fracs_a, & ! Initial relative contribution ! of each species to particle ! volume in a-bins, 0 for unused mass_fracs_b, & ! Initial relative contribution ! of each species to particle ! volume in b-bins, 0 for unused n_lognorm, & ! Number concentration for the ! log-normal modes nbin, & ! Number of size bins for ! aerosol size subranges 1 & 2 nf2a, & ! Number fraction of particles ! allocated to a-bins in ! subrange 2 b-bins will get ! 1-nf2a NH3_init, & ! Init value for ammonia nj3, & ! J3 parametrization ! 1 = condensational sink ! (Kerminen&Kulmala, 2002) ! 2 = coagulational sink ! (Lehtinen et al. 2007) ! 3 = coagS+self-coagulation ! (Anttila et al. 2010) nlcnd, & ! Condensation master switch nlcndgas, & ! Condensation of gases nlcndh2oae, & ! Condensation of H2O nlcoag, & ! Coagulation master switch nldepo, & ! Deposition master switch nldepo_vege, & ! Deposition on vegetation ! master switch nldepo_topo, & ! Deposition on topo master ! switch nldistupdate, & ! Size distribution update ! master switch nsnucl, & ! Nucleation scheme: ! 0 = off, ! 1 = binary nucleation ! 2 = activation type nucleation ! 3 = kinetic nucleation ! 4 = ternary nucleation ! 5 = nucleation with organics ! 6 = activation type of ! nucleation with H2SO4+ORG ! 7 = heteromolecular nucleation ! with H2SO4*ORG ! 8 = homomolecular nucleation ! of H2SO4 + heteromolecular ! nucleation with H2SO4*ORG ! 9 = homomolecular nucleation ! of H2SO4 and ORG + hetero- ! molecular nucleation with ! H2SO4*ORG OCNV_init, & ! Init value for non-volatile ! organic gases OCSV_init, & ! Init value for semi-volatile ! organic gases read_restart_data_salsa, & ! read restart data for ! salsa reglim, & ! Min&max diameter limits of ! size subranges salsa, & ! Master switch for SALSA salsa_source_mode,& ! 'read_from_file' or 'constant' ! or 'no_source' sigmag, & ! stdev for the initial log- ! normal modes skip_time_do_salsa, & ! Starting time of SALSA (s) van_der_waals_coagc,& ! include van der Waals forces write_binary_salsa ! Write binary for salsa line = ' ' ! !-- Try to find salsa package REWIND ( 11 ) line = ' ' DO WHILE ( INDEX( line, '&salsa_parameters' ) == 0 ) READ ( 11, '(A)', END=10 ) line ENDDO BACKSPACE ( 11 ) ! !-- Read user-defined namelist READ ( 11, salsa_parameters ) ! !-- Enable salsa (salsa switch in modules.f90) salsa = .TRUE. 10 CONTINUE END SUBROUTINE salsa_parin !------------------------------------------------------------------------------! ! Description: ! ------------ !> Check parameters routine for salsa. !------------------------------------------------------------------------------! SUBROUTINE salsa_check_parameters USE control_parameters, & ONLY: message_string IMPLICIT NONE ! !-- Checks go here (cf. check_parameters.f90). IF ( salsa .AND. .NOT. humidity ) THEN WRITE( message_string, * ) 'salsa = ', salsa, ' is ', & 'not allowed with humidity = ', humidity CALL message( 'check_parameters', 'SA0009', 1, 2, 0, 6, 0 ) ENDIF IF ( bc_salsa_b == 'dirichlet' ) THEN ibc_salsa_b = 0 ELSEIF ( bc_salsa_b == 'neumann' ) THEN ibc_salsa_b = 1 ELSE message_string = 'unknown boundary condition: bc_salsa_b = "' & // TRIM( bc_salsa_t ) // '"' CALL message( 'check_parameters', 'SA0011', 1, 2, 0, 6, 0 ) ENDIF IF ( bc_salsa_t == 'dirichlet' ) THEN ibc_salsa_t = 0 ELSEIF ( bc_salsa_t == 'neumann' ) THEN ibc_salsa_t = 1 ELSE message_string = 'unknown boundary condition: bc_salsa_t = "' & // TRIM( bc_salsa_t ) // '"' CALL message( 'check_parameters', 'SA0012', 1, 2, 0, 6, 0 ) ENDIF IF ( nj3 < 1 .OR. nj3 > 3 ) THEN message_string = 'unknown nj3 (must be 1-3)' CALL message( 'check_parameters', 'SA0044', 1, 2, 0, 6, 0 ) ENDIF END SUBROUTINE salsa_check_parameters !------------------------------------------------------------------------------! ! ! Description: ! ------------ !> Subroutine defining appropriate grid for netcdf variables. !> It is called out from subroutine netcdf. !> Same grid as for other scalars (see netcdf_interface_mod.f90) !------------------------------------------------------------------------------! SUBROUTINE salsa_define_netcdf_grid( var, found, grid_x, grid_y, grid_z ) IMPLICIT NONE CHARACTER (LEN=*), INTENT(OUT) :: grid_x !< CHARACTER (LEN=*), INTENT(OUT) :: grid_y !< CHARACTER (LEN=*), INTENT(OUT) :: grid_z !< CHARACTER (LEN=*), INTENT(IN) :: var !< LOGICAL, INTENT(OUT) :: found !< found = .TRUE. ! !-- Check for the grid IF ( var(1:2) == 'g_' ) THEN grid_x = 'x' grid_y = 'y' grid_z = 'zu' ELSEIF ( var(1:4) == 'LDSA' ) THEN grid_x = 'x' grid_y = 'y' grid_z = 'zu' ELSEIF ( var(1:5) == 'm_bin' ) THEN grid_x = 'x' grid_y = 'y' grid_z = 'zu' ELSEIF ( var(1:5) == 'N_bin' ) THEN grid_x = 'x' grid_y = 'y' grid_z = 'zu' ELSEIF ( var(1:4) == 'Ntot' ) THEN grid_x = 'x' grid_y = 'y' grid_z = 'zu' ELSEIF ( var(1:2) == 'PM' ) THEN grid_x = 'x' grid_y = 'y' grid_z = 'zu' ELSEIF ( var(1:2) == 's_' ) THEN grid_x = 'x' grid_y = 'y' grid_z = 'zu' ELSE found = .FALSE. grid_x = 'none' grid_y = 'none' grid_z = 'none' ENDIF END SUBROUTINE salsa_define_netcdf_grid !------------------------------------------------------------------------------! ! Description: ! ------------ !> Header output for new module !------------------------------------------------------------------------------! SUBROUTINE salsa_header( io ) IMPLICIT NONE INTEGER(iwp), INTENT(IN) :: io !< Unit of the output file ! !-- Write SALSA header WRITE( io, 1 ) WRITE( io, 2 ) skip_time_do_salsa WRITE( io, 3 ) dt_salsa WRITE( io, 12 ) SHAPE( aerosol_number(1)%conc ), nbins IF ( advect_particle_water ) THEN WRITE( io, 16 ) SHAPE( aerosol_mass(1)%conc ), ncc_tot*nbins, & advect_particle_water ELSE WRITE( io, 16 ) SHAPE( aerosol_mass(1)%conc ), ncc*nbins, & advect_particle_water ENDIF IF ( .NOT. salsa_gases_from_chem ) THEN WRITE( io, 17 ) SHAPE( aerosol_mass(1)%conc ), ngast, & salsa_gases_from_chem ENDIF WRITE( io, 4 ) IF ( nsnucl > 0 ) THEN WRITE( io, 5 ) nsnucl, nj3 ENDIF IF ( nlcoag ) THEN WRITE( io, 6 ) ENDIF IF ( nlcnd ) THEN WRITE( io, 7 ) nlcndgas, nlcndh2oae ENDIF IF ( nldepo ) THEN WRITE( io, 14 ) nldepo_vege, nldepo_topo ENDIF WRITE( io, 8 ) reglim, nbin, bin_low_limits WRITE( io, 15 ) nsect WRITE( io, 13 ) ncc, listspec, mass_fracs_a, mass_fracs_b IF ( .NOT. salsa_gases_from_chem ) THEN WRITE( io, 18 ) ngast, H2SO4_init, HNO3_init, NH3_init, OCNV_init, & OCSV_init ENDIF WRITE( io, 9 ) isdtyp, igctyp IF ( isdtyp == 0 ) THEN WRITE( io, 10 ) dpg, sigmag, n_lognorm ELSE WRITE( io, 11 ) ENDIF 1 FORMAT (//' SALSA information:'/ & ' ------------------------------'/) 2 FORMAT (' Starts at: skip_time_do_salsa = ', F10.2, ' s') 3 FORMAT (/' Timestep: dt_salsa = ', F6.2, ' s') 12 FORMAT (/' Array shape (z,y,x,bins):'/ & ' aerosol_number: ', 4(I3)) 16 FORMAT (/' aerosol_mass: ', 4(I3),/ & ' (advect_particle_water = ', L1, ')') 17 FORMAT (' salsa_gas: ', 4(I3),/ & ' (salsa_gases_from_chem = ', L1, ')') 4 FORMAT (/' Aerosol dynamic processes included: ') 5 FORMAT (/' nucleation (scheme = ', I1, ' and J3 parametrization = ',& I1, ')') 6 FORMAT (/' coagulation') 7 FORMAT (/' condensation (of precursor gases = ', L1, & ' and water vapour = ', L1, ')' ) 14 FORMAT (/' dry deposition (on vegetation = ', L1, & ' and on topography = ', L1, ')') 8 FORMAT (/' Aerosol bin subrange limits (in metres): ', 3(ES10.2E3), / & ' Number of size bins for each aerosol subrange: ', 2I3,/ & ' Aerosol bin limits (in metres): ', *(ES10.2E3)) 15 FORMAT (' Initial number concentration in bins at the lowest level', & ' (#/m**3):', *(ES10.2E3)) 13 FORMAT (/' Number of chemical components used: ', I1,/ & ' Species: ',7(A6),/ & ' Initial relative contribution of each species to particle', & ' volume in:',/ & ' a-bins: ', 7(F6.3),/ & ' b-bins: ', 7(F6.3)) 18 FORMAT (/' Number of gaseous tracers used: ', I1,/ & ' Initial gas concentrations:',/ & ' H2SO4: ',ES12.4E3, ' #/m**3',/ & ' HNO3: ',ES12.4E3, ' #/m**3',/ & ' NH3: ',ES12.4E3, ' #/m**3',/ & ' OCNV: ',ES12.4E3, ' #/m**3',/ & ' OCSV: ',ES12.4E3, ' #/m**3') 9 FORMAT (/' Initialising concentrations: ', / & ' Aerosol size distribution: isdtyp = ', I1,/ & ' Gas concentrations: igctyp = ', I1 ) 10 FORMAT ( ' Mode diametres: dpg(nmod) = ', 7(F7.3),/ & ' Standard deviation: sigmag(nmod) = ', 7(F7.2),/ & ' Number concentration: n_lognorm(nmod) = ', 7(ES12.4E3) ) 11 FORMAT (/' Size distribution read from a file.') END SUBROUTINE salsa_header !------------------------------------------------------------------------------! ! Description: ! ------------ !> Allocate SALSA arrays and define pointers if required !------------------------------------------------------------------------------! SUBROUTINE salsa_init_arrays USE surface_mod, & ONLY: surf_def_h, surf_def_v, surf_lsm_h, surf_lsm_v, surf_usm_h, & surf_usm_v IMPLICIT NONE INTEGER(iwp) :: gases_available !< Number of available gas components in !< the chemistry model INTEGER(iwp) :: i !< loop index for allocating INTEGER(iwp) :: l !< loop index for allocating: surfaces INTEGER(iwp) :: lsp !< loop index for chem species in the chemistry model gases_available = 0 ! !-- Allocate prognostic variables (see salsa_swap_timelevel) ! !-- Set derived indices: !-- (This does the same as the subroutine salsa_initialize in SALSA/ !-- UCLALES-SALSA) in1a = 1 ! 1st index of subrange 1a in2a = in1a + nbin(1) ! 1st index of subrange 2a fn1a = in2a - 1 ! last index of subrange 1a fn2a = fn1a + nbin(2) ! last index of subrange 2a ! !-- If the fraction of insoluble aerosols in subrange 2 is zero: do not allocate !-- arrays for them IF ( nf2a > 0.999999_wp .AND. SUM( mass_fracs_b ) < 0.00001_wp ) THEN no_insoluble = .TRUE. in2b = fn2a+1 ! 1st index of subrange 2b fn2b = fn2a ! last index of subrange 2b ELSE in2b = in2a + nbin(2) ! 1st index of subrange 2b fn2b = fn2a + nbin(2) ! last index of subrange 2b ENDIF nbins = fn2b ! total number of aerosol size bins ! !-- Create index tables for different aerosol components CALL component_index_constructor( prtcl, ncc, maxspec, listspec ) ncc_tot = ncc IF ( advect_particle_water ) ncc_tot = ncc + 1 ! Add water ! !-- Allocate: ALLOCATE( aero(nbins), bin_low_limits(nbins), nsect(nbins), massacc(nbins) ) IF ( nldepo ) ALLOCATE( sedim_vd(nzb:nzt+1,nysg:nyng,nxlg:nxrg,nbins) ) ALLOCATE( Ra_dry(nzb:nzt+1,nysg:nyng,nxlg:nxrg,nbins) ) ! !-- Aerosol number concentration ALLOCATE( aerosol_number(nbins) ) ALLOCATE( nconc_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg,nbins), & nconc_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg,nbins), & nconc_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg,nbins) ) nconc_1 = 0.0_wp nconc_2 = 0.0_wp nconc_3 = 0.0_wp DO i = 1, nbins aerosol_number(i)%conc(nzb:nzt+1,nysg:nyng,nxlg:nxrg) => nconc_1(:,:,:,i) aerosol_number(i)%conc_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg) => nconc_2(:,:,:,i) aerosol_number(i)%tconc_m(nzb:nzt+1,nysg:nyng,nxlg:nxrg) => nconc_3(:,:,:,i) ALLOCATE( aerosol_number(i)%flux_s(nzb+1:nzt,0:threads_per_task-1), & aerosol_number(i)%diss_s(nzb+1:nzt,0:threads_per_task-1), & aerosol_number(i)%flux_l(nzb+1:nzt,nys:nyn,0:threads_per_task-1),& aerosol_number(i)%diss_l(nzb+1:nzt,nys:nyn,0:threads_per_task-1),& aerosol_number(i)%init(nzb:nzt+1), & aerosol_number(i)%sums_ws_l(nzb:nzt+1,0:threads_per_task-1) ) ENDDO ! !-- Aerosol mass concentration ALLOCATE( aerosol_mass(ncc_tot*nbins) ) ALLOCATE( mconc_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg,ncc_tot*nbins), & mconc_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg,ncc_tot*nbins), & mconc_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg,ncc_tot*nbins) ) mconc_1 = 0.0_wp mconc_2 = 0.0_wp mconc_3 = 0.0_wp DO i = 1, ncc_tot*nbins aerosol_mass(i)%conc(nzb:nzt+1,nysg:nyng,nxlg:nxrg) => mconc_1(:,:,:,i) aerosol_mass(i)%conc_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg) => mconc_2(:,:,:,i) aerosol_mass(i)%tconc_m(nzb:nzt+1,nysg:nyng,nxlg:nxrg) => mconc_3(:,:,:,i) ALLOCATE( aerosol_mass(i)%flux_s(nzb+1:nzt,0:threads_per_task-1), & aerosol_mass(i)%diss_s(nzb+1:nzt,0:threads_per_task-1), & aerosol_mass(i)%flux_l(nzb+1:nzt,nys:nyn,0:threads_per_task-1),& aerosol_mass(i)%diss_l(nzb+1:nzt,nys:nyn,0:threads_per_task-1),& aerosol_mass(i)%init(nzb:nzt+1), & aerosol_mass(i)%sums_ws_l(nzb:nzt+1,0:threads_per_task-1) ) ENDDO ! !-- Surface fluxes: answs = aerosol number, amsws = aerosol mass ! !-- Horizontal surfaces: default type DO l = 0, 2 ! upward (l=0), downward (l=1) and model top (l=2) ALLOCATE( surf_def_h(l)%answs( 1:surf_def_h(l)%ns, nbins ) ) ALLOCATE( surf_def_h(l)%amsws( 1:surf_def_h(l)%ns, nbins*ncc_tot ) ) surf_def_h(l)%answs = 0.0_wp surf_def_h(l)%amsws = 0.0_wp ENDDO !-- Horizontal surfaces: natural type IF ( land_surface ) THEN ALLOCATE( surf_lsm_h%answs( 1:surf_lsm_h%ns, nbins ) ) ALLOCATE( surf_lsm_h%amsws( 1:surf_lsm_h%ns, nbins*ncc_tot ) ) surf_lsm_h%answs = 0.0_wp surf_lsm_h%amsws = 0.0_wp ENDIF !-- Horizontal surfaces: urban type IF ( urban_surface ) THEN ALLOCATE( surf_usm_h%answs( 1:surf_usm_h%ns, nbins ) ) ALLOCATE( surf_usm_h%amsws( 1:surf_usm_h%ns, nbins*ncc_tot ) ) surf_usm_h%answs = 0.0_wp surf_usm_h%amsws = 0.0_wp ENDIF ! !-- Vertical surfaces: northward (l=0), southward (l=1), eastward (l=2) and !-- westward (l=3) facing DO l = 0, 3 ALLOCATE( surf_def_v(l)%answs( 1:surf_def_v(l)%ns, nbins ) ) surf_def_v(l)%answs = 0.0_wp ALLOCATE( surf_def_v(l)%amsws( 1:surf_def_v(l)%ns, nbins*ncc_tot ) ) surf_def_v(l)%amsws = 0.0_wp IF ( land_surface) THEN ALLOCATE( surf_lsm_v(l)%answs( 1:surf_lsm_v(l)%ns, nbins ) ) surf_lsm_v(l)%answs = 0.0_wp ALLOCATE( surf_lsm_v(l)%amsws( 1:surf_lsm_v(l)%ns, nbins*ncc_tot ) ) surf_lsm_v(l)%amsws = 0.0_wp ENDIF IF ( urban_surface ) THEN ALLOCATE( surf_usm_v(l)%answs( 1:surf_usm_v(l)%ns, nbins ) ) surf_usm_v(l)%answs = 0.0_wp ALLOCATE( surf_usm_v(l)%amsws( 1:surf_usm_v(l)%ns, nbins*ncc_tot ) ) surf_usm_v(l)%amsws = 0.0_wp ENDIF ENDDO ! !-- Concentration of gaseous tracers (1. SO4, 2. HNO3, 3. NH3, 4. OCNV, 5. OCSV) !-- (number concentration (#/m3) ) ! !-- If chemistry is on, read gas phase concentrations from there. Otherwise, !-- allocate salsa_gas array. IF ( air_chemistry ) THEN DO lsp = 1, nvar IF ( TRIM( chem_species(lsp)%name ) == 'H2SO4' ) THEN gases_available = gases_available + 1 gas_index_chem(1) = lsp ELSEIF ( TRIM( chem_species(lsp)%name ) == 'HNO3' ) THEN gases_available = gases_available + 1 gas_index_chem(2) = lsp ELSEIF ( TRIM( chem_species(lsp)%name ) == 'NH3' ) THEN gases_available = gases_available + 1 gas_index_chem(3) = lsp ELSEIF ( TRIM( chem_species(lsp)%name ) == 'OCNV' ) THEN gases_available = gases_available + 1 gas_index_chem(4) = lsp ELSEIF ( TRIM( chem_species(lsp)%name ) == 'OCSV' ) THEN gases_available = gases_available + 1 gas_index_chem(5) = lsp ENDIF ENDDO IF ( gases_available == ngast ) THEN salsa_gases_from_chem = .TRUE. ELSE WRITE( message_string, * ) 'SALSA is run together with chemistry '// & 'but not all gaseous components are '// & 'provided by kpp (H2SO4, HNO3, NH3, '// & 'OCNV, OCSC)' CALL message( 'check_parameters', 'SA0024', 1, 2, 0, 6, 0 ) ENDIF ELSE ALLOCATE( salsa_gas(ngast) ) ALLOCATE( gconc_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg,ngast), & gconc_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg,ngast), & gconc_3(nzb:nzt+1,nysg:nyng,nxlg:nxrg,ngast) ) gconc_1 = 0.0_wp gconc_2 = 0.0_wp gconc_3 = 0.0_wp DO i = 1, ngast salsa_gas(i)%conc(nzb:nzt+1,nysg:nyng,nxlg:nxrg) => gconc_1(:,:,:,i) salsa_gas(i)%conc_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg) => gconc_2(:,:,:,i) salsa_gas(i)%tconc_m(nzb:nzt+1,nysg:nyng,nxlg:nxrg) => gconc_3(:,:,:,i) ALLOCATE( salsa_gas(i)%flux_s(nzb+1:nzt,0:threads_per_task-1), & salsa_gas(i)%diss_s(nzb+1:nzt,0:threads_per_task-1), & salsa_gas(i)%flux_l(nzb+1:nzt,nys:nyn,0:threads_per_task-1),& salsa_gas(i)%diss_l(nzb+1:nzt,nys:nyn,0:threads_per_task-1),& salsa_gas(i)%init(nzb:nzt+1), & salsa_gas(i)%sums_ws_l(nzb:nzt+1,0:threads_per_task-1) ) ENDDO ! !-- Surface fluxes: gtsws = gaseous tracer flux ! !-- Horizontal surfaces: default type DO l = 0, 2 ! upward (l=0), downward (l=1) and model top (l=2) ALLOCATE( surf_def_h(l)%gtsws( 1:surf_def_h(l)%ns, ngast ) ) surf_def_h(l)%gtsws = 0.0_wp ENDDO !-- Horizontal surfaces: natural type IF ( land_surface ) THEN ALLOCATE( surf_lsm_h%gtsws( 1:surf_lsm_h%ns, ngast ) ) surf_lsm_h%gtsws = 0.0_wp ENDIF !-- Horizontal surfaces: urban type IF ( urban_surface ) THEN ALLOCATE( surf_usm_h%gtsws( 1:surf_usm_h%ns, ngast ) ) surf_usm_h%gtsws = 0.0_wp ENDIF ! !-- Vertical surfaces: northward (l=0), southward (l=1), eastward (l=2) and !-- westward (l=3) facing DO l = 0, 3 ALLOCATE( surf_def_v(l)%gtsws( 1:surf_def_v(l)%ns, ngast ) ) surf_def_v(l)%gtsws = 0.0_wp IF ( land_surface ) THEN ALLOCATE( surf_lsm_v(l)%gtsws( 1:surf_lsm_v(l)%ns, ngast ) ) surf_lsm_v(l)%gtsws = 0.0_wp ENDIF IF ( urban_surface ) THEN ALLOCATE( surf_usm_v(l)%gtsws( 1:surf_usm_v(l)%ns, ngast ) ) surf_usm_v(l)%gtsws = 0.0_wp ENDIF ENDDO ENDIF END SUBROUTINE salsa_init_arrays !------------------------------------------------------------------------------! ! Description: ! ------------ !> Initialization of SALSA. Based on salsa_initialize in UCLALES-SALSA. !> Subroutines salsa_initialize, SALSAinit and DiagInitAero in UCLALES-SALSA are !> also merged here. !------------------------------------------------------------------------------! SUBROUTINE salsa_init IMPLICIT NONE INTEGER(iwp) :: b INTEGER(iwp) :: c INTEGER(iwp) :: g INTEGER(iwp) :: i INTEGER(iwp) :: j bin_low_limits = 0.0_wp nsect = 0.0_wp massacc = 1.0_wp ! !-- Indices for chemical components used (-1 = not used) i = 0 IF ( is_used( prtcl, 'SO4' ) ) THEN iso4 = get_index( prtcl,'SO4' ) i = i + 1 ENDIF IF ( is_used( prtcl,'OC' ) ) THEN ioc = get_index(prtcl, 'OC') i = i + 1 ENDIF IF ( is_used( prtcl, 'BC' ) ) THEN ibc = get_index( prtcl, 'BC' ) i = i + 1 ENDIF IF ( is_used( prtcl, 'DU' ) ) THEN idu = get_index( prtcl, 'DU' ) i = i + 1 ENDIF IF ( is_used( prtcl, 'SS' ) ) THEN iss = get_index( prtcl, 'SS' ) i = i + 1 ENDIF IF ( is_used( prtcl, 'NO' ) ) THEN ino = get_index( prtcl, 'NO' ) i = i + 1 ENDIF IF ( is_used( prtcl, 'NH' ) ) THEN inh = get_index( prtcl, 'NH' ) i = i + 1 ENDIF ! !-- All species must be known IF ( i /= ncc ) THEN message_string = 'Unknown aerosol species/component(s) given in the' // & ' initialization' CALL message( 'salsa_mod: salsa_init', 'SA0020', 1, 2, 0, 6, 0 ) ENDIF ! !-- Initialise ! !-- Aerosol size distribution (TYPE t_section) aero(:)%dwet = 1.0E-10_wp aero(:)%veqh2o = 1.0E-10_wp aero(:)%numc = nclim aero(:)%core = 1.0E-10_wp DO c = 1, maxspec+1 ! 1:SO4, 2:OC, 3:BC, 4:DU, 5:SS, 6:NO, 7:NH, 8:H2O aero(:)%volc(c) = 0.0_wp ENDDO IF ( nldepo ) sedim_vd = 0.0_wp DO b = 1, nbins IF ( .NOT. read_restart_data_salsa ) aerosol_number(b)%conc = nclim aerosol_number(b)%conc_p = 0.0_wp aerosol_number(b)%tconc_m = 0.0_wp aerosol_number(b)%flux_s = 0.0_wp aerosol_number(b)%diss_s = 0.0_wp aerosol_number(b)%flux_l = 0.0_wp aerosol_number(b)%diss_l = 0.0_wp aerosol_number(b)%init = nclim aerosol_number(b)%sums_ws_l = 0.0_wp ENDDO DO c = 1, ncc_tot*nbins IF ( .NOT. read_restart_data_salsa ) aerosol_mass(c)%conc = mclim aerosol_mass(c)%conc_p = 0.0_wp aerosol_mass(c)%tconc_m = 0.0_wp aerosol_mass(c)%flux_s = 0.0_wp aerosol_mass(c)%diss_s = 0.0_wp aerosol_mass(c)%flux_l = 0.0_wp aerosol_mass(c)%diss_l = 0.0_wp aerosol_mass(c)%init = mclim aerosol_mass(c)%sums_ws_l = 0.0_wp ENDDO IF ( .NOT. salsa_gases_from_chem ) THEN DO g = 1, ngast salsa_gas(g)%conc_p = 0.0_wp salsa_gas(g)%tconc_m = 0.0_wp salsa_gas(g)%flux_s = 0.0_wp salsa_gas(g)%diss_s = 0.0_wp salsa_gas(g)%flux_l = 0.0_wp salsa_gas(g)%diss_l = 0.0_wp salsa_gas(g)%sums_ws_l = 0.0_wp ENDDO IF ( .NOT. read_restart_data_salsa ) THEN salsa_gas(1)%conc = H2SO4_init salsa_gas(2)%conc = HNO3_init salsa_gas(3)%conc = NH3_init salsa_gas(4)%conc = OCNV_init salsa_gas(5)%conc = OCSV_init ENDIF ! !-- Set initial value for gas compound tracers and initial values salsa_gas(1)%init = H2SO4_init salsa_gas(2)%init = HNO3_init salsa_gas(3)%init = NH3_init salsa_gas(4)%init = OCNV_init salsa_gas(5)%init = OCSV_init ENDIF ! !-- Aerosol radius in each bin: dry and wet (m) Ra_dry = 1.0E-10_wp ! !-- Initialise aerosol tracers aero(:)%vhilim = 0.0_wp aero(:)%vlolim = 0.0_wp aero(:)%vratiohi = 0.0_wp aero(:)%vratiolo = 0.0_wp aero(:)%dmid = 0.0_wp ! !-- Initialise the sectional particle size distribution CALL set_sizebins() ! !-- Initialise location-dependent aerosol size distributions and !-- chemical compositions: CALL aerosol_init ! !-- Initalisation run of SALSA DO i = nxl, nxr DO j = nys, nyn CALL salsa_driver( i, j, 1 ) CALL salsa_diagnostics( i, j ) ENDDO ENDDO ! !-- Set the aerosol and gas sources IF ( salsa_source_mode == 'read_from_file' ) THEN CALL salsa_set_source ENDIF END SUBROUTINE salsa_init !------------------------------------------------------------------------------! ! Description: ! ------------ !> Initializes particle size distribution grid by calculating size bin limits !> and mid-size for *dry* particles in each bin. Called from salsa_initialize !> (only at the beginning of simulation). !> Size distribution described using: !> 1) moving center method (subranges 1 and 2) !> (Jacobson, Atmos. Env., 31, 131-144, 1997) !> 2) fixed sectional method (subrange 3) !> Size bins in each subrange are spaced logarithmically !> based on given subrange size limits and bin number. ! !> Mona changed 06/2017: Use geometric mean diameter to describe the mean !> particle diameter in a size bin, not the arithmeric mean which clearly !> overestimates the total particle volume concentration. ! !> Coded by: !> Hannele Korhonen (FMI) 2005 !> Harri Kokkola (FMI) 2006 ! !> Bug fixes for box model + updated for the new aerosol datatype: !> Juha Tonttila (FMI) 2014 !------------------------------------------------------------------------------! SUBROUTINE set_sizebins IMPLICIT NONE ! !-- Local variables INTEGER(iwp) :: cc INTEGER(iwp) :: dd REAL(wp) :: ratio_d !< ratio of the upper and lower diameter of subranges ! !-- vlolim&vhilim: min & max *dry* volumes [fxm] !-- dmid: bin mid *dry* diameter (m) !-- vratiolo&vratiohi: volume ratio between the center and low/high limit ! !-- 1) Size subrange 1: ratio_d = reglim(2) / reglim(1) ! section spacing (m) DO cc = in1a,fn1a aero(cc)%vlolim = api6 * ( reglim(1) * ratio_d ** & ( REAL( cc-1 ) / nbin(1) ) ) ** 3.0_wp aero(cc)%vhilim = api6 * ( reglim(1) * ratio_d ** & ( REAL( cc ) / nbin(1) ) ) ** 3.0_wp aero(cc)%dmid = SQRT( ( aero(cc)%vhilim / api6 ) ** ( 1.0_wp / 3.0_wp ) & * ( aero(cc)%vlolim / api6 ) ** ( 1.0_wp / 3.0_wp ) ) aero(cc)%vratiohi = aero(cc)%vhilim / ( api6 * aero(cc)%dmid ** 3.0_wp ) aero(cc)%vratiolo = aero(cc)%vlolim / ( api6 * aero(cc)%dmid ** 3.0_wp ) ENDDO ! !-- 2) Size subrange 2: !-- 2.1) Sub-subrange 2a: high hygroscopicity ratio_d = reglim(3) / reglim(2) ! section spacing DO dd = in2a, fn2a cc = dd - in2a aero(dd)%vlolim = api6 * ( reglim(2) * ratio_d ** & ( REAL( cc ) / nbin(2) ) ) ** 3.0_wp aero(dd)%vhilim = api6 * ( reglim(2) * ratio_d ** & ( REAL( cc+1 ) / nbin(2) ) ) ** 3.0_wp aero(dd)%dmid = SQRT( ( aero(dd)%vhilim / api6 ) ** ( 1.0_wp / 3.0_wp ) & * ( aero(dd)%vlolim / api6 ) ** ( 1.0_wp / 3.0_wp ) ) aero(dd)%vratiohi = aero(dd)%vhilim / ( api6 * aero(dd)%dmid ** 3.0_wp ) aero(dd)%vratiolo = aero(dd)%vlolim / ( api6 * aero(dd)%dmid ** 3.0_wp ) ENDDO ! !-- 2.2) Sub-subrange 2b: low hygroscopicity IF ( .NOT. no_insoluble ) THEN aero(in2b:fn2b)%vlolim = aero(in2a:fn2a)%vlolim aero(in2b:fn2b)%vhilim = aero(in2a:fn2a)%vhilim aero(in2b:fn2b)%dmid = aero(in2a:fn2a)%dmid aero(in2b:fn2b)%vratiohi = aero(in2a:fn2a)%vratiohi aero(in2b:fn2b)%vratiolo = aero(in2a:fn2a)%vratiolo ENDIF ! !-- Initialize the wet diameter with the bin dry diameter to avoid numerical !-- problems later aero(:)%dwet = aero(:)%dmid ! !-- Save bin limits (lower diameter) to be delivered to the host model if needed DO cc = 1, nbins bin_low_limits(cc) = ( aero(cc)%vlolim / api6 )**( 1.0_wp / 3.0_wp ) ENDDO END SUBROUTINE set_sizebins !------------------------------------------------------------------------------! ! Description: ! ------------ !> Initilize altitude-dependent aerosol size distributions and compositions. !> !> Mona added 06/2017: Correct the number and mass concentrations by normalizing !< by the given total number and mass concentration. !> !> Tomi Raatikainen, FMI, 29.2.2016 !------------------------------------------------------------------------------! SUBROUTINE aerosol_init USE arrays_3d, & ONLY: zu ! USE NETCDF USE netcdf_data_input_mod, & ONLY: get_attribute, get_variable, & netcdf_data_input_get_dimension_length, open_read_file IMPLICIT NONE INTEGER(iwp) :: b !< loop index: size bins INTEGER(iwp) :: c !< loop index: chemical components INTEGER(iwp) :: ee !< index: end INTEGER(iwp) :: g !< loop index: gases INTEGER(iwp) :: i !< loop index: x-direction INTEGER(iwp) :: id_faero !< NetCDF id of PIDS_SALSA INTEGER(iwp) :: id_fchem !< NetCDF id of PIDS_CHEM INTEGER(iwp) :: j !< loop index: y-direction INTEGER(iwp) :: k !< loop index: z-direction INTEGER(iwp) :: kk !< loop index: z-direction INTEGER(iwp) :: nz_file !< Number of grid-points in file (heights) INTEGER(iwp) :: prunmode INTEGER(iwp) :: ss !< index: start LOGICAL :: netcdf_extend = .FALSE. !< Flag indicating wether netcdf !< topography input file or not REAL(wp), DIMENSION(nbins) :: core !< size of the bin mid aerosol particle, REAL(wp) :: flag !< flag to mask topography grid points REAL(wp), DIMENSION(:,:), ALLOCATABLE :: pr_gas !< gas profiles REAL(wp), DIMENSION(:,:), ALLOCATABLE :: pr_mass_fracs_a !< mass fraction !< profiles: a REAL(wp), DIMENSION(:,:), ALLOCATABLE :: pr_mass_fracs_b !< and b REAL(wp), DIMENSION(:,:), ALLOCATABLE :: pr_nsect !< sectional size !< distribution profile REAL(wp), DIMENSION(nbins) :: nsect !< size distribution (#/m3) REAL(wp), DIMENSION(0:nz+1,nbins) :: pndist !< size dist as a function !< of height (#/m3) REAL(wp), DIMENSION(0:nz+1) :: pnf2a !< number fraction: bins 2a REAL(wp), DIMENSION(0:nz+1,maxspec) :: pvf2a !< mass distributions of !< aerosol species for a REAL(wp), DIMENSION(0:nz+1,maxspec) :: pvf2b !< and b-bins REAL(wp), DIMENSION(0:nz+1) :: pvfOC1a !< mass fraction between !< SO4 and OC in 1a REAL(wp), DIMENSION(:), ALLOCATABLE :: pr_z prunmode = 1 ! !-- Bin mean aerosol particle volume (m3) core(:) = 0.0_wp core(1:nbins) = api6 * aero(1:nbins)%dmid ** 3.0_wp ! !-- Set concentrations to zero nsect(:) = 0.0_wp pndist(:,:) = 0.0_wp pnf2a(:) = nf2a pvf2a(:,:) = 0.0_wp pvf2b(:,:) = 0.0_wp pvfOC1a(:) = 0.0_wp IF ( isdtyp == 1 ) THEN ! !-- Read input profiles from PIDS_SALSA #if defined( __netcdf ) ! !-- Location-dependent size distributions and compositions. INQUIRE( FILE='PIDS_SALSA'// TRIM( coupling_char ), EXIST=netcdf_extend ) IF ( netcdf_extend ) THEN ! !-- Open file in read-only mode CALL open_read_file( 'PIDS_SALSA' // TRIM( coupling_char ), id_faero ) ! !-- Input heights CALL netcdf_data_input_get_dimension_length( id_faero, nz_file, & "profile_z" ) ALLOCATE( pr_z(nz_file), pr_mass_fracs_a(maxspec,nz_file), & pr_mass_fracs_b(maxspec,nz_file), pr_nsect(nbins,nz_file) ) CALL get_variable( id_faero, 'profile_z', pr_z ) ! !-- Mass fracs profile: 1: H2SO4 (sulphuric acid), 2: OC (organic carbon), !-- 3: BC (black carbon), 4: DU (dust), !-- 5: SS (sea salt), 6: HNO3 (nitric acid), !-- 7: NH3 (ammonia) CALL get_variable( id_faero, "profile_mass_fracs_a", pr_mass_fracs_a,& 0, nz_file-1, 0, maxspec-1 ) CALL get_variable( id_faero, "profile_mass_fracs_b", pr_mass_fracs_b,& 0, nz_file-1, 0, maxspec-1 ) CALL get_variable( id_faero, "profile_nsect", pr_nsect, 0, nz_file-1,& 0, nbins-1 ) kk = 1 DO k = nzb, nz+1 IF ( kk < nz_file ) THEN DO WHILE ( pr_z(kk+1) <= zu(k) ) kk = kk + 1 IF ( kk == nz_file ) EXIT ENDDO ENDIF IF ( kk < nz_file ) THEN ! !-- Set initial value for gas compound tracers and initial values pvf2a(k,:) = pr_mass_fracs_a(:,kk) + ( zu(k) - pr_z(kk) ) / ( & pr_z(kk+1) - pr_z(kk) ) * ( pr_mass_fracs_a(:,kk+1)& - pr_mass_fracs_a(:,kk) ) pvf2b(k,:) = pr_mass_fracs_b(:,kk) + ( zu(k) - pr_z(kk) ) / ( & pr_z(kk+1) - pr_z(kk) ) * ( pr_mass_fracs_b(:,kk+1)& - pr_mass_fracs_b(:,kk) ) pndist(k,:) = pr_nsect(:,kk) + ( zu(k) - pr_z(kk) ) / ( & pr_z(kk+1) - pr_z(kk) ) * ( pr_nsect(:,kk+1) - & pr_nsect(:,kk) ) ELSE pvf2a(k,:) = pr_mass_fracs_a(:,kk) pvf2b(k,:) = pr_mass_fracs_b(:,kk) pndist(k,:) = pr_nsect(:,kk) ENDIF IF ( iso4 < 0 ) THEN pvf2a(k,1) = 0.0_wp pvf2b(k,1) = 0.0_wp ENDIF IF ( ioc < 0 ) THEN pvf2a(k,2) = 0.0_wp pvf2b(k,2) = 0.0_wp ENDIF IF ( ibc < 0 ) THEN pvf2a(k,3) = 0.0_wp pvf2b(k,3) = 0.0_wp ENDIF IF ( idu < 0 ) THEN pvf2a(k,4) = 0.0_wp pvf2b(k,4) = 0.0_wp ENDIF IF ( iss < 0 ) THEN pvf2a(k,5) = 0.0_wp pvf2b(k,5) = 0.0_wp ENDIF IF ( ino < 0 ) THEN pvf2a(k,6) = 0.0_wp pvf2b(k,6) = 0.0_wp ENDIF IF ( inh < 0 ) THEN pvf2a(k,7) = 0.0_wp pvf2b(k,7) = 0.0_wp ENDIF ! !-- Then normalise the mass fraction so that SUM = 1 pvf2a(k,:) = pvf2a(k,:) / SUM( pvf2a(k,:) ) IF ( SUM( pvf2b(k,:) ) > 0.0_wp ) pvf2b(k,:) = pvf2b(k,:) / & SUM( pvf2b(k,:) ) ENDDO DEALLOCATE( pr_z, pr_mass_fracs_a, pr_mass_fracs_b, pr_nsect ) ELSE message_string = 'Input file '// TRIM( 'PIDS_SALSA' ) // & TRIM( coupling_char ) // ' for SALSA missing!' CALL message( 'salsa_mod: aerosol_init', 'SA0032', 1, 2, 0, 6, 0 ) ENDIF ! netcdf_extend #endif ELSEIF ( isdtyp == 0 ) THEN ! !-- Mass fractions for species in a and b-bins IF ( iso4 > 0 ) THEN pvf2a(:,1) = mass_fracs_a(iso4) pvf2b(:,1) = mass_fracs_b(iso4) ENDIF IF ( ioc > 0 ) THEN pvf2a(:,2) = mass_fracs_a(ioc) pvf2b(:,2) = mass_fracs_b(ioc) ENDIF IF ( ibc > 0 ) THEN pvf2a(:,3) = mass_fracs_a(ibc) pvf2b(:,3) = mass_fracs_b(ibc) ENDIF IF ( idu > 0 ) THEN pvf2a(:,4) = mass_fracs_a(idu) pvf2b(:,4) = mass_fracs_b(idu) ENDIF IF ( iss > 0 ) THEN pvf2a(:,5) = mass_fracs_a(iss) pvf2b(:,5) = mass_fracs_b(iss) ENDIF IF ( ino > 0 ) THEN pvf2a(:,6) = mass_fracs_a(ino) pvf2b(:,6) = mass_fracs_b(ino) ENDIF IF ( inh > 0 ) THEN pvf2a(:,7) = mass_fracs_a(inh) pvf2b(:,7) = mass_fracs_b(inh) ENDIF DO k = nzb, nz+1 pvf2a(k,:) = pvf2a(k,:) / SUM( pvf2a(k,:) ) IF ( SUM( pvf2b(k,:) ) > 0.0_wp ) pvf2b(k,:) = pvf2b(k,:) / & SUM( pvf2b(k,:) ) ENDDO CALL size_distribution( n_lognorm, dpg, sigmag, nsect ) ! !-- Normalize by the given total number concentration nsect = nsect * SUM( n_lognorm ) * 1.0E+6_wp / SUM( nsect ) DO b = in1a, fn2b pndist(:,b) = nsect(b) ENDDO ENDIF IF ( igctyp == 1 ) THEN ! !-- Read input profiles from PIDS_CHEM #if defined( __netcdf ) ! !-- Location-dependent size distributions and compositions. INQUIRE( FILE='PIDS_CHEM' // TRIM( coupling_char ), EXIST=netcdf_extend ) IF ( netcdf_extend .AND. .NOT. salsa_gases_from_chem ) THEN ! !-- Open file in read-only mode CALL open_read_file( 'PIDS_CHEM' // TRIM( coupling_char ), id_fchem ) ! !-- Input heights CALL netcdf_data_input_get_dimension_length( id_fchem, nz_file, & "profile_z" ) ALLOCATE( pr_z(nz_file), pr_gas(ngast,nz_file) ) CALL get_variable( id_fchem, 'profile_z', pr_z ) ! !-- Gases: CALL get_variable( id_fchem, "profile_H2SO4", pr_gas(1,:) ) CALL get_variable( id_fchem, "profile_HNO3", pr_gas(2,:) ) CALL get_variable( id_fchem, "profile_NH3", pr_gas(3,:) ) CALL get_variable( id_fchem, "profile_OCNV", pr_gas(4,:) ) CALL get_variable( id_fchem, "profile_OCSV", pr_gas(5,:) ) kk = 1 DO k = nzb, nz+1 IF ( kk < nz_file ) THEN DO WHILE ( pr_z(kk+1) <= zu(k) ) kk = kk + 1 IF ( kk == nz_file ) EXIT ENDDO ENDIF IF ( kk < nz_file ) THEN ! !-- Set initial value for gas compound tracers and initial values DO g = 1, ngast salsa_gas(g)%init(k) = pr_gas(g,kk) + ( zu(k) - pr_z(kk) ) & / ( pr_z(kk+1) - pr_z(kk) ) * & ( pr_gas(g,kk+1) - pr_gas(g,kk) ) salsa_gas(g)%conc(k,:,:) = salsa_gas(g)%init(k) ENDDO ELSE DO g = 1, ngast salsa_gas(g)%init(k) = pr_gas(g,kk) salsa_gas(g)%conc(k,:,:) = salsa_gas(g)%init(k) ENDDO ENDIF ENDDO DEALLOCATE( pr_z, pr_gas ) ELSEIF ( .NOT. netcdf_extend .AND. .NOT. salsa_gases_from_chem ) THEN message_string = 'Input file '// TRIM( 'PIDS_CHEM' ) // & TRIM( coupling_char ) // ' for SALSA missing!' CALL message( 'salsa_mod: aerosol_init', 'SA0033', 1, 2, 0, 6, 0 ) ENDIF ! netcdf_extend #endif ENDIF IF ( ioc > 0 .AND. iso4 > 0 ) THEN !-- Both are there, so use the given "massDistrA" pvfOC1a(:) = pvf2a(:,2) / ( pvf2a(:,2) + pvf2a(:,1) ) ! Normalize ELSEIF ( ioc > 0 ) THEN !-- Pure organic carbon pvfOC1a(:) = 1.0_wp ELSEIF ( iso4 > 0 ) THEN !-- Pure SO4 pvfOC1a(:) = 0.0_wp ELSE message_string = 'Either OC or SO4 must be active for aerosol region 1a!' CALL message( 'salsa_mod: aerosol_init', 'SA0021', 1, 2, 0, 6, 0 ) ENDIF ! !-- Initialize concentrations DO i = nxlg, nxrg DO j = nysg, nyng DO k = nzb, nzt+1 ! !-- Predetermine flag to mask topography flag = MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_0(k,j,i), 0 ) ) ! !-- a) Number concentrations !-- Region 1: DO b = in1a, fn1a aerosol_number(b)%conc(k,j,i) = pndist(k,b) * flag IF ( prunmode == 1 ) THEN aerosol_number(b)%init = pndist(:,b) ENDIF ENDDO ! !-- Region 2: IF ( nreg > 1 ) THEN DO b = in2a, fn2a aerosol_number(b)%conc(k,j,i) = MAX( 0.0_wp, pnf2a(k) ) * & pndist(k,b) * flag IF ( prunmode == 1 ) THEN aerosol_number(b)%init = MAX( 0.0_wp, nf2a ) * pndist(:,b) ENDIF ENDDO IF ( .NOT. no_insoluble ) THEN DO b = in2b, fn2b IF ( pnf2a(k) < 1.0_wp ) THEN aerosol_number(b)%conc(k,j,i) = MAX( 0.0_wp, 1.0_wp & - pnf2a(k) ) * pndist(k,b) * flag IF ( prunmode == 1 ) THEN aerosol_number(b)%init = MAX( 0.0_wp, 1.0_wp - & nf2a ) * pndist(:,b) ENDIF ENDIF ENDDO ENDIF ENDIF ! !-- b) Aerosol mass concentrations !-- bin subrange 1: done here separately due to the SO4/OC convention !-- SO4: IF ( iso4 > 0 ) THEN ss = ( iso4 - 1 ) * nbins + in1a !< start ee = ( iso4 - 1 ) * nbins + fn1a !< end b = in1a DO c = ss, ee aerosol_mass(c)%conc(k,j,i) = MAX( 0.0_wp, 1.0_wp - & pvfOC1a(k) ) * pndist(k,b) * & core(b) * arhoh2so4 * flag IF ( prunmode == 1 ) THEN aerosol_mass(c)%init = MAX( 0.0_wp, 1.0_wp - MAXVAL( & pvfOC1a ) ) * pndist(:,b) * & core(b) * arhoh2so4 ENDIF b = b+1 ENDDO ENDIF !-- OC: IF ( ioc > 0 ) THEN ss = ( ioc - 1 ) * nbins + in1a !< start ee = ( ioc - 1 ) * nbins + fn1a !< end b = in1a DO c = ss, ee aerosol_mass(c)%conc(k,j,i) = MAX( 0.0_wp, pvfOC1a(k) ) * & pndist(k,b) * core(b) * arhooc * flag IF ( prunmode == 1 ) THEN aerosol_mass(c)%init = MAX( 0.0_wp, MAXVAL( pvfOC1a ) ) & * pndist(:,b) * core(b) * arhooc ENDIF b = b+1 ENDDO ENDIF prunmode = 3 ! Init only once ENDDO !< k ENDDO !< j ENDDO !< i ! !-- c) Aerosol mass concentrations !-- bin subrange 2: IF ( nreg > 1 ) THEN IF ( iso4 > 0 ) THEN CALL set_aero_mass( iso4, pvf2a(:,1), pvf2b(:,1), pnf2a, pndist, & core, arhoh2so4 ) ENDIF IF ( ioc > 0 ) THEN CALL set_aero_mass( ioc, pvf2a(:,2), pvf2b(:,2), pnf2a, pndist, core,& arhooc ) ENDIF IF ( ibc > 0 ) THEN CALL set_aero_mass( ibc, pvf2a(:,3), pvf2b(:,3), pnf2a, pndist, core,& arhobc ) ENDIF IF ( idu > 0 ) THEN CALL set_aero_mass( idu, pvf2a(:,4), pvf2b(:,4), pnf2a, pndist, core,& arhodu ) ENDIF IF ( iss > 0 ) THEN CALL set_aero_mass( iss, pvf2a(:,5), pvf2b(:,5), pnf2a, pndist, core,& arhoss ) ENDIF IF ( ino > 0 ) THEN CALL set_aero_mass( ino, pvf2a(:,6), pvf2b(:,6), pnf2a, pndist, core,& arhohno3 ) ENDIF IF ( inh > 0 ) THEN CALL set_aero_mass( inh, pvf2a(:,7), pvf2b(:,7), pnf2a, pndist, core,& arhonh3 ) ENDIF ENDIF END SUBROUTINE aerosol_init !------------------------------------------------------------------------------! ! Description: ! ------------ !> Create a lognormal size distribution and discretise to a sectional !> representation. !------------------------------------------------------------------------------! SUBROUTINE size_distribution( in_ntot, in_dpg, in_sigma, psd_sect ) IMPLICIT NONE !-- Log-normal size distribution: modes REAL(wp), DIMENSION(:), INTENT(in) :: in_dpg !< geometric mean diameter !< (micrometres) REAL(wp), DIMENSION(:), INTENT(in) :: in_ntot !< number conc. (#/cm3) REAL(wp), DIMENSION(:), INTENT(in) :: in_sigma !< standard deviation REAL(wp), DIMENSION(:), INTENT(inout) :: psd_sect !< sectional size !< distribution INTEGER(iwp) :: b !< running index: bin INTEGER(iwp) :: ib !< running index: iteration REAL(wp) :: d1 !< particle diameter (m, dummy) REAL(wp) :: d2 !< particle diameter (m, dummy) REAL(wp) :: delta_d !< (d2-d1)/10 REAL(wp) :: deltadp !< bin width REAL(wp) :: dmidi !< ( d1 + d2 ) / 2 DO b = in1a, fn2b !< aerosol size bins psd_sect(b) = 0.0_wp !-- Particle diameter at the low limit (largest in the bin) (m) d1 = ( aero(b)%vlolim / api6 ) ** ( 1.0_wp / 3.0_wp ) !-- Particle diameter at the high limit (smallest in the bin) (m) d2 = ( aero(b)%vhilim / api6 ) ** ( 1.0_wp / 3.0_wp ) !-- Span of particle diameter in a bin (m) delta_d = ( d2 - d1 ) / 10.0_wp !-- Iterate: DO ib = 1, 10 d1 = ( aero(b)%vlolim / api6 ) ** ( 1.0_wp / 3.0_wp ) + ( ib - 1) & * delta_d d2 = d1 + delta_d dmidi = ( d1 + d2 ) / 2.0_wp deltadp = LOG10( d2 / d1 ) !-- Size distribution !-- in_ntot = total number, total area, or total volume concentration !-- in_dpg = geometric-mean number, area, or volume diameter !-- n(k) = number, area, or volume concentration in a bin !-- n_lognorm and dpg converted to units of #/m3 and m psd_sect(b) = psd_sect(b) + SUM( in_ntot * 1.0E+6_wp * deltadp / & ( SQRT( 2.0_wp * pi ) * LOG10( in_sigma ) ) * & EXP( -LOG10( dmidi / ( 1.0E-6_wp * in_dpg ) )**2.0_wp / & ( 2.0_wp * LOG10( in_sigma ) ** 2.0_wp ) ) ) ENDDO ENDDO END SUBROUTINE size_distribution !------------------------------------------------------------------------------! ! Description: ! ------------ !> Sets the mass concentrations to aerosol arrays in 2a and 2b. !> !> Tomi Raatikainen, FMI, 29.2.2016 !------------------------------------------------------------------------------! SUBROUTINE set_aero_mass( ispec, ppvf2a, ppvf2b, ppnf2a, ppndist, pcore, prho ) IMPLICIT NONE INTEGER(iwp), INTENT(in) :: ispec !< Aerosol species index REAL(wp), INTENT(in) :: pcore(nbins) !< Aerosol bin mid core volume REAL(wp), INTENT(in) :: ppndist(0:nz+1,nbins) !< Aerosol size distribution REAL(wp), INTENT(in) :: ppnf2a(0:nz+1) !< Number fraction for 2a REAL(wp), INTENT(in) :: ppvf2a(0:nz+1) !< Mass distributions for a REAL(wp), INTENT(in) :: ppvf2b(0:nz+1) !< and b bins REAL(wp), INTENT(in) :: prho !< Aerosol density INTEGER(iwp) :: b !< loop index INTEGER(iwp) :: c !< loop index INTEGER(iwp) :: ee !< index: end INTEGER(iwp) :: i !< loop index INTEGER(iwp) :: j !< loop index INTEGER(iwp) :: k !< loop index INTEGER(iwp) :: prunmode !< 1 = initialise INTEGER(iwp) :: ss !< index: start REAL(wp) :: flag !< flag to mask topography grid points prunmode = 1 DO i = nxlg, nxrg DO j = nysg, nyng DO k = nzb, nzt+1 ! !-- Predetermine flag to mask topography flag = MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_0(k,j,i), 0 ) ) ! !-- Regime 2a: ss = ( ispec - 1 ) * nbins + in2a ee = ( ispec - 1 ) * nbins + fn2a b = in2a DO c = ss, ee aerosol_mass(c)%conc(k,j,i) = MAX( 0.0_wp, ppvf2a(k) ) * & ppnf2a(k) * ppndist(k,b) * pcore(b) * prho * flag IF ( prunmode == 1 ) THEN aerosol_mass(c)%init = MAX( 0.0_wp, MAXVAL( ppvf2a(:) ) ) * & MAXVAL( ppnf2a ) * pcore(b) * prho * & MAXVAL( ppndist(:,b) ) ENDIF b = b+1 ENDDO !-- Regime 2b: IF ( .NOT. no_insoluble ) THEN ss = ( ispec - 1 ) * nbins + in2b ee = ( ispec - 1 ) * nbins + fn2b b = in2a DO c = ss, ee aerosol_mass(c)%conc(k,j,i) = MAX( 0.0_wp, ppvf2b(k) ) * ( & 1.0_wp - ppnf2a(k) ) * ppndist(k,b) * & pcore(b) * prho * flag IF ( prunmode == 1 ) THEN aerosol_mass(c)%init = MAX( 0.0_wp, MAXVAL( ppvf2b(:) ) )& * ( 1.0_wp - MAXVAL( ppnf2a ) ) * & MAXVAL( ppndist(:,b) ) * pcore(b) * prho ENDIF b = b+1 ENDDO ENDIF prunmode = 3 ! Init only once ENDDO ENDDO ENDDO END SUBROUTINE set_aero_mass !------------------------------------------------------------------------------! ! Description: ! ------------ !> Swapping of timelevels !------------------------------------------------------------------------------! SUBROUTINE salsa_swap_timelevel( mod_count ) IMPLICIT NONE INTEGER(iwp), INTENT(IN) :: mod_count !< INTEGER(iwp) :: b !< INTEGER(iwp) :: c !< INTEGER(iwp) :: cc !< INTEGER(iwp) :: g !< IF ( simulated_time >= time_since_reference_point ) THEN SELECT CASE ( mod_count ) CASE ( 0 ) DO b = 1, nbins aerosol_number(b)%conc(nzb:nzt+1,nysg:nyng,nxlg:nxrg) => & nconc_1(:,:,:,b) aerosol_number(b)%conc_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg) => & nconc_2(:,:,:,b) DO c = 1, ncc_tot cc = ( c-1 ) * nbins + b ! required due to possible Intel18 bug aerosol_mass(cc)%conc(nzb:nzt+1,nysg:nyng,nxlg:nxrg) => & mconc_1(:,:,:,cc) aerosol_mass(cc)%conc_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg) => & mconc_2(:,:,:,cc) ENDDO ENDDO IF ( .NOT. salsa_gases_from_chem ) THEN DO g = 1, ngast salsa_gas(g)%conc(nzb:nzt+1,nysg:nyng,nxlg:nxrg) => & gconc_1(:,:,:,g) salsa_gas(g)%conc_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg) => & gconc_2(:,:,:,g) ENDDO ENDIF CASE ( 1 ) DO b = 1, nbins aerosol_number(b)%conc(nzb:nzt+1,nysg:nyng,nxlg:nxrg) => & nconc_2(:,:,:,b) aerosol_number(b)%conc_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg) => & nconc_1(:,:,:,b) DO c = 1, ncc_tot cc = ( c-1 ) * nbins + b ! required due to possible Intel18 bug aerosol_mass(cc)%conc(nzb:nzt+1,nysg:nyng,nxlg:nxrg) => & mconc_2(:,:,:,cc) aerosol_mass(cc)%conc_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg) => & mconc_1(:,:,:,cc) ENDDO ENDDO IF ( .NOT. salsa_gases_from_chem ) THEN DO g = 1, ngast salsa_gas(g)%conc(nzb:nzt+1,nysg:nyng,nxlg:nxrg) => & gconc_2(:,:,:,g) salsa_gas(g)%conc_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg) => & gconc_1(:,:,:,g) ENDDO ENDIF END SELECT ENDIF END SUBROUTINE salsa_swap_timelevel !------------------------------------------------------------------------------! ! Description: ! ------------ !> This routine reads the respective restart data. !------------------------------------------------------------------------------! SUBROUTINE salsa_rrd_local( i, k, nxlf, nxlc, nxl_on_file, nxrf, nxrc, & nxr_on_file, nynf, nync, nyn_on_file, nysf, & nysc, nys_on_file, tmp_3d, found ) IMPLICIT NONE CHARACTER (LEN=20) :: field_char !< INTEGER(iwp) :: b !< INTEGER(iwp) :: c !< INTEGER(iwp) :: g !< INTEGER(iwp) :: i !< INTEGER(iwp) :: k !< INTEGER(iwp) :: nxlc !< INTEGER(iwp) :: nxlf !< INTEGER(iwp) :: nxl_on_file !< INTEGER(iwp) :: nxrc !< INTEGER(iwp) :: nxrf !< INTEGER(iwp) :: nxr_on_file !< INTEGER(iwp) :: nync !< INTEGER(iwp) :: nynf !< INTEGER(iwp) :: nyn_on_file !< INTEGER(iwp) :: nysc !< INTEGER(iwp) :: nysf !< INTEGER(iwp) :: nys_on_file !< LOGICAL, INTENT(OUT) :: found REAL(wp), & DIMENSION(nzb:nzt+1,nys_on_file-nbgp:nyn_on_file+nbgp,nxl_on_file-nbgp:nxr_on_file+nbgp) :: tmp_3d !< found = .FALSE. IF ( read_restart_data_salsa ) THEN SELECT CASE ( restart_string(1:length) ) CASE ( 'aerosol_number' ) DO b = 1, nbins IF ( k == 1 ) READ ( 13 ) tmp_3d aerosol_number(b)%conc(:,nysc-nbgp:nync+nbgp,nxlc-nbgp:nxrc+nbgp) = & tmp_3d(:,nysf-nbgp:nynf+nbgp,nxlf-nbgp:nxrf+nbgp) found = .TRUE. ENDDO CASE ( 'aerosol_mass' ) DO c = 1, ncc_tot * nbins IF ( k == 1 ) READ ( 13 ) tmp_3d aerosol_mass(c)%conc(:,nysc-nbgp:nync+nbgp,nxlc-nbgp:nxrc+nbgp) = & tmp_3d(:,nysf-nbgp:nynf+nbgp,nxlf-nbgp:nxrf+nbgp) found = .TRUE. ENDDO CASE ( 'salsa_gas' ) DO g = 1, ngast IF ( k == 1 ) READ ( 13 ) tmp_3d salsa_gas(g)%conc(:,nysc-nbgp:nync+nbgp,nxlc-nbgp:nxrc+nbgp) = & tmp_3d(:,nysf-nbgp:nynf+nbgp,nxlf-nbgp:nxrf+nbgp) found = .TRUE. ENDDO CASE DEFAULT found = .FALSE. END SELECT ENDIF END SUBROUTINE salsa_rrd_local !------------------------------------------------------------------------------! ! Description: ! ------------ !> This routine writes the respective restart data. !> Note that the following input variables in PARIN have to be equal between !> restart runs: !> listspec, nbin, nbin2, nf2a, ncc, mass_fracs_a, mass_fracs_b !------------------------------------------------------------------------------! SUBROUTINE salsa_wrd_local IMPLICIT NONE INTEGER(iwp) :: b !< INTEGER(iwp) :: c !< INTEGER(iwp) :: g !< IF ( write_binary .AND. write_binary_salsa ) THEN CALL wrd_write_string( 'aerosol_number' ) DO b = 1, nbins WRITE ( 14 ) aerosol_number(b)%conc ENDDO CALL wrd_write_string( 'aerosol_mass' ) DO c = 1, nbins*ncc_tot WRITE ( 14 ) aerosol_mass(c)%conc ENDDO CALL wrd_write_string( 'salsa_gas' ) DO g = 1, ngast WRITE ( 14 ) salsa_gas(g)%conc ENDDO ENDIF END SUBROUTINE salsa_wrd_local !------------------------------------------------------------------------------! ! Description: ! ------------ !> Performs necessary unit and dimension conversion between the host model and !> SALSA module, and calls the main SALSA routine. !> Partially adobted form the original SALSA boxmodel version. !> Now takes masses in as kg/kg from LES!! Converted to m3/m3 for SALSA !> 05/2016 Juha: This routine is still pretty much in its original shape. !> It's dumb as a mule and twice as ugly, so implementation of !> an improved solution is necessary sooner or later. !> Juha Tonttila, FMI, 2014 !> Jaakko Ahola, FMI, 2016 !> Only aerosol processes included, Mona Kurppa, UHel, 2017 !------------------------------------------------------------------------------! SUBROUTINE salsa_driver( i, j, prunmode ) USE arrays_3d, & ONLY: pt_p, q_p, rho_air_zw, u, v, w USE plant_canopy_model_mod, & ONLY: lad_s USE surface_mod, & ONLY: surf_def_h, surf_def_v, surf_lsm_h, surf_lsm_v, surf_usm_h, & surf_usm_v IMPLICIT NONE INTEGER(iwp), INTENT(in) :: i !< loop index INTEGER(iwp), INTENT(in) :: j !< loop index INTEGER(iwp), INTENT(in) :: prunmode !< 1: Initialization call !< 2: Spinup period call !< 3: Regular runtime call !-- Local variables TYPE(t_section), DIMENSION(fn2b) :: aero_old !< helper array INTEGER(iwp) :: bb !< loop index INTEGER(iwp) :: cc !< loop index INTEGER(iwp) :: endi !< end index INTEGER(iwp) :: k_wall !< vertical index of topography top INTEGER(iwp) :: k !< loop index INTEGER(iwp) :: l !< loop index INTEGER(iwp) :: nc_h2o !< index of H2O in the prtcl index table INTEGER(iwp) :: ss !< loop index INTEGER(iwp) :: str !< start index INTEGER(iwp) :: vc !< default index in prtcl REAL(wp) :: cw_old !< previous H2O mixing ratio REAL(wp) :: flag !< flag to mask topography grid points REAL(wp), DIMENSION(nzb:nzt+1) :: in_adn !< air density (kg/m3) REAL(wp), DIMENSION(nzb:nzt+1) :: in_cs !< H2O sat. vapour conc. REAL(wp), DIMENSION(nzb:nzt+1) :: in_cw !< H2O vapour concentration REAL(wp) :: in_lad !< leaf area density (m2/m3) REAL(wp), DIMENSION(nzb:nzt+1) :: in_p !< pressure (Pa) REAL(wp) :: in_rh !< relative humidity REAL(wp), DIMENSION(nzb:nzt+1) :: in_t !< temperature (K) REAL(wp), DIMENSION(nzb:nzt+1) :: in_u !< wind magnitude (m/s) REAL(wp), DIMENSION(nzb:nzt+1) :: kvis !< kinematic viscosity of air(m2/s) REAL(wp), DIMENSION(nzb:nzt+1,fn2b) :: Sc !< particle Schmidt number REAL(wp), DIMENSION(nzb:nzt+1,fn2b) :: vd !< particle fall seed (m/s, !< sedimentation velocity) REAL(wp), DIMENSION(nzb:nzt+1) :: ppm_to_nconc !< Conversion factor !< from ppm to #/m3 REAL(wp) :: zgso4 !< SO4 REAL(wp) :: zghno3 !< HNO3 REAL(wp) :: zgnh3 !< NH3 REAL(wp) :: zgocnv !< non-volatile OC REAL(wp) :: zgocsv !< semi-volatile OC aero_old(:)%numc = 0.0_wp in_adn = 0.0_wp in_cs = 0.0_wp in_cw = 0.0_wp in_lad = 0.0_wp in_rh = 0.0_wp in_p = 0.0_wp in_t = 0.0_wp in_u = 0.0_wp kvis = 0.0_wp Sc = 0.0_wp vd = 0.0_wp ppm_to_nconc = 1.0_wp zgso4 = nclim zghno3 = nclim zgnh3 = nclim zgocnv = nclim zgocsv = nclim ! !-- Aerosol number is always set, but mass can be uninitialized DO cc = 1, nbins aero(cc)%volc = 0.0_wp aero_old(cc)%volc = 0.0_wp ENDDO ! !-- Set the salsa runtime config (How to make this more efficient?) CALL set_salsa_runtime( prunmode ) ! !-- Calculate thermodynamic quantities needed in SALSA CALL salsa_thrm_ij( i, j, p_ij=in_p, temp_ij=in_t, cw_ij=in_cw, & cs_ij=in_cs, adn_ij=in_adn ) ! !-- Magnitude of wind: needed for deposition IF ( lsdepo ) THEN in_u(nzb+1:nzt) = SQRT( & ( 0.5_wp * ( u(nzb+1:nzt,j,i) + u(nzb+1:nzt,j,i+1) ) )**2 + & ( 0.5_wp * ( v(nzb+1:nzt,j,i) + v(nzb+1:nzt,j+1,i) ) )**2 + & ( 0.5_wp * ( w(nzb:nzt-1,j,i) + w(nzb+1:nzt,j, i) ) )**2 ) ENDIF ! !-- Calculate conversion factors for gas concentrations ppm_to_nconc = for_ppm_to_nconc * in_p / in_t ! !-- Determine topography-top index on scalar grid k_wall = MAXLOC( MERGE( 1, 0, BTEST( wall_flags_0(:,j,i), 12 ) ), & DIM = 1 ) - 1 DO k = nzb+1, nzt ! !-- Predetermine flag to mask topography flag = MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_0(k,j,i), 0 ) ) ! !-- Do not run inside buildings IF ( flag == 0.0_wp ) CYCLE ! !-- Wind velocity for dry depositon on vegetation IF ( lsdepo_vege .AND. plant_canopy ) THEN in_lad = lad_s(k-k_wall,j,i) ENDIF ! !-- For initialization and spinup, limit the RH with the parameter rhlim IF ( prunmode < 3 ) THEN in_cw(k) = MIN( in_cw(k), in_cs(k) * rhlim ) ELSE in_cw(k) = in_cw(k) ENDIF cw_old = in_cw(k) !* in_adn(k) ! !-- Set volume concentrations: !-- Sulphate (SO4) or sulphuric acid H2SO4 IF ( iso4 > 0 ) THEN vc = 1 str = ( iso4-1 ) * nbins + 1 ! start index endi = iso4 * nbins ! end index cc = 1 DO ss = str, endi aero(cc)%volc(vc) = aerosol_mass(ss)%conc(k,j,i) / arhoh2so4 cc = cc+1 ENDDO aero_old(1:nbins)%volc(vc) = aero(1:nbins)%volc(vc) ENDIF !-- Organic carbon (OC) compounds IF ( ioc > 0 ) THEN vc = 2 str = ( ioc-1 ) * nbins + 1 endi = ioc * nbins cc = 1 DO ss = str, endi aero(cc)%volc(vc) = aerosol_mass(ss)%conc(k,j,i) / arhooc cc = cc+1 ENDDO aero_old(1:nbins)%volc(vc) = aero(1:nbins)%volc(vc) ENDIF !-- Black carbon (BC) IF ( ibc > 0 ) THEN vc = 3 str = ( ibc-1 ) * nbins + 1 + fn1a endi = ibc * nbins cc = 1 + fn1a DO ss = str, endi aero(cc)%volc(vc) = aerosol_mass(ss)%conc(k,j,i) / arhobc cc = cc+1 ENDDO aero_old(1:nbins)%volc(vc) = aero(1:nbins)%volc(vc) ENDIF !-- Dust (DU) IF ( idu > 0 ) THEN vc = 4 str = ( idu-1 ) * nbins + 1 + fn1a endi = idu * nbins cc = 1 + fn1a DO ss = str, endi aero(cc)%volc(vc) = aerosol_mass(ss)%conc(k,j,i) / arhodu cc = cc+1 ENDDO aero_old(1:nbins)%volc(vc) = aero(1:nbins)%volc(vc) ENDIF !-- Sea salt (SS) IF ( iss > 0 ) THEN vc = 5 str = ( iss-1 ) * nbins + 1 + fn1a endi = iss * nbins cc = 1 + fn1a DO ss = str, endi aero(cc)%volc(vc) = aerosol_mass(ss)%conc(k,j,i) / arhoss cc = cc+1 ENDDO aero_old(1:nbins)%volc(vc) = aero(1:nbins)%volc(vc) ENDIF !-- Nitrate (NO(3-)) or nitric acid HNO3 IF ( ino > 0 ) THEN vc = 6 str = ( ino-1 ) * nbins + 1 endi = ino * nbins cc = 1 DO ss = str, endi aero(cc)%volc(vc) = aerosol_mass(ss)%conc(k,j,i) / arhohno3 cc = cc+1 ENDDO aero_old(1:nbins)%volc(vc) = aero(1:nbins)%volc(vc) ENDIF !-- Ammonium (NH(4+)) or ammonia NH3 IF ( inh > 0 ) THEN vc = 7 str = ( inh-1 ) * nbins + 1 endi = inh * nbins cc = 1 DO ss = str, endi aero(cc)%volc(vc) = aerosol_mass(ss)%conc(k,j,i) / arhonh3 cc = cc+1 ENDDO aero_old(1:nbins)%volc(vc) = aero(1:nbins)%volc(vc) ENDIF !-- Water (always used) nc_h2o = get_index( prtcl,'H2O' ) vc = 8 str = ( nc_h2o-1 ) * nbins + 1 endi = nc_h2o * nbins cc = 1 IF ( advect_particle_water ) THEN DO ss = str, endi aero(cc)%volc(vc) = aerosol_mass(ss)%conc(k,j,i) / arhoh2o cc = cc+1 ENDDO ELSE aero(1:nbins)%volc(vc) = mclim ENDIF aero_old(1:nbins)%volc(vc) = aero(1:nbins)%volc(vc) ! !-- Number concentrations (numc) and particle sizes !-- (dwet = wet diameter, core = dry volume) DO bb = 1, nbins aero(bb)%numc = aerosol_number(bb)%conc(k,j,i) aero_old(bb)%numc = aero(bb)%numc IF ( aero(bb)%numc > nclim ) THEN aero(bb)%dwet = ( SUM( aero(bb)%volc(:) ) / aero(bb)%numc / api6 )& **( 1.0_wp / 3.0_wp ) aero(bb)%core = SUM( aero(bb)%volc(1:7) ) / aero(bb)%numc ELSE aero(bb)%dwet = aero(bb)%dmid aero(bb)%core = api6 * ( aero(bb)%dwet ) ** 3.0_wp ENDIF ENDDO ! !-- On EACH call of salsa_driver, calculate the ambient sizes of !-- particles by equilibrating soluble fraction of particles with water !-- using the ZSR method. in_rh = in_cw(k) / in_cs(k) IF ( prunmode==1 .OR. .NOT. advect_particle_water ) THEN CALL equilibration( in_rh, in_t(k), aero, .TRUE. ) ENDIF ! !-- Gaseous tracer concentrations in #/m3 IF ( salsa_gases_from_chem ) THEN ! !-- Convert concentrations in ppm to #/m3 zgso4 = chem_species(gas_index_chem(1))%conc(k,j,i) * ppm_to_nconc(k) zghno3 = chem_species(gas_index_chem(2))%conc(k,j,i) * ppm_to_nconc(k) zgnh3 = chem_species(gas_index_chem(3))%conc(k,j,i) * ppm_to_nconc(k) zgocnv = chem_species(gas_index_chem(4))%conc(k,j,i) * ppm_to_nconc(k) zgocsv = chem_species(gas_index_chem(5))%conc(k,j,i) * ppm_to_nconc(k) ELSE zgso4 = salsa_gas(1)%conc(k,j,i) zghno3 = salsa_gas(2)%conc(k,j,i) zgnh3 = salsa_gas(3)%conc(k,j,i) zgocnv = salsa_gas(4)%conc(k,j,i) zgocsv = salsa_gas(5)%conc(k,j,i) ENDIF ! !-- ***************************************! !-- Run SALSA ! !-- ***************************************! CALL run_salsa( in_p(k), in_cw(k), in_cs(k), in_t(k), in_u(k), & in_adn(k), in_lad, zgso4, zgocnv, zgocsv, zghno3, zgnh3,& aero, prtcl, kvis(k), Sc(k,:), vd(k,:), dt_salsa ) !-- ***************************************! IF ( lsdepo ) sedim_vd(k,j,i,:) = vd(k,:) ! !-- Calculate changes in concentrations DO bb = 1, nbins aerosol_number(bb)%conc(k,j,i) = aerosol_number(bb)%conc(k,j,i) & + ( aero(bb)%numc - aero_old(bb)%numc ) * flag ENDDO IF ( iso4 > 0 ) THEN vc = 1 str = ( iso4-1 ) * nbins + 1 endi = iso4 * nbins cc = 1 DO ss = str, endi aerosol_mass(ss)%conc(k,j,i) = aerosol_mass(ss)%conc(k,j,i) & + ( aero(cc)%volc(vc) - aero_old(cc)%volc(vc) ) & * arhoh2so4 * flag cc = cc+1 ENDDO ENDIF IF ( ioc > 0 ) THEN vc = 2 str = ( ioc-1 ) * nbins + 1 endi = ioc * nbins cc = 1 DO ss = str, endi aerosol_mass(ss)%conc(k,j,i) = aerosol_mass(ss)%conc(k,j,i) & + ( aero(cc)%volc(vc) - aero_old(cc)%volc(vc) ) & * arhooc * flag cc = cc+1 ENDDO ENDIF IF ( ibc > 0 ) THEN vc = 3 str = ( ibc-1 ) * nbins + 1 + fn1a endi = ibc * nbins cc = 1 + fn1a DO ss = str, endi aerosol_mass(ss)%conc(k,j,i) = aerosol_mass(ss)%conc(k,j,i) & + ( aero(cc)%volc(vc) - aero_old(cc)%volc(vc) ) & * arhobc * flag cc = cc+1 ENDDO ENDIF IF ( idu > 0 ) THEN vc = 4 str = ( idu-1 ) * nbins + 1 + fn1a endi = idu * nbins cc = 1 + fn1a DO ss = str, endi aerosol_mass(ss)%conc(k,j,i) = aerosol_mass(ss)%conc(k,j,i) & + ( aero(cc)%volc(vc) - aero_old(cc)%volc(vc) ) & * arhodu * flag cc = cc+1 ENDDO ENDIF IF ( iss > 0 ) THEN vc = 5 str = ( iss-1 ) * nbins + 1 + fn1a endi = iss * nbins cc = 1 + fn1a DO ss = str, endi aerosol_mass(ss)%conc(k,j,i) = aerosol_mass(ss)%conc(k,j,i) & + ( aero(cc)%volc(vc) - aero_old(cc)%volc(vc) ) & * arhoss * flag cc = cc+1 ENDDO ENDIF IF ( ino > 0 ) THEN vc = 6 str = ( ino-1 ) * nbins + 1 endi = ino * nbins cc = 1 DO ss = str, endi aerosol_mass(ss)%conc(k,j,i) = aerosol_mass(ss)%conc(k,j,i) & + ( aero(cc)%volc(vc) - aero_old(cc)%volc(vc) ) & * arhohno3 * flag cc = cc+1 ENDDO ENDIF IF ( inh > 0 ) THEN vc = 7 str = ( ino-1 ) * nbins + 1 endi = ino * nbins cc = 1 DO ss = str, endi aerosol_mass(ss)%conc(k,j,i) = aerosol_mass(ss)%conc(k,j,i) & + ( aero(cc)%volc(vc) - aero_old(cc)%volc(vc) ) & * arhonh3 * flag cc = cc+1 ENDDO ENDIF IF ( advect_particle_water ) THEN nc_h2o = get_index( prtcl,'H2O' ) vc = 8 str = ( nc_h2o-1 ) * nbins + 1 endi = nc_h2o * nbins cc = 1 DO ss = str, endi aerosol_mass(ss)%conc(k,j,i) = aerosol_mass(ss)%conc(k,j,i) & + ( aero(cc)%volc(vc) - aero_old(cc)%volc(vc) ) & * arhoh2o * flag IF ( prunmode == 1 ) THEN aerosol_mass(ss)%init(k) = MAX( aerosol_mass(ss)%init(k), & aerosol_mass(ss)%conc(k,j,i) ) ENDIF cc = cc+1 ENDDO ENDIF !-- Condensation of precursor gases IF ( lscndgas ) THEN IF ( salsa_gases_from_chem ) THEN ! !-- SO4 (or H2SO4) chem_species( gas_index_chem(1) )%conc(k,j,i) = & chem_species( gas_index_chem(1) )%conc(k,j,i) + & ( zgso4 / ppm_to_nconc(k) - & chem_species( gas_index_chem(1) )%conc(k,j,i) ) * flag ! !-- HNO3 chem_species( gas_index_chem(2) )%conc(k,j,i) = & chem_species( gas_index_chem(2) )%conc(k,j,i) + & ( zghno3 / ppm_to_nconc(k) - & chem_species( gas_index_chem(2) )%conc(k,j,i) ) * flag ! !-- NH3 chem_species( gas_index_chem(3) )%conc(k,j,i) = & chem_species( gas_index_chem(3) )%conc(k,j,i) + & ( zgnh3 / ppm_to_nconc(k) - & chem_species( gas_index_chem(3) )%conc(k,j,i) ) * flag ! !-- non-volatile OC chem_species( gas_index_chem(4) )%conc(k,j,i) = & chem_species( gas_index_chem(4) )%conc(k,j,i) + & ( zgocnv / ppm_to_nconc(k) - & chem_species( gas_index_chem(4) )%conc(k,j,i) ) * flag ! !-- semi-volatile OC chem_species( gas_index_chem(5) )%conc(k,j,i) = & chem_species( gas_index_chem(5) )%conc(k,j,i) + & ( zgocsv / ppm_to_nconc(k) - & chem_species( gas_index_chem(5) )%conc(k,j,i) ) * flag ELSE ! !-- SO4 (or H2SO4) salsa_gas(1)%conc(k,j,i) = salsa_gas(1)%conc(k,j,i) + ( zgso4 - & salsa_gas(1)%conc(k,j,i) ) * flag ! !-- HNO3 salsa_gas(2)%conc(k,j,i) = salsa_gas(2)%conc(k,j,i) + ( zghno3 - & salsa_gas(2)%conc(k,j,i) ) * flag ! !-- NH3 salsa_gas(3)%conc(k,j,i) = salsa_gas(3)%conc(k,j,i) + ( zgnh3 - & salsa_gas(3)%conc(k,j,i) ) * flag ! !-- non-volatile OC salsa_gas(4)%conc(k,j,i) = salsa_gas(4)%conc(k,j,i) + ( zgocnv - & salsa_gas(4)%conc(k,j,i) ) * flag ! !-- semi-volatile OC salsa_gas(5)%conc(k,j,i) = salsa_gas(5)%conc(k,j,i) + ( zgocsv - & salsa_gas(5)%conc(k,j,i) ) * flag ENDIF ENDIF ! !-- Tendency of water vapour mixing ratio is obtained from the !-- change in RH during SALSA run. This releases heat and changes pt. !-- Assumes no temperature change during SALSA run. !-- q = r / (1+r), Euler method for integration ! IF ( feedback_to_palm ) THEN q_p(k,j,i) = q_p(k,j,i) + 1.0_wp / ( in_cw(k) * in_adn(k) + 1.0_wp ) & ** 2.0_wp * ( in_cw(k) - cw_old ) * in_adn(k) pt_p(k,j,i) = pt_p(k,j,i) + alv / c_p * ( in_cw(k) - cw_old ) * & in_adn(k) / ( in_cw(k) / in_adn(k) + 1.0_wp ) ** 2.0_wp& * pt_p(k,j,i) / in_t(k) ENDIF ENDDO ! k ! !-- Set surfaces and wall fluxes due to deposition IF ( lsdepo_topo .AND. prunmode == 3 ) THEN IF ( .NOT. land_surface .AND. .NOT. urban_surface ) THEN CALL depo_topo( i, j, surf_def_h(0), vd, Sc, kvis, in_u, rho_air_zw ) DO l = 0, 3 CALL depo_topo( i, j, surf_def_v(l), vd, Sc, kvis, in_u, & rho_air_zw**0.0_wp ) ENDDO ELSE CALL depo_topo( i, j, surf_usm_h, vd, Sc, kvis, in_u, rho_air_zw ) DO l = 0, 3 CALL depo_topo( i, j, surf_usm_v(l), vd, Sc, kvis, in_u, & rho_air_zw**0.0_wp ) ENDDO CALL depo_topo( i, j, surf_lsm_h, vd, Sc, kvis, in_u, rho_air_zw ) DO l = 0, 3 CALL depo_topo( i, j, surf_lsm_v(l), vd, Sc, kvis, in_u, & rho_air_zw**0.0_wp ) ENDDO ENDIF ENDIF END SUBROUTINE salsa_driver !------------------------------------------------------------------------------! ! Description: ! ------------ !> The SALSA subroutine !> Modified for the new aerosol datatype, !> Juha Tonttila, FMI, 2014. !> Only aerosol processes included, Mona Kurppa, UHel, 2017 !------------------------------------------------------------------------------! SUBROUTINE run_salsa( ppres, pcw, pcs, ptemp, mag_u, adn, lad, pc_h2so4, & pc_ocnv, pc_ocsv, pc_hno3, pc_nh3, paero, prtcl, kvis, & Sc, vc, ptstep ) IMPLICIT NONE ! !-- Input parameters and variables REAL(wp), INTENT(in) :: adn !< air density (kg/m3) REAL(wp), INTENT(in) :: lad !< leaf area density (m2/m3) REAL(wp), INTENT(in) :: mag_u !< magnitude of wind (m/s) REAL(wp), INTENT(in) :: ppres !< atmospheric pressure at each grid !< point (Pa) REAL(wp), INTENT(in) :: ptemp !< temperature at each grid point (K) REAL(wp), INTENT(in) :: ptstep !< time step of salsa processes (s) TYPE(component_index), INTENT(in) :: prtcl !< part. component index table ! !-- Input variables that are changed within: REAL(wp), INTENT(inout) :: kvis !< kinematic viscosity of air (m2/s) REAL(wp), INTENT(inout) :: Sc(:) !< particle Schmidt number REAL(wp), INTENT(inout) :: vc(:) !< particle fall speed (m/s, !< sedimentation velocity) !-- Gas phase concentrations at each grid point (#/m3) REAL(wp), INTENT(inout) :: pc_h2so4 !< sulphuric acid REAL(wp), INTENT(inout) :: pc_hno3 !< nitric acid REAL(wp), INTENT(inout) :: pc_nh3 !< ammonia REAL(wp), INTENT(inout) :: pc_ocnv !< nonvolatile OC REAL(wp), INTENT(inout) :: pc_ocsv !< semivolatile OC REAL(wp), INTENT(inout) :: pcs !< Saturation concentration of water !< vapour (kg/m3) REAL(wp), INTENT(inout) :: pcw !< Water vapour concentration (kg/m3) TYPE(t_section), INTENT(inout) :: paero(fn2b) ! !-- Coagulation IF ( lscoag ) THEN CALL coagulation( paero, ptstep, ptemp, ppres ) ENDIF ! !-- Condensation IF ( lscnd ) THEN CALL condensation( paero, pc_h2so4, pc_ocnv, pc_ocsv, pc_hno3, pc_nh3, & pcw, pcs, ptemp, ppres, ptstep, prtcl ) ENDIF ! !-- Deposition IF ( lsdepo ) THEN CALL deposition( paero, ptemp, adn, mag_u, lad, kvis, Sc, vc ) ENDIF ! !-- Size distribution bin update !-- Mona: why done 3 times in SALSA-standalone? IF ( lsdistupdate ) THEN CALL distr_update( paero ) ENDIF END SUBROUTINE run_salsa !------------------------------------------------------------------------------! ! Description: ! ------------ !> Set logical switches according to the host model state and user-specified !> NAMELIST options. !> Juha Tonttila, FMI, 2014 !> Only aerosol processes included, Mona Kurppa, UHel, 2017 !------------------------------------------------------------------------------! SUBROUTINE set_salsa_runtime( prunmode ) IMPLICIT NONE INTEGER(iwp), INTENT(in) :: prunmode SELECT CASE(prunmode) CASE(1) !< Initialization lscoag = .FALSE. lscnd = .FALSE. lscndgas = .FALSE. lscndh2oae = .FALSE. lsdepo = .FALSE. lsdepo_vege = .FALSE. lsdepo_topo = .FALSE. lsdistupdate = .TRUE. CASE(2) !< Spinup period lscoag = ( .FALSE. .AND. nlcoag ) lscnd = ( .TRUE. .AND. nlcnd ) lscndgas = ( .TRUE. .AND. nlcndgas ) lscndh2oae = ( .TRUE. .AND. nlcndh2oae ) CASE(3) !< Run lscoag = nlcoag lscnd = nlcnd lscndgas = nlcndgas lscndh2oae = nlcndh2oae lsdepo = nldepo lsdepo_vege = nldepo_vege lsdepo_topo = nldepo_topo lsdistupdate = nldistupdate END SELECT END SUBROUTINE set_salsa_runtime !------------------------------------------------------------------------------! ! Description: ! ------------ !> Calculates the absolute temperature (using hydrostatic pressure), saturation !> vapour pressure and mixing ratio over water, relative humidity and air !> density needed in the SALSA model. !> NOTE, no saturation adjustment takes place -> the resulting water vapour !> mixing ratio can be supersaturated, allowing the microphysical calculations !> in SALSA. ! !> Juha Tonttila, FMI, 2014 (original SALSAthrm) !> Mona Kurppa, UHel, 2017 (adjustment for PALM and only aerosol processes) !------------------------------------------------------------------------------! SUBROUTINE salsa_thrm_ij( i, j, p_ij, temp_ij, cw_ij, cs_ij, adn_ij ) USE arrays_3d, & ONLY: p, pt, q, zu USE basic_constants_and_equations_mod, & ONLY: barometric_formula, exner_function, ideal_gas_law_rho, magnus USE control_parameters, & ONLY: pt_surface, surface_pressure IMPLICIT NONE INTEGER(iwp), INTENT(in) :: i INTEGER(iwp), INTENT(in) :: j REAL(wp), DIMENSION(:), INTENT(inout) :: adn_ij REAL(wp), DIMENSION(:), INTENT(inout) :: p_ij REAL(wp), DIMENSION(:), INTENT(inout) :: temp_ij REAL(wp), DIMENSION(:), INTENT(inout), OPTIONAL :: cw_ij REAL(wp), DIMENSION(:), INTENT(inout), OPTIONAL :: cs_ij REAL(wp), DIMENSION(nzb:nzt+1) :: e_s !< saturation vapour pressure !< over water (Pa) REAL(wp) :: t_surface !< absolute surface temperature (K) ! !-- Pressure p_ijk (Pa) = hydrostatic pressure + perturbation pressure (p) t_surface = pt_surface * exner_function( surface_pressure * 100.0_wp ) p_ij(:) = barometric_formula( zu, t_surface, surface_pressure * 100.0_wp ) & + p(:,j,i) ! !-- Absolute ambient temperature (K) temp_ij(:) = pt(:,j,i) * exner_function( p_ij(:) ) ! !-- Air density adn_ij(:) = ideal_gas_law_rho( p_ij(:), temp_ij(:) ) ! !-- Water vapour concentration r_v (kg/m3) IF ( PRESENT( cw_ij ) ) THEN cw_ij(:) = ( q(:,j,i) / ( 1.0_wp - q(:,j,i) ) ) * adn_ij(:) ENDIF ! !-- Saturation mixing ratio r_s (kg/kg) from vapour pressure at temp (Pa) IF ( PRESENT( cs_ij ) ) THEN e_s(:) = 611.0_wp * EXP( alv_d_rv * ( 3.6609E-3_wp - 1.0_wp / & temp_ij(:) ) )! magnus( temp_ij(:) ) cs_ij(:) = ( 0.622_wp * e_s / ( p_ij(:) - e_s(:) ) ) * adn_ij(:) ENDIF END SUBROUTINE salsa_thrm_ij !------------------------------------------------------------------------------! ! Description: ! ------------ !> Calculates ambient sizes of particles by equilibrating soluble fraction of !> particles with water using the ZSR method (Stokes and Robinson, 1966). !> Method: !> Following chemical components are assumed water-soluble !> - (ammonium) sulphate (100%) !> - sea salt (100 %) !> - organic carbon (epsoc * 100%) !> Exact thermodynamic considerations neglected. !> - If particles contain no sea salt, calculation according to sulphate !> properties !> - If contain sea salt but no sulphate, calculation according to sea salt !> properties !> - If contain both sulphate and sea salt -> the molar fraction of these !> compounds determines which one of them is used as the basis of calculation. !> If sulphate and sea salt coexist in a particle, it is assumed that the Cl is !> replaced by sulphate; thus only either sulphate + organics or sea salt + !> organics is included in the calculation of soluble fraction. !> Molality parameterizations taken from Table 1 of Tang: Thermodynamic and !> optical properties of mixed-salt aerosols of atmospheric importance, !> J. Geophys. Res., 102 (D2), 1883-1893 (1997) ! !> Coded by: !> Hannele Korhonen (FMI) 2005 !> Harri Kokkola (FMI) 2006 !> Matti Niskanen(FMI) 2012 !> Anton Laakso (FMI) 2013 !> Modified for the new aerosol datatype, Juha Tonttila (FMI) 2014 ! !> fxm: should sea salt form a solid particle when prh is very low (even though !> it could be mixed with e.g. sulphate)? !> fxm: crashes if no sulphate or sea salt !> fxm: do we really need to consider Kelvin effect for subrange 2 !------------------------------------------------------------------------------! SUBROUTINE equilibration( prh, ptemp, paero, init ) IMPLICIT NONE ! !-- Input variables LOGICAL, INTENT(in) :: init !< TRUE: Initialization call !< FALSE: Normal runtime: update water !< content only for 1a REAL(wp), INTENT(in) :: prh !< relative humidity [0-1] REAL(wp), INTENT(in) :: ptemp !< temperature (K) ! !-- Output variables TYPE(t_section), INTENT(inout) :: paero(fn2b) ! !-- Local INTEGER(iwp) :: b !< loop index INTEGER(iwp) :: counti !< loop index REAL(wp) :: zaw !< water activity [0-1] REAL(wp) :: zbinmol(7) !< binary molality of each components (mol/kg) REAL(wp) :: zcore !< Volume of dry particle REAL(wp) :: zdold !< Old diameter REAL(wp) :: zdwet !< Wet diameter or mean droplet diameter REAL(wp) :: zke !< Kelvin term in the Köhler equation REAL(wp) :: zlwc !< liquid water content [kg/m3-air] REAL(wp) :: zrh !< Relative humidity REAL(wp) :: zvpart(7) !< volume of chem. compounds in one particle zaw = 0.0_wp zbinmol = 0.0_wp zcore = 0.0_wp zdold = 0.0_wp zdwet = 0.0_wp zlwc = 0.0_wp zrh = 0.0_wp ! !-- Relative humidity: zrh = prh zrh = MAX( zrh, 0.05_wp ) zrh = MIN( zrh, 0.98_wp) ! !-- 1) Regime 1: sulphate and partly water-soluble OC. Done for every CALL DO b = in1a, fn1a ! size bin zbinmol = 0.0_wp zdold = 1.0_wp zke = 1.02_wp IF ( paero(b)%numc > nclim ) THEN ! !-- Volume in one particle zvpart = 0.0_wp zvpart(1:2) = paero(b)%volc(1:2) / paero(b)%numc zvpart(6:7) = paero(b)%volc(6:7) / paero(b)%numc ! !-- Total volume and wet diameter of one dry particle zcore = SUM( zvpart(1:2) ) zdwet = paero(b)%dwet counti = 0 DO WHILE ( ABS( zdwet / zdold - 1.0_wp ) > 1.0E-2_wp ) zdold = MAX( zdwet, 1.0E-20_wp ) zaw = MAX( 1.0E-3_wp, zrh / zke ) ! To avoid underflow ! !-- Binary molalities (mol/kg): !-- Sulphate zbinmol(1) = 1.1065495E+2_wp - 3.6759197E+2_wp * zaw & + 5.0462934E+2_wp * zaw**2.0_wp & - 3.1543839E+2_wp * zaw**3.0_wp & + 6.770824E+1_wp * zaw**4.0_wp !-- Organic carbon zbinmol(2) = 1.0_wp / ( zaw * amh2o ) - 1.0_wp / amh2o !-- Nitric acid zbinmol(6) = 2.306844303E+1_wp - 3.563608869E+1_wp * zaw & - 6.210577919E+1_wp * zaw**2.0_wp & + 5.510176187E+2_wp * zaw**3.0_wp & - 1.460055286E+3_wp * zaw**4.0_wp & + 1.894467542E+3_wp * zaw**5.0_wp & - 1.220611402E+3_wp * zaw**6.0_wp & + 3.098597737E+2_wp * zaw**7.0_wp ! !-- Calculate the liquid water content (kg/m3-air) using ZSR (see e.g. !-- Eq. 10.98 in Seinfeld and Pandis (2006)) zlwc = ( paero(b)%volc(1) * ( arhoh2so4 / amh2so4 ) ) / & zbinmol(1) + epsoc * paero(b)%volc(2) * ( arhooc / amoc ) & / zbinmol(2) + ( paero(b)%volc(6) * ( arhohno3/amhno3 ) ) & / zbinmol(6) ! !-- Particle wet diameter (m) zdwet = ( zlwc / paero(b)%numc / arhoh2o / api6 + & ( SUM( zvpart(6:7) ) / api6 ) + & zcore / api6 )**( 1.0_wp / 3.0_wp ) ! !-- Kelvin effect (Eq. 10.85 in in Seinfeld and Pandis (2006)). Avoid !-- overflow. zke = EXP( MIN( 50.0_wp, & 4.0_wp * surfw0 * amvh2so4 / ( abo * ptemp * zdwet ) ) ) counti = counti + 1 IF ( counti > 1000 ) THEN message_string = 'Subrange 1: no convergence!' CALL message( 'salsa_mod: equilibration', 'SA0042', & 1, 2, 0, 6, 0 ) ENDIF ENDDO ! !-- Instead of lwc, use the volume concentration of water from now on !-- (easy to convert...) paero(b)%volc(8) = zlwc / arhoh2o ! !-- If this is initialization, update the core and wet diameter IF ( init ) THEN paero(b)%dwet = zdwet paero(b)%core = zcore ENDIF ELSE !-- If initialization !-- 1.2) empty bins given bin average values IF ( init ) THEN paero(b)%dwet = paero(b)%dmid paero(b)%core = api6 * paero(b)%dmid ** 3.0_wp ENDIF ENDIF ENDDO !< b ! !-- 2) Regime 2a: sulphate, OC, BC and sea salt !-- This is done only for initialization call, otherwise the water contents !-- are computed via condensation IF ( init ) THEN DO b = in2a, fn2b !-- Initialize zke = 1.02_wp zbinmol = 0.0_wp zdold = 1.0_wp ! !-- 1) Particle properties calculated for non-empty bins IF ( paero(b)%numc > nclim ) THEN ! !-- Volume in one particle [fxm] zvpart = 0.0_wp zvpart(1:7) = paero(b)%volc(1:7) / paero(b)%numc ! !-- Total volume and wet diameter of one dry particle [fxm] zcore = SUM( zvpart(1:5) ) zdwet = paero(b)%dwet counti = 0 DO WHILE ( ABS( zdwet / zdold - 1.0_wp ) > 1.0E-12_wp ) zdold = MAX( zdwet, 1.0E-20_wp ) zaw = zrh / zke ! !-- Binary molalities (mol/kg): !-- Sulphate zbinmol(1) = 1.1065495E+2_wp - 3.6759197E+2_wp * zaw & + 5.0462934E+2_wp * zaw**2 - 3.1543839E+2_wp * zaw**3 & + 6.770824E+1_wp * zaw**4 !-- Organic carbon zbinmol(2) = 1.0_wp / ( zaw * amh2o ) - 1.0_wp / amh2o !-- Nitric acid zbinmol(6) = 2.306844303E+1_wp - 3.563608869E+1_wp * zaw & - 6.210577919E+1_wp * zaw**2 + 5.510176187E+2_wp * zaw**3 & - 1.460055286E+3_wp * zaw**4 + 1.894467542E+3_wp * zaw**5 & - 1.220611402E+3_wp * zaw**6 + 3.098597737E+2_wp * zaw**7 !-- Sea salt (natrium chloride) zbinmol(5) = 5.875248E+1_wp - 1.8781997E+2_wp * zaw & + 2.7211377E+2_wp * zaw**2 - 1.8458287E+2_wp * zaw**3 & + 4.153689E+1_wp * zaw**4 ! !-- Calculate the liquid water content (kg/m3-air) zlwc = ( paero(b)%volc(1) * ( arhoh2so4 / amh2so4 ) ) / & zbinmol(1) + epsoc * ( paero(b)%volc(2) * ( arhooc / & amoc ) ) / zbinmol(2) + ( paero(b)%volc(6) * ( arhohno3 & / amhno3 ) ) / zbinmol(6) + ( paero(b)%volc(5) * & ( arhoss / amss ) ) / zbinmol(5) !-- Particle wet radius (m) zdwet = ( zlwc / paero(b)%numc / arhoh2o / api6 + & ( SUM( zvpart(6:7) ) / api6 ) + & zcore / api6 ) ** ( 1.0_wp / 3.0_wp ) ! !-- Kelvin effect (Eq. 10.85 in Seinfeld and Pandis (2006)) zke = EXP( MIN( 50.0_wp, & 4.0_wp * surfw0 * amvh2so4 / ( abo * zdwet * ptemp ) ) ) counti = counti + 1 IF ( counti > 1000 ) THEN message_string = 'Subrange 2: no convergence!' CALL message( 'salsa_mod: equilibration', 'SA0043', & 1, 2, 0, 6, 0 ) ENDIF ENDDO ! !-- Liquid water content; instead of LWC use the volume concentration paero(b)%volc(8) = zlwc / arhoh2o paero(b)%dwet = zdwet paero(b)%core = zcore ELSE !-- 2.2) empty bins given bin average values paero(b)%dwet = paero(b)%dmid paero(b)%core = api6 * paero(b)%dmid ** 3.0_wp ENDIF ENDDO ! b ENDIF END SUBROUTINE equilibration !------------------------------------------------------------------------------! !> Description: !> ------------ !> Calculation of the settling velocity vc (m/s) per aerosol size bin and !> deposition on plant canopy (lsdepo_vege). ! !> Deposition is based on either the scheme presented in: !> Zhang et al. (2001), Atmos. Environ. 35, 549-560 (includes collection due to !> Brownian diffusion, impaction, interception and sedimentation) !> OR !> Petroff & Zhang (2010), Geosci. Model Dev. 3, 753-769 (includes also !> collection due to turbulent impaction) ! !> Equation numbers refer to equation in Jacobson (2005): Fundamentals of !> Atmospheric Modeling, 2nd Edition. ! !> Subroutine follows closely sedim_SALSA in UCLALES-SALSA written by Juha !> Tonttila (KIT/FMI) and Zubair Maalick (UEF). !> Rewritten to PALM by Mona Kurppa (UH), 2017. ! !> Call for grid point i,j,k !------------------------------------------------------------------------------! SUBROUTINE deposition( paero, tk, adn, mag_u, lad, kvis, Sc, vc ) USE plant_canopy_model_mod, & ONLY: cdc IMPLICIT NONE REAL(wp), INTENT(in) :: adn !< air density (kg/m3) REAL(wp), INTENT(out) :: kvis !< kinematic viscosity of air (m2/s) REAL(wp), INTENT(in) :: lad !< leaf area density (m2/m3) REAL(wp), INTENT(in) :: mag_u !< wind velocity (m/s) REAL(wp), INTENT(out) :: Sc(:) !< particle Schmidt number REAL(wp), INTENT(in) :: tk !< abs.temperature (K) REAL(wp), INTENT(out) :: vc(:) !< critical fall speed i.e. settling !< velocity of an aerosol particle (m/s) TYPE(t_section), INTENT(inout) :: paero(fn2b) INTEGER(iwp) :: b !< loop index INTEGER(iwp) :: c !< loop index REAL(wp) :: avis !< molecular viscocity of air (kg/(m*s)) REAL(wp), PARAMETER :: c_A = 1.249_wp !< Constants A, B and C for REAL(wp), PARAMETER :: c_B = 0.42_wp !< calculating the Cunningham REAL(wp), PARAMETER :: c_C = 0.87_wp !< slip-flow correction (Cc) !< according to Jacobson (2005), !< Eq. 15.30 REAL(wp) :: Cc !< Cunningham slip-flow correction factor REAL(wp) :: Kn !< Knudsen number REAL(wp) :: lambda !< molecular mean free path (m) REAL(wp) :: mdiff !< particle diffusivity coefficient REAL(wp) :: pdn !< particle density (kg/m3) REAL(wp) :: ustar !< friction velocity (m/s) REAL(wp) :: va !< thermal speed of an air molecule (m/s) REAL(wp) :: zdwet !< wet diameter (m) ! !-- Initialise Cc = 0.0_wp Kn = 0.0_wp mdiff = 0.0_wp pdn = 1500.0_wp ! default value ustar = 0.0_wp ! !-- Molecular viscosity of air (Eq. 4.54) avis = 1.8325E-5_wp * ( 416.16_wp / ( tk + 120.0_wp ) ) * ( tk / & 296.16_wp )**1.5_wp ! !-- Kinematic viscosity (Eq. 4.55) kvis = avis / adn ! !-- Thermal velocity of an air molecule (Eq. 15.32) va = SQRT( 8.0_wp * abo * tk / ( pi * am_airmol ) ) ! !-- Mean free path (m) (Eq. 15.24) lambda = 2.0_wp * avis / ( adn * va ) DO b = 1, nbins IF ( paero(b)%numc < nclim ) CYCLE zdwet = paero(b)%dwet ! !-- Knudsen number (Eq. 15.23) Kn = MAX( 1.0E-2_wp, lambda / ( zdwet * 0.5_wp ) ) ! To avoid underflow ! !-- Cunningham slip-flow correction (Eq. 15.30) Cc = 1.0_wp + Kn * ( c_A + c_B * EXP( -c_C / Kn ) ) !-- Particle diffusivity coefficient (Eq. 15.29) mdiff = ( abo * tk * Cc ) / ( 3.0_wp * pi * avis * zdwet ) ! !-- Particle Schmidt number (Eq. 15.36) Sc(b) = kvis / mdiff ! !-- Critical fall speed i.e. settling velocity (Eq. 20.4) vc(b) = MIN( 1.0_wp, terminal_vel( 0.5_wp * zdwet, pdn, adn, avis, Cc) ) IF ( lsdepo_vege .AND. plant_canopy .AND. lad > 0.0_wp ) THEN ! !-- Friction velocity calculated following Prandtl (1925): ustar = SQRT( cdc ) * mag_u CALL depo_vege( paero, b, vc(b), mag_u, ustar, kvis, Sc(b), lad ) ENDIF ENDDO END SUBROUTINE deposition !------------------------------------------------------------------------------! ! Description: ! ------------ !> Calculate change in number and volume concentrations due to deposition on !> plant canopy. !------------------------------------------------------------------------------! SUBROUTINE depo_vege( paero, b, vc, mag_u, ustar, kvis_a, Sc, lad ) IMPLICIT NONE INTEGER(iwp), INTENT(in) :: b !< loop index REAL(wp), INTENT(in) :: kvis_a !< kinematic viscosity of air (m2/s) REAL(wp), INTENT(in) :: lad !< leaf area density (m2/m3) REAL(wp), INTENT(in) :: mag_u !< wind velocity (m/s) REAL(wp), INTENT(in) :: Sc !< particle Schmidt number REAL(wp), INTENT(in) :: ustar !< friction velocity (m/s) REAL(wp), INTENT(in) :: vc !< terminal velocity (m/s) TYPE(t_section), INTENT(inout) :: paero(fn2b) INTEGER(iwp) :: c !< loop index REAL(wp), PARAMETER :: c_A = 1.249_wp !< Constants A, B and C for REAL(wp), PARAMETER :: c_B = 0.42_wp !< calculating the Cunningham REAL(wp), PARAMETER :: c_C = 0.87_wp !< slip-flow correction (Cc) !< according to Jacobson (2005), !< Eq. 15.30 REAL(wp) :: alpha !< parameter, Table 3 in Zhang et al. (2001) REAL(wp) :: depo !< deposition efficiency REAL(wp) :: C_Br !< coefficient for Brownian diffusion REAL(wp) :: C_IM !< coefficient for inertial impaction REAL(wp) :: C_IN !< coefficient for interception REAL(wp) :: C_IT !< coefficient for turbulent impaction REAL(wp) :: gamma !< parameter, Table 3 in Zhang et al. (2001) REAL(wp) :: par_A !< parameter A for the characteristic radius of !< collectors, Table 3 in Zhang et al. (2001) REAL(wp) :: rt !< the overall quasi-laminar resistance for !< particles REAL(wp) :: St !< Stokes number for smooth surfaces or bluff !< surface elements REAL(wp) :: tau_plus !< dimensionless particle relaxation time REAL(wp) :: v_bd !< deposition velocity due to Brownian diffusion REAL(wp) :: v_im !< deposition velocity due to impaction REAL(wp) :: v_in !< deposition velocity due to interception REAL(wp) :: v_it !< deposition velocity due to turbulent impaction ! !-- Initialise depo = 0.0_wp rt = 0.0_wp St = 0.0_wp tau_plus = 0.0_wp v_bd = 0.0_wp v_im = 0.0_wp v_in = 0.0_wp v_it = 0.0_wp IF ( depo_vege_type == 'zhang2001' ) THEN ! !-- Parameters for the land use category 'deciduous broadleaf trees'(Table 3) par_A = 5.0E-3_wp alpha = 0.8_wp gamma = 0.56_wp ! !-- Stokes number for vegetated surfaces (Seinfeld & Pandis (2006): Eq.19.24) St = vc * ustar / ( g * par_A ) ! !-- The overall quasi-laminar resistance for particles (Zhang et al., Eq. 5) rt = MAX( EPSILON( 1.0_wp ), ( 3.0_wp * ustar * EXP( -St**0.5_wp ) * & ( Sc**( -gamma ) + ( St / ( alpha + St ) )**2.0_wp + & 0.5_wp * ( paero(b)%dwet / par_A )**2.0_wp ) ) ) depo = ( rt + vc ) * lad paero(b)%numc = paero(b)%numc - depo * paero(b)%numc * dt_salsa DO c = 1, maxspec+1 paero(b)%volc(c) = paero(b)%volc(c) - depo * paero(b)%volc(c) * & dt_salsa ENDDO ELSEIF ( depo_vege_type == 'petroff2010' ) THEN ! !-- vd = v_BD + v_IN + v_IM + v_IT + vc !-- Deposition efficiencies from Table 1. Constants from Table 2. C_Br = 1.262_wp C_IM = 0.130_wp C_IN = 0.216_wp C_IT = 0.056_wp par_A = 0.03_wp ! Here: leaf width (m) ! !-- Stokes number for vegetated surfaces (Seinfeld & Pandis (2006): Eq.19.24) St = vc * ustar / ( g * par_A ) ! !-- Non-dimensional relexation time of the particle on top of canopy tau_plus = vc * ustar**2.0_wp / ( kvis_a * g ) ! !-- Brownian diffusion v_bd = mag_u * C_Br * Sc**( -2.0_wp / 3.0_wp ) * & ( mag_u * par_A / kvis_a )**( -0.5_wp ) ! !-- Interception v_in = mag_u * C_IN * paero(b)%dwet / par_A * ( 2.0_wp + LOG( 2.0_wp * & par_A / paero(b)%dwet ) ) ! !-- Impaction: Petroff (2009) Eq. 18 v_im = mag_u * C_IM * ( St / ( St + 0.47_wp ) )**2.0_wp IF ( tau_plus < 20.0_wp ) THEN v_it = 2.5E-3_wp * C_IT * tau_plus**2.0_wp ELSE v_it = C_IT ENDIF depo = ( v_bd + v_in + v_im + v_it + vc ) * lad paero(b)%numc = paero(b)%numc - depo * paero(b)%numc * dt_salsa DO c = 1, maxspec+1 paero(b)%volc(c) = paero(b)%volc(c) - depo * paero(b)%volc(c) * & dt_salsa ENDDO ENDIF END SUBROUTINE depo_vege !------------------------------------------------------------------------------! ! Description: ! ------------ !> Calculate deposition on horizontal and vertical surfaces. Implement as !> surface flux. !------------------------------------------------------------------------------! SUBROUTINE depo_topo( i, j, surf, vc, Sc, kvis, mag_u, norm ) USE surface_mod, & ONLY: surf_type IMPLICIT NONE INTEGER(iwp), INTENT(in) :: i !< loop index INTEGER(iwp), INTENT(in) :: j !< loop index REAL(wp), INTENT(in) :: kvis(:) !< kinematic viscosity of air (m2/s) REAL(wp), INTENT(in) :: mag_u(:) !< wind velocity (m/s) REAL(wp), INTENT(in) :: norm(:) !< normalisation (usually air density) REAL(wp), INTENT(in) :: Sc(:,:) !< particle Schmidt number REAL(wp), INTENT(in) :: vc(:,:) !< terminal velocity (m/s) TYPE(surf_type), INTENT(inout) :: surf !< respective surface type INTEGER(iwp) :: b !< loop index INTEGER(iwp) :: c !< loop index INTEGER(iwp) :: k !< loop index INTEGER(iwp) :: m !< loop index INTEGER(iwp) :: surf_e !< End index of surface elements at (j,i)-gridpoint INTEGER(iwp) :: surf_s !< Start index of surface elements at (j,i)-gridpoint REAL(wp) :: alpha !< parameter, Table 3 in Zhang et al. (2001) REAL(wp) :: C_Br !< coefficient for Brownian diffusion REAL(wp) :: C_IM !< coefficient for inertial impaction REAL(wp) :: C_IN !< coefficient for interception REAL(wp) :: C_IT !< coefficient for turbulent impaction REAL(wp) :: depo !< deposition efficiency REAL(wp) :: gamma !< parameter, Table 3 in Zhang et al. (2001) REAL(wp) :: par_A !< parameter A for the characteristic radius of !< collectors, Table 3 in Zhang et al. (2001) REAL(wp) :: rt !< the overall quasi-laminar resistance for !< particles REAL(wp) :: St !< Stokes number for bluff surface elements REAL(wp) :: tau_plus !< dimensionless particle relaxation time REAL(wp) :: v_bd !< deposition velocity due to Brownian diffusion REAL(wp) :: v_im !< deposition velocity due to impaction REAL(wp) :: v_in !< deposition velocity due to interception REAL(wp) :: v_it !< deposition velocity due to turbulent impaction ! !-- Initialise rt = 0.0_wp St = 0.0_wp tau_plus = 0.0_wp v_bd = 0.0_wp v_im = 0.0_wp v_in = 0.0_wp v_it = 0.0_wp surf_s = surf%start_index(j,i) surf_e = surf%end_index(j,i) DO m = surf_s, surf_e k = surf%k(m) DO b = 1, nbins IF ( aerosol_number(b)%conc(k,j,i) <= nclim .OR. & Sc(k+1,b) < 1.0_wp ) CYCLE IF ( depo_topo_type == 'zhang2001' ) THEN ! !-- Parameters for the land use category 'urban' in Table 3 alpha = 1.5_wp gamma = 0.56_wp par_A = 10.0E-3_wp ! !-- Stokes number for smooth surfaces or surfaces with bluff roughness !-- elements (Seinfeld and Pandis, 2nd edition (2006): Eq. 19.23) St = MAX( 0.01_wp, vc(k+1,b) * surf%us(m) ** 2.0_wp / & ( g * kvis(k+1) ) ) ! !-- The overall quasi-laminar resistance for particles (Eq. 5) rt = MAX( EPSILON( 1.0_wp ), ( 3.0_wp * surf%us(m) * ( & Sc(k+1,b)**( -gamma ) + ( St / ( alpha + St ) )**2.0_wp & + 0.5_wp * ( Ra_dry(k,j,i,b) / par_A )**2.0_wp ) * & EXP( -St**0.5_wp ) ) ) depo = vc(k+1,b) + rt ELSEIF ( depo_topo_type == 'petroff2010' ) THEN ! !-- vd = v_BD + v_IN + v_IM + v_IT + vc !-- Deposition efficiencies from Table 1. Constants from Table 2. C_Br = 1.262_wp C_IM = 0.130_wp C_IN = 0.216_wp C_IT = 0.056_wp par_A = 0.03_wp ! Here: leaf width (m) ! !-- Stokes number for smooth surfaces or surfaces with bluff roughness !-- elements (Seinfeld and Pandis, 2nd edition (2006): Eq. 19.23) St = MAX( 0.01_wp, vc(k+1,b) * surf%us(m) ** 2.0_wp / & ( g * kvis(k+1) ) ) ! !-- Non-dimensional relexation time of the particle on top of canopy tau_plus = vc(k+1,b) * surf%us(m)**2.0_wp / ( kvis(k+1) * g ) ! !-- Brownian diffusion v_bd = mag_u(k+1) * C_Br * Sc(k+1,b)**( -2.0_wp / 3.0_wp ) * & ( mag_u(k+1) * par_A / kvis(k+1) )**( -0.5_wp ) ! !-- Interception v_in = mag_u(k+1) * C_IN * Ra_dry(k,j,i,b)/ par_A * ( 2.0_wp + & LOG( 2.0_wp * par_A / Ra_dry(k,j,i,b) ) ) ! !-- Impaction: Petroff (2009) Eq. 18 v_im = mag_u(k+1) * C_IM * ( St / ( St + 0.47_wp ) )**2.0_wp IF ( tau_plus < 20.0_wp ) THEN v_it = 2.5E-3_wp * C_IT * tau_plus**2.0_wp ELSE v_it = C_IT ENDIF depo = v_bd + v_in + v_im + v_it + vc(k+1,b) ENDIF IF ( lod_aero == 3 .OR. salsa_source_mode == 'no_source' ) THEN surf%answs(m,b) = -depo * norm(k) * aerosol_number(b)%conc(k,j,i) DO c = 1, ncc_tot surf%amsws(m,(c-1)*nbins+b) = -depo * norm(k) * & aerosol_mass((c-1)*nbins+b)%conc(k,j,i) ENDDO ! c ELSE surf%answs(m,b) = SUM( aerosol_number(b)%source(:,j,i) ) - & MAX( 0.0_wp, depo * norm(k) * & aerosol_number(b)%conc(k,j,i) ) DO c = 1, ncc_tot surf%amsws(m,(c-1)*nbins+b) = SUM( & aerosol_mass((c-1)*nbins+b)%source(:,j,i) ) - & MAX( 0.0_wp, depo * norm(k) * & aerosol_mass((c-1)*nbins+b)%conc(k,j,i) ) ENDDO ENDIF ENDDO ! b ENDDO ! m END SUBROUTINE depo_topo !------------------------------------------------------------------------------! ! Description: ! ------------ ! Function for calculating terminal velocities for different particles sizes. !------------------------------------------------------------------------------! REAL(wp) FUNCTION terminal_vel( radius, rhop, rhoa, visc, beta ) IMPLICIT NONE REAL(wp), INTENT(in) :: beta !< Cunningham correction factor REAL(wp), INTENT(in) :: radius !< particle radius (m) REAL(wp), INTENT(in) :: rhop !< particle density (kg/m3) REAL(wp), INTENT(in) :: rhoa !< air density (kg/m3) REAL(wp), INTENT(in) :: visc !< molecular viscosity of air (kg/(m*s)) REAL(wp), PARAMETER :: rhoa_ref = 1.225_wp ! reference air density (kg/m3) ! !-- Stokes law with Cunningham slip correction factor terminal_vel = ( 4.0_wp * radius**2.0_wp ) * ( rhop - rhoa ) * g * beta / & ( 18.0_wp * visc ) ! (m/s) END FUNCTION terminal_vel !------------------------------------------------------------------------------! ! Description: ! ------------ !> Calculates particle loss and change in size distribution due to (Brownian) !> coagulation. Only for particles with dwet < 30 micrometres. ! !> Method: !> Semi-implicit, non-iterative method: (Jacobson, 1994) !> Volume concentrations of the smaller colliding particles added to the bin of !> the larger colliding particles. Start from first bin and use the updated !> number and volume for calculation of following bins. NB! Our bin numbering !> does not follow particle size in subrange 2. ! !> Schematic for bin numbers in different subranges: !> 1 2 !> +-------------------------------------------+ !> a | 1 | 2 | 3 || 4 | 5 | 6 | 7 | 8 | 9 | 10|| !> b | ||11 |12 |13 |14 | 15 | 16 | 17|| !> +-------------------------------------------+ ! !> Exact coagulation coefficients for each pressure level are scaled according !> to current particle wet size (linear scaling). !> Bins are organized in terms of the dry size of the condensation nucleus, !> while coagulation kernell is calculated with the actual hydrometeor !> size. ! !> Called from salsa_driver !> fxm: Process selection should be made smarter - now just lots of IFs inside !> loops ! !> Coded by: !> Hannele Korhonen (FMI) 2005 !> Harri Kokkola (FMI) 2006 !> Tommi Bergman (FMI) 2012 !> Matti Niskanen(FMI) 2012 !> Anton Laakso (FMI) 2013 !> Juha Tonttila (FMI) 2014 !------------------------------------------------------------------------------! SUBROUTINE coagulation( paero, ptstep, ptemp, ppres ) IMPLICIT NONE !-- Input and output variables TYPE(t_section), INTENT(inout) :: paero(fn2b) !< Aerosol properties REAL(wp), INTENT(in) :: ppres !< ambient pressure (Pa) REAL(wp), INTENT(in) :: ptemp !< ambient temperature (K) REAL(wp), INTENT(in) :: ptstep !< time step (s) !-- Local variables INTEGER(iwp) :: index_2a !< corresponding bin in subrange 2a INTEGER(iwp) :: index_2b !< corresponding bin in subrange 2b INTEGER(iwp) :: b !< loop index INTEGER(iwp) :: ll !< loop index INTEGER(iwp) :: mm !< loop index INTEGER(iwp) :: nn !< loop index REAL(wp) :: pressi !< pressure REAL(wp) :: temppi !< temperature REAL(wp) :: zcc(fn2b,fn2b) !< updated coagulation coefficients (m3/s) REAL(wp) :: zdpart_mm !< diameter of particle (m) REAL(wp) :: zdpart_nn !< diameter of particle (m) REAL(wp) :: zminusterm !< coagulation loss in a bin (1/s) REAL(wp) :: zplusterm(8) !< coagulation gain in a bin (fxm/s) !< (for each chemical compound) REAL(wp) :: zmpart(fn2b) !< approximate mass of particles (kg) zcc = 0.0_wp zmpart = 0.0_wp zdpart_mm = 0.0_wp zdpart_nn = 0.0_wp ! !-- 1) Coagulation to coarse mode calculated in a simplified way: !-- CoagSink ~ Dp in continuum subrange, thus we calculate 'effective' !-- number concentration of coarse particles !-- 2) Updating coagulation coefficients ! !-- Aerosol mass (kg). Density of 1500 kg/m3 assumed zmpart(1:fn2b) = api6 * ( MIN( paero(1:fn2b)%dwet, 30.0E-6_wp )**3.0_wp ) & * 1500.0_wp temppi = ptemp pressi = ppres zcc = 0.0_wp ! !-- Aero-aero coagulation DO mm = 1, fn2b ! smaller colliding particle IF ( paero(mm)%numc < nclim ) CYCLE DO nn = mm, fn2b ! larger colliding particle IF ( paero(nn)%numc < nclim ) CYCLE zdpart_mm = MIN( paero(mm)%dwet, 30.0E-6_wp ) ! Limit to 30 um zdpart_nn = MIN( paero(nn)%dwet, 30.0E-6_wp ) ! Limit to 30 um ! !-- Coagulation coefficient of particles (m3/s) zcc(mm,nn) = coagc( zdpart_mm, zdpart_nn, zmpart(mm), zmpart(nn), & temppi, pressi ) zcc(nn,mm) = zcc(mm,nn) ENDDO ENDDO ! !-- 3) New particle and volume concentrations after coagulation: !-- Calculated according to Jacobson (2005) eq. 15.9 ! !-- Aerosols in subrange 1a: DO b = in1a, fn1a IF ( paero(b)%numc < nclim ) CYCLE zminusterm = 0.0_wp zplusterm(:) = 0.0_wp ! !-- Particles lost by coagulation with larger aerosols DO ll = b+1, fn2b zminusterm = zminusterm + zcc(b,ll) * paero(ll)%numc ENDDO ! !-- Coagulation gain in a bin: change in volume conc. (cm3/cm3): DO ll = in1a, b-1 zplusterm(1:2) = zplusterm(1:2) + zcc(ll,b) * paero(ll)%volc(1:2) zplusterm(6:7) = zplusterm(6:7) + zcc(ll,b) * paero(ll)%volc(6:7) zplusterm(8) = zplusterm(8) + zcc(ll,b) * paero(ll)%volc(8) ENDDO ! !-- Volume and number concentrations after coagulation update [fxm] paero(b)%volc(1:2) = ( paero(b)%volc(1:2) + ptstep * zplusterm(1:2) * & paero(b)%numc ) / ( 1.0_wp + ptstep * zminusterm ) paero(b)%volc(6:7) = ( paero(b)%volc(6:7) + ptstep * zplusterm(6:7) * & paero(b)%numc ) / ( 1.0_wp + ptstep * zminusterm ) paero(b)%volc(8) = ( paero(b)%volc(8) + ptstep * zplusterm(8) * & paero(b)%numc ) / ( 1.0_wp + ptstep * zminusterm ) paero(b)%numc = paero(b)%numc / ( 1.0_wp + ptstep * zminusterm + & 0.5_wp * ptstep * zcc(b,b) * paero(b)%numc ) ENDDO ! !-- Aerosols in subrange 2a: DO b = in2a, fn2a IF ( paero(b)%numc < nclim ) CYCLE zminusterm = 0.0_wp zplusterm(:) = 0.0_wp ! !-- Find corresponding size bin in subrange 2b index_2b = b - in2a + in2b ! !-- Particles lost by larger particles in 2a DO ll = b+1, fn2a zminusterm = zminusterm + zcc(b,ll) * paero(ll)%numc ENDDO ! !-- Particles lost by larger particles in 2b IF ( .NOT. no_insoluble ) THEN DO ll = index_2b+1, fn2b zminusterm = zminusterm + zcc(b,ll) * paero(ll)%numc ENDDO ENDIF ! !-- Particle volume gained from smaller particles in subranges 1, 2a and 2b DO ll = in1a, b-1 zplusterm(1:2) = zplusterm(1:2) + zcc(ll,b) * paero(ll)%volc(1:2) zplusterm(6:7) = zplusterm(6:7) + zcc(ll,b) * paero(ll)%volc(6:7) zplusterm(8) = zplusterm(8) + zcc(ll,b) * paero(ll)%volc(8) ENDDO ! !-- Particle volume gained from smaller particles in 2a !-- (Note, for components not included in the previous loop!) DO ll = in2a, b-1 zplusterm(3:5) = zplusterm(3:5) + zcc(ll,b)*paero(ll)%volc(3:5) ENDDO ! !-- Particle volume gained from smaller (and equal) particles in 2b IF ( .NOT. no_insoluble ) THEN DO ll = in2b, index_2b zplusterm(1:8) = zplusterm(1:8) + zcc(ll,b) * paero(ll)%volc(1:8) ENDDO ENDIF ! !-- Volume and number concentrations after coagulation update [fxm] paero(b)%volc(1:8) = ( paero(b)%volc(1:8) + ptstep * zplusterm(1:8) * & paero(b)%numc ) / ( 1.0_wp + ptstep * zminusterm ) paero(b)%numc = paero(b)%numc / ( 1.0_wp + ptstep * zminusterm + & 0.5_wp * ptstep * zcc(b,b) * paero(b)%numc ) ENDDO ! !-- Aerosols in subrange 2b: IF ( .NOT. no_insoluble ) THEN DO b = in2b, fn2b IF ( paero(b)%numc < nclim ) CYCLE zminusterm = 0.0_wp zplusterm(:) = 0.0_wp ! !-- Find corresponding size bin in subsubrange 2a index_2a = b - in2b + in2a ! !-- Particles lost to larger particles in subranges 2b DO ll = b+1, fn2b zminusterm = zminusterm + zcc(b,ll) * paero(ll)%numc ENDDO ! !-- Particles lost to larger and equal particles in 2a DO ll = index_2a, fn2a zminusterm = zminusterm + zcc(b,ll) * paero(ll)%numc ENDDO ! !-- Particle volume gained from smaller particles in subranges 1 & 2a DO ll = in1a, index_2a-1 zplusterm(1:8) = zplusterm(1:8) + zcc(ll,b) * paero(ll)%volc(1:8) ENDDO ! !-- Particle volume gained from smaller particles in 2b DO ll = in2b, b-1 zplusterm(1:8) = zplusterm(1:8) + zcc(ll,b) * paero(ll)%volc(1:8) ENDDO ! !-- Volume and number concentrations after coagulation update [fxm] paero(b)%volc(1:8) = ( paero(b)%volc(1:8) + ptstep * zplusterm(1:8)& * paero(b)%numc ) / ( 1.0_wp + ptstep * zminusterm ) paero(b)%numc = paero(b)%numc / ( 1.0_wp + ptstep * zminusterm + & 0.5_wp * ptstep * zcc(b,b) * paero(b)%numc ) ENDDO ENDIF END SUBROUTINE coagulation !------------------------------------------------------------------------------! ! Description: ! ------------ !> Calculation of coagulation coefficients. Extended version of the function !> originally found in mo_salsa_init. ! !> J. Tonttila, FMI, 05/2014 !------------------------------------------------------------------------------! REAL(wp) FUNCTION coagc( diam1, diam2, mass1, mass2, temp, pres ) IMPLICIT NONE ! !-- Input and output variables REAL(wp), INTENT(in) :: diam1 !< diameter of colliding particle 1 (m) REAL(wp), INTENT(in) :: diam2 !< diameter of colliding particle 2 (m) REAL(wp), INTENT(in) :: mass1 !< mass of colliding particle 1 (kg) REAL(wp), INTENT(in) :: mass2 !< mass of colliding particle 2 (kg) REAL(wp), INTENT(in) :: pres !< ambient pressure (Pa?) [fxm] REAL(wp), INTENT(in) :: temp !< ambient temperature (K) ! !-- Local variables REAL(wp) :: fmdist !< distance of flux matching (m) REAL(wp) :: knud_p !< particle Knudsen number REAL(wp) :: mdiam !< mean diameter of colliding particles (m) REAL(wp) :: mfp !< mean free path of air molecules (m) REAL(wp) :: visc !< viscosity of air (kg/(m s)) REAL(wp), DIMENSION (2) :: beta !< Cunningham correction factor REAL(wp), DIMENSION (2) :: dfpart !< particle diffusion coefficient !< (m2/s) REAL(wp), DIMENSION (2) :: diam !< diameters of particles (m) REAL(wp), DIMENSION (2) :: flux !< flux in continuum and free molec. !< regime (m/s) REAL(wp), DIMENSION (2) :: knud !< particle Knudsen number REAL(wp), DIMENSION (2) :: mpart !< masses of particles (kg) REAL(wp), DIMENSION (2) :: mtvel !< particle mean thermal velocity (m/s) REAL(wp), DIMENSION (2) :: omega !< particle mean free path REAL(wp), DIMENSION (2) :: tva !< temporary variable (m) ! !-- Initialisation coagc = 0.0_wp ! !-- 1) Initializing particle and ambient air variables diam = (/ diam1, diam2 /) !< particle diameters (m) mpart = (/ mass1, mass2 /) !< particle masses (kg) !-- Viscosity of air (kg/(m s)) visc = ( 7.44523E-3_wp * temp ** 1.5_wp ) / & ( 5093.0_wp * ( temp + 110.4_wp ) ) !-- Mean free path of air (m) mfp = ( 1.656E-10_wp * temp + 1.828E-8_wp ) * ( p_0 + 1325.0_wp ) / pres ! !-- 2) Slip correction factor for small particles knud = 2.0_wp * EXP( LOG(mfp) - LOG(diam) )! Knudsen number for air (15.23) !-- Cunningham correction factor (Allen and Raabe, Aerosol Sci. Tech. 4, 269) beta = 1.0_wp + knud * ( 1.142_wp + 0.558_wp * EXP( -0.999_wp / knud ) ) ! !-- 3) Particle properties !-- Diffusion coefficient (m2/s) (Jacobson (2005) eq. 15.29) dfpart = beta * abo * temp / ( 3.0_wp * pi * visc * diam ) !-- Mean thermal velocity (m/s) (Jacobson (2005) eq. 15.32) mtvel = SQRT( ( 8.0_wp * abo * temp ) / ( pi * mpart ) ) !-- Particle mean free path (m) (Jacobson (2005) eq. 15.34 ) omega = 8.0_wp * dfpart / ( pi * mtvel ) !-- Mean diameter (m) mdiam = 0.5_wp * ( diam(1) + diam(2) ) ! !-- 4) Calculation of fluxes (Brownian collision kernels) and flux matching !-- following Jacobson (2005): !-- Flux in continuum regime (m3/s) (eq. 15.28) flux(1) = 4.0_wp * pi * mdiam * ( dfpart(1) + dfpart(2) ) !-- Flux in free molec. regime (m3/s) (eq. 15.31) flux(2) = pi * SQRT( ( mtvel(1)**2.0_wp ) + ( mtvel(2)**2.0_wp ) ) * & ( mdiam**2.0_wp ) !-- temporary variables (m) to calculate flux matching distance (m) tva(1) = ( ( mdiam + omega(1) )**3.0_wp - ( mdiam**2.0_wp + & omega(1)**2.0_wp ) * SQRT( ( mdiam**2.0_wp + omega(1)**2.0_wp ) & ) ) / ( 3.0_wp * mdiam * omega(1) ) - mdiam tva(2) = ( ( mdiam + omega(2) )**3.0_wp - ( mdiam**2.0_wp + & omega(2)**2.0_wp ) * SQRT( ( mdiam**2 + omega(2)**2 ) ) ) / & ( 3.0_wp * mdiam * omega(2) ) - mdiam !-- Flux matching distance (m) i.e. the mean distance from the centre of a !-- sphere reached by particles leaving sphere's surface and travelling a !-- distance of particle mean free path mfp (eq. 15 34) fmdist = SQRT( tva(1)**2 + tva(2)**2.0_wp) ! !-- 5) Coagulation coefficient (m3/s) (eq. 15.33). Here assumed !-- coalescence efficiency 1!! coagc = flux(1) / ( mdiam / ( mdiam + fmdist) + flux(1) / flux(2) ) !-- coagulation coefficient = coalescence efficiency * collision kernel ! !-- Corrected collision kernel following Karl et al., 2016 (ACP): !-- Inclusion of van der Waals and viscous forces IF ( van_der_waals_coagc ) THEN knud_p = SQRT( omega(1)**2 + omega(2)**2 ) / mdiam IF ( knud_p >= 0.1_wp .AND. knud_p <= 10.0_wp ) THEN coagc = coagc * ( 2.0_wp + 0.4_wp * LOG( knud_p ) ) ELSE coagc = coagc * 3.0_wp ENDIF ENDIF END FUNCTION coagc !------------------------------------------------------------------------------! ! Description: ! ------------ !> Calculates the change in particle volume and gas phase !> concentrations due to nucleation, condensation and dissolutional growth. ! !> Sulphuric acid and organic vapour: only condensation and no evaporation. ! !> New gas and aerosol phase concentrations calculated according to Jacobson !> (1997): Numerical techniques to solve condensational and dissolutional growth !> equations when growth is coupled to reversible reactions, Aerosol Sci. Tech., !> 27, pp 491-498. ! !> Following parameterization has been used: !> Molecular diffusion coefficient of condensing vapour (m2/s) !> (Reid et al. (1987): Properties of gases and liquids, McGraw-Hill, New York.) !> D = {1.d-7*sqrt(1/M_air + 1/M_gas)*T^1.75} / & ! {p_atm/p_stand * (d_air^(1/3) + d_gas^(1/3))^2 } ! M_air = 28.965 : molar mass of air (g/mol) ! d_air = 19.70 : diffusion volume of air ! M_h2so4 = 98.08 : molar mass of h2so4 (g/mol) ! d_h2so4 = 51.96 : diffusion volume of h2so4 ! !> Called from main aerosol model ! !> fxm: calculated for empty bins too !> fxm: same diffusion coefficients and mean free paths used for sulphuric acid !> and organic vapours (average values? 'real' values for each?) !> fxm: one should really couple with vapour production and loss terms as well !> should nucleation be coupled here as well???? ! ! Coded by: ! Hannele Korhonen (FMI) 2005 ! Harri Kokkola (FMI) 2006 ! Juha Tonttila (FMI) 2014 ! Rewritten to PALM by Mona Kurppa (UHel) 2017 !------------------------------------------------------------------------------! SUBROUTINE condensation( paero, pcsa, pcocnv, pcocsv, pchno3, pcnh3, pcw, pcs,& ptemp, ppres, ptstep, prtcl ) IMPLICIT NONE !-- Input and output variables REAL(wp), INTENT(IN) :: ppres !< ambient pressure (Pa) REAL(wp), INTENT(IN) :: pcs !< Water vapour saturation concentration !< (kg/m3) REAL(wp), INTENT(IN) :: ptemp !< ambient temperature (K) REAL(wp), INTENT(IN) :: ptstep !< timestep (s) TYPE(component_index), INTENT(in) :: prtcl !< Keeps track which substances !< are used REAL(wp), INTENT(INOUT) :: pchno3 !< Gas concentrations (#/m3): !< nitric acid HNO3 REAL(wp), INTENT(INOUT) :: pcnh3 !< ammonia NH3 REAL(wp), INTENT(INOUT) :: pcocnv !< non-volatile organics REAL(wp), INTENT(INOUT) :: pcocsv !< semi-volatile organics REAL(wp), INTENT(INOUT) :: pcsa !< sulphuric acid H2SO4 REAL(wp), INTENT(INOUT) :: pcw !< Water vapor concentration (kg/m3) TYPE(t_section), INTENT(inout) :: paero(fn2b) !< Aerosol properties !-- Local variables REAL(wp) :: zbeta(fn2b) !< transitional correction factor for aerosols REAL(wp) :: zcolrate(fn2b) !< collision rate of molecules to particles !< (1/s) REAL(wp) :: zcolrate_ocnv(fn2b) !< collision rate of organic molecules !< to particles (1/s) REAL(wp) :: zcs_ocnv !< condensation sink of nonvolatile organics (1/s) REAL(wp) :: zcs_ocsv !< condensation sink of semivolatile organics (1/s) REAL(wp) :: zcs_su !< condensation sink of sulfate (1/s) REAL(wp) :: zcs_tot!< total condensation sink (1/s) (gases) !-- vapour concentration after time step (#/m3) REAL(wp) :: zcvap_new1 !< sulphuric acid REAL(wp) :: zcvap_new2 !< nonvolatile organics REAL(wp) :: zcvap_new3 !< semivolatile organics REAL(wp) :: zdfpart(in1a+1) !< particle diffusion coefficient (m2/s) REAL(wp) :: zdfvap !< air diffusion coefficient (m2/s) !-- change in vapour concentration (#/m3) REAL(wp) :: zdvap1 !< sulphuric acid REAL(wp) :: zdvap2 !< nonvolatile organics REAL(wp) :: zdvap3 !< semivolatile organics REAL(wp) :: zdvoloc(fn2b) !< change of organics volume in each bin [fxm] REAL(wp) :: zdvolsa(fn2b) !< change of sulphate volume in each bin [fxm] REAL(wp) :: zj3n3(2) !< Formation massrate of molecules in !< nucleation, (molec/m3s). 1: H2SO4 !< and 2: organic vapor REAL(wp) :: zknud(fn2b) !< particle Knudsen number REAL(wp) :: zmfp !< mean free path of condensing vapour (m) REAL(wp) :: zrh !< Relative humidity [0-1] REAL(wp) :: zvisc !< viscosity of air (kg/(m s)) REAL(wp) :: zn_vs_c !< ratio of nucleation of all mass transfer in the !< smallest bin REAL(wp) :: zxocnv !< ratio of organic vapour in 3nm particles REAL(wp) :: zxsa !< Ratio in 3nm particles: sulphuric acid zj3n3 = 0.0_wp zrh = pcw / pcs zxocnv = 0.0_wp zxsa = 0.0_wp ! !-- Nucleation IF ( nsnucl > 0 ) THEN CALL nucleation( paero, ptemp, zrh, ppres, pcsa, pcocnv, pcnh3, ptstep, & zj3n3, zxsa, zxocnv ) ENDIF ! !-- Condensation on pre-existing particles IF ( lscndgas ) THEN ! !-- Initialise: zdvolsa = 0.0_wp zdvoloc = 0.0_wp zcolrate = 0.0_wp ! !-- 1) Properties of air and condensing gases: !-- Viscosity of air (kg/(m s)) (Eq. 4.54 in Jabonson (2005)) zvisc = ( 7.44523E-3_wp * ptemp ** 1.5_wp ) / ( 5093.0_wp * & ( ptemp + 110.4_wp ) ) !-- Diffusion coefficient of air (m2/s) zdfvap = 5.1111E-10_wp * ptemp ** 1.75_wp * ( p_0 + 1325.0_wp ) / ppres !-- Mean free path (m): same for H2SO4 and organic compounds zmfp = 3.0_wp * zdfvap * SQRT( pi * amh2so4 / ( 8.0_wp * argas * ptemp ) ) ! !-- 2) Transition regime correction factor zbeta for particles: !-- Fuchs and Sutugin (1971), In: Hidy et al. (ed.) Topics in current !-- aerosol research, Pergamon. Size of condensing molecule considered !-- only for nucleation mode (3 - 20 nm) ! !-- Particle Knudsen number: condensation of gases on aerosols zknud(in1a:in1a+1) = 2.0_wp * zmfp / ( paero(in1a:in1a+1)%dwet + d_sa ) zknud(in1a+2:fn2b) = 2.0_wp * zmfp / paero(in1a+2:fn2b)%dwet ! !-- Transitional correction factor: aerosol + gas (the semi-empirical Fuchs- !-- Sutugin interpolation function (Fuchs and Sutugin, 1971)) zbeta = ( zknud + 1.0_wp ) / ( 0.377_wp * zknud + 1.0_wp + 4.0_wp / & ( 3.0_wp * massacc ) * ( zknud + zknud ** 2.0_wp ) ) ! !-- 3) Collision rate of molecules to particles !-- Particle diffusion coefficient considered only for nucleation mode !-- (3 - 20 nm) ! !-- Particle diffusion coefficient (m2/s) (e.g. Eq. 15.29 in Jacobson (2005)) zdfpart = abo * ptemp * zbeta(in1a:in1a+1) / ( 3.0_wp * pi * zvisc * & paero(in1a:in1a+1)%dwet ) ! !-- Collision rate (mass-transfer coefficient): gases on aerosols (1/s) !-- (Eq. 16.64 in Jacobson (2005)) zcolrate(in1a:in1a+1) = MERGE( 2.0_wp * pi * & ( paero(in1a:in1a+1)%dwet + d_sa ) * & ( zdfvap + zdfpart ) * zbeta(in1a:in1a+1)& * paero(in1a:in1a+1)%numc, 0.0_wp, & paero(in1a:in1a+1)%numc > nclim ) zcolrate(in1a+2:fn2b) = MERGE( 2.0_wp * pi * paero(in1a+2:fn2b)%dwet * & zdfvap * zbeta(in1a+2:fn2b) * & paero(in1a+2:fn2b)%numc, 0.0_wp, & paero(in1a+2:fn2b)%numc > nclim ) ! !-- 4) Condensation sink (1/s) zcs_tot = SUM( zcolrate ) ! total sink ! !-- 5) Changes in gas-phase concentrations and particle volume ! !-- 5.1) Organic vapours ! !-- 5.1.1) Non-volatile organic compound: condenses onto all bins IF ( pcocnv > 1.0E+10_wp .AND. zcs_tot > 1.0E-30_wp .AND. & is_used( prtcl,'OC' ) ) & THEN !-- Ratio of nucleation vs. condensation rates in the smallest bin zn_vs_c = 0.0_wp IF ( zj3n3(2) > 1.0_wp ) THEN zn_vs_c = ( zj3n3(2) ) / ( zj3n3(2) + pcocnv * zcolrate(in1a) ) ENDIF ! !-- Collision rate in the smallest bin, including nucleation and !-- condensation(see Jacobson, Fundamentals of Atmospheric Modeling, 2nd !-- Edition (2005), equation (16.73) ) zcolrate_ocnv = zcolrate zcolrate_ocnv(in1a) = zcolrate_ocnv(in1a) + zj3n3(2) / pcocnv ! !-- Total sink for organic vapor zcs_ocnv = zcs_tot + zj3n3(2) / pcocnv ! !-- New gas phase concentration (#/m3) zcvap_new2 = pcocnv / ( 1.0_wp + ptstep * zcs_ocnv ) ! !-- Change in gas concentration (#/m3) zdvap2 = pcocnv - zcvap_new2 ! !-- Updated vapour concentration (#/m3) pcocnv = zcvap_new2 ! !-- Volume change of particles (m3(OC)/m3(air)) zdvoloc = zcolrate_ocnv(in1a:fn2b) / zcs_ocnv * amvoc * zdvap2 ! !-- Change of volume due to condensation in 1a-2b paero(in1a:fn2b)%volc(2) = paero(in1a:fn2b)%volc(2) + zdvoloc ! !-- Change of number concentration in the smallest bin caused by !-- nucleation (Jacobson (2005), equation (16.75)). If zxocnv = 0, then !-- the chosen nucleation mechanism doesn't take into account the non- !-- volatile organic vapors and thus the paero doesn't have to be updated. IF ( zxocnv > 0.0_wp ) THEN paero(in1a)%numc = paero(in1a)%numc + zn_vs_c * zdvoloc(in1a) / & amvoc / ( n3 * zxocnv ) ENDIF ENDIF ! !-- 5.1.2) Semivolatile organic compound: all bins except subrange 1 zcs_ocsv = SUM( zcolrate(in2a:fn2b) ) !< sink for semi-volatile organics IF ( pcocsv > 1.0E+10_wp .AND. zcs_ocsv > 1.0E-30 .AND. & is_used( prtcl,'OC') ) & THEN ! !-- New gas phase concentration (#/m3) zcvap_new3 = pcocsv / ( 1.0_wp + ptstep * zcs_ocsv ) ! !-- Change in gas concentration (#/m3) zdvap3 = pcocsv - zcvap_new3 ! !-- Updated gas concentration (#/m3) pcocsv = zcvap_new3 ! !-- Volume change of particles (m3(OC)/m3(air)) zdvoloc(in2a:fn2b) = zdvoloc(in2a:fn2b) + zcolrate(in2a:fn2b) / & zcs_ocsv * amvoc * zdvap3 ! !-- Change of volume due to condensation in 1a-2b paero(in1a:fn2b)%volc(2) = paero(in1a:fn2b)%volc(2) + zdvoloc ENDIF ! !-- 5.2) Sulphate: condensed on all bins IF ( pcsa > 1.0E+10_wp .AND. zcs_tot > 1.0E-30_wp .AND. & is_used( prtcl,'SO4' ) ) & THEN ! !-- Ratio of mass transfer between nucleation and condensation zn_vs_c = 0.0_wp IF ( zj3n3(1) > 1.0_wp ) THEN zn_vs_c = ( zj3n3(1) ) / ( zj3n3(1) + pcsa * zcolrate(in1a) ) ENDIF ! !-- Collision rate in the smallest bin, including nucleation and !-- condensation (see Jacobson, Fundamentals of Atmospheric Modeling, 2nd !-- Edition (2005), equation (16.73)) zcolrate(in1a) = zcolrate(in1a) + zj3n3(1) / pcsa ! !-- Total sink for sulfate (1/s) zcs_su = zcs_tot + zj3n3(1) / pcsa ! !-- Sulphuric acid: !-- New gas phase concentration (#/m3) zcvap_new1 = pcsa / ( 1.0_wp + ptstep * zcs_su ) ! !-- Change in gas concentration (#/m3) zdvap1 = pcsa - zcvap_new1 ! !-- Updating vapour concentration (#/m3) pcsa = zcvap_new1 ! !-- Volume change of particles (m3(SO4)/m3(air)) by condensation zdvolsa = zcolrate(in1a:fn2b) / zcs_su * amvh2so4 * zdvap1 !-- For validation: zdvolsa = 5.5 mum3/cm3 per 12 h ! zdvolsa = zdvolsa / SUM( zdvolsa ) * 5.5E-12_wp * dt_salsa / 43200.0_wp !0.3E-12_wp, 0.6E-12_wp, 11.0E-12_wp, 4.6E-12_wp, 9.2E-12_wp ! !-- Change of volume concentration of sulphate in aerosol [fxm] paero(in1a:fn2b)%volc(1) = paero(in1a:fn2b)%volc(1) + zdvolsa ! !-- Change of number concentration in the smallest bin caused by nucleation !-- (Jacobson (2005), equation (16.75)) IF ( zxsa > 0.0_wp ) THEN paero(in1a)%numc = paero(in1a)%numc + zn_vs_c * zdvolsa(in1a) / & amvh2so4 / ( n3 * zxsa ) ENDIF ENDIF ENDIF ! ! !-- Condensation of water vapour IF ( lscndh2oae ) THEN CALL gpparth2o( paero, ptemp, ppres, pcs, pcw, ptstep ) ENDIF ! ! !-- Partitioning of H2O, HNO3, and NH3: Dissolutional growth IF ( lscndgas .AND. ino > 0 .AND. inh > 0 .AND. & ( pchno3 > 1.0E+10_wp .OR. pcnh3 > 1.0E+10_wp ) ) & THEN CALL gpparthno3( ppres, ptemp, paero, pchno3, pcnh3, pcw, pcs, zbeta, & ptstep ) ENDIF END SUBROUTINE condensation !------------------------------------------------------------------------------! ! Description: ! ------------ !> Calculates the particle number and volume increase, and gas-phase !> concentration decrease due to nucleation subsequent growth to detectable size !> of 3 nm. ! !> Method: !> When the formed clusters grow by condensation (possibly also by self- !> coagulation), their number is reduced due to scavenging to pre-existing !> particles. Thus, the apparent nucleation rate at 3 nm is significantly lower !> than the real nucleation rate (at ~1 nm). ! !> Calculation of the formation rate of detectable particles at 3 nm (i.e. J3): !> nj3 = 1: Kerminen, V.-M. and Kulmala, M. (2002), J. Aerosol Sci.,33, 609-622. !> nj3 = 2: Lehtinen et al. (2007), J. Aerosol Sci., 38(9), 988-994. !> nj3 = 3: Anttila et al. (2010), J. Aerosol Sci., 41(7), 621-636. ! !> Called from subroutine condensation (in module salsa_dynamics_mod.f90) ! !> Calls one of the following subroutines: !> - binnucl !> - ternucl !> - kinnucl !> - actnucl ! !> fxm: currently only sulphuric acid grows particles from 1 to 3 nm !> (if asked from Markku, this is terribly wrong!!!) ! !> Coded by: !> Hannele Korhonen (FMI) 2005 !> Harri Kokkola (FMI) 2006 !> Matti Niskanen(FMI) 2012 !> Anton Laakso (FMI) 2013 !------------------------------------------------------------------------------! SUBROUTINE nucleation( paero, ptemp, prh, ppres, pcsa, pcocnv, pcnh3, ptstep, & pj3n3, pxsa, pxocnv ) IMPLICIT NONE ! !-- Input and output variables REAL(wp), INTENT(in) :: pcnh3 !< ammonia concentration (#/m3) REAL(wp), INTENT(in) :: pcocnv !< conc. of non-volatile OC (#/m3) REAL(wp), INTENT(in) :: pcsa !< sulphuric acid conc. (#/m3) REAL(wp), INTENT(in) :: ppres !< ambient air pressure (Pa) REAL(wp), INTENT(in) :: prh !< ambient rel. humidity [0-1] REAL(wp), INTENT(in) :: ptemp !< ambient temperature (K) REAL(wp), INTENT(in) :: ptstep !< time step (s) of SALSA TYPE(t_section), INTENT(inout) :: paero(fn2b) !< aerosol properties REAL(wp), INTENT(inout) :: pj3n3(2) !< formation mass rate of molecules !< (molec/m3s) for 1: H2SO4 and !< 2: organic vapour REAL(wp), INTENT(out) :: pxocnv !< ratio of non-volatile organic vapours in !< 3nm aerosol particles REAL(wp), INTENT(out) :: pxsa !< ratio of H2SO4 in 3nm aerosol particles !-- Local variables INTEGER(iwp) :: iteration REAL(wp) :: zbeta(fn2b) !< transitional correction factor REAL(wp) :: zc_h2so4 !< H2SO4 conc. (#/cm3) !UNITS! REAL(wp) :: zc_org !< organic vapour conc. (#/cm3) REAL(wp) :: zCoagStot !< total losses due to coagulation, including !< condensation and self-coagulation REAL(wp) :: zcocnv_local !< organic vapour conc. (#/m3) REAL(wp) :: zcsink !< condensational sink (#/m2) REAL(wp) :: zcsa_local !< H2SO4 conc. (#/m3) REAL(wp) :: zdcrit !< diameter of critical cluster (m) REAL(wp) :: zdelta_vap !< change of H2SO4 and organic vapour !< concentration (#/m3) REAL(wp) :: zdfvap !< air diffusion coefficient (m2/s) REAL(wp) :: zdmean !< mean diameter of existing particles (m) REAL(wp) :: zeta !< constant: proportional to ratio of CS/GR (m) !< (condensation sink / growth rate) REAL(wp) :: zgamma !< proportionality factor ((nm2*m2)/h) REAL(wp) :: zGRclust !< growth rate of formed clusters (nm/h) REAL(wp) :: zGRtot !< total growth rate REAL(wp) :: zj3 !< number conc. of formed 3nm particles (#/m3) REAL(wp) :: zjnuc !< nucleation rate at ~1nm (#/m3s) REAL(wp) :: zKeff !< effective cogulation coefficient between !< freshly nucleated particles REAL(wp) :: zknud(fn2b) !< particle Knudsen number REAL(wp) :: zkocnv !< lever: zkocnv=1 --> organic compounds involved !< in nucleation REAL(wp) :: zksa !< lever: zksa=1 --> H2SO4 involved in nucleation REAL(wp) :: zlambda !< parameter for adjusting the growth rate due to !< self-coagulation REAL(wp) :: zmfp !< mean free path of condesing vapour(m) REAL(wp) :: zmixnh3 !< ammonia mixing ratio (ppt) REAL(wp) :: zNnuc !< number of clusters/particles at the size range !< d1-dx (#/m3) REAL(wp) :: znoc !< number of organic molecules in critical cluster REAL(wp) :: znsa !< number of H2SO4 molecules in critical cluster ! !-- Variable determined for the m-parameter REAL(wp) :: zCc_2(fn2b) !< REAL(wp) :: zCc_c !< REAL(wp) :: zCc_x !< REAL(wp) :: zCoagS_c !< REAL(wp) :: zCoagS_x !< REAL(wp) :: zcv_2(fn2b) !< REAL(wp) :: zcv_c !< REAL(wp) :: zcv_c2(fn2b) !< REAL(wp) :: zcv_x !< REAL(wp) :: zcv_x2(fn2b) !< REAL(wp) :: zDc_2(fn2b) !< REAL(wp) :: zDc_c(fn2b) !< REAL(wp) :: zDc_c2(fn2b) !< REAL(wp) :: zDc_x(fn2b) !< REAL(wp) :: zDc_x2(fn2b) !< REAL(wp) :: zgammaF_2(fn2b) !< REAL(wp) :: zgammaF_c(fn2b) !< REAL(wp) :: zgammaF_x(fn2b) !< REAL(wp) :: zK_c2(fn2b) !< REAL(wp) :: zK_x2(fn2b) !< REAL(wp) :: zknud_2(fn2b) !< REAL(wp) :: zknud_c !< REAL(wp) :: zknud_x !< REAL(wp) :: zm_2(fn2b) !< REAL(wp) :: zm_c !< REAL(wp) :: zm_para !< REAL(wp) :: zm_x !< REAL(wp) :: zmyy !< REAL(wp) :: zomega_2c(fn2b) !< REAL(wp) :: zomega_2x(fn2b) !< REAL(wp) :: zomega_c(fn2b) !< REAL(wp) :: zomega_x(fn2b) !< REAL(wp) :: zRc2(fn2b) !< REAL(wp) :: zRx2(fn2b) !< REAL(wp) :: zsigma_c2(fn2b) !< REAL(wp) :: zsigma_x2(fn2b) !< ! !-- 1) Nucleation rate (zjnuc) and diameter of critical cluster (zdcrit) zjnuc = 0.0_wp znsa = 0.0_wp znoc = 0.0_wp zdcrit = 0.0_wp zksa = 0.0_wp zkocnv = 0.0_wp SELECT CASE ( nsnucl ) CASE(1) ! Binary H2SO4-H2O nucleation zc_h2so4 = pcsa * 1.0E-6_wp ! sulphuric acid conc. to #/cm3 CALL binnucl( zc_h2so4, ptemp, prh, zjnuc, znsa, znoc, zdcrit, zksa, & zkocnv ) CASE(2) ! Activation type nucleation zc_h2so4 = pcsa * 1.0E-6_wp ! sulphuric acid conc. to #/cm3 CALL binnucl( zc_h2so4, ptemp, prh, zjnuc, znsa, znoc, zdcrit, zksa, & zkocnv ) CALL actnucl( pcsa, zjnuc, zdcrit, znsa, znoc, zksa, zkocnv, act_coeff ) CASE(3) ! Kinetically limited nucleation of (NH4)HSO4 clusters zc_h2so4 = pcsa * 1.0E-6_wp ! sulphuric acid conc. to #/cm3 CALL binnucl( zc_h2so4, ptemp, prh, zjnuc, znsa, znoc, zdcrit, zksa, & zkocnv ) CALL kinnucl( zc_h2so4, zjnuc, zdcrit, znsa, znoc, zksa, zkocnv ) CASE(4) ! Ternary H2SO4-H2O-NH3 nucleation zmixnh3 = pcnh3 * ptemp * argas / ( ppres * avo ) zc_h2so4 = pcsa * 1.0E-6_wp ! sulphuric acid conc. to #/cm3 CALL ternucl( zc_h2so4, zmixnh3, ptemp, prh, zjnuc, znsa, znoc, zdcrit, & zksa, zkocnv ) CASE(5) ! Organic nucleation, J~[ORG] or J~[ORG]**2 zc_org = pcocnv * 1.0E-6_wp ! conc. of non-volatile OC to #/cm3 zc_h2so4 = pcsa * 1.0E-6_wp ! sulphuric acid conc. to #/cm3 CALL binnucl( zc_h2so4, ptemp, prh, zjnuc, znsa, znoc, zdcrit, zksa, & zkocnv ) CALL orgnucl( pcocnv, zjnuc, zdcrit, znsa, znoc, zksa, zkocnv ) CASE(6) ! Sum of H2SO4 and organic activation type nucleation, ! J~[H2SO4]+[ORG] zc_h2so4 = pcsa * 1.0E-6_wp ! sulphuric acid conc. to #/cm3 CALL binnucl( zc_h2so4, ptemp, prh, zjnuc, znsa, znoc, zdcrit, zksa, & zkocnv ) CALL sumnucl( pcsa, pcocnv, zjnuc, zdcrit, znsa, znoc, zksa, zkocnv ) CASE(7) ! Heteromolecular nucleation, J~[H2SO4]*[ORG] zc_h2so4 = pcsa * 1.0E-6_wp ! sulphuric acid conc. to #/cm3 zc_org = pcocnv * 1.0E-6_wp ! conc. of non-volatile OC to #/cm3 CALL binnucl( zc_h2so4, ptemp, prh, zjnuc, znsa, znoc, zdcrit, zksa, & zkocnv ) CALL hetnucl( zc_h2so4, zc_org, zjnuc, zdcrit, znsa, znoc, zksa, zkocnv ) CASE(8) ! Homomolecular nucleation of H2SO4 and heteromolecular ! nucleation of H2SO4 and organic vapour, ! J~[H2SO4]**2 + [H2SO4]*[ORG] (EUCAARI project) zc_h2so4 = pcsa * 1.0E-6_wp ! sulphuric acid conc. to #/cm3 zc_org = pcocnv * 1.0E-6_wp ! conc. of non-volatile OC to #/cm3 CALL binnucl( zc_h2so4, ptemp, prh, zjnuc, znsa, znoc, zdcrit, zksa, & zkocnv ) CALL SAnucl( zc_h2so4, zc_org, zjnuc, zdcrit, znsa, znoc, zksa, zkocnv ) CASE(9) ! Homomolecular nucleation of H2SO4 and organic vapour and ! heteromolecular nucleation of H2SO4 and organic vapour, ! J~[H2SO4]**2 + [H2SO4]*[ORG]+[ORG]**2 (EUCAARI project) zc_h2so4 = pcsa * 1.0E-6_wp ! sulphuric acid conc. to #/cm3 zc_org = pcocnv * 1.0E-6_wp ! conc. of non-volatile OC to #/cm3 CALL binnucl( zc_h2so4, ptemp, prh, zjnuc, znsa, znoc, zdcrit, zksa, & zkocnv ) CALL SAORGnucl( zc_h2so4, zc_org, zjnuc, zdcrit, znsa, znoc, zksa, & zkocnv ) END SELECT zcsa_local = pcsa zcocnv_local = pcocnv ! !-- 2) Change of particle and gas concentrations due to nucleation ! !-- 2.1) Check that there is enough H2SO4 and organic vapour to produce the !-- nucleation IF ( nsnucl <= 4 ) THEN !-- If the chosen nucleation scheme is 1-4, nucleation occurs only due to !-- H2SO4. All of the total vapour concentration that is taking part to the !-- nucleation is there for sulphuric acid (sa = H2SO4) and non-volatile !-- organic vapour is zero. pxsa = 1.0_wp ! ratio of sulphuric acid in 3nm particles pxocnv = 0.0_wp ! ratio of non-volatile origanic vapour ! in 3nm particles ELSEIF ( nsnucl > 4 ) THEN !-- If the chosen nucleation scheme is 5-9, nucleation occurs due to organic !-- vapour or the combination of organic vapour and H2SO4. The number of !-- needed molecules depends on the chosen nucleation type and it has an !-- effect also on the minimum ratio of the molecules present. IF ( pcsa * znsa + pcocnv * znoc < 1.E-14_wp ) THEN pxsa = 0.0_wp pxocnv = 0.0_wp ELSE pxsa = pcsa * znsa / ( pcsa * znsa + pcocnv * znoc ) pxocnv = pcocnv * znoc / ( pcsa * znsa + pcocnv * znoc ) ENDIF ENDIF ! !-- The change in total vapour concentration is the sum of the concentrations !-- of the vapours taking part to the nucleation (depends on the chosen !-- nucleation scheme) zdelta_vap = MIN( zjnuc * ( znoc + znsa ), ( pcocnv * zkocnv + pcsa * & zksa ) / ptstep ) ! !-- Nucleation rate J at ~1nm (#/m3s) zjnuc = zdelta_vap / ( znoc + znsa ) ! !-- H2SO4 concentration after nucleation in #/m3 zcsa_local = MAX( 1.0_wp, pcsa - zdelta_vap * pxsa ) ! !-- Non-volative organic vapour concentration after nucleation (#/m3) zcocnv_local = MAX( 1.0_wp, pcocnv - zdelta_vap * pxocnv ) ! !-- 2.2) Formation rate of 3 nm particles (Kerminen & Kulmala, 2002) ! !-- 2.2.1) Growth rate of clusters formed by H2SO4 ! !-- GR = 3.0e-15 / dens_clus * sum( molecspeed * molarmass * conc ) ! !-- dens_clus = density of the clusters (here 1830 kg/m3) !-- molarmass = molar mass of condensing species (here 98.08 g/mol) !-- conc = concentration of condensing species [#/m3] !-- molecspeed = molecular speed of condensing species [m/s] !-- = sqrt( 8.0 * R * ptemp / ( pi * molarmass ) ) !-- (Seinfeld & Pandis, 1998) ! !-- Growth rate by H2SO4 and organic vapour in nm/h (Eq. 21) zGRclust = 2.3623E-15_wp * SQRT( ptemp ) * ( zcsa_local + zcocnv_local ) ! !-- 2.2.2) Condensational sink of pre-existing particle population ! !-- Diffusion coefficient (m2/s) zdfvap = 5.1111E-10_wp * ptemp ** 1.75_wp * ( p_0 + 1325.0_wp ) / ppres !-- Mean free path of condensing vapour (m) (Jacobson (2005), Eq. 15.25 and !-- 16.29) zmfp = 3.0_wp * zdfvap * SQRT( pi * amh2so4 / ( 8.0_wp * argas * ptemp ) ) !-- Knudsen number zknud = 2.0_wp * zmfp / ( paero(:)%dwet + d_sa ) !-- Transitional regime correction factor (zbeta) according to Fuchs and !-- Sutugin (1971), In: Hidy et al. (ed.), Topics in current aerosol research, !-- Pergamon. (Eq. 4 in Kerminen and Kulmala, 2002) zbeta = ( zknud + 1.0_wp) / ( 0.377_wp * zknud + 1.0_wp + 4.0_wp / & ( 3.0_wp * massacc ) * ( zknud + zknud ** 2 ) ) !-- Condensational sink (#/m2) (Eq. 3) zcsink = SUM( paero(:)%dwet * zbeta * paero(:)%numc ) ! !-- Parameterised formation rate of detectable 3 nm particles (i.e. J3) IF ( nj3 == 1 ) THEN ! Kerminen and Kulmala (2002) !-- 2.2.3) Parameterised formation rate of detectable 3 nm particles !-- Constants needed for the parameterisation: !-- dapp = 3 nm and dens_nuc = 1830 kg/m3 IF ( zcsink < 1.0E-30_wp ) THEN zeta = 0._dp ELSE !-- Mean diameter of backgroud population (nm) zdmean = 1.0_wp / SUM( paero(:)%numc ) * SUM( paero(:)%numc * & paero(:)%dwet ) * 1.0E+9_wp !-- Proportionality factor (nm2*m2/h) (Eq. 22) zgamma = 0.23_wp * ( zdcrit * 1.0E+9_wp ) ** 0.2_wp * ( zdmean / & 150.0_wp ) ** 0.048_wp * ( ptemp / 293.0_wp ) ** ( -0.75_wp ) & * ( arhoh2so4 / 1000.0_wp ) ** ( -0.33_wp ) !-- Factor eta (nm) (Eq. 11) zeta = MIN( zgamma * zcsink / zGRclust, zdcrit * 1.0E11_wp ) ENDIF ! !-- Number conc. of clusters surviving to 3 nm in a time step (#/m3) (Eq.14) zj3 = zjnuc * EXP( MIN( 0.0_wp, zeta / 3.0_wp - zeta / & ( zdcrit * 1.0E9_wp ) ) ) ELSEIF ( nj3 > 1 ) THEN !-- Defining the value for zm_para. The growth is investigated between !-- [d1,reglim(1)] = [zdcrit,3nm] !-- m = LOG( CoagS_dx / CoagX_zdcrit ) / LOG( reglim / zdcrit ) !-- (Lehtinen et al. 2007, Eq. 5) !-- The steps for the coagulation sink for reglim = 3nm and zdcrit ~= 1nm are !-- explained in article of Kulmala et al. (2001). The particles of diameter !-- zdcrit ~1.14 nm and reglim = 3nm are both in turn the "number 1" !-- variables (Kulmala et al. 2001). !-- c = critical (1nm), x = 3nm, 2 = wet or mean droplet !-- Sum of the radii, R12 = R1 + zR2 (m) of two particles 1 and 2 zRc2 = zdcrit / 2.0_wp + paero(:)%dwet / 2.0_wp zRx2 = reglim(1) / 2.0_wp + paero(:)%dwet / 2.0_wp ! !-- The mass of particle (kg) (comes only from H2SO4) zm_c = 4.0_wp / 3.0_wp * pi * ( zdcrit / 2.0_wp ) ** 3.0_wp * arhoh2so4 zm_x = 4.0_wp / 3.0_wp * pi * ( reglim(1) / 2.0_wp ) ** 3.0_wp * & arhoh2so4 zm_2 = 4.0_wp / 3.0_wp * pi * ( paero(:)%dwet / 2.0_wp )** 3.0_wp * & arhoh2so4 ! !-- Mean relative thermal velocity between the particles (m/s) zcv_c = SQRT( 8.0_wp * abo * ptemp / ( pi * zm_c ) ) zcv_x = SQRT( 8.0_wp * abo * ptemp / ( pi * zm_x ) ) zcv_2 = SQRT( 8.0_wp * abo * ptemp / ( pi * zm_2 ) ) ! !-- Average velocity after coagulation zcv_c2 = SQRT( zcv_c ** 2.0_wp + zcv_2 ** 2.0_wp ) zcv_x2 = SQRT( zcv_x ** 2.0_wp + zcv_2 ** 2.0_wp ) ! !-- Knudsen number (zmfp = mean free path of condensing vapour) zknud_c = 2.0_wp * zmfp / zdcrit zknud_x = 2.0_wp * zmfp / reglim(1) zknud_2 = MAX( 0.0_wp, 2.0_wp * zmfp / paero(:)%dwet ) ! !-- Cunningham correction factor zCc_c = 1.0_wp + zknud_c * ( 1.142_wp + 0.558_wp * & EXP( -0.999_wp / zknud_c ) ) zCc_x = 1.0_wp + zknud_x * ( 1.142_wp + 0.558_wp * & EXP( -0.999_wp / zknud_x ) ) zCc_2 = 1.0_wp + zknud_2 * ( 1.142_wp + 0.558_wp * & EXP( -0.999_wp / zknud_2 ) ) ! !-- Gas dynamic viscosity (N*s/m2). !-- Viscocity(air @20C) = 1.81e-5_dp N/m2 *s (Hinds, p. 25) zmyy = 1.81E-5_wp * ( ptemp / 293.0_wp) ** ( 0.74_wp ) ! !-- Particle diffusion coefficient (m2/s) zDc_c = abo * ptemp * zCc_c / ( 3.0_wp * pi * zmyy * zdcrit ) zDc_x = abo * ptemp * zCc_x / ( 3.0_wp * pi * zmyy * reglim(1) ) zDc_2 = abo * ptemp * zCc_2 / ( 3.0_wp * pi * zmyy * paero(:)%dwet ) ! !-- D12 = D1+D2 (Seinfield and Pandis, 2nd ed. Eq. 13.38) zDc_c2 = zDc_c + zDc_2 zDc_x2 = zDc_x + zDc_2 ! !-- zgammaF = 8*D/pi/zcv (m) for calculating zomega zgammaF_c = 8.0_wp * zDc_c / pi / zcv_c zgammaF_x = 8.0_wp * zDc_x / pi / zcv_x zgammaF_2 = 8.0_wp * zDc_2 / pi / zcv_2 ! !-- zomega (m) for calculating zsigma zomega_c = ( ( zRc2 + zgammaF_c ) ** 3 - ( zRc2 ** 2 + & zgammaF_c ) ** ( 3.0_wp / 2.0_wp ) ) / ( 3.0_wp * & zRc2 * zgammaF_c ) - zRc2 zomega_x = ( ( zRx2 + zgammaF_x ) ** 3.0_wp - ( zRx2 ** 2.0_wp + & zgammaF_x ) ** ( 3.0_wp / 2.0_wp ) ) / ( 3.0_wp * & zRx2 * zgammaF_x ) - zRx2 zomega_2c = ( ( zRc2 + zgammaF_2 ) ** 3.0_wp - ( zRc2 ** 2.0_wp + & zgammaF_2 ) ** ( 3.0_wp / 2.0_wp ) ) / ( 3.0_wp * & zRc2 * zgammaF_2 ) - zRc2 zomega_2x = ( ( zRx2 + zgammaF_2 ) ** 3.0_wp - ( zRx2 ** 2.0_wp + & zgammaF_2 ) ** ( 3.0_wp / 2.0_wp ) ) / ( 3.0_wp * & zRx2 * zgammaF_2 ) - zRx2 ! !-- The distance (m) at which the two fluxes are matched (condensation and !-- coagulation sinks?) zsigma_c2 = SQRT( zomega_c ** 2.0_wp + zomega_2c ** 2.0_wp ) zsigma_x2 = SQRT( zomega_x ** 2.0_wp + zomega_2x ** 2.0_wp ) ! !-- Coagulation coefficient in the continuum regime (m*m2/s) zK_c2 = 4.0_wp * pi * zRc2 * zDc_c2 / ( zRc2 / ( zRc2 + zsigma_c2 ) + & 4.0_wp * zDc_c2 / ( zcv_c2 * zRc2 ) ) zK_x2 = 4.0_wp * pi * zRx2 * zDc_x2 / ( zRx2 / ( zRx2 + zsigma_x2 ) + & 4.0_wp * zDc_x2 / ( zcv_x2 * zRx2 ) ) ! !-- Coagulation sink (1/s) zCoagS_c = MAX( 1.0E-20_wp, SUM( zK_c2 * paero(:)%numc ) ) zCoagS_x = MAX( 1.0E-20_wp, SUM( zK_x2 * paero(:)%numc ) ) ! !-- Parameter m for calculating the coagulation sink onto background !-- particles (Eq. 5&6 in Lehtinen et al. 2007) zm_para = LOG( zCoagS_x / zCoagS_c ) / LOG( reglim(1) / zdcrit ) ! !-- Parameter gamma for calculating the formation rate J of particles having !-- a diameter zdcrit < d < reglim(1) (Anttila et al. 2010, eq. 5) zgamma = ( ( ( reglim(1) / zdcrit ) ** ( zm_para + 1.0_wp ) ) - 1.0_wp )& / ( zm_para + 1.0_wp ) IF ( nj3 == 2 ) THEN ! Coagulation sink ! !-- Formation rate J before iteration (#/m3s) zj3 = zjnuc * EXP( MIN( 0.0_wp, -zgamma * zdcrit * zCoagS_c / & ( zGRclust * 1.0E-9_wp / ( 60.0_wp ** 2.0_wp ) ) ) ) ELSEIF ( nj3 == 3 ) THEN ! Coagulation sink and self-coag. !-- IF polluted air... then the self-coagulation becomes important. !-- Self-coagulation of small particles < 3 nm. ! !-- "Effective" coagulation coefficient between freshly-nucleated !-- particles: zKeff = 5.0E-16_wp ! cm3/s ! !-- zlambda parameter for "adjusting" the growth rate due to the !-- self-coagulation zlambda = 6.0_wp IF ( reglim(1) >= 10.0E-9_wp ) THEN ! for particles >10 nm: zKeff = 5.0E-17_wp zlambda = 3.0_wp ENDIF ! !-- Initial values for coagulation sink and growth rate (m/s) zCoagStot = zCoagS_c zGRtot = zGRclust * 1.0E-9_wp / 60.0_wp ** 2.0_wp ! !-- Number of clusters/particles at the size range [d1,dx] (#/m3): zNnuc = zjnuc / zCoagStot !< Initial guess ! !-- Coagulation sink and growth rate due to self-coagulation: DO iteration = 1, 5 zCoagStot = zCoagS_c + zKeff * zNnuc * 1.0E-6_wp ! (1/s) zGRtot = zGRclust * 1.0E-9_wp / ( 3600.0_wp ) + 1.5708E-6_wp * & zlambda * zdcrit ** 3.0_wp * ( zNnuc * 1.0E-6_wp ) * & zcv_c * avo * 1.0E-9_wp / 3600.0_wp zeta = - zCoagStot / ( ( zm_para + 1.0_wp ) * zGRtot * ( zdcrit **& zm_para ) ) ! Eq. 7b (Anttila) zNnuc = zNnuc_tayl( zdcrit, reglim(1), zm_para, zjnuc, zeta, & zGRtot ) ENDDO ! !-- Calculate the final values with new zNnuc: zCoagStot = zCoagS_c + zKeff * zNnuc * 1.0E-6_wp ! (1/s) zGRtot = zGRclust * 1.0E-9_wp / 3600.0_wp + 1.5708E-6_wp * zlambda & * zdcrit ** 3.0_wp * ( zNnuc * 1.0E-6_wp ) * zcv_c * avo * & 1.0E-9_wp / 3600.0_wp !< (m/s) zj3 = zjnuc * EXP( MIN( 0.0_wp, -zgamma * zdcrit * zCoagStot / & zGRtot ) ) ! (Eq. 5a) (#/m3s) ENDIF ENDIF !-- If J3 very small (< 1 #/cm3), neglect particle formation. In real atmosphere !-- this would mean that clusters form but coagulate to pre-existing particles !-- who gain sulphate. Since CoagS ~ CS (4piD*CS'), we do *not* update H2SO4 !-- concentration here but let condensation take care of it. !-- Formation mass rate of molecules (molec/m3s) for 1: H2SO4 and 2: organic !-- vapour pj3n3(1) = zj3 * n3 * pxsa pj3n3(2) = zj3 * n3 * pxocnv END SUBROUTINE nucleation !------------------------------------------------------------------------------! ! Description: ! ------------ !> Calculate the nucleation rate and the size of critical clusters assuming !> binary nucleation. !> Parametrisation according to Vehkamaki et al. (2002), J. Geophys. Res., !> 107(D22), 4622. Called from subroutine nucleation. !------------------------------------------------------------------------------! SUBROUTINE binnucl( pc_sa, ptemp, prh, pnuc_rate, pn_crit_sa, pn_crit_ocnv, & pd_crit, pk_sa, pk_ocnv ) IMPLICIT NONE ! !-- Input and output variables REAL(wp), INTENT(in) :: pc_sa !< H2SO4 conc. (#/cm3) REAL(wp), INTENT(in) :: prh !< relative humidity [0-1] REAL(wp), INTENT(in) :: ptemp !< ambient temperature (K) REAL(wp), INTENT(out) :: pnuc_rate !< nucleation rate (#/(m3 s)) REAL(wp), INTENT(out) :: pn_crit_sa !< number of H2SO4 molecules in !< cluster (#) REAL(wp), INTENT(out) :: pn_crit_ocnv !< number of organic molecules in !< cluster (#) REAL(wp), INTENT(out) :: pd_crit !< diameter of critical cluster (m) REAL(wp), INTENT(out) :: pk_sa !< Lever: if pk_sa = 1, H2SO4 is !< involved in nucleation. REAL(wp), INTENT(out) :: pk_ocnv !< Lever: if pk_ocnv = 1, organic !< compounds are involved in !< nucleation. !-- Local variables REAL(wp) :: zx !< mole fraction of sulphate in critical cluster REAL(wp) :: zntot !< number of molecules in critical cluster REAL(wp) :: zt !< temperature REAL(wp) :: zpcsa !< sulfuric acid concentration REAL(wp) :: zrh !< relative humidity REAL(wp) :: zma !< REAL(wp) :: zmw !< REAL(wp) :: zxmass!< REAL(wp) :: za !< REAL(wp) :: zb !< REAL(wp) :: zc !< REAL(wp) :: zroo !< REAL(wp) :: zm1 !< REAL(wp) :: zm2 !< REAL(wp) :: zv1 !< REAL(wp) :: zv2 !< REAL(wp) :: zcoll !< pnuc_rate = 0.0_wp pd_crit = 1.0E-9_wp ! !-- 1) Checking that we are in the validity range of the parameterization zt = MAX( ptemp, 190.15_wp ) zt = MIN( zt, 300.15_wp ) zpcsa = MAX( pc_sa, 1.0E4_wp ) zpcsa = MIN( zpcsa, 1.0E11_wp ) zrh = MAX( prh, 0.0001_wp ) zrh = MIN( zrh, 1.0_wp ) ! !-- 2) Mole fraction of sulphate in a critical cluster (Eq. 11) zx = 0.7409967177282139_wp & - 0.002663785665140117_wp * zt & + 0.002010478847383187_wp * LOG( zrh ) & - 0.0001832894131464668_wp* zt * LOG( zrh ) & + 0.001574072538464286_wp * LOG( zrh ) ** 2 & - 0.00001790589121766952_wp * zt * LOG( zrh ) ** 2 & + 0.0001844027436573778_wp * LOG( zrh ) ** 3 & - 1.503452308794887E-6_wp * zt * LOG( zrh ) ** 3 & - 0.003499978417957668_wp * LOG( zpcsa ) & + 0.0000504021689382576_wp * zt * LOG( zpcsa ) ! !-- 3) Nucleation rate (Eq. 12) pnuc_rate = 0.1430901615568665_wp & + 2.219563673425199_wp * zt & - 0.02739106114964264_wp * zt ** 2 & + 0.00007228107239317088_wp * zt ** 3 & + 5.91822263375044_wp / zx & + 0.1174886643003278_wp * LOG( zrh ) & + 0.4625315047693772_wp * zt * LOG( zrh ) & - 0.01180591129059253_wp * zt ** 2 * LOG( zrh ) & + 0.0000404196487152575_wp * zt ** 3 * LOG( zrh ) & + ( 15.79628615047088_wp * LOG( zrh ) ) / zx & - 0.215553951893509_wp * LOG( zrh ) ** 2 & - 0.0810269192332194_wp * zt * LOG( zrh ) ** 2 & + 0.001435808434184642_wp * zt ** 2 * LOG( zrh ) ** 2 & - 4.775796947178588E-6_wp * zt ** 3 * LOG( zrh ) ** 2 & - (2.912974063702185_wp * LOG( zrh ) ** 2 ) / zx & - 3.588557942822751_wp * LOG( zrh ) ** 3 & + 0.04950795302831703_wp * zt * LOG( zrh ) ** 3 & - 0.0002138195118737068_wp * zt ** 2 * LOG( zrh ) ** 3 & + 3.108005107949533E-7_wp * zt ** 3 * LOG( zrh ) ** 3 & - ( 0.02933332747098296_wp * LOG( zrh ) ** 3 ) / zx & + 1.145983818561277_wp * LOG( zpcsa ) & - 0.6007956227856778_wp * zt * LOG( zpcsa ) & + 0.00864244733283759_wp * zt ** 2 * LOG( zpcsa ) & - 0.00002289467254710888_wp * zt ** 3 * LOG( zpcsa ) & - ( 8.44984513869014_wp * LOG( zpcsa ) ) / zx & + 2.158548369286559_wp * LOG( zrh ) * LOG( zpcsa ) & + 0.0808121412840917_wp * zt * LOG( zrh ) * LOG( zpcsa ) & - 0.0004073815255395214_wp * zt ** 2 * LOG( zrh ) * LOG( zpcsa ) & - 4.019572560156515E-7_wp * zt ** 3 * LOG( zrh ) * LOG( zpcsa ) & + ( 0.7213255852557236_wp * LOG( zrh ) * LOG( zpcsa ) ) / zx & + 1.62409850488771_wp * LOG( zrh ) ** 2 * LOG( zpcsa ) & - 0.01601062035325362_wp * zt * LOG( zrh ) ** 2 * LOG( zpcsa ) & + 0.00003771238979714162_wp*zt**2* LOG( zrh )**2 * LOG( zpcsa ) & + 3.217942606371182E-8_wp * zt**3 * LOG( zrh )**2 * LOG( zpcsa ) & - (0.01132550810022116_wp * LOG( zrh )**2 * LOG( zpcsa ) ) / zx & + 9.71681713056504_wp * LOG( zpcsa ) ** 2 & - 0.1150478558347306_wp * zt * LOG( zpcsa ) ** 2 & + 0.0001570982486038294_wp * zt ** 2 * LOG( zpcsa ) ** 2 & + 4.009144680125015E-7_wp * zt ** 3 * LOG( zpcsa ) ** 2 & + ( 0.7118597859976135_wp * LOG( zpcsa ) ** 2 ) / zx & - 1.056105824379897_wp * LOG( zrh ) * LOG( zpcsa ) ** 2 & + 0.00903377584628419_wp * zt * LOG( zrh ) * LOG( zpcsa )**2 & - 0.00001984167387090606_wp*zt**2*LOG( zrh )*LOG( zpcsa )**2 & + 2.460478196482179E-8_wp * zt**3 * LOG( zrh ) * LOG( zpcsa )**2 & - ( 0.05790872906645181_wp * LOG( zrh ) * LOG( zpcsa )**2 ) / zx & - 0.1487119673397459_wp * LOG( zpcsa ) ** 3 & + 0.002835082097822667_wp * zt * LOG( zpcsa ) ** 3 & - 9.24618825471694E-6_wp * zt ** 2 * LOG( zpcsa ) ** 3 & + 5.004267665960894E-9_wp * zt ** 3 * LOG( zpcsa ) ** 3 & - ( 0.01270805101481648_wp * LOG( zpcsa ) ** 3 ) / zx ! !-- Nucleation rate in #/(cm3 s) pnuc_rate = EXP( pnuc_rate ) ! !-- Check the validity of parameterization IF ( pnuc_rate < 1.0E-7_wp ) THEN pnuc_rate = 0.0_wp pd_crit = 1.0E-9_wp ENDIF ! !-- 4) Total number of molecules in the critical cluster (Eq. 13) zntot = - 0.002954125078716302_wp & - 0.0976834264241286_wp * zt & + 0.001024847927067835_wp * zt ** 2 & - 2.186459697726116E-6_wp * zt ** 3 & - 0.1017165718716887_wp / zx & - 0.002050640345231486_wp * LOG( zrh ) & - 0.007585041382707174_wp * zt * LOG( zrh ) & + 0.0001926539658089536_wp * zt ** 2 * LOG( zrh ) & - 6.70429719683894E-7_wp * zt ** 3 * LOG( zrh ) & - ( 0.2557744774673163_wp * LOG( zrh ) ) / zx & + 0.003223076552477191_wp * LOG( zrh ) ** 2 & + 0.000852636632240633_wp * zt * LOG( zrh ) ** 2 & - 0.00001547571354871789_wp * zt ** 2 * LOG( zrh ) ** 2 & + 5.666608424980593E-8_wp * zt ** 3 * LOG( zrh ) ** 2 & + ( 0.03384437400744206_wp * LOG( zrh ) ** 2 ) / zx & + 0.04743226764572505_wp * LOG( zrh ) ** 3 & - 0.0006251042204583412_wp * zt * LOG( zrh ) ** 3 & + 2.650663328519478E-6_wp * zt ** 2 * LOG( zrh ) ** 3 & - 3.674710848763778E-9_wp * zt ** 3 * LOG( zrh ) ** 3 & - ( 0.0002672510825259393_wp * LOG( zrh ) ** 3 ) / zx & - 0.01252108546759328_wp * LOG( zpcsa ) & + 0.005806550506277202_wp * zt * LOG( zpcsa ) & - 0.0001016735312443444_wp * zt ** 2 * LOG( zpcsa ) & + 2.881946187214505E-7_wp * zt ** 3 * LOG( zpcsa ) & + ( 0.0942243379396279_wp * LOG( zpcsa ) ) / zx & - 0.0385459592773097_wp * LOG( zrh ) * LOG( zpcsa ) & - 0.0006723156277391984_wp * zt * LOG( zrh ) * LOG( zpcsa ) & + 2.602884877659698E-6_wp * zt ** 2 * LOG( zrh ) * LOG( zpcsa ) & + 1.194163699688297E-8_wp * zt ** 3 * LOG( zrh ) * LOG( zpcsa ) & - ( 0.00851515345806281_wp * LOG( zrh ) * LOG( zpcsa ) ) / zx & - 0.01837488495738111_wp * LOG( zrh ) ** 2 * LOG( zpcsa ) & + 0.0001720723574407498_wp * zt * LOG( zrh ) ** 2 * LOG( zpcsa ) & - 3.717657974086814E-7_wp * zt**2 * LOG( zrh )**2 * LOG( zpcsa ) & - 5.148746022615196E-10_wp * zt**3 * LOG( zrh )**2 * LOG( zpcsa ) & + ( 0.0002686602132926594_wp * LOG(zrh)**2 * LOG(zpcsa) ) / zx & - 0.06199739728812199_wp * LOG( zpcsa ) ** 2 & + 0.000906958053583576_wp * zt * LOG( zpcsa ) ** 2 & - 9.11727926129757E-7_wp * zt ** 2 * LOG( zpcsa ) ** 2 & - 5.367963396508457E-9_wp * zt ** 3 * LOG( zpcsa ) ** 2 & - ( 0.007742343393937707_wp * LOG( zpcsa ) ** 2 ) / zx & + 0.0121827103101659_wp * LOG( zrh ) * LOG( zpcsa ) ** 2 & - 0.0001066499571188091_wp * zt * LOG( zrh ) * LOG( zpcsa ) ** 2 & + 2.534598655067518E-7_wp * zt**2 * LOG( zrh ) * LOG( zpcsa )**2 & - 3.635186504599571E-10_wp * zt**3 * LOG( zrh ) * LOG( zpcsa )**2 & + ( 0.0006100650851863252_wp * LOG( zrh ) * LOG( zpcsa ) **2 )/ zx & + 0.0003201836700403512_wp * LOG( zpcsa ) ** 3 & - 0.0000174761713262546_wp * zt * LOG( zpcsa ) ** 3 & + 6.065037668052182E-8_wp * zt ** 2 * LOG( zpcsa ) ** 3 & - 1.421771723004557E-11_wp * zt ** 3 * LOG( zpcsa ) ** 3 & + ( 0.0001357509859501723_wp * LOG( zpcsa ) ** 3 ) / zx zntot = EXP( zntot ) ! in # ! !-- 5) Size of the critical cluster pd_crit (m) (diameter) (Eq. 14) pn_crit_sa = zx * zntot pd_crit = 2.0E-9_wp * EXP( -1.6524245_wp + 0.42316402_wp * zx + & 0.33466487_wp * LOG( zntot ) ) ! !-- 6) Organic compounds not involved when binary nucleation is assumed pn_crit_ocnv = 0.0_wp ! number of organic molecules pk_sa = 1.0_wp ! if = 1, H2SO4 involved in nucleation pk_ocnv = 0.0_wp ! if = 1, organic compounds involved ! !-- Set nucleation rate to collision rate IF ( pn_crit_sa < 4.0_wp ) THEN ! !-- Volumes of the colliding objects zma = 96.0_wp ! molar mass of SO4 in g/mol zmw = 18.0_wp ! molar mass of water in g/mol zxmass = 1.0_wp ! mass fraction of H2SO4 za = 0.7681724_wp + zxmass * ( 2.1847140_wp + zxmass * ( & 7.1630022_wp + zxmass * ( -44.31447_wp + zxmass * ( & 88.75606 + zxmass * ( -75.73729_wp + zxmass * & 23.43228_wp ) ) ) ) ) zb = 1.808225E-3_wp + zxmass * ( -9.294656E-3_wp + zxmass * & ( -0.03742148_wp + zxmass * ( 0.2565321_wp + zxmass * & ( -0.5362872_wp + zxmass * ( 0.4857736 - zxmass * & 0.1629592_wp ) ) ) ) ) zc = - 3.478524E-6_wp + zxmass * ( 1.335867E-5_wp + zxmass * & ( 5.195706E-5_wp + zxmass * ( -3.717636E-4_wp + zxmass * & ( 7.990811E-4_wp + zxmass * ( -7.458060E-4_wp + zxmass * & 2.58139E-4_wp ) ) ) ) ) ! !-- Density for the sulphuric acid solution (Eq. 10 in Vehkamaki) zroo = za + zt * ( zb + zc * zt ) ! g/cm^3 zroo = zroo * 1.0E+3_wp ! kg/m^3 zm1 = 0.098_wp ! molar mass of H2SO4 in kg/mol zm2 = zm1 zv1 = zm1 / avo / zroo ! volume zv2 = zv1 ! !-- Collision rate zcoll = zpcsa * zpcsa * ( 3.0_wp * pi / 4.0_wp ) ** ( 1.0_wp / 6.0_wp )& * SQRT( 6.0_wp * argas * zt / zm1 + 6.0_wp * argas * zt / zm2 )& * ( zv1 ** ( 1.0_wp / 3.0_wp ) + zv2 ** ( 1.0_wp /3.0_wp ) ) **& 2.0_wp * 1.0E+6_wp ! m3 -> cm3 zcoll = MIN( zcoll, 1.0E+10_wp ) pnuc_rate = zcoll ! (#/(cm3 s)) ELSE pnuc_rate = MIN( pnuc_rate, 1.0E+10_wp ) ENDIF pnuc_rate = pnuc_rate * 1.0E+6_wp ! (#/(m3 s)) END SUBROUTINE binnucl !------------------------------------------------------------------------------! ! Description: ! ------------ !> Calculate the nucleation rate and the size of critical clusters assuming !> ternary nucleation. Parametrisation according to: !> Napari et al. (2002), J. Chem. Phys., 116, 4221-4227 and !> Napari et al. (2002), J. Geophys. Res., 107(D19), AAC 6-1-ACC 6-6. !> Called from subroutine nucleation. !------------------------------------------------------------------------------! SUBROUTINE ternucl( pc_sa, pc_nh3, ptemp, prh, pnuc_rate, pn_crit_sa, & pn_crit_ocnv, pd_crit, pk_sa, pk_ocnv ) IMPLICIT NONE !-- Input and output variables REAL(wp), INTENT(in) :: pc_nh3 !< ammonia mixing ratio (ppt) REAL(wp), INTENT(in) :: pc_sa !< H2SO4 conc. (#/cm3) REAL(wp), INTENT(in) :: prh !< relative humidity [0-1] REAL(wp), INTENT(in) :: ptemp !< ambient temperature (K) REAL(wp), INTENT(out) :: pd_crit !< diameter of critical !< cluster (m) REAL(wp), INTENT(out) :: pk_ocnv !< Lever: if pk_ocnv = 1,organic compounds !< are involved in nucleation REAL(wp), INTENT(out) :: pk_sa !< Lever: if pk_sa = 1, H2SO4 is involved !< in nucleation REAL(wp), INTENT(out) :: pn_crit_ocnv !< number of organic molecules in !< cluster (#) REAL(wp), INTENT(out) :: pn_crit_sa !< number of H2SO4 molecules in !< cluster (#) REAL(wp), INTENT(out) :: pnuc_rate !< nucleation rate (#/(m3 s)) !-- Local variables REAL(wp) :: zlnj !< logarithm of nucleation rate !-- 1) Checking that we are in the validity range of the parameterization. !-- Validity of parameterization : DO NOT REMOVE! IF ( ptemp < 240.0_wp .OR. ptemp > 300.0_wp ) THEN message_string = 'Invalid input value: ptemp' CALL message( 'salsa_mod: ternucl', 'SA0045', 1, 2, 0, 6, 0 ) ENDIF IF ( prh < 0.05_wp .OR. prh > 0.95_wp ) THEN message_string = 'Invalid input value: prh' CALL message( 'salsa_mod: ternucl', 'SA0046', 1, 2, 0, 6, 0 ) ENDIF IF ( pc_sa < 1.0E+4_wp .OR. pc_sa > 1.0E+9_wp ) THEN message_string = 'Invalid input value: pc_sa' CALL message( 'salsa_mod: ternucl', 'SA0047', 1, 2, 0, 6, 0 ) ENDIF IF ( pc_nh3 < 0.1_wp .OR. pc_nh3 > 100.0_wp ) THEN message_string = 'Invalid input value: pc_nh3' CALL message( 'salsa_mod: ternucl', 'SA0048', 1, 2, 0, 6, 0 ) ENDIF ! !-- 2) Nucleation rate (Eq. 7 in Napari et al., 2002: Parameterization of !-- ternary nucleation of sulfuric acid - ammonia - water. zlnj = - 84.7551114741543_wp & + 0.3117595133628944_wp * prh & + 1.640089605712946_wp * prh * ptemp & - 0.003438516933381083_wp * prh * ptemp ** 2.0_wp & - 0.00001097530402419113_wp * prh * ptemp ** 3.0_wp & - 0.3552967070274677_wp / LOG( pc_sa ) & - ( 0.06651397829765026_wp * prh ) / LOG( pc_sa ) & - ( 33.84493989762471_wp * ptemp ) / LOG( pc_sa ) & - ( 7.823815852128623_wp * prh * ptemp ) / LOG( pc_sa) & + ( 0.3453602302090915_wp * ptemp ** 2.0_wp ) / LOG( pc_sa ) & + ( 0.01229375748100015_wp * prh * ptemp ** 2.0_wp ) / LOG( pc_sa ) & - ( 0.000824007160514956_wp *ptemp ** 3.0_wp ) / LOG( pc_sa ) & + ( 0.00006185539100670249_wp * prh * ptemp ** 3.0_wp ) & / LOG( pc_sa ) & + 3.137345238574998_wp * LOG( pc_sa ) & + 3.680240980277051_wp * prh * LOG( pc_sa ) & - 0.7728606202085936_wp * ptemp * LOG( pc_sa ) & - 0.204098217156962_wp * prh * ptemp * LOG( pc_sa ) & + 0.005612037586790018_wp * ptemp ** 2.0_wp * LOG( pc_sa ) & + 0.001062588391907444_wp * prh * ptemp ** 2.0_wp * LOG( pc_sa ) & - 9.74575691760229E-6_wp * ptemp ** 3.0_wp * LOG( pc_sa ) & - 1.265595265137352E-6_wp * prh * ptemp ** 3.0_wp * LOG( pc_sa ) & + 19.03593713032114_wp * LOG( pc_sa ) ** 2.0_wp & - 0.1709570721236754_wp * ptemp * LOG( pc_sa ) ** 2.0_wp & + 0.000479808018162089_wp * ptemp ** 2.0_wp * LOG( pc_sa ) ** 2.0_wp& - 4.146989369117246E-7_wp * ptemp ** 3.0_wp * LOG( pc_sa ) ** 2.0_wp& + 1.076046750412183_wp * LOG( pc_nh3 ) & + 0.6587399318567337_wp * prh * LOG( pc_nh3 ) & + 1.48932164750748_wp * ptemp * LOG( pc_nh3 ) & + 0.1905424394695381_wp * prh * ptemp * LOG( pc_nh3 ) & - 0.007960522921316015_wp * ptemp ** 2.0_wp * LOG( pc_nh3 ) & - 0.001657184248661241_wp * prh * ptemp ** 2.0_wp * LOG( pc_nh3 ) & + 7.612287245047392E-6_wp * ptemp ** 3.0_wp * LOG( pc_nh3 ) & + 3.417436525881869E-6_wp * prh * ptemp ** 3.0_wp * LOG( pc_nh3 ) & + ( 0.1655358260404061_wp * LOG( pc_nh3 ) ) / LOG( pc_sa) & + ( 0.05301667612522116_wp * prh * LOG( pc_nh3 ) ) / LOG( pc_sa ) & + ( 3.26622914116752_wp * ptemp * LOG( pc_nh3 ) ) / LOG( pc_sa ) & - ( 1.988145079742164_wp * prh * ptemp * LOG( pc_nh3 ) ) & / LOG( pc_sa ) & - ( 0.04897027401984064_wp * ptemp ** 2.0_wp * LOG( pc_nh3) ) & / LOG( pc_sa ) & + ( 0.01578269253599732_wp * prh * ptemp ** 2.0_wp * LOG( pc_nh3 ) & ) / LOG( pc_sa ) & + ( 0.0001469672236351303_wp * ptemp ** 3.0_wp * LOG( pc_nh3 ) ) & / LOG( pc_sa ) & - ( 0.00002935642836387197_wp * prh * ptemp ** 3.0_wp *LOG( pc_nh3 )& ) / LOG( pc_sa ) & + 6.526451177887659_wp * LOG( pc_sa ) * LOG( pc_nh3 ) & - 0.2580021816722099_wp * ptemp * LOG( pc_sa ) * LOG( pc_nh3 ) & + 0.001434563104474292_wp * ptemp ** 2.0_wp * LOG( pc_sa ) & * LOG( pc_nh3 ) & - 2.020361939304473E-6_wp * ptemp ** 3.0_wp * LOG( pc_sa ) & * LOG( pc_nh3 ) & - 0.160335824596627_wp * LOG( pc_sa ) ** 2.0_wp * LOG( pc_nh3 ) & + 0.00889880721460806_wp * ptemp * LOG( pc_sa ) ** 2.0_wp & * LOG( pc_nh3 ) & - 0.00005395139051155007_wp * ptemp ** 2.0_wp & * LOG( pc_sa) ** 2.0_wp * LOG( pc_nh3 ) & + 8.39521718689596E-8_wp * ptemp ** 3.0_wp * LOG( pc_sa ) ** 2.0_wp& * LOG( pc_nh3 ) & + 6.091597586754857_wp * LOG( pc_nh3 ) ** 2.0_wp & + 8.5786763679309_wp * prh * LOG( pc_nh3 ) ** 2.0_wp & - 1.253783854872055_wp * ptemp * LOG( pc_nh3 ) ** 2.0_wp & - 0.1123577232346848_wp * prh * ptemp * LOG( pc_nh3 ) ** 2.0_wp & + 0.00939835595219825_wp * ptemp ** 2.0_wp * LOG( pc_nh3 ) ** 2.0_wp& + 0.0004726256283031513_wp * prh * ptemp ** 2.0_wp & * LOG( pc_nh3) ** 2.0_wp & - 0.00001749269360523252_wp * ptemp ** 3.0_wp & * LOG( pc_nh3 ) ** 2.0_wp & - 6.483647863710339E-7_wp * prh * ptemp ** 3.0_wp & * LOG( pc_nh3 ) ** 2.0_wp & + ( 0.7284285726576598_wp * LOG( pc_nh3 ) ** 2.0_wp ) / LOG( pc_sa )& + ( 3.647355600846383_wp * ptemp * LOG( pc_nh3 ) ** 2.0_wp ) & / LOG( pc_sa ) & - ( 0.02742195276078021_wp * ptemp ** 2.0_wp & * LOG( pc_nh3) ** 2.0_wp ) / LOG( pc_sa ) & + ( 0.00004934777934047135_wp * ptemp ** 3.0_wp & * LOG( pc_nh3 ) ** 2.0_wp ) / LOG( pc_sa ) & + 41.30162491567873_wp * LOG( pc_sa ) * LOG( pc_nh3 ) ** 2.0_wp & - 0.357520416800604_wp * ptemp * LOG( pc_sa ) & * LOG( pc_nh3 ) ** 2.0_wp & + 0.000904383005178356_wp * ptemp ** 2.0_wp * LOG( pc_sa ) & * LOG( pc_nh3 ) ** 2.0_wp & - 5.737876676408978E-7_wp * ptemp ** 3.0_wp * LOG( pc_sa ) & * LOG( pc_nh3 ) ** 2.0_wp & - 2.327363918851818_wp * LOG( pc_sa ) ** 2.0_wp & * LOG( pc_nh3 ) ** 2.0_wp & + 0.02346464261919324_wp * ptemp * LOG( pc_sa ) ** 2.0_wp & * LOG( pc_nh3 ) ** 2.0_wp & - 0.000076518969516405_wp * ptemp ** 2.0_wp & * LOG( pc_sa ) ** 2.0_wp * LOG( pc_nh3 ) ** 2.0_wp & + 8.04589834836395E-8_wp * ptemp ** 3.0_wp * LOG( pc_sa ) ** 2.0_wp & * LOG( pc_nh3 ) ** 2.0_wp & - 0.02007379204248076_wp * LOG( prh ) & - 0.7521152446208771_wp * ptemp * LOG( prh ) & + 0.005258130151226247_wp * ptemp ** 2.0_wp * LOG( prh ) & - 8.98037634284419E-6_wp * ptemp ** 3.0_wp * LOG( prh ) & + ( 0.05993213079516759_wp * LOG( prh ) ) / LOG( pc_sa ) & + ( 5.964746463184173_wp * ptemp * LOG( prh ) ) / LOG( pc_sa ) & - ( 0.03624322255690942_wp * ptemp ** 2.0_wp * LOG( prh ) ) & / LOG( pc_sa ) & + ( 0.00004933369382462509_wp * ptemp ** 3.0_wp * LOG( prh ) ) & / LOG( pc_sa ) & - 0.7327310805365114_wp * LOG( pc_nh3 ) * LOG( prh ) & - 0.01841792282958795_wp * ptemp * LOG( pc_nh3 ) * LOG( prh ) & + 0.0001471855981005184_wp * ptemp ** 2.0_wp * LOG( pc_nh3 ) & * LOG( prh ) & - 2.377113195631848E-7_wp * ptemp ** 3.0_wp * LOG( pc_nh3 ) & * LOG( prh ) pnuc_rate = EXP( zlnj ) ! (#/(cm3 s)) ! !-- Check validity of parametrization IF ( pnuc_rate < 1.0E-5_wp ) THEN pnuc_rate = 0.0_wp pd_crit = 1.0E-9_wp ELSEIF ( pnuc_rate > 1.0E6_wp ) THEN message_string = 'Invalid output value: nucleation rate > 10^6 1/cm3s' CALL message( 'salsa_mod: ternucl', 'SA0049', 1, 2, 0, 6, 0 ) ENDIF pnuc_rate = pnuc_rate * 1.0E6_wp ! (#/(m3 s)) ! !-- 3) Number of H2SO4 molecules in a critical cluster (Eq. 9) pn_crit_sa = 38.16448247950508_wp + 0.7741058259731187_wp * zlnj + & 0.002988789927230632_wp * zlnj ** 2.0_wp - & 0.3576046920535017_wp * ptemp - & 0.003663583011953248_wp * zlnj * ptemp + & 0.000855300153372776_wp * ptemp ** 2.0_wp !-- Kinetic limit: at least 2 H2SO4 molecules in a cluster pn_crit_sa = MAX( pn_crit_sa, 2.0E0_wp ) ! !-- 4) Size of the critical cluster in nm (Eq. 12) pd_crit = 0.1410271086638381_wp - 0.001226253898894878_wp * zlnj - & 7.822111731550752E-6_wp * zlnj ** 2.0_wp - & 0.001567273351921166_wp * ptemp - & 0.00003075996088273962_wp * zlnj * ptemp + & 0.00001083754117202233_wp * ptemp ** 2.0_wp pd_crit = pd_crit * 2.0E-9_wp ! Diameter in m ! !-- 5) Organic compounds not involved when ternary nucleation assumed pn_crit_ocnv = 0.0_wp pk_sa = 1.0_wp pk_ocnv = 0.0_wp END SUBROUTINE ternucl !------------------------------------------------------------------------------! ! Description: ! ------------ !> Calculate the nucleation rate and the size of critical clusters assuming !> kinetic nucleation. Each sulphuric acid molecule forms an (NH4)HSO4 molecule !> in the atmosphere and two colliding (NH4)HSO4 molecules form a stable !> cluster. See Sihto et al. (2006), Atmos. Chem. Phys., 6(12), 4079-4091. !> !> Below the following assumption have been made: !> nucrate = coagcoeff*zpcsa**2 !> coagcoeff = 8*sqrt(3*boltz*ptemp*r_abs/dens_abs) !> r_abs = 0.315d-9 radius of bisulphate molecule [m] !> dens_abs = 1465 density of - " - [kg/m3] !------------------------------------------------------------------------------! SUBROUTINE kinnucl( pc_sa, pnuc_rate, pd_crit, pn_crit_sa, pn_crit_ocnv, & pk_sa, pk_ocnv ) IMPLICIT NONE !-- Input and output variables REAL(wp), INTENT(in) :: pc_sa !< H2SO4 conc. (#/m3) REAL(wp), INTENT(out) :: pd_crit !< critical diameter of clusters (m) REAL(wp), INTENT(out) :: pk_ocnv !< Lever: if pk_ocnv = 1, organic !< compounds are involved in nucleation REAL(wp), INTENT(out) :: pk_sa !< Lever: if pk_sa = 1, H2SO4 is involved !< in nucleation REAL(wp), INTENT(out) :: pn_crit_ocnv !< number of organic molecules in !< cluster (#) REAL(wp), INTENT(out) :: pn_crit_sa !< number of H2SO4 molecules in !< cluster (#) REAL(wp), INTENT(out) :: pnuc_rate !< nucl. rate (#/(m3 s)) !-- Nucleation rate (#/(m3 s)) pnuc_rate = 5.0E-13_wp * pc_sa ** 2.0_wp * 1.0E+6_wp !-- Organic compounds not involved when kinetic nucleation is assumed. pn_crit_sa = 2.0_wp pn_crit_ocnv = 0.0_wp pk_sa = 1.0_wp pk_ocnv = 0.0_wp pd_crit = 7.9375E-10_wp ! (m) END SUBROUTINE kinnucl !------------------------------------------------------------------------------! ! Description: ! ------------ !> Calculate the nucleation rate and the size of critical clusters assuming !> activation type nucleation. !> See Riipinen et al. (2007), Atmos. Chem. Phys., 7(8), 1899-1914. !------------------------------------------------------------------------------! SUBROUTINE actnucl( psa_conc, pnuc_rate, pd_crit, pn_crit_sa, pn_crit_ocnv, & pk_sa, pk_ocnv, activ ) IMPLICIT NONE !-- Input and output variables REAL(wp), INTENT(in) :: psa_conc !< H2SO4 conc. (#/m3) REAL(wp), INTENT(in) :: activ !< REAL(wp), INTENT(out) :: pd_crit !< critical diameter of clusters (m) REAL(wp), INTENT(out) :: pk_ocnv !< Lever: if pk_ocnv = 1, organic !< compounds are involved in nucleation REAL(wp), INTENT(out) :: pk_sa !< Lever: if pk_sa = 1, H2SO4 is involved !< in nucleation REAL(wp), INTENT(out) :: pn_crit_ocnv !< number of organic molecules in !< cluster (#) REAL(wp), INTENT(out) :: pn_crit_sa !< number of H2SO4 molecules in !< cluster (#) REAL(wp), INTENT(out) :: pnuc_rate !< nucl. rate (#/(m3 s)) !-- act_coeff 1e-7 by default pnuc_rate = activ * psa_conc ! (#/(m3 s)) !-- Organic compounds not involved when kinetic nucleation is assumed. pn_crit_sa = 2.0_wp pn_crit_ocnv = 0.0_wp pk_sa = 1.0_wp pk_ocnv = 0.0_wp pd_crit = 7.9375E-10_wp ! (m) END SUBROUTINE actnucl !------------------------------------------------------------------------------! ! Description: ! ------------ !> Conciders only the organic matter in nucleation. Paasonen et al. (2010) !> determined particle formation rates for 2 nm particles, J2, from different !> kind of combinations of sulphuric acid and organic matter concentration. !> See Paasonen et al. (2010), Atmos. Chem. Phys., 10, 11223-11242. !------------------------------------------------------------------------------! SUBROUTINE orgnucl( pc_org, pnuc_rate, pd_crit, pn_crit_sa, pn_crit_ocnv, & pk_sa, pk_ocnv ) IMPLICIT NONE !-- Input and output variables REAL(wp), INTENT(in) :: pc_org !< organic vapour concentration (#/m3) REAL(wp), INTENT(out) :: pd_crit !< critical diameter of clusters (m) REAL(wp), INTENT(out) :: pk_ocnv !< Lever: if pk_ocnv = 1, organic !< compounds are involved in nucleation REAL(wp), INTENT(out) :: pk_sa !< Lever: if pk_sa = 1, H2SO4 is involved !< in nucleation REAL(wp), INTENT(out) :: pn_crit_ocnv !< number of organic molecules in !< cluster (#) REAL(wp), INTENT(out) :: pn_crit_sa !< number of H2SO4 molecules in !< cluster (#) REAL(wp), INTENT(out) :: pnuc_rate !< nucl. rate (#/(m3 s)) !-- Local variables REAL(wp) :: Aorg = 1.3E-7_wp !< (1/s) (Paasonen et al. Table 4: median) !-- Homomolecular nuleation - which one? pnuc_rate = Aorg * pc_org !-- H2SO4 not involved when pure organic nucleation is assumed. pn_crit_sa = 0.0_wp pn_crit_ocnv = 1.0_wp pk_sa = 0.0_wp pk_ocnv = 1.0_wp pd_crit = 1.5E-9_wp ! (m) END SUBROUTINE orgnucl !------------------------------------------------------------------------------! ! Description: ! ------------ !> Conciders both the organic vapor and H2SO4 in nucleation - activation type !> of nucleation. !> See Paasonen et al. (2010), Atmos. Chem. Phys., 10, 11223-11242. !------------------------------------------------------------------------------! SUBROUTINE sumnucl( pc_sa, pc_org, pnuc_rate, pd_crit, pn_crit_sa, & pn_crit_ocnv, pk_sa, pk_ocnv ) IMPLICIT NONE !-- Input and output variables REAL(wp), INTENT(in) :: pc_org !< organic vapour concentration (#/m3) REAL(wp), INTENT(in) :: pc_sa !< H2SO4 conc. (#/m3) REAL(wp), INTENT(out) :: pd_crit !< critical diameter of clusters (m) REAL(wp), INTENT(out) :: pk_ocnv !< Lever: if pk_ocnv = 1, organic !< compounds are involved in nucleation REAL(wp), INTENT(out) :: pk_sa !< Lever: if pk_sa = 1, H2SO4 is involved !< in nucleation REAL(wp), INTENT(out) :: pn_crit_ocnv !< number of organic molecules in !< cluster (#) REAL(wp), INTENT(out) :: pn_crit_sa !< number of H2SO4 molecules in !< cluster (#) REAL(wp), INTENT(out) :: pnuc_rate !< nucl. rate (#/(m3 s)) !-- Local variables REAL(wp) :: As1 = 6.1E-7_wp !< (1/s) REAL(wp) :: As2 = 0.39E-7_wp !< (1/s) (Paasonen et al. Table 3.) !-- Nucleation rate (#/m3/s) pnuc_rate = As1 * pc_sa + As2 * pc_org !-- Both Organic compounds and H2SO4 are involved when SUMnucleation is assumed. pn_crit_sa = 1.0_wp pn_crit_ocnv = 1.0_wp pk_sa = 1.0_wp pk_ocnv = 1.0_wp pd_crit = 1.5E-9_wp ! (m) END SUBROUTINE sumnucl !------------------------------------------------------------------------------! ! Description: ! ------------ !> Conciders both the organic vapor and H2SO4 in nucleation - heteromolecular !> nucleation. !> See Paasonen et al. (2010), Atmos. Chem. Phys., 10, 11223-11242. !------------------------------------------------------------------------------! SUBROUTINE hetnucl( pc_sa, pc_org, pnuc_rate, pd_crit, pn_crit_sa, & pn_crit_ocnv, pk_sa, pk_ocnv ) IMPLICIT NONE !-- Input and output variables REAL(wp), INTENT(in) :: pc_org !< organic vapour concentration (#/m3) REAL(wp), INTENT(in) :: pc_sa !< H2SO4 conc. (#/m3) REAL(wp), INTENT(out) :: pd_crit !< critical diameter of clusters (m) REAL(wp), INTENT(out) :: pk_ocnv !< Lever: if pk_ocnv = 1, organic !< compounds are involved in nucleation REAL(wp), INTENT(out) :: pk_sa !< Lever: if pk_sa = 1, H2SO4 is involved !< in nucleation REAL(wp), INTENT(out) :: pn_crit_ocnv !< number of organic molecules in !< cluster (#) REAL(wp), INTENT(out) :: pn_crit_sa !< number of H2SO4 molecules in !< cluster (#) REAL(wp), INTENT(out) :: pnuc_rate !< nucl. rate (#/(m3 s)) !-- Local variables REAL(wp) :: zKhet = 4.1E-14_wp !< (cm3/s) (Paasonen et al. Table 4: median) !-- Nucleation rate (#/m3/s) pnuc_rate = zKhet * pc_sa * pc_org * 1.0E6_wp !-- Both Organic compounds and H2SO4 are involved when heteromolecular !-- nucleation is assumed. pn_crit_sa = 1.0_wp pn_crit_ocnv = 1.0_wp pk_sa = 1.0_wp pk_ocnv = 1.0_wp pd_crit = 1.5E-9_wp ! (m) END SUBROUTINE hetnucl !------------------------------------------------------------------------------! ! Description: ! ------------ !> Takes into account the homomolecular nucleation of sulphuric acid H2SO4 with !> both of the available vapours. !> See Paasonen et al. (2010), Atmos. Chem. Phys., 10, 11223-11242. !------------------------------------------------------------------------------! SUBROUTINE SAnucl( pc_sa, pc_org, pnuc_rate, pd_crit, pn_crit_sa, & pn_crit_ocnv, pk_sa, pk_ocnv ) IMPLICIT NONE !-- Input and output variables REAL(wp), INTENT(in) :: pc_org !< organic vapour concentration (#/m3) REAL(wp), INTENT(in) :: pc_sa !< H2SO4 conc. (#/m3) REAL(wp), INTENT(out) :: pd_crit !< critical diameter of clusters (m) REAL(wp), INTENT(out) :: pk_ocnv !< Lever: if pk_ocnv = 1, organic !< compounds are involved in nucleation REAL(wp), INTENT(out) :: pk_sa !< Lever: if pk_sa = 1, H2SO4 is involved !< in nucleation REAL(wp), INTENT(out) :: pn_crit_ocnv !< number of organic molecules in !< cluster (#) REAL(wp), INTENT(out) :: pn_crit_sa !< number of H2SO4 molecules in !< cluster (#) REAL(wp), INTENT(out) :: pnuc_rate !< nucleation rate (#/(m3 s)) !-- Local variables REAL(wp) :: zKsa1 = 1.1E-14_wp !< (cm3/s) REAL(wp) :: zKsa2 = 3.2E-14_wp !< (cm3/s) (Paasonen et al. Table 3.) !-- Nucleation rate (#/m3/s) pnuc_rate = ( zKsa1 * pc_sa ** 2.0_wp + zKsa2 * pc_sa * pc_org ) * 1.0E+6_wp !-- Both Organic compounds and H2SO4 are involved when SAnucleation is assumed. pn_crit_sa = 3.0_wp pn_crit_ocnv = 1.0_wp pk_sa = 1.0_wp pk_ocnv = 1.0_wp pd_crit = 1.5E-9_wp ! (m) END SUBROUTINE SAnucl !------------------------------------------------------------------------------! ! Description: ! ------------ !> Takes into account the homomolecular nucleation of both sulphuric acid and !> Lorganic with heteromolecular nucleation. !> See Paasonen et al. (2010), Atmos. Chem. Phys., 10, 11223-11242. !------------------------------------------------------------------------------! SUBROUTINE SAORGnucl( pc_sa, pc_org, pnuc_rate, pd_crit, pn_crit_sa, & pn_crit_ocnv, pk_sa, pk_ocnv ) IMPLICIT NONE !-- Input and output variables REAL(wp), INTENT(in) :: pc_org !< organic vapour concentration (#/m3) REAL(wp), INTENT(in) :: pc_sa !< H2SO4 conc. (#/m3) REAL(wp), INTENT(out) :: pd_crit !< critical diameter of clusters (m) REAL(wp), INTENT(out) :: pk_ocnv !< Lever: if pk_ocnv = 1, organic !< compounds are involved in nucleation REAL(wp), INTENT(out) :: pk_sa !< Lever: if pk_sa = 1, H2SO4 is involved !< in nucleation REAL(wp), INTENT(out) :: pn_crit_ocnv !< number of organic molecules in !< cluster (#) REAL(wp), INTENT(out) :: pn_crit_sa !< number of H2SO4 molecules in !< cluster (#) REAL(wp), INTENT(out) :: pnuc_rate !< nucl. rate (#/(m3 s)) !-- Local variables REAL(wp) :: zKs1 = 1.4E-14_wp !< (cm3/s]) REAL(wp) :: zKs2 = 2.6E-14_wp !< (cm3/s]) REAL(wp) :: zKs3 = 0.037E-14_wp !< (cm3/s]) (Paasonen et al. Table 3.) !-- Nucleation rate (#/m3/s) pnuc_rate = ( zKs1 * pc_sa **2 + zKs2 * pc_sa * pc_org + zKs3 * & pc_org ** 2.0_wp ) * 1.0E+6_wp !-- Organic compounds not involved when kinetic nucleation is assumed. pn_crit_sa = 3.0_wp pn_crit_ocnv = 3.0_wp pk_sa = 1.0_wp pk_ocnv = 1.0_wp pd_crit = 1.5E-9_wp ! (m) END SUBROUTINE SAORGnucl !------------------------------------------------------------------------------! ! Description: ! ------------ !> Function zNnuc_tayl is connected to the calculation of self-coagualtion of !> small particles. It calculates number of the particles in the size range !> [zdcrit,dx] using Taylor-expansion (please note that the expansion is not !> valid for certain rational numbers, e.g. -4/3 and -3/2) !------------------------------------------------------------------------------! FUNCTION zNnuc_tayl( d1, dx, zm_para, zjnuc_t, zeta, zGRtot ) IMPLICIT NONE INTEGER(iwp) :: i REAL(wp) :: d1 REAL(wp) :: dx REAL(wp) :: zjnuc_t REAL(wp) :: zeta REAL(wp) :: term1 REAL(wp) :: term2 REAL(wp) :: term3 REAL(wp) :: term4 REAL(wp) :: term5 REAL(wp) :: zNnuc_tayl REAL(wp) :: zGRtot REAL(wp) :: zm_para zNnuc_tayl = 0.0_wp DO i = 0, 29 IF ( i == 0 .OR. i == 1 ) THEN term1 = 1.0_wp ELSE term1 = term1 * REAL( i, SELECTED_REAL_KIND(12,307) ) END IF term2 = ( REAL( i, SELECTED_REAL_KIND(12,307) ) * ( zm_para + 1.0_wp & ) + 1.0_wp ) * term1 term3 = zeta ** i term4 = term3 / term2 term5 = REAL( i, SELECTED_REAL_KIND(12,307) ) * ( zm_para + 1.0_wp ) & + 1.0_wp zNnuc_tayl = zNnuc_tayl + term4 * ( dx ** term5 - d1 ** term5 ) ENDDO zNnuc_tayl = zNnuc_tayl * zjnuc_t * EXP( -zeta * & ( d1 ** ( zm_para + 1 ) ) ) / zGRtot END FUNCTION zNnuc_tayl !------------------------------------------------------------------------------! ! Description: ! ------------ !> Calculates the condensation of water vapour on aerosol particles. Follows the !> analytical predictor method by Jacobson (2005). !> For equations, see Jacobson (2005), Fundamentals of atmospheric modelling !> (2nd edition). !------------------------------------------------------------------------------! SUBROUTINE gpparth2o( paero, ptemp, ppres, pcs, pcw, ptstep ) IMPLICIT NONE ! !-- Input and output variables REAL(wp), INTENT(in) :: ppres !< Air pressure (Pa) REAL(wp), INTENT(in) :: pcs !< Water vapour saturation !< concentration (kg/m3) REAL(wp), INTENT(in) :: ptemp !< Ambient temperature (K) REAL(wp), INTENT(in) :: ptstep !< timestep (s) REAL(wp), INTENT(inout) :: pcw !< Water vapour concentration !< (kg/m3) TYPE(t_section), INTENT(inout) :: paero(nbins) !< Aerosol properties !-- Local variables INTEGER(iwp) :: b !< loop index INTEGER(iwp) :: nstr REAL(wp) :: adt !< internal timestep in this subroutine REAL(wp) :: adtc(nbins) REAL(wp) :: rhoair REAL(wp) :: ttot REAL(wp) :: zact !< Water activity REAL(wp) :: zaelwc1 !< Current aerosol water content REAL(wp) :: zaelwc2 !< New aerosol water content after !< equilibrium calculation REAL(wp) :: zbeta !< Transitional correction factor REAL(wp) :: zcwc !< Current water vapour mole concentration REAL(wp) :: zcwcae(nbins) !< Current water mole concentrations !< in aerosols REAL(wp) :: zcwint !< Current and new water vapour mole concentrations REAL(wp) :: zcwintae(nbins) !< Current and new water mole concentrations !< in aerosols REAL(wp) :: zcwn !< New water vapour mole concentration REAL(wp) :: zcwnae(nbins) !< New water mole concentration in aerosols REAL(wp) :: zcwsurfae(nbins) !< Surface mole concentration REAL(wp) :: zcwtot !< Total water mole concentration REAL(wp) :: zdfh2o REAL(wp) :: zhlp1 REAL(wp) :: zhlp2 REAL(wp) :: zhlp3 REAL(wp) :: zka(nbins) !< Activity coefficient REAL(wp) :: zkelvin(nbins) !< Kelvin effect REAL(wp) :: zknud REAL(wp) :: zmfph2o !< mean free path of H2O gas molecule REAL(wp) :: zmtae(nbins) !< Mass transfer coefficients REAL(wp) :: zrh !< Relative humidity [0-1] REAL(wp) :: zthcond REAL(wp) :: zwsatae(nbins) !< Water saturation ratio above aerosols ! !-- Relative humidity [0-1] zrh = pcw / pcs !-- Calculate the condensation only for 2a/2b aerosol bins nstr = in2a !-- Save the current aerosol water content, 8 in paero is H2O zaelwc1 = SUM( paero(in1a:fn2b)%volc(8) ) * arhoh2o ! !-- Equilibration: IF ( advect_particle_water ) THEN IF ( zrh < 0.98_wp .OR. .NOT. lscndh2oae ) THEN CALL equilibration( zrh, ptemp, paero, .TRUE. ) ELSE CALL equilibration( zrh, ptemp, paero, .FALSE. ) ENDIF ENDIF ! !-- The new aerosol water content after equilibrium calculation zaelwc2 = SUM( paero(in1a:fn2b)%volc(8) ) * arhoh2o !-- New water vapour mixing ratio (kg/m3) pcw = pcw - ( zaelwc2 - zaelwc1 ) * ppres * amdair / ( argas * ptemp ) ! !-- Initialise variables adtc(:) = 0.0_wp zcwc = 0.0_wp zcwcae = 0.0_wp zcwint = 0.0_wp zcwintae = 0.0_wp zcwn = 0.0_wp zcwnae = 0.0_wp zhlp1 = 0.0_wp zwsatae = 0.0_wp ! !-- Air: !-- Density (kg/m3) rhoair = amdair * ppres / ( argas * ptemp ) !-- Thermal conductivity of air zthcond = 0.023807_wp + 7.1128E-5_wp * ( ptemp - 273.16_wp ) ! !-- Water vapour: ! !-- Molecular diffusion coefficient (cm2/s) (eq.16.17) zdfh2o = ( 5.0_wp / ( 16.0_wp * avo * rhoair * 1.0E-3_wp * & ( 3.11E-8_wp ) ** 2.0_wp ) ) * SQRT( argas * 1.0E+7_wp * ptemp * & amdair * 1.0E+3_wp * ( amh2o + amdair ) * 1.0E+3_wp / ( 2.0_wp * & pi * amh2o * 1.0E+3_wp ) ) zdfh2o = zdfh2o * 1.0E-4 ! Unit change to m^2/s ! !-- Mean free path (eq. 15.25 & 16.29) zmfph2o = 3.0_wp * zdfh2o * SQRT( pi * amh2o / ( 8.0_wp * argas * ptemp ) ) zka = 1.0_wp ! Assume activity coefficients as 1 for now. ! !-- Kelvin effect (eq. 16.33) zkelvin = 1.0_wp zkelvin(1:nbins) = EXP( 4.0_wp * surfw0 * amh2o / ( argas * ptemp * & arhoh2o * paero(1:nbins)%dwet) ) ! ! --Aerosols: zmtae(:) = 0.0_wp ! mass transfer coefficient zcwsurfae(:) = 0.0_wp ! surface mole concentrations DO b = 1, nbins IF ( paero(b)%numc > nclim .AND. zrh > 0.98_wp ) THEN ! !-- Water activity zact = acth2o( paero(b) ) ! !-- Saturation mole concentration over flat surface. Limit the super- !-- saturation to max 1.01 for the mass transfer. Experimental! zcwsurfae(b) = MAX( pcs, pcw / 1.01_wp ) * rhoair / amh2o ! !-- Equilibrium saturation ratio zwsatae(b) = zact * zkelvin(b) ! !-- Knudsen number (eq. 16.20) zknud = 2.0_wp * zmfph2o / paero(b)%dwet ! !-- Transitional correction factor (Fuks & Sutugin, 1971) zbeta = ( zknud + 1.0_wp ) / ( 0.377_wp * zknud + 1.0_wp + 4.0_wp / & ( 3.0_wp * massacc(b) ) * ( zknud + zknud ** 2.0_wp ) ) ! !-- Mass transfer of H2O: Eq. 16.64 but here D^eff = zdfh2o * zbeta zhlp1 = paero(b)%numc * 2.0_wp * pi * paero(b)%dwet * zdfh2o * & zbeta !-- 1st term on the left side of the denominator in eq. 16.55 zhlp2 = amh2o * zdfh2o * alv * zwsatae(b) * zcwsurfae(b) / & ( zthcond * ptemp ) !-- 2nd term on the left side of the denominator in eq. 16.55 zhlp3 = ( (alv * amh2o ) / ( argas * ptemp ) ) - 1.0_wp !-- Full eq. 16.64: Mass transfer coefficient (1/s) zmtae(b) = zhlp1 / ( zhlp2 * zhlp3 + 1.0_wp ) ENDIF ENDDO ! !-- Current mole concentrations of water zcwc = pcw * rhoair / amh2o ! as vapour zcwcae(1:nbins) = paero(1:nbins)%volc(8) * arhoh2o / amh2o ! in aerosols zcwtot = zcwc + SUM( zcwcae ) ! total water concentration ttot = 0.0_wp adtc = 0.0_wp zcwintae = zcwcae ! !-- Substepping loop zcwint = 0.0_wp DO WHILE ( ttot < ptstep ) adt = 2.0E-2_wp ! internal timestep ! !-- New vapour concentration: (eq. 16.71) zhlp1 = zcwc + adt * ( SUM( zmtae(nstr:nbins) * zwsatae(nstr:nbins) * & zcwsurfae(nstr:nbins) ) ) ! numerator zhlp2 = 1.0_wp + adt * ( SUM( zmtae(nstr:nbins) ) ) ! denomin. zcwint = zhlp1 / zhlp2 ! new vapour concentration zcwint = MIN( zcwint, zcwtot ) IF ( ANY( paero(:)%numc > nclim ) .AND. zrh > 0.98_wp ) THEN DO b = nstr, nbins zcwintae(b) = zcwcae(b) + MIN( MAX( adt * zmtae(b) * & ( zcwint - zwsatae(b) * zcwsurfae(b) ), & -0.02_wp * zcwcae(b) ), 0.05_wp * zcwcae(b) ) zwsatae(b) = acth2o( paero(b), zcwintae(b) ) * zkelvin(b) ENDDO ENDIF zcwintae(nstr:nbins) = MAX( zcwintae(nstr:nbins), 0.0_wp ) ! !-- Update vapour concentration for consistency zcwint = zcwtot - SUM( zcwintae(1:nbins) ) !-- Update "old" values for next cycle zcwcae = zcwintae ttot = ttot + adt ENDDO ! ADT zcwn = zcwint zcwnae = zcwintae pcw = zcwn * amh2o / rhoair paero(1:nbins)%volc(8) = MAX( 0.0_wp, zcwnae(1:nbins) * amh2o / arhoh2o ) END SUBROUTINE gpparth2o !------------------------------------------------------------------------------! ! Description: ! ------------ !> Calculates the activity coefficient of liquid water !------------------------------------------------------------------------------! REAL(wp) FUNCTION acth2o( ppart, pcw ) IMPLICIT NONE TYPE(t_section), INTENT(in) :: ppart !< Aerosol properties of a bin REAL(wp), INTENT(in), OPTIONAL :: pcw !< molar concentration of water !< (mol/m3) REAL(wp) :: zns !< molar concentration of solutes (mol/m3) REAL(wp) :: znw !< molar concentration of water (mol/m3) zns = ( 3.0_wp * ( ppart%volc(1) * arhoh2so4 / amh2so4 ) + & ( ppart%volc(2) * arhooc / amoc ) + & 2.0_wp * ( ppart%volc(5) * arhoss / amss ) + & ( ppart%volc(6) * arhohno3 / amhno3 ) + & ( ppart%volc(7) * arhonh3 / amnh3 ) ) IF ( PRESENT(pcw) ) THEN znw = pcw ELSE znw = ppart%volc(8) * arhoh2o / amh2o ENDIF !-- Activity = partial pressure of water vapour / !-- sat. vapour pressure of water over a bulk liquid surface !-- = molality * activity coefficient (Jacobson, 2005: eq. 17.20-21) !-- Assume activity coefficient of 1 for water acth2o = MAX( 0.1_wp, znw / MAX( EPSILON( 1.0_wp ),( znw + zns ) ) ) END FUNCTION acth2o !------------------------------------------------------------------------------! ! Description: ! ------------ !> Calculates the dissolutional growth of particles (i.e. gas transfers to a !> particle surface and dissolves in liquid water on the surface). Treated here !> as a non-equilibrium (time-dependent) process. Gases: HNO3 and NH3 !> (Chapter 17.14 in Jacobson, 2005). ! !> Called from subroutine condensation. !> Coded by: !> Harri Kokkola (FMI) !------------------------------------------------------------------------------! SUBROUTINE gpparthno3( ppres, ptemp, paero, pghno3, pgnh3, pcw, pcs, pbeta, & ptstep ) IMPLICIT NONE ! !-- Input and output variables REAL(wp), INTENT(in) :: pbeta(nbins) !< transitional correction factor for !< aerosols REAL(wp), INTENT(in) :: ppres !< ambient pressure (Pa) REAL(wp), INTENT(in) :: pcs !< water vapour saturation !< concentration (kg/m3) REAL(wp), INTENT(in) :: ptemp !< ambient temperature (K) REAL(wp), INTENT(in) :: ptstep !< time step (s) REAL(wp), INTENT(inout) :: pghno3 !< nitric acid concentration (#/m3) REAL(wp), INTENT(inout) :: pgnh3 !< ammonia conc. (#/m3) REAL(wp), INTENT(inout) :: pcw !< water vapour concentration (kg/m3) TYPE(t_section), INTENT(inout) :: paero(nbins) !< Aerosol properties ! !-- Local variables INTEGER(iwp) :: b !< loop index REAL(wp) :: adt !< timestep REAL(wp) :: zachhso4ae(nbins) !< Activity coefficients for HHSO4 REAL(wp) :: zacnh3ae(nbins) !< Activity coefficients for NH3 REAL(wp) :: zacnh4hso2ae(nbins)!< Activity coefficients for NH4HSO2 REAL(wp) :: zacno3ae(nbins) !< Activity coefficients for HNO3 REAL(wp) :: zcgnh3eqae(nbins) !< Equilibrium gas concentration: NH3 REAL(wp) :: zcgno3eqae(nbins) !< Equilibrium gas concentration: HNO3 REAL(wp) :: zcgwaeqae(nbins) !< Equilibrium gas concentration: H2O REAL(wp) :: zcnh3c !< Current NH3 gas concentration REAL(wp) :: zcnh3int !< Intermediate NH3 gas concentration REAL(wp) :: zcnh3intae(nbins) !< Intermediate NH3 aerosol concentration REAL(wp) :: zcnh3n !< New NH3 gas concentration REAL(wp) :: zcnh3cae(nbins) !< Current NH3 in aerosols REAL(wp) :: zcnh3nae(nbins) !< New NH3 in aerosols REAL(wp) :: zcnh3tot !< Total NH3 concentration REAL(wp) :: zcno3c !< Current HNO3 gas concentration REAL(wp) :: zcno3int !< Intermediate HNO3 gas concentration REAL(wp) :: zcno3intae(nbins) !< Intermediate HNO3 aerosol concentration REAL(wp) :: zcno3n !< New HNO3 gas concentration REAL(wp) :: zcno3cae(nbins) !< Current HNO3 in aerosols REAL(wp) :: zcno3nae(nbins) !< New HNO3 in aerosols REAL(wp) :: zcno3tot !< Total HNO3 concentration REAL(wp) :: zdfvap !< Diffusion coefficient for vapors REAL(wp) :: zhlp1 !< helping variable REAL(wp) :: zhlp2 !< helping variable REAL(wp) :: zkelnh3ae(nbins) !< Kelvin effects for NH3 REAL(wp) :: zkelno3ae(nbins) !< Kelvin effect for HNO3 REAL(wp) :: zmolsae(nbins,7) !< Ion molalities from pdfite REAL(wp) :: zmtnh3ae(nbins) !< Mass transfer coefficients for NH3 REAL(wp) :: zmtno3ae(nbins) !< Mass transfer coefficients for HNO3 REAL(wp) :: zrh !< relative humidity REAL(wp) :: zsathno3ae(nbins) !< HNO3 saturation ratio REAL(wp) :: zsatnh3ae(nbins) !< NH3 saturation ratio = the partial !< pressure of a gas divided by its !< saturation vapor pressure over a surface ! !-- Initialise: adt = ptstep zachhso4ae = 0.0_wp zacnh3ae = 0.0_wp zacnh4hso2ae = 0.0_wp zacno3ae = 0.0_wp zcgnh3eqae = 0.0_wp zcgno3eqae = 0.0_wp zcnh3c = 0.0_wp zcnh3cae = 0.0_wp zcnh3int = 0.0_wp zcnh3intae = 0.0_wp zcnh3n = 0.0_wp zcnh3nae = 0.0_wp zcnh3tot = 0.0_wp zcno3c = 0.0_wp zcno3cae = 0.0_wp zcno3int = 0.0_wp zcno3intae = 0.0_wp zcno3n = 0.0_wp zcno3nae = 0.0_wp zcno3tot = 0.0_wp zhlp1 = 0.0_wp zhlp2 = 0.0_wp zkelno3ae = 1.0_wp zkelnh3ae = 1.0_wp zmolsae = 0.0_wp zmtno3ae = 0.0_wp zmtnh3ae = 0.0_wp zrh = 0.0_wp zsatnh3ae = 1.0_wp zsathno3ae = 1.0_wp ! !-- Diffusion coefficient (m2/s) zdfvap = 5.1111E-10_wp * ptemp ** 1.75_wp * ( p_0 + 1325.0_wp ) / ppres ! !-- Kelvin effects (Jacobson (2005), eq. 16.33) zkelno3ae(1:nbins) = EXP( 4.0_wp * surfw0 * amvhno3 / ( abo * ptemp * & paero(1:nbins)%dwet ) ) zkelnh3ae(1:nbins) = EXP( 4.0_wp * surfw0 * amvnh3 / ( abo * ptemp * & paero(1:nbins)%dwet ) ) ! !-- Current vapour mole concentrations (mol/m3) zcno3c = pghno3 / avo ! HNO3 zcnh3c = pgnh3 / avo ! NH3 ! !-- Current particle mole concentrations (mol/m3) zcno3cae(1:nbins) = paero(1:nbins)%volc(6) * arhohno3 / amhno3 zcnh3cae(1:nbins) = paero(1:nbins)%volc(7) * arhonh3 / amnh3 ! !-- Total mole concentrations: gas and particle phase zcno3tot = zcno3c + SUM( zcno3cae(1:nbins) ) zcnh3tot = zcnh3c + SUM( zcnh3cae(1:nbins) ) ! !-- Relative humidity [0-1] zrh = pcw / pcs ! !-- Mass transfer coefficients (Jacobson, Eq. 16.64) zmtno3ae(1:nbins) = 2.0_wp * pi * paero(1:nbins)%dwet * zdfvap * & paero(1:nbins)%numc * pbeta(1:nbins) zmtnh3ae(1:nbins) = 2.0_wp * pi * paero(1:nbins)%dwet * zdfvap * & paero(1:nbins)%numc * pbeta(1:nbins) ! !-- Get the equilibrium concentrations above aerosols CALL NONHEquil( zrh, ptemp, paero, zcgno3eqae, zcgnh3eqae, zacno3ae, & zacnh3ae, zacnh4hso2ae, zachhso4ae, zmolsae ) ! !-- NH4/HNO3 saturation ratios for aerosols CALL SVsat( ptemp, paero, zacno3ae, zacnh3ae, zacnh4hso2ae, zachhso4ae, & zcgno3eqae, zcno3cae, zcnh3cae, zkelno3ae, zkelnh3ae, & zsathno3ae, zsatnh3ae, zmolsae ) ! !-- Intermediate concentrations zhlp1 = SUM( zcno3cae(1:nbins) / ( 1.0_wp + adt * zmtno3ae(1:nbins) * & zsathno3ae(1:nbins) ) ) zhlp2 = SUM( zmtno3ae(1:nbins) / ( 1.0_wp + adt * zmtno3ae(1:nbins) * & zsathno3ae(1:nbins) ) ) zcno3int = ( zcno3tot - zhlp1 ) / ( 1.0_wp + adt * zhlp2 ) zhlp1 = SUM( zcnh3cae(1:nbins) / ( 1.0_wp + adt * zmtnh3ae(1:nbins) * & zsatnh3ae(1:nbins) ) ) zhlp2 = SUM( zmtnh3ae(1:nbins) / ( 1.0_wp + adt * zmtnh3ae(1:nbins) * & zsatnh3ae(1:nbins) ) ) zcnh3int = ( zcnh3tot - zhlp1 )/( 1.0_wp + adt * zhlp2 ) zcno3int = MIN(zcno3int, zcno3tot) zcnh3int = MIN(zcnh3int, zcnh3tot) ! !-- Calculate the new particle concentrations zcno3intae = zcno3cae zcnh3intae = zcnh3cae DO b = 1, nbins zcno3intae(b) = ( zcno3cae(b) + adt * zmtno3ae(b) * zcno3int ) / & ( 1.0_wp + adt * zmtno3ae(b) * zsathno3ae(b) ) zcnh3intae(b) = ( zcnh3cae(b) + adt * zmtnh3ae(b) * zcnh3int ) / & ( 1.0_wp + adt * zmtnh3ae(b) * zsatnh3ae(b) ) ENDDO zcno3intae(1:nbins) = MAX( zcno3intae(1:nbins), 0.0_wp ) zcnh3intae(1:nbins) = MAX( zcnh3intae(1:nbins), 0.0_wp ) zcno3n = zcno3int ! Final molar gas concentration of HNO3 zcno3nae = zcno3intae ! Final molar particle concentration of HNO3 zcnh3n = zcnh3int ! Final molar gas concentration of NH3 zcnh3nae = zcnh3intae ! Final molar particle concentration of NH3 ! !-- Model timestep reached - update the new arrays pghno3 = zcno3n * avo pgnh3 = zcnh3n * avo DO b = in1a, fn2b paero(b)%volc(6) = zcno3nae(b) * amhno3 / arhohno3 paero(b)%volc(7) = zcnh3nae(b) * amnh3 / arhonh3 ENDDO END SUBROUTINE gpparthno3 !------------------------------------------------------------------------------! ! Description: ! ------------ !> Calculate the equilibrium concentrations above aerosols (reference?) !------------------------------------------------------------------------------! SUBROUTINE NONHEquil( prh, ptemp, ppart, pcgno3eq, pcgnh3eq, pgammano, & pgammanh, pgammanh4hso2, pgammahhso4, pmols ) IMPLICIT NONE REAL(wp), INTENT(in) :: prh !< relative humidity REAL(wp), INTENT(in) :: ptemp !< ambient temperature (K) TYPE(t_section), INTENT(inout) :: ppart(nbins) !< Aerosol properties !-- Equilibrium molar concentration above aerosols: REAL(wp), INTENT(inout) :: pcgnh3eq(nbins) !< of NH3 REAL(wp), INTENT(inout) :: pcgno3eq(nbins) !< of HNO3 !< Activity coefficients: REAL(wp), INTENT(inout) :: pgammahhso4(nbins) !< HHSO4 REAL(wp), INTENT(inout) :: pgammanh(nbins) !< NH3 REAL(wp), INTENT(inout) :: pgammanh4hso2(nbins) !< NH4HSO2 REAL(wp), INTENT(inout) :: pgammano(nbins) !< HNO3 REAL(wp), INTENT(inout) :: pmols(nbins,7) !< Ion molalities INTEGER(iwp) :: b REAL(wp) :: zgammas(7) !< Activity coefficients REAL(wp) :: zhlp !< Dummy variable REAL(wp) :: zions(7) !< molar concentration of ion (mol/m3) REAL(wp) :: zphcl !< Equilibrium vapor pressures (Pa??) REAL(wp) :: zphno3 !< Equilibrium vapor pressures (Pa??) REAL(wp) :: zpnh3 !< Equilibrium vapor pressures (Pa??) REAL(wp) :: zwatertotal !< Total water in particles (mol/m3) ??? zgammas = 0.0_wp zhlp = 0.0_wp zions = 0.0_wp zphcl = 0.0_wp zphno3 = 0.0_wp zpnh3 = 0.0_wp zwatertotal = 0.0_wp DO b = 1, nbins IF ( ppart(b)%numc < nclim ) CYCLE ! !-- 2*H2SO4 + CL + NO3 - Na - NH4 zhlp = 2.0_wp * ppart(b)%volc(1) * arhoh2so4 / amh2so4 + & ppart(b)%volc(5) * arhoss / amss + & ppart(b)%volc(6) * arhohno3 / amhno3 - & ppart(b)%volc(5) * arhoss / amss - & ppart(b)%volc(7) * arhonh3 / amnh3 zhlp = MAX( zhlp, 1.0E-30_wp ) zions(1) = zhlp ! H+ zions(2) = ppart(b)%volc(7) * arhonh3 / amnh3 ! NH4+ zions(3) = ppart(b)%volc(5) * arhoss / amss ! Na+ zions(4) = ppart(b)%volc(1) * arhoh2so4 / amh2so4 ! SO4(2-) zions(5) = 0.0_wp ! HSO4- zions(6) = ppart(b)%volc(6) * arhohno3 / amhno3 ! NO3- zions(7) = ppart(b)%volc(5) * arhoss / amss ! Cl- zwatertotal = ppart(b)%volc(8) * arhoh2o / amh2o IF ( zwatertotal > 1.0E-30_wp ) THEN CALL inorganic_pdfite( prh, ptemp, zions, zwatertotal, zphno3, zphcl,& zpnh3, zgammas, pmols(b,:) ) ENDIF ! !-- Activity coefficients pgammano(b) = zgammas(1) ! HNO3 pgammanh(b) = zgammas(3) ! NH3 pgammanh4hso2(b) = zgammas(6) ! NH4HSO2 pgammahhso4(b) = zgammas(7) ! HHSO4 ! !-- Equilibrium molar concentrations (mol/m3) from equlibrium pressures (Pa) pcgno3eq(b) = zphno3 / ( argas * ptemp ) pcgnh3eq(b) = zpnh3 / ( argas * ptemp ) ENDDO END SUBROUTINE NONHEquil !------------------------------------------------------------------------------! ! Description: ! ------------ !> Calculate saturation ratios of NH4 and HNO3 for aerosols !------------------------------------------------------------------------------! SUBROUTINE SVsat( ptemp, ppart, pachno3, pacnh3, pacnh4hso2, pachhso4, & pchno3eq, pchno3, pcnh3, pkelhno3, pkelnh3, psathno3, & psatnh3, pmols ) IMPLICIT NONE REAL(wp), INTENT(in) :: ptemp !< ambient temperature (K) TYPE(t_section), INTENT(inout) :: ppart(nbins) !< Aerosol properties !-- Activity coefficients REAL(wp), INTENT(in) :: pachhso4(nbins) !< REAL(wp), INTENT(in) :: pacnh3(nbins) !< REAL(wp), INTENT(in) :: pacnh4hso2(nbins) !< REAL(wp), INTENT(in) :: pachno3(nbins) !< REAL(wp), INTENT(in) :: pchno3eq(nbins) !< Equilibrium surface concentration !< of HNO3 REAL(wp), INTENT(in) :: pchno3(nbins) !< Current particle mole !< concentration of HNO3 (mol/m3) REAL(wp), INTENT(in) :: pcnh3(nbins) !< Current particle mole !< concentration of NH3 (mol/m3) REAL(wp), INTENT(in) :: pkelhno3(nbins) !< Kelvin effect for HNO3 REAL(wp), INTENT(in) :: pkelnh3(nbins) !< Kelvin effect for NH3 REAL(wp), INTENT(in) :: pmols(nbins,7) !-- Saturation ratios REAL(wp), INTENT(out) :: psathno3(nbins) !< REAL(wp), INTENT(out) :: psatnh3(nbins) !< INTEGER :: b !< running index for aerosol bins !-- Constants for calculating equilibrium constants: REAL(wp), PARAMETER :: a1 = -22.52_wp !< REAL(wp), PARAMETER :: a2 = -1.50_wp !< REAL(wp), PARAMETER :: a3 = 13.79_wp !< REAL(wp), PARAMETER :: a4 = 29.17_wp !< REAL(wp), PARAMETER :: b1 = 26.92_wp !< REAL(wp), PARAMETER :: b2 = 26.92_wp !< REAL(wp), PARAMETER :: b3 = -5.39_wp !< REAL(wp), PARAMETER :: b4 = 16.84_wp !< REAL(wp), PARAMETER :: K01 = 1.01E-14_wp !< REAL(wp), PARAMETER :: K02 = 1.81E-5_wp !< REAL(wp), PARAMETER :: K03 = 57.64_wp !< REAL(wp), PARAMETER :: K04 = 2.51E+6_wp !< !-- Equilibrium constants of equilibrium reactions REAL(wp) :: KllH2O !< H2O(aq) <--> H+ + OH- (mol/kg) REAL(wp) :: KllNH3 !< NH3(aq) + H2O(aq) <--> NH4+ + OH- (mol/kg) REAL(wp) :: KglNH3 !< NH3(g) <--> NH3(aq) (mol/kg/atm) REAL(wp) :: KglHNO3 !< HNO3(g) <--> H+ + NO3- (mol2/kg2/atm) REAL(wp) :: zmolno3 !< molality of NO3- (mol/kg) REAL(wp) :: zmolhp !< molality of H+ (mol/kg) REAL(wp) :: zmolso4 !< molality of SO4(2-) (mol/kg) REAL(wp) :: zmolcl !< molality of Cl (mol/kg) REAL(wp) :: zmolnh4 !< Molality of NH4 (mol/kg) REAL(wp) :: zmolna !< Molality of Na (mol/kg) REAL(wp) :: zhlp1 !< REAL(wp) :: zhlp2 !< REAL(wp) :: zhlp3 !< REAL(wp) :: zxi !< REAL(wp) :: zt0 !< Reference temp zhlp1 = 0.0_wp zhlp2 = 0.0_wp zhlp3 = 0.0_wp zmolcl = 0.0_wp zmolhp = 0.0_wp zmolna = 0.0_wp zmolnh4 = 0.0_wp zmolno3 = 0.0_wp zmolso4 = 0.0_wp zt0 = 298.15_wp zxi = 0.0_wp ! !-- Calculates equlibrium rate constants based on Table B.7 in Jacobson (2005) !-- K^ll_H20, K^ll_NH3, K^gl_NH3, K^gl_HNO3 zhlp1 = zt0 / ptemp zhlp2 = zhlp1 - 1.0_wp zhlp3 = 1.0_wp + LOG( zhlp1 ) - zhlp1 KllH2O = K01 * EXP( a1 * zhlp2 + b1 * zhlp3 ) KllNH3 = K02 * EXP( a2 * zhlp2 + b2 * zhlp3 ) KglNH3 = K03 * EXP( a3 * zhlp2 + b3 * zhlp3 ) KglHNO3 = K04 * EXP( a4 * zhlp2 + b4 * zhlp3 ) DO b = 1, nbins IF ( ppart(b)%numc > nclim .AND. ppart(b)%volc(8) > 1.0E-30_wp ) THEN ! !-- Molality of H+ and NO3- zhlp1 = pcnh3(b) * amnh3 + ppart(b)%volc(1) * arhoh2so4 + & ppart(b)%volc(2) * arhooc + ppart(b)%volc(5) * arhoss + & ppart(b)%volc(8) * arhoh2o zmolno3 = pchno3(b) / zhlp1 !< mol/kg ! !-- Particle mole concentration ratio: (NH3+SS)/H2SO4 zxi = ( pcnh3(b) + ppart(b)%volc(5) * arhoss / amss ) / & ( ppart(b)%volc(1) * arhoh2so4 / amh2so4 ) IF ( zxi <= 2.0_wp ) THEN ! !-- Molality of SO4(2-) zhlp1 = pcnh3(b) * amnh3 + pchno3(b) * amhno3 + & ppart(b)%volc(2) * arhooc + ppart(b)%volc(5) * arhoss + & ppart(b)%volc(8) * arhoh2o zmolso4 = ( ppart(b)%volc(1) * arhoh2so4 / amh2so4 ) / zhlp1 ! !-- Molality of Cl- zhlp1 = pcnh3(b) * amnh3 + pchno3(b) * amhno3 + & ppart(b)%volc(2) * arhooc + ppart(b)%volc(1) * arhoh2so4 & + ppart(b)%volc(8) * arhoh2o zmolcl = ( ppart(b)%volc(5) * arhoss / amss ) / zhlp1 ! !-- Molality of NH4+ zhlp1 = pchno3(b) * amhno3 + ppart(b)%volc(1) * arhoh2so4 + & ppart(b)%volc(2) * arhooc + ppart(b)%volc(5) * arhoss + & ppart(b)%volc(8) * arhoh2o zmolnh4 = pcnh3(b) / zhlp1 ! !-- Molality of Na+ zmolna = zmolcl ! !-- Molality of H+ zmolhp = 2.0_wp * zmolso4 + zmolno3 + zmolcl - ( zmolnh4 + zmolna ) ELSE zhlp2 = pkelhno3(b) * zmolno3 * pachno3(b) ** 2.0_wp ! !-- Mona debugging IF ( zhlp2 > 1.0E-30_wp ) THEN zmolhp = KglHNO3 * pchno3eq(b) / zhlp2 ! Eq. 17.38 ELSE zmolhp = 0.0_wp ENDIF ENDIF zhlp1 = ppart(b)%volc(8) * arhoh2o * argas * ptemp * KglHNO3 ! !-- Saturation ratio for NH3 and for HNO3 IF ( zmolhp > 0.0_wp ) THEN zhlp2 = pkelnh3(b) / ( zhlp1 * zmolhp ) zhlp3 = KllH2O / ( KllNH3 + KglNH3 ) psatnh3(b) = zhlp2 * ( ( pacnh4hso2(b) / pachhso4(b) ) **2.0_wp ) & * zhlp3 psathno3(b) = ( pkelhno3(b) * zmolhp * pachno3(b)**2.0_wp ) / zhlp1 ELSE psatnh3(b) = 1.0_wp psathno3(b) = 1.0_wp ENDIF ELSE psatnh3(b) = 1.0_wp psathno3(b) = 1.0_wp ENDIF ENDDO END SUBROUTINE SVsat !------------------------------------------------------------------------------! ! Description: ! ------------ !> Prototype module for calculating the water content of a mixed inorganic/ !> organic particle + equilibrium water vapour pressure above the solution !> (HNO3, HCL, NH3 and representative organic compounds. Efficient calculation !> of the partitioning of species between gas and aerosol. Based in a chamber !> study. ! !> Written by Dave Topping. Pure organic component properties predicted by Mark !> Barley based on VOCs predicted in MCM simulations performed by Mike Jenkin. !> Delivered by Gordon McFiggans as Deliverable D22 from WP1.4 in the EU FP6 !> EUCAARI Integrated Project. ! !> Queries concerning the use of this code through Gordon McFiggans, !> g.mcfiggans@manchester.ac.uk, !> Ownership: D. Topping, Centre for Atmospheric Sciences, University of !> Manchester, 2007 ! !> Rewritten to PALM by Mona Kurppa, UHel, 2017 !------------------------------------------------------------------------------! SUBROUTINE inorganic_pdfite( RH, temp, ions, water_total, Press_HNO3, & Press_HCL, Press_NH3, gamma_out, mols_out ) IMPLICIT NONE REAL(wp), DIMENSION(:) :: gamma_out !< Activity coefficient for calculating !< the non-ideal dissociation constants !< 1: HNO3, 2: HCL, 3: NH4+/H+ (NH3) !< 4: HHSO4**2/H2SO4, !< 5: H2SO4**3/HHSO4**2 !< 6: NH4HSO2, 7: HHSO4 REAL(wp), DIMENSION(:) :: ions !< ion molarities (mol/m3) !< 1: H+, 2: NH4+, 3: Na+, 4: SO4(2-), !< 5: HSO4-, 6: NO3-, 7: Cl- REAL(wp), DIMENSION(7) :: ions_mol !< ion molalities (mol/kg) !< 1: H+, 2: NH4+, 3: Na+, 4: SO4(2-), !< 5: HSO4-, 6: NO3-, 7: Cl- REAL(wp), DIMENSION(:) :: mols_out !< ion molality output (mol/kg) !< 1: H+, 2: NH4+, 3: Na+, 4: SO4(2-), !< 5: HSO4-, 6: NO3-, 7: Cl- REAL(wp) :: act_product !< ionic activity coef. product: !< = (gamma_h2so4**3d0) / !< (gamma_hhso4**2d0) REAL(wp) :: ammonium_chloride !< REAL(wp) :: ammonium_chloride_eq_frac !< REAL(wp) :: ammonium_nitrate !< REAL(wp) :: ammonium_nitrate_eq_frac !< REAL(wp) :: ammonium_sulphate !< REAL(wp) :: ammonium_sulphate_eq_frac !< REAL(wp) :: binary_h2so4 !< binary H2SO4 activity coeff. REAL(wp) :: binary_hcl !< binary HCL activity coeff. REAL(wp) :: binary_hhso4 !< binary HHSO4 activity coeff. REAL(wp) :: binary_hno3 !< binary HNO3 activity coeff. REAL(wp) :: binary_nh4hso4 !< binary NH4HSO4 activity coeff. REAL(wp) :: charge_sum !< sum of ionic charges REAL(wp) :: gamma_h2so4 !< activity coefficient REAL(wp) :: gamma_hcl !< activity coefficient REAL(wp) :: gamma_hhso4 !< activity coeffient REAL(wp) :: gamma_hno3 !< activity coefficient REAL(wp) :: gamma_nh3 !< activity coefficient REAL(wp) :: gamma_nh4hso4 !< activity coefficient REAL(wp) :: h_out !< REAL(wp) :: h_real !< new hydrogen ion conc. REAL(wp) :: H2SO4_hcl !< contribution of H2SO4 REAL(wp) :: H2SO4_hno3 !< contribution of H2SO4 REAL(wp) :: H2SO4_nh3 !< contribution of H2SO4 REAL(wp) :: H2SO4_nh4hso4 !< contribution of H2SO4 REAL(wp) :: HCL_h2so4 !< contribution of HCL REAL(wp) :: HCL_hhso4 !< contribution of HCL REAL(wp) :: HCL_hno3 !< contribution of HCL REAL(wp) :: HCL_nh3 !< contribution of HCL REAL(wp) :: HCL_nh4hso4 !< contribution of HCL REAL(wp) :: henrys_temp_dep !< temperature dependence of !< Henry's Law REAL(wp) :: HNO3_h2so4 !< contribution of HNO3 REAL(wp) :: HNO3_hcl !< contribution of HNO3 REAL(wp) :: HNO3_hhso4 !< contribution of HNO3 REAL(wp) :: HNO3_nh3 !< contribution of HNO3 REAL(wp) :: HNO3_nh4hso4 !< contribution of HNO3 REAL(wp) :: hso4_out !< REAL(wp) :: hso4_real !< new bisulphate ion conc. REAL(wp) :: hydrochloric_acid !< REAL(wp) :: hydrochloric_acid_eq_frac !< REAL(wp) :: Kh !< equilibrium constant for H+ REAL(wp) :: K_hcl !< equilibrium constant of HCL REAL(wp) :: K_hno3 !< equilibrium constant of HNO3 REAL(wp) :: Knh4 !< equilibrium constant for NH4+ REAL(wp) :: Kw !< equil. const. for water_surface REAL(wp) :: Ln_h2so4_act !< gamma_h2so4 = EXP(Ln_h2so4_act) REAL(wp) :: Ln_HCL_act !< gamma_hcl = EXP( Ln_HCL_act ) REAL(wp) :: Ln_hhso4_act !< gamma_hhso4 = EXP(Ln_hhso4_act) REAL(wp) :: Ln_HNO3_act !< gamma_hno3 = EXP( Ln_HNO3_act ) REAL(wp) :: Ln_NH4HSO4_act !< gamma_nh4hso4 = !< EXP( Ln_NH4HSO4_act ) REAL(wp) :: molality_ratio_nh3 !< molality ratio of NH3 !< (NH4+ and H+) REAL(wp) :: Na2SO4_h2so4 !< contribution of Na2SO4 REAL(wp) :: Na2SO4_hcl !< contribution of Na2SO4 REAL(wp) :: Na2SO4_hhso4 !< contribution of Na2SO4 REAL(wp) :: Na2SO4_hno3 !< contribution of Na2SO4 REAL(wp) :: Na2SO4_nh3 !< contribution of Na2SO4 REAL(wp) :: Na2SO4_nh4hso4 !< contribution of Na2SO4 REAL(wp) :: NaCl_h2so4 !< contribution of NaCl REAL(wp) :: NaCl_hcl !< contribution of NaCl REAL(wp) :: NaCl_hhso4 !< contribution of NaCl REAL(wp) :: NaCl_hno3 !< contribution of NaCl REAL(wp) :: NaCl_nh3 !< contribution of NaCl REAL(wp) :: NaCl_nh4hso4 !< contribution of NaCl REAL(wp) :: NaNO3_h2so4 !< contribution of NaNO3 REAL(wp) :: NaNO3_hcl !< contribution of NaNO3 REAL(wp) :: NaNO3_hhso4 !< contribution of NaNO3 REAL(wp) :: NaNO3_hno3 !< contribution of NaNO3 REAL(wp) :: NaNO3_nh3 !< contribution of NaNO3 REAL(wp) :: NaNO3_nh4hso4 !< contribution of NaNO3 REAL(wp) :: NH42SO4_h2so4 !< contribution of NH42SO4 REAL(wp) :: NH42SO4_hcl !< contribution of NH42SO4 REAL(wp) :: NH42SO4_hhso4 !< contribution of NH42SO4 REAL(wp) :: NH42SO4_hno3 !< contribution of NH42SO4 REAL(wp) :: NH42SO4_nh3 !< contribution of NH42SO4 REAL(wp) :: NH42SO4_nh4hso4 !< contribution of NH42SO4 REAL(wp) :: NH4Cl_h2so4 !< contribution of NH4Cl REAL(wp) :: NH4Cl_hcl !< contribution of NH4Cl REAL(wp) :: NH4Cl_hhso4 !< contribution of NH4Cl REAL(wp) :: NH4Cl_hno3 !< contribution of NH4Cl REAL(wp) :: NH4Cl_nh3 !< contribution of NH4Cl REAL(wp) :: NH4Cl_nh4hso4 !< contribution of NH4Cl REAL(wp) :: NH4NO3_h2so4 !< contribution of NH4NO3 REAL(wp) :: NH4NO3_hcl !< contribution of NH4NO3 REAL(wp) :: NH4NO3_hhso4 !< contribution of NH4NO3 REAL(wp) :: NH4NO3_hno3 !< contribution of NH4NO3 REAL(wp) :: NH4NO3_nh3 !< contribution of NH4NO3 REAL(wp) :: NH4NO3_nh4hso4 !< contribution of NH4NO3 REAL(wp) :: nitric_acid !< REAL(wp) :: nitric_acid_eq_frac !< Equivalent fractions REAL(wp) :: Press_HCL !< partial pressure of HCL REAL(wp) :: Press_HNO3 !< partial pressure of HNO3 REAL(wp) :: Press_NH3 !< partial pressure of NH3 REAL(wp) :: RH !< relative humidity [0-1] REAL(wp) :: temp !< temperature REAL(wp) :: so4_out !< REAL(wp) :: so4_real !< new sulpate ion concentration REAL(wp) :: sodium_chloride !< REAL(wp) :: sodium_chloride_eq_frac !< REAL(wp) :: sodium_nitrate !< REAL(wp) :: sodium_nitrate_eq_frac !< REAL(wp) :: sodium_sulphate !< REAL(wp) :: sodium_sulphate_eq_frac !< REAL(wp) :: solutes !< REAL(wp) :: sulphuric_acid !< REAL(wp) :: sulphuric_acid_eq_frac !< REAL(wp) :: water_total !< REAL(wp) :: a !< auxiliary variable REAL(wp) :: b !< auxiliary variable REAL(wp) :: c !< auxiliary variable REAL(wp) :: root1 !< auxiliary variable REAL(wp) :: root2 !< auxiliary variable INTEGER(iwp) :: binary_case INTEGER(iwp) :: full_complexity ! !-- Value initialisation binary_h2so4 = 0.0_wp binary_hcl = 0.0_wp binary_hhso4 = 0.0_wp binary_hno3 = 0.0_wp binary_nh4hso4 = 0.0_wp henrys_temp_dep = ( 1.0_wp / temp - 1.0_wp / 298.0_wp ) HCL_hno3 = 1.0_wp H2SO4_hno3 = 1.0_wp NH42SO4_hno3 = 1.0_wp NH4NO3_hno3 = 1.0_wp NH4Cl_hno3 = 1.0_wp Na2SO4_hno3 = 1.0_wp NaNO3_hno3 = 1.0_wp NaCl_hno3 = 1.0_wp HNO3_hcl = 1.0_wp H2SO4_hcl = 1.0_wp NH42SO4_hcl = 1.0_wp NH4NO3_hcl = 1.0_wp NH4Cl_hcl = 1.0_wp Na2SO4_hcl = 1.0_wp NaNO3_hcl = 1.0_wp NaCl_hcl = 1.0_wp HNO3_nh3 = 1.0_wp HCL_nh3 = 1.0_wp H2SO4_nh3 = 1.0_wp NH42SO4_nh3 = 1.0_wp NH4NO3_nh3 = 1.0_wp NH4Cl_nh3 = 1.0_wp Na2SO4_nh3 = 1.0_wp NaNO3_nh3 = 1.0_wp NaCl_nh3 = 1.0_wp HNO3_hhso4 = 1.0_wp HCL_hhso4 = 1.0_wp NH42SO4_hhso4 = 1.0_wp NH4NO3_hhso4 = 1.0_wp NH4Cl_hhso4 = 1.0_wp Na2SO4_hhso4 = 1.0_wp NaNO3_hhso4 = 1.0_wp NaCl_hhso4 = 1.0_wp HNO3_h2so4 = 1.0_wp HCL_h2so4 = 1.0_wp NH42SO4_h2so4 = 1.0_wp NH4NO3_h2so4 = 1.0_wp NH4Cl_h2so4 = 1.0_wp Na2SO4_h2so4 = 1.0_wp NaNO3_h2so4 = 1.0_wp NaCl_h2so4 = 1.0_wp !-- New NH3 variables HNO3_nh4hso4 = 1.0_wp HCL_nh4hso4 = 1.0_wp H2SO4_nh4hso4 = 1.0_wp NH42SO4_nh4hso4 = 1.0_wp NH4NO3_nh4hso4 = 1.0_wp NH4Cl_nh4hso4 = 1.0_wp Na2SO4_nh4hso4 = 1.0_wp NaNO3_nh4hso4 = 1.0_wp NaCl_nh4hso4 = 1.0_wp ! !-- Juha Tonttila added mols_out = 0.0_wp Press_HNO3 = 0.0_wp Press_HCL = 0.0_wp Press_NH3 = 0.0_wp !< Initialising vapour pressure over the !< multicomponent particle gamma_out = 1.0_wp !< i.e. don't alter the ideal mixing ratios if !< there's nothing there. ! !-- 1) - COMPOSITION DEFINITIONS ! !-- a) Inorganic ion pairing: !-- In order to calculate the water content, which is also used in !-- calculating vapour pressures, one needs to pair the anions and cations !-- for use in the ZSR mixing rule. The equation provided by Clegg et al. !-- (2001) is used for ion pairing. The solutes chosen comprise of 9 !-- inorganic salts and acids which provide a pairing between each anion and !-- cation: (NH4)2SO4, NH4NO3, NH4Cl, Na2SO4, NaNO3, NaCl, H2SO4, HNO3, HCL. !-- The organic compound is treated as a seperate solute. !-- Ions: 1: H+, 2: NH4+, 3: Na+, 4: SO4(2-), 5: HSO4-, 6: NO3-, 7: Cl- ! charge_sum = ions(1) + ions(2) + ions(3) + 2.0_wp * ions(4) + ions(5) + & ions(6) + ions(7) nitric_acid = 0.0_wp ! HNO3 nitric_acid = ( 2.0_wp * ions(1) * ions(6) * & ( ( 1.0_wp / 1.0_wp ) ** 0.5_wp ) ) / ( charge_sum ) hydrochloric_acid = 0.0_wp ! HCL hydrochloric_acid = ( 2.0_wp * ions(1) * ions(7) * & ( ( 1.0_wp / 1.0_wp ) ** 0.5_wp ) ) / ( charge_sum ) sulphuric_acid = 0.0_wp ! H2SO4 sulphuric_acid = ( 2.0_wp * ions(1) * ions(4) * & ( ( 2.0_wp / 2.0_wp ) ** 0.5_wp ) ) / ( charge_sum ) ammonium_sulphate = 0.0_wp ! (NH4)2SO4 ammonium_sulphate = ( 2.0_wp * ions(2) * ions(4) * & ( ( 2.0_wp / 2.0_wp ) ** 0.5_wp ) ) / ( charge_sum ) ammonium_nitrate = 0.0_wp ! NH4NO3 ammonium_nitrate = ( 2.0_wp * ions(2) * ions(6) * & ( ( 1.0_wp / 1.0_wp ) ** 0.5_wp ) ) / ( charge_sum ) ammonium_chloride = 0.0_wp ! NH4Cl ammonium_chloride = ( 2.0_wp * ions(2) * ions(7) * & ( ( 1.0_wp / 1.0_wp ) ** 0.5_wp ) ) / ( charge_sum ) sodium_sulphate = 0.0_wp ! Na2SO4 sodium_sulphate = ( 2.0_wp * ions(3) * ions(4) * & ( ( 2.0_wp / 2.0_wp ) ** 0.5_wp ) ) / ( charge_sum ) sodium_nitrate = 0.0_wp ! NaNO3 sodium_nitrate = ( 2.0_wp * ions(3) *ions(6) * & ( ( 1.0_wp / 1.0_wp ) ** 0.5_wp ) ) / ( charge_sum ) sodium_chloride = 0.0_wp ! NaCl sodium_chloride = ( 2.0_wp * ions(3) * ions(7) * & ( ( 1.0_wp / 1.0_wp ) ** 0.5_wp ) ) / ( charge_sum ) solutes = 0.0_wp solutes = 3.0_wp * sulphuric_acid + 2.0_wp * hydrochloric_acid + & 2.0_wp * nitric_acid + 3.0_wp * ammonium_sulphate + & 2.0_wp * ammonium_nitrate + 2.0_wp * ammonium_chloride + & 3.0_wp * sodium_sulphate + 2.0_wp * sodium_nitrate + & 2.0_wp * sodium_chloride ! !-- b) Inorganic equivalent fractions: !-- These values are calculated so that activity coefficients can be !-- expressed by a linear additive rule, thus allowing more efficient !-- calculations and future expansion (see more detailed description below) nitric_acid_eq_frac = 2.0_wp * nitric_acid / ( solutes ) hydrochloric_acid_eq_frac = 2.0_wp * hydrochloric_acid / ( solutes ) sulphuric_acid_eq_frac = 3.0_wp * sulphuric_acid / ( solutes ) ammonium_sulphate_eq_frac = 3.0_wp * ammonium_sulphate / ( solutes ) ammonium_nitrate_eq_frac = 2.0_wp * ammonium_nitrate / ( solutes ) ammonium_chloride_eq_frac = 2.0_wp * ammonium_chloride / ( solutes ) sodium_sulphate_eq_frac = 3.0_wp * sodium_sulphate / ( solutes ) sodium_nitrate_eq_frac = 2.0_wp * sodium_nitrate / ( solutes ) sodium_chloride_eq_frac = 2.0_wp * sodium_chloride / ( solutes ) ! !-- Inorganic ion molalities ions_mol(:) = 0.0_wp ions_mol(1) = ions(1) / ( water_total * 18.01528E-3_wp ) ! H+ ions_mol(2) = ions(2) / ( water_total * 18.01528E-3_wp ) ! NH4+ ions_mol(3) = ions(3) / ( water_total * 18.01528E-3_wp ) ! Na+ ions_mol(4) = ions(4) / ( water_total * 18.01528E-3_wp ) ! SO4(2-) ions_mol(5) = ions(5) / ( water_total * 18.01528E-3_wp ) ! HSO4(2-) ions_mol(6) = ions(6) / ( water_total * 18.01528E-3_wp ) ! NO3- ions_mol(7) = ions(7) / ( water_total * 18.01528E-3_wp ) ! Cl- !-- *** !-- At this point we may need to introduce a method for prescribing H+ when !-- there is no 'real' value for H+..i.e. in the sulphate poor domain !-- This will give a value for solve quadratic proposed by Zaveri et al. 2005 ! !-- 2) - WATER CALCULATION ! !-- a) The water content is calculated using the ZSR rule with solute !-- concentrations calculated using 1a above. Whilst the usual approximation of !-- ZSR relies on binary data consisting of 5th or higher order polynomials, in !-- this code 4 different RH regimes are used, each housing cubic equations for !-- the water associated with each solute listed above. Binary water contents !-- for inorganic components were calculated using AIM online (Clegg et al !-- 1998). The water associated with the organic compound is calculated assuming !-- ideality and that aw = RH. ! !-- b) Molality of each inorganic ion and organic solute (initial input) is !-- calculated for use in vapour pressure calculation. ! !-- 3) - BISULPHATE ION DISSOCIATION CALCULATION ! !-- The dissociation of the bisulphate ion is calculated explicitly. A solution !-- to the equilibrium equation between the bisulphate ion, hydrogen ion and !-- sulphate ion is found using tabulated equilibrium constants (referenced). It !-- is necessary to calculate the activity coefficients of HHSO4 and H2SO4 in a !-- non-iterative manner. These are calculated using the same format as !-- described in 4) below, where both activity coefficients were fit to the !-- output from ADDEM (Topping et al 2005a,b) covering an extensive composition !-- space, providing the activity coefficients and bisulphate ion dissociation !-- as a function of equivalent mole fractions and relative humidity. ! !-- NOTE: the flags "binary_case" and "full_complexity" are not used in this !-- prototype. They are used for simplification of the fit expressions when !-- using limited composition regions. This section of code calculates the !-- bisulphate ion concentration ! IF ( ions(1) > 0.0_wp .AND. ions(4) > 0.0_wp ) THEN ! !-- HHSO4: binary_case = 1 IF ( RH > 0.1_wp .AND. RH < 0.9_wp ) THEN binary_hhso4 = - 4.9521_wp * ( RH**3 ) + 9.2881_wp * ( RH**2 ) - & 10.777_wp * RH + 6.0534_wp ELSEIF ( RH >= 0.9_wp .AND. RH < 0.955_wp ) THEN binary_hhso4 = - 6.3777_wp * RH + 5.962_wp ELSEIF ( RH >= 0.955_wp .AND. RH < 0.99_wp ) THEN binary_hhso4 = 2367.2_wp * ( RH**3 ) - 6849.7_wp * ( RH**2 ) + & 6600.9_wp * RH - 2118.7_wp ELSEIF ( RH >= 0.99_wp .AND. RH < 0.9999_wp ) THEN binary_hhso4 = 3E-7_wp * ( RH**5 ) - 2E-5_wp * ( RH**4 ) + & 0.0004_wp * ( RH**3 ) - 0.0035_wp * ( RH**2 ) + & 0.0123_wp * RH - 0.3025_wp ENDIF IF ( nitric_acid > 0.0_wp ) THEN HNO3_hhso4 = - 4.2204_wp * ( RH**4 ) + 12.193_wp * ( RH**3 ) - & 12.481_wp * ( RH**2 ) + 6.459_wp * RH - 1.9004_wp ENDIF IF ( hydrochloric_acid > 0.0_wp ) THEN HCL_hhso4 = - 54.845_wp * ( RH**7 ) + 209.54_wp * ( RH**6 ) - & 336.59_wp * ( RH**5 ) + 294.21_wp * ( RH**4 ) - & 150.07_wp * ( RH**3 ) + 43.767_wp * ( RH**2 ) - & 6.5495_wp * RH + 0.60048_wp ENDIF IF ( ammonium_sulphate > 0.0_wp ) THEN NH42SO4_hhso4 = 16.768_wp * ( RH**3 ) - 28.75_wp * ( RH**2 ) + & 20.011_wp * RH - 8.3206_wp ENDIF IF ( ammonium_nitrate > 0.0_wp ) THEN NH4NO3_hhso4 = - 17.184_wp * ( RH**4 ) + 56.834_wp * ( RH**3 ) - & 65.765_wp * ( RH**2 ) + 35.321_wp * RH - 9.252_wp ENDIF IF (ammonium_chloride > 0.0_wp ) THEN IF ( RH < 0.2_wp .AND. RH >= 0.1_wp ) THEN NH4Cl_hhso4 = 3.2809_wp * RH - 2.0637_wp ELSEIF ( RH >= 0.2_wp .AND. RH < 0.99_wp ) THEN NH4Cl_hhso4 = - 1.2981_wp * ( RH**3 ) + 4.7461_wp * ( RH**2 ) - & 2.3269_wp * RH - 1.1259_wp ENDIF ENDIF IF ( sodium_sulphate > 0.0_wp ) THEN Na2SO4_hhso4 = 118.87_wp * ( RH**6 ) - 358.63_wp * ( RH**5 ) + & 435.85_wp * ( RH**4 ) - 272.88_wp * ( RH**3 ) + & 94.411_wp * ( RH**2 ) - 18.21_wp * RH + 0.45935_wp ENDIF IF ( sodium_nitrate > 0.0_wp ) THEN IF ( RH < 0.2_wp .AND. RH >= 0.1_wp ) THEN NaNO3_hhso4 = 4.8456_wp * RH - 2.5773_wp ELSEIF ( RH >= 0.2_wp .AND. RH < 0.99_wp ) THEN NaNO3_hhso4 = 0.5964_wp * ( RH**3 ) - 0.38967_wp * ( RH**2 ) + & 1.7918_wp * RH - 1.9691_wp ENDIF ENDIF IF ( sodium_chloride > 0.0_wp ) THEN IF ( RH < 0.2_wp ) THEN NaCl_hhso4 = 0.51995_wp * RH - 1.3981_wp ELSEIF ( RH >= 0.2_wp .AND. RH < 0.99_wp ) THEN NaCl_hhso4 = 1.6539_wp * RH - 1.6101_wp ENDIF ENDIF Ln_hhso4_act = binary_hhso4 + & nitric_acid_eq_frac * HNO3_hhso4 + & hydrochloric_acid_eq_frac * HCL_hhso4 + & ammonium_sulphate_eq_frac * NH42SO4_hhso4 + & ammonium_nitrate_eq_frac * NH4NO3_hhso4 + & ammonium_chloride_eq_frac * NH4Cl_hhso4 + & sodium_sulphate_eq_frac * Na2SO4_hhso4 + & sodium_nitrate_eq_frac * NaNO3_hhso4 + & sodium_chloride_eq_frac * NaCl_hhso4 gamma_hhso4 = EXP( Ln_hhso4_act ) ! molal activity coefficient of HHSO4 !-- H2SO4 (sulphuric acid): IF ( RH >= 0.1_wp .AND. RH < 0.9_wp ) THEN binary_h2so4 = 2.4493_wp * ( RH**2 ) - 6.2326_wp * RH + 2.1763_wp ELSEIF ( RH >= 0.9_wp .AND. RH < 0.98 ) THEN binary_h2so4 = 914.68_wp * ( RH**3 ) - 2502.3_wp * ( RH**2 ) + & 2281.9_wp * RH - 695.11_wp ELSEIF ( RH >= 0.98 .AND. RH < 0.9999 ) THEN binary_h2so4 = 3E-8_wp * ( RH**4 ) - 5E-6_wp * ( RH**3 ) + & 0.0003_wp * ( RH**2 ) - 0.0022_wp * RH - 1.1305_wp ENDIF IF ( nitric_acid > 0.0_wp ) THEN HNO3_h2so4 = - 16.382_wp * ( RH**5 ) + 46.677_wp * ( RH**4 ) - & 54.149_wp * ( RH**3 ) + 34.36_wp * ( RH**2 ) - & 12.54_wp * RH + 2.1368_wp ENDIF IF ( hydrochloric_acid > 0.0_wp ) THEN HCL_h2so4 = - 14.409_wp * ( RH**5 ) + 42.804_wp * ( RH**4 ) - & 47.24_wp * ( RH**3 ) + 24.668_wp * ( RH**2 ) - & 5.8015_wp * RH + 0.084627_wp ENDIF IF ( ammonium_sulphate > 0.0_wp ) THEN NH42SO4_h2so4 = 66.71_wp * ( RH**5 ) - 187.5_wp * ( RH**4 ) + & 210.57_wp * ( RH**3 ) - 121.04_wp * ( RH**2 ) + & 39.182_wp * RH - 8.0606_wp ENDIF IF ( ammonium_nitrate > 0.0_wp ) THEN NH4NO3_h2so4 = - 22.532_wp * ( RH**4 ) + 66.615_wp * ( RH**3 ) - & 74.647_wp * ( RH**2 ) + 37.638_wp * RH - 6.9711_wp ENDIF IF ( ammonium_chloride > 0.0_wp ) THEN IF ( RH >= 0.1_wp .AND. RH < 0.2_wp ) THEN NH4Cl_h2so4 = - 0.32089_wp * RH + 0.57738_wp ELSEIF ( RH >= 0.2_wp .AND. RH < 0.9_wp ) THEN NH4Cl_h2so4 = 18.089_wp * ( RH**5 ) - 51.083_wp * ( RH**4 ) + & 50.32_wp * ( RH**3 ) - 17.012_wp * ( RH**2 ) - & 0.93435_wp * RH + 1.0548_wp ELSEIF ( RH >= 0.9_wp .AND. RH < 0.99_wp ) THEN NH4Cl_h2so4 = - 1.5749_wp * RH + 1.7002_wp ENDIF ENDIF IF ( sodium_sulphate > 0.0_wp ) THEN Na2SO4_h2so4 = 29.843_wp * ( RH**4 ) - 69.417_wp * ( RH**3 ) + & 61.507_wp * ( RH**2 ) - 29.874_wp * RH + 7.7556_wp ENDIF IF ( sodium_nitrate > 0.0_wp ) THEN NaNO3_h2so4 = - 122.37_wp * ( RH**6 ) + 427.43_wp * ( RH**5 ) - & 604.68_wp * ( RH**4 ) + 443.08_wp * ( RH**3 ) - & 178.61_wp * ( RH**2 ) + 37.242_wp * RH - 1.9564_wp ENDIF IF ( sodium_chloride > 0.0_wp ) THEN NaCl_h2so4 = - 40.288_wp * ( RH**5 ) + 115.61_wp * ( RH**4 ) - & 129.99_wp * ( RH**3 ) + 72.652_wp * ( RH**2 ) - & 22.124_wp * RH + 4.2676_wp ENDIF Ln_h2so4_act = binary_h2so4 + & nitric_acid_eq_frac * HNO3_h2so4 + & hydrochloric_acid_eq_frac * HCL_h2so4 + & ammonium_sulphate_eq_frac * NH42SO4_h2so4 + & ammonium_nitrate_eq_frac * NH4NO3_h2so4 + & ammonium_chloride_eq_frac * NH4Cl_h2so4 + & sodium_sulphate_eq_frac * Na2SO4_h2so4 + & sodium_nitrate_eq_frac * NaNO3_h2so4 + & sodium_chloride_eq_frac * NaCl_h2so4 gamma_h2so4 = EXP( Ln_h2so4_act ) ! molal activity coefficient ! !-- Export activity coefficients IF ( gamma_h2so4 > 1.0E-10_wp ) THEN gamma_out(4) = ( gamma_hhso4**2.0_wp ) / gamma_h2so4 ENDIF IF ( gamma_hhso4 > 1.0E-10_wp ) THEN gamma_out(5) = ( gamma_h2so4**3.0_wp ) / ( gamma_hhso4**2.0_wp ) ENDIF ! !-- Ionic activity coefficient product act_product = ( gamma_h2so4**3.0_wp ) / ( gamma_hhso4**2.0_wp ) ! !-- Solve the quadratic equation (i.e. x in ax**2 + bx + c = 0) a = 1.0_wp b = - 1.0_wp * ( ions(4) + ions(1) + ( ( water_total * 18.0E-3_wp ) / & ( 99.0_wp * act_product ) ) ) c = ions(4) * ions(1) root1 = ( ( -1.0_wp * b ) + ( ( ( b**2 ) - 4.0_wp * a * c )**0.5_wp & ) ) / ( 2 * a ) root2 = ( ( -1.0_wp * b ) - ( ( ( b**2 ) - 4.0_wp * a * c) **0.5_wp & ) ) / ( 2 * a ) IF ( root1 > ions(1) .OR. root1 < 0.0_wp ) THEN root1 = 0.0_wp ENDIF IF ( root2 > ions(1) .OR. root2 < 0.0_wp ) THEN root2 = 0.0_wp ENDIF ! !-- Calculate the new hydrogen ion, bisulphate ion and sulphate ion !-- concentration hso4_real = 0.0_wp h_real = ions(1) so4_real = ions(4) IF ( root1 == 0.0_wp ) THEN hso4_real = root2 ELSEIF ( root2 == 0.0_wp ) THEN hso4_real = root1 ENDIF h_real = ions(1) - hso4_real so4_real = ions(4) - hso4_real ! !-- Recalculate ion molalities ions_mol(1) = h_real / ( water_total * 18.01528E-3_wp ) ! H+ ions_mol(4) = so4_real / ( water_total * 18.01528E-3_wp ) ! SO4(2-) ions_mol(5) = hso4_real / ( water_total * 18.01528E-3_wp ) ! HSO4(2-) h_out = h_real hso4_out = hso4_real so4_out = so4_real ELSEIF ( ions(1) == 0.0_wp .OR. ions(4) == 0.0_wp ) THEN h_out = ions(1) hso4_out = 0.0_wp so4_out = ions(4) ENDIF ! !-- 4) ACTIVITY COEFFICIENTS -for vapour pressures of HNO3,HCL and NH3 ! !-- This section evaluates activity coefficients and vapour pressures using the !-- water content calculated above) for each inorganic condensing species: !-- a - HNO3, b - NH3, c - HCL. !-- The following procedure is used: !-- Zaveri et al (2005) found that one could express the variation of activity !-- coefficients linearly in log-space if equivalent mole fractions were used. !-- So, by a taylor series expansion LOG( activity coefficient ) = !-- LOG( binary activity coefficient at a given RH ) + !-- (equivalent mole fraction compound A) * !-- ('interaction' parameter between A and condensing species) + !-- equivalent mole fraction compound B) * !-- ('interaction' parameter between B and condensing species). !-- Here, the interaction parameters have been fit to ADDEM by searching the !-- whole compositon space and fit usign the Levenberg-Marquardt non-linear !-- least squares algorithm. ! !-- They are given as a function of RH and vary with complexity ranging from !-- linear to 5th order polynomial expressions, the binary activity coefficients !-- were calculated using AIM online. !-- NOTE: for NH3, no binary activity coefficient was used and the data were fit !-- to the ratio of the activity coefficients for the ammonium and hydrogen !-- ions. Once the activity coefficients are obtained the vapour pressure can be !-- easily calculated using tabulated equilibrium constants (referenced). This !-- procedure differs from that of Zaveri et al (2005) in that it is not assumed !-- one can carry behaviour from binary mixtures in multicomponent systems. To !-- this end we have fit the 'interaction' parameters explicitly to a general !-- inorganic equilibrium model (ADDEM - Topping et al. 2005a,b). Such !-- parameters take into account bisulphate ion dissociation and water content. !-- This also allows us to consider one regime for all composition space, rather !-- than defining sulphate rich and sulphate poor regimes !-- NOTE: The flags "binary_case" and "full_complexity" are not used in this !-- prototype. They are used for simplification of the fit expressions when !-- using limited composition regions. ! !-- a) - ACTIVITY COEFF/VAPOUR PRESSURE - HNO3 IF ( ions(1) > 0.0_wp .AND. ions(6) > 0.0_wp ) THEN binary_case = 1 IF ( RH > 0.1_wp .AND. RH < 0.98_wp ) THEN IF ( binary_case == 1 ) THEN binary_hno3 = 1.8514_wp * ( RH**3 ) - 4.6991_wp * ( RH**2 ) + & 1.5514_wp * RH + 0.90236_wp ELSEIF ( binary_case == 2 ) THEN binary_hno3 = - 1.1751_wp * ( RH**2 ) - 0.53794_wp * RH + & 1.2808_wp ENDIF ELSEIF ( RH >= 0.98_wp .AND. RH < 0.9999_wp ) THEN binary_hno3 = 1244.69635941351_wp * ( RH**3 ) - & 2613.93941099991_wp * ( RH**2 ) + & 1525.0684974546_wp * RH -155.946764059316_wp ENDIF ! !-- Contributions from other solutes full_complexity = 1 IF ( hydrochloric_acid > 0.0_wp ) THEN ! HCL IF ( full_complexity == 1 .OR. RH < 0.4_wp ) THEN HCL_hno3 = 16.051_wp * ( RH**4 ) - 44.357_wp * ( RH**3 ) + & 45.141_wp * ( RH**2 ) - 21.638_wp * RH + 4.8182_wp ELSEIF ( full_complexity == 0 .AND. RH > 0.4_wp ) THEN HCL_hno3 = - 1.5833_wp * RH + 1.5569_wp ENDIF ENDIF IF ( sulphuric_acid > 0.0_wp ) THEN ! H2SO4 IF ( full_complexity == 1 .OR. RH < 0.4_wp ) THEN H2SO4_hno3 = - 3.0849_wp * ( RH**3 ) + 5.9609_wp * ( RH**2 ) - & 4.468_wp * RH + 1.5658_wp ELSEIF ( full_complexity == 0 .AND. RH > 0.4_wp ) THEN H2SO4_hno3 = - 0.93473_wp * RH + 0.9363_wp ENDIF ENDIF IF ( ammonium_sulphate > 0.0_wp ) THEN ! NH42SO4 NH42SO4_hno3 = 16.821_wp * ( RH**3 ) - 28.391_wp * ( RH**2 ) + & 18.133_wp * RH - 6.7356_wp ENDIF IF ( ammonium_nitrate > 0.0_wp ) THEN ! NH4NO3 NH4NO3_hno3 = 11.01_wp * ( RH**3 ) - 21.578_wp * ( RH**2 ) + & 14.808_wp * RH - 4.2593_wp ENDIF IF ( ammonium_chloride > 0.0_wp ) THEN ! NH4Cl IF ( full_complexity == 1 .OR. RH <= 0.4_wp ) THEN NH4Cl_hno3 = - 1.176_wp * ( RH**3 ) + 5.0828_wp * ( RH**2 ) - & 3.8792_wp * RH - 0.05518_wp ELSEIF ( full_complexity == 0 .AND. RH > 0.4_wp ) THEN NH4Cl_hno3 = 2.6219_wp * ( RH**2 ) - 2.2609_wp * RH - 0.38436_wp ENDIF ENDIF IF ( sodium_sulphate > 0.0_wp ) THEN ! Na2SO4 Na2SO4_hno3 = 35.504_wp * ( RH**4 ) - 80.101_wp * ( RH**3 ) + & 67.326_wp * ( RH**2 ) - 28.461_wp * RH + 5.6016_wp ENDIF IF ( sodium_nitrate > 0.0_wp ) THEN ! NaNO3 IF ( full_complexity == 1 .OR. RH <= 0.4_wp ) THEN NaNO3_hno3 = 23.659_wp * ( RH**5 ) - 66.917_wp * ( RH**4 ) + & 74.686_wp * ( RH**3 ) - 40.795_wp * ( RH**2 ) + & 10.831_wp * RH - 1.4701_wp ELSEIF ( full_complexity == 0 .AND. RH > 0.4_wp ) THEN NaNO3_hno3 = 14.749_wp * ( RH**4 ) - 35.237_wp * ( RH**3 ) + & 31.196_wp * ( RH**2 ) - 12.076_wp * RH + 1.3605_wp ENDIF ENDIF IF ( sodium_chloride > 0.0_wp ) THEN ! NaCl IF ( full_complexity == 1 .OR. RH <= 0.4_wp ) THEN NaCl_hno3 = 13.682_wp * ( RH**4 ) - 35.122_wp * ( RH**3 ) + & 33.397_wp * ( RH**2 ) - 14.586_wp * RH + 2.6276_wp ELSEIF ( full_complexity == 0 .AND. RH > 0.4_wp ) THEN NaCl_hno3 = 1.1882_wp * ( RH**3 ) - 1.1037_wp * ( RH**2 ) - & 0.7642_wp * RH + 0.6671_wp ENDIF ENDIF Ln_HNO3_act = binary_hno3 + & hydrochloric_acid_eq_frac * HCL_hno3 + & sulphuric_acid_eq_frac * H2SO4_hno3 + & ammonium_sulphate_eq_frac * NH42SO4_hno3 + & ammonium_nitrate_eq_frac * NH4NO3_hno3 + & ammonium_chloride_eq_frac * NH4Cl_hno3 + & sodium_sulphate_eq_frac * Na2SO4_hno3 + & sodium_nitrate_eq_frac * NaNO3_hno3 + & sodium_chloride_eq_frac * NaCl_hno3 gamma_hno3 = EXP( Ln_HNO3_act ) ! Molal activity coefficient of HNO3 gamma_out(1) = gamma_hno3 ! !-- Partial pressure calculation !-- K_hno3 = 2.51 * ( 10**6 ) !-- K_hno3 = 2.628145923d6 !< calculated by AIM online (Clegg et al 1998) !-- after Chameides (1984) (and NIST database) K_hno3 = 2.6E6_wp * EXP( 8700.0_wp * henrys_temp_dep) Press_HNO3 = ( ions_mol(1) * ions_mol(6) * ( gamma_hno3**2 ) ) / & K_hno3 ENDIF ! !-- b) - ACTIVITY COEFF/VAPOUR PRESSURE - NH3 !-- Follow the two solute approach of Zaveri et al. (2005) IF ( ions(2) > 0.0_wp .AND. ions_mol(1) > 0.0_wp ) THEN !-- NH4HSO4: binary_nh4hso4 = 56.907_wp * ( RH**6 ) - 155.32_wp * ( RH**5 ) + & 142.94_wp * ( RH**4 ) - 32.298_wp * ( RH**3 ) - & 27.936_wp * ( RH**2 ) + 19.502_wp * RH - 4.2618_wp IF ( nitric_acid > 0.0_wp) THEN ! HNO3 HNO3_nh4hso4 = 104.8369_wp * ( RH**8 ) - 288.8923_wp * ( RH**7 ) + & 129.3445_wp * ( RH**6 ) + 373.0471_wp * ( RH**5 ) - & 571.0385_wp * ( RH**4 ) + 326.3528_wp * ( RH**3 ) - & 74.169_wp * ( RH**2 ) - 2.4999_wp * RH + 3.17_wp ENDIF IF ( hydrochloric_acid > 0.0_wp) THEN ! HCL HCL_nh4hso4 = - 7.9133_wp * ( RH**8 ) + 126.6648_wp * ( RH**7 ) - & 460.7425_wp * ( RH**6 ) + 731.606_wp * ( RH**5 ) - & 582.7467_wp * ( RH**4 ) + 216.7197_wp * ( RH**3 ) - & 11.3934_wp * ( RH**2 ) - 17.7728_wp * RH + 5.75_wp ENDIF IF ( sulphuric_acid > 0.0_wp) THEN ! H2SO4 H2SO4_nh4hso4 = 195.981_wp * ( RH**8 ) - 779.2067_wp * ( RH**7 ) + & 1226.3647_wp * ( RH**6 ) - 964.0261_wp * ( RH**5 ) + & 391.7911_wp * ( RH**4 ) - 84.1409_wp * ( RH**3 ) + & 20.0602_wp * ( RH**2 ) - 10.2663_wp * RH + 3.5817_wp ENDIF IF ( ammonium_sulphate > 0.0_wp) THEN ! NH42SO4 NH42SO4_nh4hso4 = 617.777_wp * ( RH**8 ) - 2547.427_wp * ( RH**7 ) & + 4361.6009_wp * ( RH**6 ) - 4003.162_wp * ( RH**5 ) & + 2117.8281_wp * ( RH**4 ) - 640.0678_wp * ( RH**3 ) & + 98.0902_wp * ( RH**2 ) - 2.2615_wp * RH - 2.3811_wp ENDIF IF ( ammonium_nitrate > 0.0_wp) THEN ! NH4NO3 NH4NO3_nh4hso4 = - 104.4504_wp * ( RH**8 ) + 539.5921_wp * & ( RH**7 ) - 1157.0498_wp * ( RH**6 ) + 1322.4507_wp * & ( RH**5 ) - 852.2475_wp * ( RH**4 ) + 298.3734_wp * & ( RH**3 ) - 47.0309_wp * ( RH**2 ) + 1.297_wp * RH - & 0.8029_wp ENDIF IF ( ammonium_chloride > 0.0_wp) THEN ! NH4Cl NH4Cl_nh4hso4 = 258.1792_wp * ( RH**8 ) - 1019.3777_wp * & ( RH**7 ) + 1592.8918_wp * ( RH**6 ) - 1221.0726_wp * & ( RH**5 ) + 442.2548_wp * ( RH**4 ) - 43.6278_wp * & ( RH**3 ) - 7.5282_wp * ( RH**2 ) - 3.8459_wp * RH + 2.2728_wp ENDIF IF ( sodium_sulphate > 0.0_wp) THEN ! Na2SO4 Na2SO4_nh4hso4 = 225.4238_wp * ( RH**8 ) - 732.4113_wp * & ( RH**7 ) + 843.7291_wp * ( RH**6 ) - 322.7328_wp * & ( RH**5 ) - 88.6252_wp * ( RH**4 ) + 72.4434_wp * & ( RH**3 ) + 22.9252_wp * ( RH**2 ) - 25.3954_wp * RH + & 4.6971_wp ENDIF IF ( sodium_nitrate > 0.0_wp) THEN ! NaNO3 NaNO3_nh4hso4 = 96.1348_wp * ( RH**8 ) - 341.6738_wp * ( RH**7 ) + & 406.5314_wp * ( RH**6 ) - 98.5777_wp * ( RH**5 ) - & 172.8286_wp * ( RH**4 ) + 149.3151_wp * ( RH**3 ) - & 38.9998_wp * ( RH**2 ) - 0.2251 * RH + 0.4953_wp ENDIF IF ( sodium_chloride > 0.0_wp) THEN ! NaCl NaCl_nh4hso4 = 91.7856_wp * ( RH**8 ) - 316.6773_wp * ( RH**7 ) + & 358.2703_wp * ( RH**6 ) - 68.9142 * ( RH**5 ) - & 156.5031_wp * ( RH**4 ) + 116.9592_wp * ( RH**3 ) - & 22.5271_wp * ( RH**2 ) - 3.7716_wp * RH + 1.56_wp ENDIF Ln_NH4HSO4_act = binary_nh4hso4 + & nitric_acid_eq_frac * HNO3_nh4hso4 + & hydrochloric_acid_eq_frac * HCL_nh4hso4 + & sulphuric_acid_eq_frac * H2SO4_nh4hso4 + & ammonium_sulphate_eq_frac * NH42SO4_nh4hso4 + & ammonium_nitrate_eq_frac * NH4NO3_nh4hso4 + & ammonium_chloride_eq_frac * NH4Cl_nh4hso4 + & sodium_sulphate_eq_frac * Na2SO4_nh4hso4 + & sodium_nitrate_eq_frac * NaNO3_nh4hso4 + & sodium_chloride_eq_frac * NaCl_nh4hso4 gamma_nh4hso4 = EXP( Ln_NH4HSO4_act ) ! molal act. coefficient of NH4HSO4 !-- Molal activity coefficient of NO3- gamma_out(6) = gamma_nh4hso4 !-- Molal activity coefficient of NH4+ gamma_nh3 = ( gamma_nh4hso4**2 ) / ( gamma_hhso4**2 ) gamma_out(3) = gamma_nh3 ! !-- This actually represents the ratio of the ammonium to hydrogen ion !-- activity coefficients (see Zaveri paper) - multiply this by the ratio !-- of the ammonium to hydrogen ion molality and the ratio of appropriate !-- equilibrium constants ! !-- Equilibrium constants !-- Kh = 57.64d0 ! Zaveri et al. (2005) Kh = 5.8E1_wp * EXP( 4085.0_wp * henrys_temp_dep ) ! after Chameides ! ! (1984) (and NIST database) !-- Knh4 = 1.81E-5_wp ! Zaveri et al. (2005) Knh4 = 1.7E-5_wp * EXP( -4325.0_wp * henrys_temp_dep ) ! Chameides ! (1984) !-- Kw = 1.01E-14_wp ! Zaveri et al (2005) Kw = 1.E-14_wp * EXP( -6716.0_wp * henrys_temp_dep ) ! Chameides ! (1984) ! molality_ratio_nh3 = ions_mol(2) / ions_mol(1) !-- Partial pressure calculation Press_NH3 = molality_ratio_nh3 * gamma_nh3 * ( Kw / ( Kh * Knh4 ) ) ENDIF ! !-- c) - ACTIVITY COEFF/VAPOUR PRESSURE - HCL IF ( ions(1) > 0.0_wp .AND. ions(7) > 0.0_wp ) THEN binary_case = 1 IF ( RH > 0.1_wp .AND. RH < 0.98 ) THEN IF ( binary_case == 1 ) THEN binary_hcl = - 5.0179_wp * ( RH**3 ) + 9.8816_wp * ( RH**2 ) - & 10.789_wp * RH + 5.4737_wp ELSEIF ( binary_case == 2 ) THEN binary_hcl = - 4.6221_wp * RH + 4.2633_wp ENDIF ELSEIF ( RH >= 0.98_wp .AND. RH < 0.9999_wp ) THEN binary_hcl = 775.6111008626_wp * ( RH**3 ) - 2146.01320888771_wp * & ( RH**2 ) + 1969.01979670259_wp * RH - 598.878230033926_wp ENDIF ENDIF IF ( nitric_acid > 0.0_wp ) THEN ! HNO3 IF ( full_complexity == 1 .OR. RH <= 0.4_wp ) THEN HNO3_hcl = 9.6256_wp * ( RH**4 ) - 26.507_wp * ( RH**3 ) + & 27.622_wp * ( RH**2 ) - 12.958_wp * RH + 2.2193_wp ELSEIF ( full_complexity == 0 .AND. RH > 0.4_wp ) THEN HNO3_hcl = 1.3242_wp * ( RH**2 ) - 1.8827_wp * RH + 0.55706_wp ENDIF ENDIF IF ( sulphuric_acid > 0.0_wp ) THEN ! H2SO4 IF ( full_complexity == 1 .OR. RH <= 0.4 ) THEN H2SO4_hcl = 1.4406_wp * ( RH**3 ) - 2.7132_wp * ( RH**2 ) + & 1.014_wp * RH + 0.25226_wp ELSEIF ( full_complexity == 0 .AND. RH > 0.4_wp ) THEN H2SO4_hcl = 0.30993_wp * ( RH**2 ) - 0.99171_wp * RH + 0.66913_wp ENDIF ENDIF IF ( ammonium_sulphate > 0.0_wp ) THEN ! NH42SO4 NH42SO4_hcl = 22.071_wp * ( RH**3 ) - 40.678_wp * ( RH**2 ) + & 27.893_wp * RH - 9.4338_wp ENDIF IF ( ammonium_nitrate > 0.0_wp ) THEN ! NH4NO3 NH4NO3_hcl = 19.935_wp * ( RH**3 ) - 42.335_wp * ( RH**2 ) + & 31.275_wp * RH - 8.8675_wp ENDIF IF ( ammonium_chloride > 0.0_wp ) THEN ! NH4Cl IF ( full_complexity == 1 .OR. RH <= 0.4_wp ) THEN NH4Cl_hcl = 2.8048_wp * ( RH**3 ) - 4.3182_wp * ( RH**2 ) + & 3.1971_wp * RH - 1.6824_wp ELSEIF ( full_complexity == 0 .AND. RH > 0.4_wp ) THEN NH4Cl_hcl = 1.2304_wp * ( RH**2 ) - 0.18262_wp * RH - 1.0643_wp ENDIF ENDIF IF ( sodium_sulphate > 0.0_wp ) THEN ! Na2SO4 Na2SO4_hcl = 36.104_wp * ( RH**4 ) - 78.658_wp * ( RH**3 ) + & 63.441_wp * ( RH**2 ) - 26.727_wp * RH + 5.7007_wp ENDIF IF ( sodium_nitrate > 0.0_wp ) THEN ! NaNO3 IF ( full_complexity == 1 .OR. RH <= 0.4_wp ) THEN NaNO3_hcl = 54.471_wp * ( RH**5 ) - 159.42_wp * ( RH**4 ) + & 180.25_wp * ( RH**3 ) - 98.176_wp * ( RH**2 ) + & 25.309_wp * RH - 2.4275_wp ELSEIF ( full_complexity == 0 .AND. RH > 0.4_wp ) THEN NaNO3_hcl = 21.632_wp * ( RH**4 ) - 53.088_wp * ( RH**3 ) + & 47.285_wp * ( RH**2 ) - 18.519_wp * RH + 2.6846_wp ENDIF ENDIF IF ( sodium_chloride > 0.0_wp ) THEN ! NaCl IF ( full_complexity == 1 .OR. RH <= 0.4_wp ) THEN NaCl_hcl = 5.4138_wp * ( RH**4 ) - 12.079_wp * ( RH**3 ) + & 9.627_wp * ( RH**2 ) - 3.3164_wp * RH + 0.35224_wp ELSEIF ( full_complexity == 0 .AND. RH > 0.4_wp ) THEN NaCl_hcl = 2.432_wp * ( RH**3 ) - 4.3453_wp * ( RH**2 ) + & 2.3834_wp * RH - 0.4762_wp ENDIF ENDIF Ln_HCL_act = binary_hcl + & nitric_acid_eq_frac * HNO3_hcl + & sulphuric_acid_eq_frac * H2SO4_hcl + & ammonium_sulphate_eq_frac * NH42SO4_hcl + & ammonium_nitrate_eq_frac * NH4NO3_hcl + & ammonium_chloride_eq_frac * NH4Cl_hcl + & sodium_sulphate_eq_frac * Na2SO4_hcl + & sodium_nitrate_eq_frac * NaNO3_hcl + & sodium_chloride_eq_frac * NaCl_hcl gamma_hcl = EXP( Ln_HCL_act ) ! Molal activity coefficient gamma_out(2) = gamma_hcl ! !-- Equilibrium constant after Wagman et al. (1982) (and NIST database) K_hcl = 2E6_wp * EXP( 9000.0_wp * henrys_temp_dep ) Press_HCL = ( ions_mol(1) * ions_mol(7) * ( gamma_hcl**2 ) ) / K_hcl ! !-- 5) Ion molility output mols_out = ions_mol ! !-- REFERENCES !-- Clegg et al. (1998) A Thermodynamic Model of the System !-- H+-NH4+-Na+-SO42- -NO3--Cl--H2O at 298.15 K, J. Phys. Chem., 102A, !-- 2155-2171. !-- Clegg et al. (2001) Thermodynamic modelling of aqueous aerosols containing !-- electrolytes and dissolved organic compounds. Journal of Aerosol Science !-- 2001;32(6):713-738. !-- Topping et al. (2005a) A curved multi-component aerosol hygroscopicity model !-- framework: Part 1 - Inorganic compounds. Atmospheric Chemistry and !-- Physics 2005;5:1205-1222. !-- Topping et al. (2005b) A curved multi-component aerosol hygroscopicity model !-- framework: Part 2 - Including organic compounds. Atmospheric Chemistry !-- and Physics 2005;5:1223-1242. !-- Wagman et al. (1982). The NBS tables of chemical thermodynamic properties: !-- selected values for inorganic and C₁ and C₂ organic substances in SI !-- units (book) !-- Zaveri et al. (2005). A new method for multicomponent activity coefficients !-- of electrolytes in aqueous atmospheric aerosols, JGR, 110, D02201, 2005. END SUBROUTINE inorganic_pdfite !------------------------------------------------------------------------------! ! Description: ! ------------ !> Update the particle size distribution. Put particles into corrects bins. !> !> Moving-centre method assumed, i.e. particles are allowed to grow to their !> exact size as long as they are not crossing the fixed diameter bin limits. !> If the particles in a size bin cross the lower or upper diameter limit, they !> are all moved to the adjacent diameter bin and their volume is averaged with !> the particles in the new bin, which then get a new diameter. ! !> Moving-centre method minimises numerical diffusion. !------------------------------------------------------------------------------! SUBROUTINE distr_update( paero ) IMPLICIT NONE !-- Input and output variables TYPE(t_section), INTENT(inout) :: paero(fn2b) !< Aerosols particle !< size distribution and properties !-- Local variables INTEGER(iwp) :: b !< loop index INTEGER(iwp) :: mm !< loop index INTEGER(iwp) :: counti LOGICAL :: within_bins !< logical (particle belongs to the bin?) REAL(wp) :: znfrac !< number fraction to be moved to the larger bin REAL(wp) :: zvfrac !< volume fraction to be moved to the larger bin REAL(wp) :: zVexc !< Volume in the grown bin which exceeds the bin !< upper limit REAL(wp) :: zVihi !< particle volume at the high end of the bin REAL(wp) :: zVilo !< particle volume at the low end of the bin REAL(wp) :: zvpart !< particle volume (m3) REAL(wp) :: zVrat !< volume ratio of a size bin zvpart = 0.0_wp zvfrac = 0.0_wp within_bins = .FALSE. ! !-- Check if the volume of the bin is within bin limits after update counti = 0 DO WHILE ( .NOT. within_bins ) within_bins = .TRUE. DO b = fn2b-1, in1a, -1 mm = 0 IF ( paero(b)%numc > nclim ) THEN zvpart = 0.0_wp zvfrac = 0.0_wp IF ( b == fn2a ) CYCLE ! !-- Dry volume zvpart = SUM( paero(b)%volc(1:7) ) / paero(b)%numc ! !-- Smallest bin cannot decrease IF ( paero(b)%vlolim > zvpart .AND. b == in1a ) CYCLE ! !-- Decreasing bins IF ( paero(b)%vlolim > zvpart ) THEN mm = b - 1 IF ( b == in2b ) mm = fn1a ! 2b goes to 1a paero(mm)%numc = paero(mm)%numc + paero(b)%numc paero(b)%numc = 0.0_wp paero(mm)%volc(:) = paero(mm)%volc(:) + paero(b)%volc(:) paero(b)%volc(:) = 0.0_wp CYCLE ENDIF ! !-- If size bin has not grown, cycle !-- Changed by Mona: compare to the arithmetic mean volume, as done !-- originally. Now particle volume is derived from the geometric mean !-- diameter, not arithmetic (see SUBROUTINE set_sizebins). IF ( zvpart <= api6 * ( ( aero(b)%vhilim + aero(b)%vlolim ) / & ( 2.0_wp * api6 ) ) ) CYCLE IF ( ABS( zvpart - api6 * paero(b)%dmid ** 3.0_wp ) < & 1.0E-35_wp ) CYCLE ! Mona: to avoid precision problems ! !-- Volume ratio of the size bin zVrat = paero(b)%vhilim / paero(b)%vlolim !-- Particle volume at the low end of the bin zVilo = 2.0_wp * zvpart / ( 1.0_wp + zVrat ) !-- Particle volume at the high end of the bin zVihi = zVrat * zVilo !-- Volume in the grown bin which exceeds the bin upper limit zVexc = 0.5_wp * ( zVihi + paero(b)%vhilim ) !-- Number fraction to be moved to the larger bin znfrac = MIN( 1.0_wp, ( zVihi - paero(b)%vhilim) / & ( zVihi - zVilo ) ) !-- Volume fraction to be moved to the larger bin zvfrac = MIN( 0.99_wp, znfrac * zVexc / zvpart ) IF ( zvfrac < 0.0_wp ) THEN message_string = 'Error: zvfrac < 0' CALL message( 'salsa_mod: distr_update', 'SA0050', & 1, 2, 0, 6, 0 ) ENDIF ! !-- Update bin mm = b + 1 !-- Volume (cm3/cm3) paero(mm)%volc(:) = paero(mm)%volc(:) + znfrac * paero(b)%numc * & zVexc * paero(b)%volc(:) / & SUM( paero(b)%volc(1:7) ) paero(b)%volc(:) = paero(b)%volc(:) - znfrac * paero(b)%numc * & zVexc * paero(b)%volc(:) / & SUM( paero(b)%volc(1:7) ) !-- Number concentration (#/m3) paero(mm)%numc = paero(mm)%numc + znfrac * paero(b)%numc paero(b)%numc = paero(b)%numc * ( 1.0_wp - znfrac ) ENDIF ! nclim IF ( paero(b)%numc > nclim ) THEN zvpart = SUM( paero(b)%volc(1:7) ) / paero(b)%numc within_bins = ( paero(b)%vlolim < zvpart .AND. & zvpart < paero(b)%vhilim ) ENDIF ENDDO ! - b counti = counti + 1 IF ( counti > 100 ) THEN message_string = 'Error: Aerosol bin update not converged' CALL message( 'salsa_mod: distr_update', 'SA0051', 1, 2, 0, 6, 0 ) ENDIF ENDDO ! - within bins END SUBROUTINE distr_update !------------------------------------------------------------------------------! ! Description: ! ------------ !> salsa_diagnostics: Update properties for the current timestep: !> !> Juha Tonttila, FMI, 2014 !> Tomi Raatikainen, FMI, 2016 !------------------------------------------------------------------------------! SUBROUTINE salsa_diagnostics( i, j ) USE arrays_3d, & ONLY: p, pt, zu USE basic_constants_and_equations_mod, & ONLY: g USE control_parameters, & ONLY: pt_surface, surface_pressure USE cpulog, & ONLY: cpu_log, log_point_s IMPLICIT NONE INTEGER(iwp), INTENT(in) :: i !< INTEGER(iwp), INTENT(in) :: j !< INTEGER(iwp) :: b !< INTEGER(iwp) :: c !< INTEGER(iwp) :: gt !< INTEGER(iwp) :: k !< INTEGER(iwp) :: nc !< REAL(wp), DIMENSION(nzb:nzt+1) :: flag !< flag to mask topography REAL(wp), DIMENSION(nzb:nzt+1) :: flag_zddry !< flag to mask zddry REAL(wp), DIMENSION(nzb:nzt+1) :: in_adn !< air density (kg/m3) REAL(wp), DIMENSION(nzb:nzt+1) :: in_p !< pressure REAL(wp), DIMENSION(nzb:nzt+1) :: in_t !< temperature (K) REAL(wp), DIMENSION(nzb:nzt+1) :: mcsum !< sum of mass concentration REAL(wp), DIMENSION(nzb:nzt+1) :: ppm_to_nconc !< Conversion factor !< from ppm to #/m3 REAL(wp), DIMENSION(nzb:nzt+1) :: zddry !< REAL(wp), DIMENSION(nzb:nzt+1) :: zvol !< flag_zddry = 0.0_wp in_adn = 0.0_wp in_p = 0.0_wp in_t = 0.0_wp ppm_to_nconc = 1.0_wp zddry = 0.0_wp zvol = 0.0_wp CALL cpu_log( log_point_s(94), 'salsa diagnostics ', 'start' ) ! !-- Calculate thermodynamic quantities needed in SALSA CALL salsa_thrm_ij( i, j, p_ij=in_p, temp_ij=in_t, adn_ij=in_adn ) ! !-- Calculate conversion factors for gas concentrations ppm_to_nconc = for_ppm_to_nconc * in_p / in_t ! !-- Predetermine flag to mask topography flag = MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_0(:,j,i), 0 ) ) DO b = 1, nbins ! aerosol size bins ! !-- Remove negative values aerosol_number(b)%conc(:,j,i) = MAX( nclim, & aerosol_number(b)%conc(:,j,i) ) * flag mcsum = 0.0_wp ! total mass concentration DO c = 1, ncc_tot ! !-- Remove negative concentrations aerosol_mass((c-1)*nbins+b)%conc(:,j,i) = MAX( mclim, & aerosol_mass((c-1)*nbins+b)%conc(:,j,i) ) & * flag mcsum = mcsum + aerosol_mass((c-1)*nbins+b)%conc(:,j,i) * flag ENDDO ! !-- Check that number and mass concentration match qualitatively IF ( ANY ( aerosol_number(b)%conc(:,j,i) > nclim .AND. & mcsum <= 0.0_wp ) ) & THEN DO k = nzb+1, nzt IF ( aerosol_number(b)%conc(k,j,i) > nclim .AND. & mcsum(k) <= 0.0_wp ) & THEN aerosol_number(b)%conc(k,j,i) = nclim * flag(k) DO c = 1, ncc_tot aerosol_mass((c-1)*nbins+b)%conc(k,j,i) = mclim * flag(k) ENDDO ENDIF ENDDO ENDIF ! !-- Update aerosol particle radius CALL bin_mixrat( 'dry', b, i, j, zvol ) zvol = zvol / arhoh2so4 ! Why on sulphate? ! !-- Particles smaller then 0.1 nm diameter are set to zero zddry = ( zvol / MAX( nclim, aerosol_number(b)%conc(:,j,i) ) / api6 )** & ( 1.0_wp / 3.0_wp ) flag_zddry = MERGE( 1.0_wp, 0.0_wp, ( zddry < 1.0E-10_wp .AND. & aerosol_number(b)%conc(:,j,i) > nclim ) ) ! !-- Volatile species to the gas phase IF ( is_used( prtcl, 'SO4' ) .AND. lscndgas ) THEN nc = get_index( prtcl, 'SO4' ) c = ( nc - 1 ) * nbins + b IF ( salsa_gases_from_chem ) THEN chem_species( gas_index_chem(1) )%conc(:,j,i) = & chem_species( gas_index_chem(1) )%conc(:,j,i) + & aerosol_mass(c)%conc(:,j,i) * avo * flag * & flag_zddry / ( amh2so4 * ppm_to_nconc ) ELSE salsa_gas(1)%conc(:,j,i) = salsa_gas(1)%conc(:,j,i) + & aerosol_mass(c)%conc(:,j,i) / amh2so4 *& avo * flag * flag_zddry ENDIF ENDIF IF ( is_used( prtcl, 'OC' ) .AND. lscndgas ) THEN nc = get_index( prtcl, 'OC' ) c = ( nc - 1 ) * nbins + b IF ( salsa_gases_from_chem ) THEN chem_species( gas_index_chem(5) )%conc(:,j,i) = & chem_species( gas_index_chem(5) )%conc(:,j,i) + & aerosol_mass(c)%conc(:,j,i) * avo * flag * & flag_zddry / ( amoc * ppm_to_nconc ) ELSE salsa_gas(5)%conc(:,j,i) = salsa_gas(5)%conc(:,j,i) + & aerosol_mass(c)%conc(:,j,i) / amoc * & avo * flag * flag_zddry ENDIF ENDIF IF ( is_used( prtcl, 'NO' ) .AND. lscndgas ) THEN nc = get_index( prtcl, 'NO' ) c = ( nc - 1 ) * nbins + b IF ( salsa_gases_from_chem ) THEN chem_species( gas_index_chem(2) )%conc(:,j,i) = & chem_species( gas_index_chem(2) )%conc(:,j,i) + & aerosol_mass(c)%conc(:,j,i) * avo * flag * & flag_zddry / ( amhno3 * ppm_to_nconc ) ELSE salsa_gas(2)%conc(:,j,i) = salsa_gas(2)%conc(:,j,i) + & aerosol_mass(c)%conc(:,j,i) / amhno3 * & avo * flag * flag_zddry ENDIF ENDIF IF ( is_used( prtcl, 'NH' ) .AND. lscndgas ) THEN nc = get_index( prtcl, 'NH' ) c = ( nc - 1 ) * nbins + b IF ( salsa_gases_from_chem ) THEN chem_species( gas_index_chem(3) )%conc(:,j,i) = & chem_species( gas_index_chem(3) )%conc(:,j,i) + & aerosol_mass(c)%conc(:,j,i) * avo * flag * & flag_zddry / ( amnh3 * ppm_to_nconc ) ELSE salsa_gas(3)%conc(:,j,i) = salsa_gas(3)%conc(:,j,i) + & aerosol_mass(c)%conc(:,j,i) / amnh3 * & avo * flag * flag_zddry ENDIF ENDIF ! !-- Mass and number to zero (insoluble species and water are lost) DO c = 1, ncc_tot aerosol_mass((c-1)*nbins+b)%conc(:,j,i) = MERGE( mclim * flag, & aerosol_mass((c-1)*nbins+b)%conc(:,j,i), & flag_zddry > 0.0_wp ) ENDDO aerosol_number(b)%conc(:,j,i) = MERGE( nclim * flag, & aerosol_number(b)%conc(:,j,i), & flag_zddry > 0.0_wp ) Ra_dry(:,j,i,b) = MAX( 1.0E-10_wp, 0.5_wp * zddry ) ENDDO IF ( .NOT. salsa_gases_from_chem ) THEN DO gt = 1, ngast salsa_gas(gt)%conc(:,j,i) = MAX( nclim, salsa_gas(gt)%conc(:,j,i) ) & * flag ENDDO ENDIF CALL cpu_log( log_point_s(94), 'salsa diagnostics ', 'stop' ) END SUBROUTINE salsa_diagnostics ! !------------------------------------------------------------------------------! ! Description: ! ------------ !> Calculate the tendencies for aerosol number and mass concentrations. !> Cache-optimized. !------------------------------------------------------------------------------! SUBROUTINE salsa_tendency_ij( id, rs_p, rs, trs_m, i, j, i_omp_start, tn, b, & c, flux_s, diss_s, flux_l, diss_l, rs_init ) USE advec_ws, & ONLY: advec_s_ws USE advec_s_pw_mod, & ONLY: advec_s_pw USE advec_s_up_mod, & ONLY: advec_s_up USE arrays_3d, & ONLY: ddzu, hyp, pt, rdf_sc, tend USE diffusion_s_mod, & ONLY: diffusion_s USE indices, & ONLY: wall_flags_0 USE pegrid, & ONLY: threads_per_task, myid USE surface_mod, & ONLY : surf_def_h, surf_def_v, surf_lsm_h, surf_lsm_v, surf_usm_h, & surf_usm_v IMPLICIT NONE CHARACTER (LEN = *) :: id INTEGER(iwp) :: b !< bin index in derived type aerosol_size_bin INTEGER(iwp) :: c !< bin index in derived type aerosol_size_bin INTEGER(iwp) :: i !< INTEGER(iwp) :: i_omp_start !< INTEGER(iwp) :: j !< INTEGER(iwp) :: k !< INTEGER(iwp) :: nc !< (c-1)*nbins+b INTEGER(iwp) :: tn !< REAL(wp), DIMENSION(nzb+1:nzt,nys:nyn,0:threads_per_task-1) :: diss_l !< REAL(wp), DIMENSION(nzb+1:nzt,0:threads_per_task-1) :: diss_s !< REAL(wp), DIMENSION(nzb+1:nzt,nys:nyn,0:threads_per_task-1) :: flux_l !< REAL(wp), DIMENSION(nzb+1:nzt,0:threads_per_task-1) :: flux_s !< REAL(wp), DIMENSION(nzb:nzt+1) :: rs_init !< REAL(wp), DIMENSION(:,:,:), POINTER :: rs_p !< REAL(wp), DIMENSION(:,:,:), POINTER :: rs !< REAL(wp), DIMENSION(:,:,:), POINTER :: trs_m !< nc = (c-1)*nbins+b ! !-- Tendency-terms for reactive scalar tend(:,j,i) = 0.0_wp IF ( id == 'aerosol_number' .AND. lod_aero == 3 ) THEN tend(:,j,i) = tend(:,j,i) + aerosol_number(b)%source(:,j,i) ELSEIF ( id == 'aerosol_mass' .AND. lod_aero == 3 ) THEN tend(:,j,i) = tend(:,j,i) + aerosol_mass(nc)%source(:,j,i) ENDIF ! !-- Advection terms IF ( timestep_scheme(1:5) == 'runge' ) THEN IF ( ws_scheme_sca ) THEN CALL advec_s_ws( i, j, rs, id, flux_s, diss_s, flux_l, diss_l, & i_omp_start, tn ) ELSE CALL advec_s_pw( i, j, rs ) ENDIF ELSE CALL advec_s_up( i, j, rs ) ENDIF ! !-- Diffusion terms IF ( id == 'aerosol_number' ) THEN CALL diffusion_s( i, j, rs, surf_def_h(0)%answs(:,b), & surf_def_h(1)%answs(:,b), surf_def_h(2)%answs(:,b), & surf_lsm_h%answs(:,b), surf_usm_h%answs(:,b), & surf_def_v(0)%answs(:,b), surf_def_v(1)%answs(:,b), & surf_def_v(2)%answs(:,b), surf_def_v(3)%answs(:,b), & surf_lsm_v(0)%answs(:,b), surf_lsm_v(1)%answs(:,b), & surf_lsm_v(2)%answs(:,b), surf_lsm_v(3)%answs(:,b), & surf_usm_v(0)%answs(:,b), surf_usm_v(1)%answs(:,b), & surf_usm_v(2)%answs(:,b), surf_usm_v(3)%answs(:,b) ) ! !-- Sedimentation for aerosol number and mass IF ( lsdepo ) THEN tend(nzb+1:nzt,j,i) = tend(nzb+1:nzt,j,i) - MAX( 0.0_wp, & ( rs(nzb+2:nzt+1,j,i) * sedim_vd(nzb+2:nzt+1,j,i,b) - & rs(nzb+1:nzt,j,i) * sedim_vd(nzb+1:nzt,j,i,b) ) * & ddzu(nzb+1:nzt) ) * MERGE( 1.0_wp, 0.0_wp, & BTEST( wall_flags_0(nzb:nzt-1,j,i), 0 ) ) ENDIF ELSEIF ( id == 'aerosol_mass' ) THEN CALL diffusion_s( i, j, rs, surf_def_h(0)%amsws(:,nc), & surf_def_h(1)%amsws(:,nc), surf_def_h(2)%amsws(:,nc), & surf_lsm_h%amsws(:,nc), surf_usm_h%amsws(:,nc), & surf_def_v(0)%amsws(:,nc), surf_def_v(1)%amsws(:,nc), & surf_def_v(2)%amsws(:,nc), surf_def_v(3)%amsws(:,nc), & surf_lsm_v(0)%amsws(:,nc), surf_lsm_v(1)%amsws(:,nc), & surf_lsm_v(2)%amsws(:,nc), surf_lsm_v(3)%amsws(:,nc), & surf_usm_v(0)%amsws(:,nc), surf_usm_v(1)%amsws(:,nc), & surf_usm_v(2)%amsws(:,nc), surf_usm_v(3)%amsws(:,nc) ) ! !-- Sedimentation for aerosol number and mass IF ( lsdepo ) THEN tend(nzb+1:nzt,j,i) = tend(nzb+1:nzt,j,i) - MAX( 0.0_wp, & ( rs(nzb+2:nzt+1,j,i) * sedim_vd(nzb+2:nzt+1,j,i,b) - & rs(nzb+1:nzt,j,i) * sedim_vd(nzb+1:nzt,j,i,b) ) * & ddzu(nzb+1:nzt) ) * MERGE( 1.0_wp, 0.0_wp, & BTEST( wall_flags_0(nzb:nzt-1,j,i), 0 ) ) ENDIF ELSEIF ( id == 'salsa_gas' ) THEN CALL diffusion_s( i, j, rs, surf_def_h(0)%gtsws(:,b), & surf_def_h(1)%gtsws(:,b), surf_def_h(2)%gtsws(:,b), & surf_lsm_h%gtsws(:,b), surf_usm_h%gtsws(:,b), & surf_def_v(0)%gtsws(:,b), surf_def_v(1)%gtsws(:,b), & surf_def_v(2)%gtsws(:,b), surf_def_v(3)%gtsws(:,b), & surf_lsm_v(0)%gtsws(:,b), surf_lsm_v(1)%gtsws(:,b), & surf_lsm_v(2)%gtsws(:,b), surf_lsm_v(3)%gtsws(:,b), & surf_usm_v(0)%gtsws(:,b), surf_usm_v(1)%gtsws(:,b), & surf_usm_v(2)%gtsws(:,b), surf_usm_v(3)%gtsws(:,b) ) ENDIF ! !-- Prognostic equation for a scalar DO k = nzb+1, nzt rs_p(k,j,i) = rs(k,j,i) + ( dt_3d * ( tsc(2) * tend(k,j,i) + & tsc(3) * trs_m(k,j,i) ) & - tsc(5) * rdf_sc(k) & * ( rs(k,j,i) - rs_init(k) ) ) & * MERGE( 1.0_wp, 0.0_wp, & BTEST( wall_flags_0(k,j,i), 0 ) ) IF ( rs_p(k,j,i) < 0.0_wp ) rs_p(k,j,i) = 0.1_wp * rs(k,j,i) ENDDO ! !-- Calculate tendencies for the next Runge-Kutta step IF ( timestep_scheme(1:5) == 'runge' ) THEN IF ( intermediate_timestep_count == 1 ) THEN DO k = nzb+1, nzt trs_m(k,j,i) = tend(k,j,i) ENDDO ELSEIF ( intermediate_timestep_count < & intermediate_timestep_count_max ) THEN DO k = nzb+1, nzt trs_m(k,j,i) = -9.5625_wp * tend(k,j,i) + 5.3125_wp * trs_m(k,j,i) ENDDO ENDIF ENDIF END SUBROUTINE salsa_tendency_ij ! !------------------------------------------------------------------------------! ! Description: ! ------------ !> Calculate the tendencies for aerosol number and mass concentrations. !> Vector-optimized. !------------------------------------------------------------------------------! SUBROUTINE salsa_tendency( id, rs_p, rs, trs_m, b, c, rs_init ) USE advec_ws, & ONLY: advec_s_ws USE advec_s_pw_mod, & ONLY: advec_s_pw USE advec_s_up_mod, & ONLY: advec_s_up USE arrays_3d, & ONLY: ddzu, hyp, pt, rdf_sc, tend USE diffusion_s_mod, & ONLY: diffusion_s USE indices, & ONLY: wall_flags_0 USE surface_mod, & ONLY : surf_def_h, surf_def_v, surf_lsm_h, surf_lsm_v, surf_usm_h, & surf_usm_v IMPLICIT NONE CHARACTER (LEN = *) :: id INTEGER(iwp) :: b !< bin index in derived type aerosol_size_bin INTEGER(iwp) :: c !< bin index in derived type aerosol_size_bin INTEGER(iwp) :: i !< INTEGER(iwp) :: j !< INTEGER(iwp) :: k !< INTEGER(iwp) :: nc !< (c-1)*nbins+b REAL(wp), DIMENSION(nzb:nzt+1) :: rs_init !< REAL(wp), DIMENSION(:,:,:), POINTER :: rs_p !< REAL(wp), DIMENSION(:,:,:), POINTER :: rs !< REAL(wp), DIMENSION(:,:,:), POINTER :: trs_m !< nc = (c-1)*nbins+b ! !-- Tendency-terms for reactive scalar tend = 0.0_wp IF ( id == 'aerosol_number' .AND. lod_aero == 3 ) THEN tend = tend + aerosol_number(b)%source ELSEIF ( id == 'aerosol_mass' .AND. lod_aero == 3 ) THEN tend = tend + aerosol_mass(nc)%source ENDIF ! !-- Advection terms IF ( timestep_scheme(1:5) == 'runge' ) THEN IF ( ws_scheme_sca ) THEN CALL advec_s_ws( rs, id ) ELSE CALL advec_s_pw( rs ) ENDIF ELSE CALL advec_s_up( rs ) ENDIF ! !-- Diffusion terms IF ( id == 'aerosol_number' ) THEN CALL diffusion_s( rs, surf_def_h(0)%answs(:,b), & surf_def_h(1)%answs(:,b), surf_def_h(2)%answs(:,b), & surf_lsm_h%answs(:,b), surf_usm_h%answs(:,b), & surf_def_v(0)%answs(:,b), surf_def_v(1)%answs(:,b), & surf_def_v(2)%answs(:,b), surf_def_v(3)%answs(:,b), & surf_lsm_v(0)%answs(:,b), surf_lsm_v(1)%answs(:,b), & surf_lsm_v(2)%answs(:,b), surf_lsm_v(3)%answs(:,b), & surf_usm_v(0)%answs(:,b), surf_usm_v(1)%answs(:,b), & surf_usm_v(2)%answs(:,b), surf_usm_v(3)%answs(:,b) ) ELSEIF ( id == 'aerosol_mass' ) THEN CALL diffusion_s( rs, surf_def_h(0)%amsws(:,nc), & surf_def_h(1)%amsws(:,nc), surf_def_h(2)%amsws(:,nc), & surf_lsm_h%amsws(:,nc), surf_usm_h%amsws(:,nc), & surf_def_v(0)%amsws(:,nc), surf_def_v(1)%amsws(:,nc), & surf_def_v(2)%amsws(:,nc), surf_def_v(3)%amsws(:,nc), & surf_lsm_v(0)%amsws(:,nc), surf_lsm_v(1)%amsws(:,nc), & surf_lsm_v(2)%amsws(:,nc), surf_lsm_v(3)%amsws(:,nc), & surf_usm_v(0)%amsws(:,nc), surf_usm_v(1)%amsws(:,nc), & surf_usm_v(2)%amsws(:,nc), surf_usm_v(3)%amsws(:,nc) ) ELSEIF ( id == 'salsa_gas' ) THEN CALL diffusion_s( rs, surf_def_h(0)%gtsws(:,b), & surf_def_h(1)%gtsws(:,b), surf_def_h(2)%gtsws(:,b), & surf_lsm_h%gtsws(:,b), surf_usm_h%gtsws(:,b), & surf_def_v(0)%gtsws(:,b), surf_def_v(1)%gtsws(:,b), & surf_def_v(2)%gtsws(:,b), surf_def_v(3)%gtsws(:,b), & surf_lsm_v(0)%gtsws(:,b), surf_lsm_v(1)%gtsws(:,b), & surf_lsm_v(2)%gtsws(:,b), surf_lsm_v(3)%gtsws(:,b), & surf_usm_v(0)%gtsws(:,b), surf_usm_v(1)%gtsws(:,b), & surf_usm_v(2)%gtsws(:,b), surf_usm_v(3)%gtsws(:,b) ) ENDIF ! !-- Prognostic equation for a scalar DO i = nxl, nxr DO j = nys, nyn IF ( id == 'salsa_gas' .AND. lod_gases == 3 ) THEN tend(:,j,i) = tend(:,j,i) + salsa_gas(b)%source(:,j,i) * & for_ppm_to_nconc * hyp(:) / pt(:,j,i) * ( hyp(:) / & 100000.0_wp )**0.286_wp ! ppm to #/m3 ELSEIF ( id == 'aerosol_mass' .OR. id == 'aerosol_number') THEN ! !-- Sedimentation for aerosol number and mass IF ( lsdepo ) THEN tend(nzb+1:nzt,j,i) = tend(nzb+1:nzt,j,i) - MAX( 0.0_wp, & ( rs(nzb+2:nzt+1,j,i) * sedim_vd(nzb+2:nzt+1,j,i,b) - & rs(nzb+1:nzt,j,i) * sedim_vd(nzb+1:nzt,j,i,b) ) * & ddzu(nzb+1:nzt) ) * MERGE( 1.0_wp, 0.0_wp, & BTEST( wall_flags_0(nzb:nzt-1,j,i), 0 ) ) ENDIF ENDIF DO k = nzb+1, nzt rs_p(k,j,i) = rs(k,j,i) + ( dt_3d * ( tsc(2) * tend(k,j,i) + & tsc(3) * trs_m(k,j,i) ) & - tsc(5) * rdf_sc(k) & * ( rs(k,j,i) - rs_init(k) ) )& * MERGE( 1.0_wp, 0.0_wp, & BTEST( wall_flags_0(k,j,i), 0 ) ) IF ( rs_p(k,j,i) < 0.0_wp ) rs_p(k,j,i) = 0.1_wp * rs(k,j,i) ENDDO ENDDO ENDDO ! !-- Calculate tendencies for the next Runge-Kutta step IF ( timestep_scheme(1:5) == 'runge' ) THEN IF ( intermediate_timestep_count == 1 ) THEN DO i = nxl, nxr DO j = nys, nyn DO k = nzb+1, nzt trs_m(k,j,i) = tend(k,j,i) ENDDO ENDDO ENDDO ELSEIF ( intermediate_timestep_count < & intermediate_timestep_count_max ) THEN DO i = nxl, nxr DO j = nys, nyn DO k = nzb+1, nzt trs_m(k,j,i) = -9.5625_wp * tend(k,j,i) & + 5.3125_wp * trs_m(k,j,i) ENDDO ENDDO ENDDO ENDIF ENDIF END SUBROUTINE salsa_tendency !------------------------------------------------------------------------------! ! Description: ! ------------ !> Boundary conditions for prognostic variables in SALSA !------------------------------------------------------------------------------! SUBROUTINE salsa_boundary_conds USE surface_mod, & ONLY : bc_h IMPLICIT NONE INTEGER(iwp) :: b !< index for aerosol size bins INTEGER(iwp) :: c !< index for chemical compounds in aerosols INTEGER(iwp) :: g !< idex for gaseous compounds INTEGER(iwp) :: i !< grid index x direction INTEGER(iwp) :: j !< grid index y direction INTEGER(iwp) :: k !< grid index y direction INTEGER(iwp) :: kb !< variable to set respective boundary value, depends on !< facing. INTEGER(iwp) :: l !< running index boundary type, for up- and downward- !< facing walls INTEGER(iwp) :: m !< running index surface elements ! !-- Surface conditions: IF ( ibc_salsa_b == 0 ) THEN ! Dirichlet ! !-- Run loop over all non-natural and natural walls. Note, in wall-datatype !-- the k coordinate belongs to the atmospheric grid point, therefore, set !-- s_p at k-1 DO l = 0, 1 ! !-- Set kb, for upward-facing surfaces value at topography top (k-1) is !-- set, for downward-facing surfaces at topography bottom (k+1) kb = MERGE ( -1, 1, l == 0 ) !$OMP PARALLEL PRIVATE( b, c, g, i, j, k ) !$OMP DO DO m = 1, bc_h(l)%ns i = bc_h(l)%i(m) j = bc_h(l)%j(m) k = bc_h(l)%k(m) DO b = 1, nbins aerosol_number(b)%conc_p(k+kb,j,i) = & aerosol_number(b)%conc(k+kb,j,i) DO c = 1, ncc_tot aerosol_mass((c-1)*nbins+b)%conc_p(k+kb,j,i) = & aerosol_mass((c-1)*nbins+b)%conc(k+kb,j,i) ENDDO ENDDO IF ( .NOT. salsa_gases_from_chem ) THEN DO g = 1, ngast salsa_gas(g)%conc_p(k+kb,j,i) = salsa_gas(g)%conc(k+kb,j,i) ENDDO ENDIF ENDDO !$OMP END PARALLEL ENDDO ELSE ! Neumann DO l = 0, 1 ! !-- Set kb, for upward-facing surfaces value at topography top (k-1) is !-- set, for downward-facing surfaces at topography bottom (k+1) kb = MERGE( -1, 1, l == 0 ) !$OMP PARALLEL PRIVATE( b, c, g, i, j, k ) !$OMP DO DO m = 1, bc_h(l)%ns i = bc_h(l)%i(m) j = bc_h(l)%j(m) k = bc_h(l)%k(m) DO b = 1, nbins aerosol_number(b)%conc_p(k+kb,j,i) = & aerosol_number(b)%conc_p(k,j,i) DO c = 1, ncc_tot aerosol_mass((c-1)*nbins+b)%conc_p(k+kb,j,i) = & aerosol_mass((c-1)*nbins+b)%conc_p(k,j,i) ENDDO ENDDO IF ( .NOT. salsa_gases_from_chem ) THEN DO g = 1, ngast salsa_gas(g)%conc_p(k+kb,j,i) = salsa_gas(g)%conc_p(k,j,i) ENDDO ENDIF ENDDO !$OMP END PARALLEL ENDDO ENDIF ! !--Top boundary conditions: IF ( ibc_salsa_t == 0 ) THEN ! Dirichlet DO b = 1, nbins aerosol_number(b)%conc_p(nzt+1,:,:) = & aerosol_number(b)%conc(nzt+1,:,:) DO c = 1, ncc_tot aerosol_mass((c-1)*nbins+b)%conc_p(nzt+1,:,:) = & aerosol_mass((c-1)*nbins+b)%conc(nzt+1,:,:) ENDDO ENDDO IF ( .NOT. salsa_gases_from_chem ) THEN DO g = 1, ngast salsa_gas(g)%conc_p(nzt+1,:,:) = salsa_gas(g)%conc(nzt+1,:,:) ENDDO ENDIF ELSEIF ( ibc_salsa_t == 1 ) THEN ! Neumann DO b = 1, nbins aerosol_number(b)%conc_p(nzt+1,:,:) = & aerosol_number(b)%conc_p(nzt,:,:) DO c = 1, ncc_tot aerosol_mass((c-1)*nbins+b)%conc_p(nzt+1,:,:) = & aerosol_mass((c-1)*nbins+b)%conc_p(nzt,:,:) ENDDO ENDDO IF ( .NOT. salsa_gases_from_chem ) THEN DO g = 1, ngast salsa_gas(g)%conc_p(nzt+1,:,:) = salsa_gas(g)%conc_p(nzt,:,:) ENDDO ENDIF ENDIF ! !-- Lateral boundary conditions at the outflow IF ( bc_radiation_s ) THEN DO b = 1, nbins aerosol_number(b)%conc_p(:,nys-1,:) = aerosol_number(b)%conc_p(:,nys,:) DO c = 1, ncc_tot aerosol_mass((c-1)*nbins+b)%conc_p(:,nys-1,:) = & aerosol_mass((c-1)*nbins+b)%conc_p(:,nys,:) ENDDO ENDDO ELSEIF ( bc_radiation_n ) THEN DO b = 1, nbins aerosol_number(b)%conc_p(:,nyn+1,:) = aerosol_number(b)%conc_p(:,nyn,:) DO c = 1, ncc_tot aerosol_mass((c-1)*nbins+b)%conc_p(:,nyn+1,:) = & aerosol_mass((c-1)*nbins+b)%conc_p(:,nyn,:) ENDDO ENDDO ELSEIF ( bc_radiation_l ) THEN DO b = 1, nbins aerosol_number(b)%conc_p(:,nxl-1,:) = aerosol_number(b)%conc_p(:,nxl,:) DO c = 1, ncc_tot aerosol_mass((c-1)*nbins+b)%conc_p(:,nxl-1,:) = & aerosol_mass((c-1)*nbins+b)%conc_p(:,nxl,:) ENDDO ENDDO ELSEIF ( bc_radiation_r ) THEN DO b = 1, nbins aerosol_number(b)%conc_p(:,nxr+1,:) = aerosol_number(b)%conc_p(:,nxr,:) DO c = 1, ncc_tot aerosol_mass((c-1)*nbins+b)%conc_p(:,nxr+1,:) = & aerosol_mass((c-1)*nbins+b)%conc_p(:,nxr,:) ENDDO ENDDO ENDIF END SUBROUTINE salsa_boundary_conds !------------------------------------------------------------------------------! ! Description: ! ------------ ! Undoing of the previously done cyclic boundary conditions. !------------------------------------------------------------------------------! SUBROUTINE salsa_boundary_conds_decycle ( sq, sq_init ) IMPLICIT NONE INTEGER(iwp) :: boundary !< INTEGER(iwp) :: ee !< INTEGER(iwp) :: copied !< INTEGER(iwp) :: i !< INTEGER(iwp) :: j !< INTEGER(iwp) :: k !< INTEGER(iwp) :: ss !< REAL(wp), DIMENSION(nzb:nzt+1) :: sq_init REAL(wp), DIMENSION(nzb:nzt+1,nysg:nyng,nxlg:nxrg) :: sq REAL(wp) :: flag !< flag to mask topography grid points flag = 0.0_wp ! !-- Left and right boundaries IF ( decycle_lr .AND. ( bc_lr_cyc .OR. bc_lr == 'nested' ) ) THEN DO boundary = 1, 2 IF ( decycle_method(boundary) == 'dirichlet' ) THEN ! !-- Initial profile is copied to ghost and first three layers ss = 1 ee = 0 IF ( boundary == 1 .AND. nxl == 0 ) THEN ss = nxlg ee = nxl+2 ELSEIF ( boundary == 2 .AND. nxr == nx ) THEN ss = nxr-2 ee = nxrg ENDIF DO i = ss, ee DO j = nysg, nyng DO k = nzb+1, nzt flag = MERGE( 1.0_wp, 0.0_wp, & BTEST( wall_flags_0(k,j,i), 0 ) ) sq(k,j,i) = sq_init(k) * flag ENDDO ENDDO ENDDO ELSEIF ( decycle_method(boundary) == 'neumann' ) THEN ! !-- The value at the boundary is copied to the ghost layers to simulate !-- an outlet with zero gradient ss = 1 ee = 0 IF ( boundary == 1 .AND. nxl == 0 ) THEN ss = nxlg ee = nxl-1 copied = nxl ELSEIF ( boundary == 2 .AND. nxr == nx ) THEN ss = nxr+1 ee = nxrg copied = nxr ENDIF DO i = ss, ee DO j = nysg, nyng DO k = nzb+1, nzt flag = MERGE( 1.0_wp, 0.0_wp, & BTEST( wall_flags_0(k,j,i), 0 ) ) sq(k,j,i) = sq(k,j,copied) * flag ENDDO ENDDO ENDDO ELSE WRITE(message_string,*) & 'unknown decycling method: decycle_method (', & boundary, ') ="' // TRIM( decycle_method(boundary) ) // '"' CALL message( 'salsa_boundary_conds_decycle', 'SA0029', & 1, 2, 0, 6, 0 ) ENDIF ENDDO ENDIF ! !-- South and north boundaries IF ( decycle_ns .AND. ( bc_ns_cyc .OR. bc_ns == 'nested' ) ) THEN DO boundary = 3, 4 IF ( decycle_method(boundary) == 'dirichlet' ) THEN ! !-- Initial profile is copied to ghost and first three layers ss = 1 ee = 0 IF ( boundary == 3 .AND. nys == 0 ) THEN ss = nysg ee = nys+2 ELSEIF ( boundary == 4 .AND. nyn == ny ) THEN ss = nyn-2 ee = nyng ENDIF DO i = nxlg, nxrg DO j = ss, ee DO k = nzb+1, nzt flag = MERGE( 1.0_wp, 0.0_wp, & BTEST( wall_flags_0(k,j,i), 0 ) ) sq(k,j,i) = sq_init(k) * flag ENDDO ENDDO ENDDO ELSEIF ( decycle_method(boundary) == 'neumann' ) THEN ! !-- The value at the boundary is copied to the ghost layers to simulate !-- an outlet with zero gradient ss = 1 ee = 0 IF ( boundary == 3 .AND. nys == 0 ) THEN ss = nysg ee = nys-1 copied = nys ELSEIF ( boundary == 4 .AND. nyn == ny ) THEN ss = nyn+1 ee = nyng copied = nyn ENDIF DO i = nxlg, nxrg DO j = ss, ee DO k = nzb+1, nzt flag = MERGE( 1.0_wp, 0.0_wp, & BTEST( wall_flags_0(k,j,i), 0 ) ) sq(k,j,i) = sq(k,copied,i) * flag ENDDO ENDDO ENDDO ELSE WRITE(message_string,*) & 'unknown decycling method: decycle_method (', & boundary, ') ="' // TRIM( decycle_method(boundary) ) // '"' CALL message( 'salsa_boundary_conds_decycle', 'SA0030', & 1, 2, 0, 6, 0 ) ENDIF ENDDO ENDIF END SUBROUTINE salsa_boundary_conds_decycle !------------------------------------------------------------------------------! ! Description: ! ------------ !> Calculates the total dry or wet mass concentration for individual bins !> Juha Tonttila (FMI) 2015 !> Tomi Raatikainen (FMI) 2016 !------------------------------------------------------------------------------! SUBROUTINE bin_mixrat( itype, ibin, i, j, mconc ) IMPLICIT NONE CHARACTER(len=*), INTENT(in) :: itype !< 'dry' or 'wet' INTEGER(iwp), INTENT(in) :: ibin !< index of the chemical component INTEGER(iwp), INTENT(in) :: i !< loop index for x-direction INTEGER(iwp), INTENT(in) :: j !< loop index for y-direction REAL(wp), DIMENSION(:), INTENT(out) :: mconc !< total dry or wet mass !< concentration INTEGER(iwp) :: c !< loop index for mass bin number INTEGER(iwp) :: iend !< end index: include water or not !-- Number of components IF ( itype == 'dry' ) THEN iend = get_n_comp( prtcl ) - 1 ELSE IF ( itype == 'wet' ) THEN iend = get_n_comp( prtcl ) ELSE STOP 1 ! "INFO for Developer: please use the message routine to pass the output string" bin_mixrat: Error in itype ENDIF mconc = 0.0_wp DO c = ibin, iend*nbins+ibin, nbins !< every nbins'th element mconc = mconc + aerosol_mass(c)%conc(:,j,i) ENDDO END SUBROUTINE bin_mixrat !------------------------------------------------------------------------------! !> Description: !> ------------ !> Define aerosol fluxes: constant or read from a from file !------------------------------------------------------------------------------! SUBROUTINE salsa_set_source ! USE date_and_time_mod, & ! ONLY: index_dd, index_hh, index_mm #if defined( __netcdf ) USE NETCDF USE netcdf_data_input_mod, & ONLY: get_attribute, get_variable, & netcdf_data_input_get_dimension_length, open_read_file USE surface_mod, & ONLY: surf_def_h, surf_lsm_h, surf_usm_h IMPLICIT NONE INTEGER(iwp), PARAMETER :: ndm = 3 !< number of default modes INTEGER(iwp), PARAMETER :: ndc = 4 !< number of default categories CHARACTER (LEN=10) :: unita !< Unit of aerosol fluxes CHARACTER (LEN=10) :: unitg !< Unit of gaseous fluxes INTEGER(iwp) :: b !< loop index: aerosol number bins INTEGER(iwp) :: c !< loop index: aerosol chemical components INTEGER(iwp) :: ee !< loop index: end INTEGER(iwp), ALLOCATABLE, DIMENSION(:) :: eci !< emission category index INTEGER(iwp) :: g !< loop index: gaseous tracers INTEGER(iwp) :: i !< loop index: x-direction INTEGER(iwp) :: id_faero !< NetCDF id of aerosol source input file INTEGER(iwp) :: id_fchem !< NetCDF id of aerosol source input file INTEGER(iwp) :: id_sa !< NetCDF id of variable: source INTEGER(iwp) :: j !< loop index: y-direction INTEGER(iwp) :: k !< loop index: z-direction INTEGER(iwp) :: kg !< loop index: z-direction (gases) INTEGER(iwp) :: n_dt !< number of time steps in the emission file INTEGER(iwp) :: nc_stat !< local variable for storing the result of !< netCDF calls for error message handling INTEGER(iwp) :: nb_file !< Number of grid-points in file (bins) INTEGER(iwp) :: ncat !< Number of emission categories INTEGER(iwp) :: ng_file !< Number of grid-points in file (gases) INTEGER(iwp) :: num_vars !< number of variables in input file INTEGER(iwp) :: nz_file !< number of grid-points in file INTEGER(iwp) :: n !< loop index INTEGER(iwp) :: ni !< loop index INTEGER(iwp) :: ss !< loop index LOGICAL :: netcdf_extend = .FALSE. !< Flag indicating wether netcdf !< topography input file or not REAL(wp), ALLOCATABLE, DIMENSION(:,:,:,:) :: dum_var_4d !< variable for !< temporary data REAL(wp) :: fillval !< fill value REAL(wp) :: flag !< flag to mask topography grid points REAL(wp), DIMENSION(nbins) :: nsect_emission !< sectional emission (lod1) REAL(wp), ALLOCATABLE, DIMENSION(:,:,:) :: pm_emission !< aerosol mass !< emission (lod1) REAL(wp), DIMENSION(nbins) :: source_ijka !< aerosol source at (k,j,i) ! !-- The default size distribution and mass composition per emission category: !-- 1 = traffic, 2 = road dust, 3 = wood combustion, 4 = other !-- Mass fractions: H2SO4, OC, BC, DU, SS, HNO3, NH3 CHARACTER(LEN=15), DIMENSION(ndc) :: cat_name_table = &!< emission category (/'road traffic ','road dust ',& 'wood combustion','other '/) REAL(wp), DIMENSION(ndc) :: avg_density !< average density REAL(wp), DIMENSION(ndc) :: conversion_factor !< unit conversion factor !< for aerosol emissions REAL(wp), DIMENSION(ndm), PARAMETER :: dpg_table = & !< mean diameter (mum) (/ 13.5E-3_wp, 1.4_wp, 5.4E-2_wp/) REAL(wp), DIMENSION(ndm) :: ntot_table REAL(wp), DIMENSION(maxspec,ndc), PARAMETER :: mass_fraction_table = & RESHAPE( (/ 0.04_wp, 0.48_wp, 0.48_wp, 0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp, & 0.0_wp, 0.05_wp, 0.0_wp, 0.95_wp, 0.0_wp, 0.0_wp, 0.0_wp, & 0.0_wp, 0.5_wp, 0.5_wp, 0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp, & 0.0_wp, 0.5_wp, 0.5_wp, 0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp & /), (/maxspec,ndc/) ) REAL(wp), DIMENSION(ndm,ndc), PARAMETER :: PMfrac_table = & !< rel. mass RESHAPE( (/ 0.016_wp, 0.000_wp, 0.984_wp, & 0.000_wp, 1.000_wp, 0.000_wp, & 0.000_wp, 0.000_wp, 1.000_wp, & 1.000_wp, 0.000_wp, 1.000_wp & /), (/ndm,ndc/) ) REAL(wp), DIMENSION(ndm), PARAMETER :: sigmag_table = & !< mode std (/1.6_wp, 1.4_wp, 1.7_wp/) avg_density = 1.0_wp nb_file = 0 ng_file = 0 nsect_emission = 0.0_wp nz_file = 0 source_ijka = 0.0_wp ! !-- First gases, if needed: IF ( .NOT. salsa_gases_from_chem ) THEN ! !-- Read sources from PIDS_CHEM INQUIRE( FILE='PIDS_CHEM' // TRIM( coupling_char ), EXIST=netcdf_extend ) IF ( .NOT. netcdf_extend ) THEN message_string = 'Input file '// TRIM( 'PIDS_CHEM' ) // & TRIM( coupling_char ) // ' for SALSA missing!' CALL message( 'salsa_mod: salsa_set_source', 'SA0027', 1, 2, 0, 6, 0 ) ENDIF ! netcdf_extend CALL location_message( ' salsa_set_source: NOTE! Gaseous emissions'//& ' should be provided with following emission indices:'// & ' 1=H2SO4, 2=HNO3, 3=NH3, 4=OCNV, 5=OCSV', .TRUE. ) CALL location_message( ' salsa_set_source: No time dependency for '//& 'gaseous emissions. Use emission_values '// & 'directly.', .TRUE. ) ! !-- Open PIDS_CHEM in read-only mode CALL open_read_file( 'PIDS_CHEM' // TRIM( coupling_char ), id_fchem ) ! !-- Inquire the level of detail (lod) CALL get_attribute( id_fchem, 'lod', lod_gases, .FALSE., & "emission_values" ) IF ( lod_gases == 2 ) THEN ! !-- Index of gaseous compounds CALL netcdf_data_input_get_dimension_length( id_fchem, ng_file, & "nspecies" ) IF ( ng_file < 5 ) THEN message_string = 'Some gaseous emissions missing.' CALL message( 'salsa_mod: salsa_set_source', 'SA0041', & 1, 2, 0, 6, 0 ) ENDIF ! !-- Get number of emission categories CALL netcdf_data_input_get_dimension_length( id_fchem, ncat, "ncat" ) ! !-- Inquire the unit of gaseous fluxes CALL get_attribute( id_fchem, 'units', unitg, .FALSE., & "emission_values") ! !-- Inquire the fill value CALL get_attribute( id_fchem, '_FillValue', fillval, .FALSE., & "emission_values" ) ! !-- Read surface emission data (x,y) PE-wise ALLOCATE( dum_var_4d(ng_file,ncat,nys:nyn,nxl:nxr) ) CALL get_variable( id_fchem, 'emission_values', dum_var_4d, nxl, nxr,& nys, nyn, 0, ncat-1, 0, ng_file-1 ) DO g = 1, ngast ALLOCATE( salsa_gas(g)%source(ncat,nys:nyn,nxl:nxr) ) salsa_gas(g)%source = 0.0_wp salsa_gas(g)%source = salsa_gas(g)%source + dum_var_4d(g,:,:,:) ENDDO ! !-- Set surface fluxes of gaseous compounds on horizontal surfaces. !-- Set fluxes only for either default, land or urban surface. IF ( .NOT. land_surface .AND. .NOT. urban_surface ) THEN CALL set_gas_flux( surf_def_h(0), ncat, unitg ) ELSE CALL set_gas_flux( surf_lsm_h, ncat, unitg ) CALL set_gas_flux( surf_usm_h, ncat, unitg ) ENDIF DEALLOCATE( dum_var_4d ) DO g = 1, ngast DEALLOCATE( salsa_gas(g)%source ) ENDDO ELSE message_string = 'Input file PIDS_CHEM needs to have lod = 2 when '//& 'SALSA is applied but not the chemistry module!' CALL message( 'salsa_mod: salsa_set_source', 'SA0039', 1, 2, 0, 6, 0 ) ENDIF ENDIF ! !-- Read sources from PIDS_SALSA INQUIRE( FILE='PIDS_SALSA' // TRIM( coupling_char ), EXIST=netcdf_extend ) IF ( .NOT. netcdf_extend ) THEN message_string = 'Input file '// TRIM( 'PIDS_SALSA' ) // & TRIM( coupling_char ) // ' for SALSA missing!' CALL message( 'salsa_mod: salsa_set_source', 'SA0034', 1, 2, 0, 6, 0 ) ENDIF ! netcdf_extend ! !-- Open file in read-only mode CALL open_read_file( 'PIDS_SALSA' // TRIM( coupling_char ), id_faero ) ! !-- Get number of emission categories and their indices CALL netcdf_data_input_get_dimension_length( id_faero, ncat, "ncat" ) ! !-- Get emission category indices ALLOCATE( eci(1:ncat) ) CALL get_variable( id_faero, 'emission_category_index', eci ) ! !-- Inquire the level of detail (lod) CALL get_attribute( id_faero, 'lod', lod_aero, .FALSE., & "aerosol_emission_values" ) IF ( lod_aero < 3 .AND. ibc_salsa_b == 0 ) THEN message_string = 'lod1/2 for aerosol emissions requires '// & 'bc_salsa_b = "Neumann"' CALL message( 'salsa_mod: salsa_set_source','SA0025', 1, 2, 0, 6, 0 ) ENDIF ! !-- Inquire the fill value CALL get_attribute( id_faero, '_FillValue', fillval, .FALSE., & "aerosol_emission_values" ) ! !-- Aerosol chemical composition: ALLOCATE( emission_mass_fracs(1:ncat,1:maxspec) ) emission_mass_fracs = 0.0_wp !-- Chemical composition: 1: H2SO4 (sulphuric acid), 2: OC (organic carbon), !-- 3: BC (black carbon), 4: DU (dust), !-- 5: SS (sea salt), 6: HNO3 (nitric acid), !-- 7: NH3 (ammonia) DO n = 1, ncat IF ( lod_aero < 2 ) THEN emission_mass_fracs(n,:) = mass_fraction_table(:,n) ELSE CALL get_variable( id_faero, "emission_mass_fracs", & emission_mass_fracs(n,:) ) ENDIF ! !-- If the chemical component is not activated, set its mass fraction to 0 !-- to avoid inbalance between number and mass flux IF ( iso4 < 0 ) emission_mass_fracs(n,1) = 0.0_wp IF ( ioc < 0 ) emission_mass_fracs(n,2) = 0.0_wp IF ( ibc < 0 ) emission_mass_fracs(n,3) = 0.0_wp IF ( idu < 0 ) emission_mass_fracs(n,4) = 0.0_wp IF ( iss < 0 ) emission_mass_fracs(n,5) = 0.0_wp IF ( ino < 0 ) emission_mass_fracs(n,6) = 0.0_wp IF ( inh < 0 ) emission_mass_fracs(n,7) = 0.0_wp !-- Then normalise the mass fraction so that SUM = 1 emission_mass_fracs(n,:) = emission_mass_fracs(n,:) / & SUM( emission_mass_fracs(n,:) ) ENDDO IF ( lod_aero > 1 ) THEN ! !-- Aerosol geometric mean diameter CALL netcdf_data_input_get_dimension_length( id_faero, nb_file, 'Dmid' ) IF ( nb_file /= nbins ) THEN message_string = 'The number of size bins in aerosol input data '// & 'does not correspond to the model set-up' CALL message( 'salsa_mod: salsa_set_source','SA0040', 1, 2, 0, 6, 0 ) ENDIF ENDIF IF ( lod_aero < 3 ) THEN CALL location_message( ' salsa_set_source: No time dependency for '//& 'aerosol emissions. Use aerosol_emission_values'//& ' directly.', .TRUE. ) ! !-- Allocate source arrays DO b = 1, nbins ALLOCATE( aerosol_number(b)%source(1:ncat,nys:nyn,nxl:nxr) ) aerosol_number(b)%source = 0.0_wp ENDDO DO c = 1, ncc_tot*nbins ALLOCATE( aerosol_mass(c)%source(1:ncat,nys:nyn,nxl:nxr) ) aerosol_mass(c)%source = 0.0_wp ENDDO IF ( lod_aero == 1 ) THEN DO n = 1, ncat avg_density(n) = emission_mass_fracs(n,1) * arhoh2so4 + & emission_mass_fracs(n,2) * arhooc + & emission_mass_fracs(n,3) * arhobc + & emission_mass_fracs(n,4) * arhodu + & emission_mass_fracs(n,5) * arhoss + & emission_mass_fracs(n,6) * arhohno3 + & emission_mass_fracs(n,7) * arhonh3 ENDDO ! !-- Emission unit CALL get_attribute( id_faero, 'units', unita, .FALSE., & "aerosol_emission_values") conversion_factor = 1.0_wp IF ( unita == 'kg/m2/yr' ) THEN conversion_factor = 3.170979e-8_wp / avg_density ELSEIF ( unita == 'g/m2/yr' ) THEN conversion_factor = 3.170979e-8_wp * 1.0E-3_wp / avg_density ELSEIF ( unita == 'kg/m2/s' ) THEN conversion_factor = 1.0_wp / avg_density ELSEIF ( unita == 'g/m2/s' ) THEN conversion_factor = 1.0E-3_wp / avg_density ELSE message_string = 'unknown unit for aerosol emissions: ' & // TRIM( unita ) // ' (lod1)' CALL message( 'salsa_mod: salsa_set_source','SA0035', & 1, 2, 0, 6, 0 ) ENDIF ! !-- Read surface emission data (x,y) PE-wise ALLOCATE( pm_emission(ncat,nys:nyn,nxl:nxr) ) CALL get_variable( id_faero, 'aerosol_emission_values', pm_emission, & nxl, nxr, nys, nyn, 0, ncat-1 ) DO ni = 1, SIZE( eci ) n = eci(ni) ! !-- Calculate the number concentration of a log-normal size !-- distribution following Jacobson (2005): Eq 13.25. ntot_table = 6.0_wp * PMfrac_table(:,n) / ( pi * dpg_table**3 * & EXP( 4.5_wp * LOG( sigmag_table )**2 ) ) * 1.0E+12_wp ! !-- Sectional size distibution from a log-normal one CALL size_distribution( ntot_table, dpg_table, sigmag_table, & nsect_emission ) DO b = 1, nbins aerosol_number(b)%source(ni,:,:) = & aerosol_number(b)%source(ni,:,:) + & pm_emission(ni,:,:) * conversion_factor(n) & * nsect_emission(b) ENDDO ENDDO ELSEIF ( lod_aero == 2 ) THEN ! !-- Read surface emission data (x,y) PE-wise ALLOCATE( dum_var_4d(nb_file,ncat,nys:nyn,nxl:nxr) ) CALL get_variable( id_faero, 'aerosol_emission_values', dum_var_4d, & nxl, nxr, nys, nyn, 0, ncat-1, 0, nb_file-1 ) DO b = 1, nbins aerosol_number(b)%source = dum_var_4d(b,:,:,:) ENDDO DEALLOCATE( dum_var_4d ) ENDIF ! !-- Set surface fluxes of aerosol number and mass on horizontal surfaces. !-- Set fluxes only for either default, land or urban surface. IF ( .NOT. land_surface .AND. .NOT. urban_surface ) THEN CALL set_flux( surf_def_h(0), ncat ) ELSE CALL set_flux( surf_usm_h, ncat ) CALL set_flux( surf_lsm_h, ncat ) ENDIF ELSEIF ( lod_aero == 3 ) THEN ! !-- Inquire aerosol emission rate per bin (#/(m3s)) nc_stat = NF90_INQ_VARID( id_faero, "aerosol_emission_values", id_sa ) ! !-- Emission time step CALL netcdf_data_input_get_dimension_length( id_faero, n_dt, & 'dt_emission' ) IF ( n_dt > 1 ) THEN CALL location_message( ' salsa_set_source: hourly emission data'//& ' provided but currently the value of the '// & ' first hour is applied.', .TRUE. ) ENDIF ! !-- Allocate source arrays DO b = 1, nbins ALLOCATE( aerosol_number(b)%source(nzb:nzt+1,nys:nyn,nxl:nxr) ) aerosol_number(b)%source = 0.0_wp ENDDO DO c = 1, ncc_tot*nbins ALLOCATE( aerosol_mass(c)%source(nzb:nzt+1,nys:nyn,nxl:nxr) ) aerosol_mass(c)%source = 0.0_wp ENDDO ! !-- Get dimension of z-axis: CALL netcdf_data_input_get_dimension_length( id_faero, nz_file, 'z' ) ! !-- Read surface emission data (x,y) PE-wise DO i = nxl, nxr DO j = nys, nyn DO k = 0, nz_file-1 ! !-- Predetermine flag to mask topography flag = MERGE( 1.0_wp, 0.0_wp, BTEST(wall_flags_0(k,j,i), 0 )) ! !-- No sources inside buildings ! IF ( flag == 0.0_wp ) CYCLE ! !-- Read volume source: nc_stat = NF90_GET_VAR( id_faero, id_sa, source_ijka, & start = (/ i+1, j+1, k+1, 1, 1 /), & count = (/ 1, 1, 1, 1, nb_file /) ) IF ( nc_stat /= NF90_NOERR ) THEN message_string = 'error in aerosol emissions: lod3' CALL message( 'salsa_mod: salsa_set_source','SA0038', 1, 2, & 0, 6, 0 ) ENDIF ! !-- Set mass fluxes. First bins include only SO4 and/or OC. Call !-- subroutine set_mass_source for larger bins. ! !-- Sulphate and organic carbon IF ( iso4 > 0 .AND. ioc > 0 ) THEN !-- First sulphate: ss = ( iso4 - 1 ) * nbins + in1a ! start ee = ( iso4 - 1 ) * nbins + fn1a ! end b = in1a DO c = ss, ee IF ( source_ijka(b) /= fillval ) & aerosol_mass(c)%source(k,j,i) = & aerosol_mass(c)%source(k,j,i) + & emission_mass_fracs(1,1) / ( emission_mass_fracs(1,1) & + emission_mass_fracs(1,2) ) * source_ijka(b) * & aero(b)%core * arhoh2so4 b = b+1 ENDDO !-- Then organic carbon: ss = ( ioc - 1 ) * nbins + in1a ! start ee = ( ioc - 1 ) * nbins + fn1a ! end b = in1a DO c = ss, ee IF ( source_ijka(b) /= fillval ) & aerosol_mass(c)%source(k,j,i) = & aerosol_mass(c)%source(k,j,i) + & emission_mass_fracs(1,2) / ( emission_mass_fracs(1,1) & + emission_mass_fracs(1,2) ) * source_ijka(b) * & aero(b)%core * arhooc b = b+1 ENDDO CALL set_mass_source( k, j, i, iso4, & emission_mass_fracs(1,1), arhoh2so4, & source_ijka, fillval ) CALL set_mass_source( k, j, i, ioc, emission_mass_fracs(1,2),& arhooc, source_ijka, fillval ) !-- Only sulphate: ELSEIF ( iso4 > 0 .AND. ioc < 0 ) THEN ss = ( iso4 - 1 ) * nbins + in1a ! start ee = ( iso4 - 1 ) * nbins + fn1a ! end b = in1a DO c = ss, ee IF ( source_ijka(b) /= fillval ) & aerosol_mass(c)%source(k,j,i) = & aerosol_mass(c)%source(k,j,i) + source_ijka(b) * & aero(b)%core * arhoh2so4 b = b+1 ENDDO CALL set_mass_source( k, j, i, iso4, & emission_mass_fracs(1,1), arhoh2so4, & source_ijka, fillval ) !-- Only organic carbon: ELSEIF ( iso4 < 0 .AND. ioc > 0 ) THEN ss = ( ioc - 1 ) * nbins + in1a ! start ee = ( ioc - 1 ) * nbins + fn1a ! end b = in1a DO c = ss, ee IF ( source_ijka(b) /= fillval ) & aerosol_mass(c)%source(k,j,i) = & aerosol_mass(c)%source(k,j,i) + source_ijka(b) * & aero(b)%core * arhooc b = b+1 ENDDO CALL set_mass_source( k, j, i, ioc, emission_mass_fracs(1,2),& arhooc, source_ijka, fillval ) ENDIF !-- Black carbon IF ( ibc > 0 ) THEN CALL set_mass_source( k, j, i, ibc, emission_mass_fracs(1,3),& arhobc, source_ijka, fillval ) ENDIF !-- Dust IF ( idu > 0 ) THEN CALL set_mass_source( k, j, i, idu, emission_mass_fracs(1,4),& arhodu, source_ijka, fillval ) ENDIF !-- Sea salt IF ( iss > 0 ) THEN CALL set_mass_source( k, j, i, iss, emission_mass_fracs(1,5),& arhoss, source_ijka, fillval ) ENDIF !-- Nitric acid IF ( ino > 0 ) THEN CALL set_mass_source( k, j, i, ino, emission_mass_fracs(1,6),& arhohno3, source_ijka, fillval ) ENDIF !-- Ammonia IF ( inh > 0 ) THEN CALL set_mass_source( k, j, i, inh, emission_mass_fracs(1,7),& arhonh3, source_ijka, fillval ) ENDIF ! !-- Save aerosol number sources in the end DO b = 1, nbins IF ( source_ijka(b) /= fillval ) & aerosol_number(b)%source(k,j,i) = & aerosol_number(b)%source(k,j,i) + source_ijka(b) ENDDO ENDDO ! k ENDDO ! j ENDDO ! i ELSE message_string = 'NetCDF attribute lod is not set properly.' CALL message( 'salsa_mod: salsa_set_source','SA0026', 1, 2, 0, 6, 0 ) ENDIF #endif END SUBROUTINE salsa_set_source !------------------------------------------------------------------------------! ! Description: ! ------------ !> Sets the gaseous fluxes !------------------------------------------------------------------------------! SUBROUTINE set_gas_flux( surface, ncat_emission, unit ) USE arrays_3d, & ONLY: dzw, hyp, pt, rho_air_zw USE grid_variables, & ONLY: dx, dy USE surface_mod, & ONLY: surf_type IMPLICIT NONE CHARACTER(LEN=*) :: unit !< flux unit in the input file INTEGER(iwp) :: ncat_emission !< number of emission categories TYPE(surf_type), INTENT(inout) :: surface !< respective surface type INTEGER(iwp) :: g !< loop index INTEGER(iwp) :: i !< loop index INTEGER(iwp) :: j !< loop index INTEGER(iwp) :: k !< loop index INTEGER(iwp) :: m !< running index for surface elements INTEGER(iwp) :: n !< running index for emission categories REAL(wp), DIMENSION(ngast) :: conversion_factor conversion_factor = 1.0_wp DO m = 1, surface%ns ! !-- Get indices of respective grid point i = surface%i(m) j = surface%j(m) k = surface%k(m) IF ( unit == '#/m2/s' ) THEN conversion_factor = 1.0_wp ELSEIF ( unit == 'g/m2/s' ) THEN conversion_factor(1) = avo / ( amh2so4 * 1000.0_wp ) conversion_factor(2) = avo / ( amhno3 * 1000.0_wp ) conversion_factor(3) = avo / ( amnh3 * 1000.0_wp ) conversion_factor(4) = avo / ( amoc * 1000.0_wp ) conversion_factor(5) = avo / ( amoc * 1000.0_wp ) ELSEIF ( unit == 'ppm/m2/s' ) THEN conversion_factor = for_ppm_to_nconc * hyp(k) / pt(k,j,i) * ( hyp(k) & / 100000.0_wp )**0.286_wp * dx * dy * dzw(k) ELSEIF ( unit == 'mumol/m2/s' ) THEN conversion_factor = 1.0E-6_wp * avo ELSE message_string = 'Unknown unit for gaseous emissions!' CALL message( 'salsa_mod: set_gas_flux', 'SA0031', 1, 2, 0, 6, 0 ) ENDIF DO n = 1, ncat_emission DO g = 1, ngast IF ( salsa_gas(g)%source(n,j,i) < 0.0_wp ) THEN salsa_gas(g)%source(n,j,i) = 0.0_wp CYCLE ENDIF surface%gtsws(m,g) = surface%gtsws(m,g) + & salsa_gas(g)%source(n,j,i) * rho_air_zw(k-1) & * conversion_factor(g) ENDDO ENDDO ENDDO END SUBROUTINE set_gas_flux !------------------------------------------------------------------------------! ! Description: ! ------------ !> Sets the aerosol flux to aerosol arrays in 2a and 2b. !------------------------------------------------------------------------------! SUBROUTINE set_flux( surface, ncat_emission ) USE arrays_3d, & ONLY: hyp, pt, rho_air_zw USE surface_mod, & ONLY: surf_type IMPLICIT NONE INTEGER(iwp) :: ncat_emission !< number of emission categories TYPE(surf_type), INTENT(inout) :: surface !< respective surface type INTEGER(iwp) :: b !< loop index INTEGER(iwp) :: ee !< loop index INTEGER(iwp) :: g !< loop index INTEGER(iwp) :: i !< loop index INTEGER(iwp) :: j !< loop index INTEGER(iwp) :: k !< loop index INTEGER(iwp) :: m !< running index for surface elements INTEGER(iwp) :: n !< loop index for emission categories INTEGER(iwp) :: c !< loop index INTEGER(iwp) :: ss !< loop index DO m = 1, surface%ns ! !-- Get indices of respective grid point i = surface%i(m) j = surface%j(m) k = surface%k(m) DO n = 1, ncat_emission DO b = 1, nbins IF ( aerosol_number(b)%source(n,j,i) < 0.0_wp ) THEN aerosol_number(b)%source(n,j,i) = 0.0_wp CYCLE ENDIF ! !-- Set mass fluxes. First bins include only SO4 and/or OC. IF ( b <= fn1a ) THEN ! !-- Both sulphate and organic carbon IF ( iso4 > 0 .AND. ioc > 0 ) THEN c = ( iso4 - 1 ) * nbins + b surface%amsws(m,c) = surface%amsws(m,c) + & emission_mass_fracs(n,1) / & ( emission_mass_fracs(n,1) + & emission_mass_fracs(n,2) ) * & aerosol_number(b)%source(n,j,i) * & api6 * aero(b)%dmid**3.0_wp * & arhoh2so4 * rho_air_zw(k-1) aerosol_mass(c)%source(n,j,i) = & aerosol_mass(c)%source(n,j,i) + surface%amsws(m,c) c = ( ioc - 1 ) * nbins + b surface%amsws(m,c) = surface%amsws(m,c) + & emission_mass_fracs(n,2) / & ( emission_mass_fracs(n,1) + & emission_mass_fracs(n,2) ) * & aerosol_number(b)%source(n,j,i) * & api6 * aero(b)%dmid**3.0_wp * arhooc & * rho_air_zw(k-1) aerosol_mass(c)%source(n,j,i) = & aerosol_mass(c)%source(n,j,i) + surface%amsws(m,c) ! !-- Only sulphates ELSEIF ( iso4 > 0 .AND. ioc < 0 ) THEN c = ( iso4 - 1 ) * nbins + b surface%amsws(m,c) = surface%amsws(m,c) + & aerosol_number(b)%source(n,j,i) * api6 & * aero(b)%dmid**3.0_wp * arhoh2so4 & * rho_air_zw(k-1) aerosol_mass(c)%source(n,j,i) = & aerosol_mass(c)%source(n,j,i) + surface%amsws(m,c) ! !-- Only organic carbon ELSEIF ( iso4 < 0 .AND. ioc > 0 ) THEN c = ( ioc - 1 ) * nbins + b surface%amsws(m,c) = surface%amsws(m,c) + & aerosol_number(b)%source(n,j,i) * api6 & * aero(b)%dmid**3.0_wp * arhooc & * rho_air_zw(k-1) aerosol_mass(c)%source(n,j,i) = & aerosol_mass(c)%source(n,j,i) + surface%amsws(m,c) ENDIF ELSEIF ( b > fn1a ) THEN ! !-- Sulphate IF ( iso4 > 0 ) THEN CALL set_mass_flux( surface, m, b, iso4, n, & emission_mass_fracs(n,1), arhoh2so4, & aerosol_number(b)%source(n,j,i) ) ENDIF ! !-- Organic carbon IF ( ioc > 0 ) THEN CALL set_mass_flux( surface, m, b, ioc, n, & emission_mass_fracs(n,2), arhooc, & aerosol_number(b)%source(n,j,i) ) ENDIF ! !-- Black carbon IF ( ibc > 0 ) THEN CALL set_mass_flux( surface, m, b, ibc, n, & emission_mass_fracs(n,3), arhobc, & aerosol_number(b)%source(n,j,i) ) ENDIF ! !-- Dust IF ( idu > 0 ) THEN CALL set_mass_flux( surface, m, b, idu, n, & emission_mass_fracs(n,4), arhodu, & aerosol_number(b)%source(n,j,i) ) ENDIF ! !-- Sea salt IF ( iss > 0 ) THEN CALL set_mass_flux( surface, m, b, iss, n, & emission_mass_fracs(n,5), arhoss, & aerosol_number(b)%source(n,j,i) ) ENDIF ! !-- Nitric acid IF ( ino > 0 ) THEN CALL set_mass_flux( surface, m, b, ino, n, & emission_mass_fracs(n,6), arhohno3, & aerosol_number(b)%source(n,j,i) ) ENDIF ! !-- Ammonia IF ( inh > 0 ) THEN CALL set_mass_flux( surface, m, b, inh, n, & emission_mass_fracs(n,7), arhonh3, & aerosol_number(b)%source(n,j,i) ) ENDIF ENDIF ! !-- Save number fluxes in the end surface%answs(m,b) = surface%answs(m,b) + & aerosol_number(b)%source(n,j,i) * rho_air_zw(k-1) aerosol_number(b)%source(n,j,i) = surface%answs(m,b) ENDDO ENDDO ENDDO END SUBROUTINE set_flux !------------------------------------------------------------------------------! ! Description: ! ------------ !> Sets the mass emissions to aerosol arrays in 2a and 2b. !------------------------------------------------------------------------------! SUBROUTINE set_mass_flux( surface, surf_num, b, ispec, n, mass_frac, prho, & nsource ) USE arrays_3d, & ONLY: rho_air_zw USE surface_mod, & ONLY: surf_type IMPLICIT NONE INTEGER(iwp), INTENT(in) :: b !< Aerosol size bin index INTEGER(iwp), INTENT(in) :: ispec !< Aerosol species index INTEGER(iwp), INTENT(in) :: n !< emission category number INTEGER(iwp), INTENT(in) :: surf_num !< index surface elements REAL(wp), INTENT(in) :: mass_frac !< mass fraction of a chemical !< compound in all bins REAL(wp), INTENT(in) :: nsource !< number source (#/m2/s) REAL(wp), INTENT(in) :: prho !< Aerosol density TYPE(surf_type), INTENT(inout) :: surface !< respective surface type INTEGER(iwp) :: ee !< index: end INTEGER(iwp) :: i !< loop index INTEGER(iwp) :: j !< loop index INTEGER(iwp) :: k !< loop index INTEGER(iwp) :: c !< loop index INTEGER(iwp) :: ss ! Sets the mass sources to aerosol arrays in 2a and 2b. !------------------------------------------------------------------------------! SUBROUTINE set_mass_source( k, j, i, ispec, mass_frac, prho, nsource, fillval ) USE surface_mod, & ONLY: surf_type IMPLICIT NONE INTEGER(iwp), INTENT(in) :: ispec !< Aerosol species index REAL(wp), INTENT(in) :: fillval !< _FillValue in the NetCDF file REAL(wp), INTENT(in) :: mass_frac !< mass fraction of a chemical !< compound in all bins REAL(wp), INTENT(in), DIMENSION(:) :: nsource !< number source REAL(wp), INTENT(in) :: prho !< Aerosol density INTEGER(iwp) :: b !< loop index INTEGER(iwp) :: ee !< index: end INTEGER(iwp) :: i !< loop index INTEGER(iwp) :: j !< loop index INTEGER(iwp) :: k !< loop index INTEGER(iwp) :: c !< loop index INTEGER(iwp) :: ss ! Check data output for salsa. !------------------------------------------------------------------------------! SUBROUTINE salsa_check_data_output( var, unit ) USE control_parameters, & ONLY: message_string IMPLICIT NONE CHARACTER (LEN=*) :: unit !< CHARACTER (LEN=*) :: var !< SELECT CASE ( TRIM( var ) ) CASE ( 'g_H2SO4', 'g_HNO3', 'g_NH3', 'g_OCNV', 'g_OCSV', & 'N_bin1', 'N_bin2', 'N_bin3', 'N_bin4', 'N_bin5', 'N_bin6', & 'N_bin7', 'N_bin8', 'N_bin9', 'N_bin10', 'N_bin11', 'N_bin12', & 'Ntot' ) IF ( .NOT. salsa ) THEN message_string = 'output of "' // TRIM( var ) // '" requi' // & 'res salsa = .TRUE.' CALL message( 'check_parameters', 'SA0006', 1, 2, 0, 6, 0 ) ENDIF unit = '#/m3' CASE ( 'LDSA' ) IF ( .NOT. salsa ) THEN message_string = 'output of "' // TRIM( var ) // '" requi' // & 'res salsa = .TRUE.' CALL message( 'check_parameters', 'SA0003', 1, 2, 0, 6, 0 ) ENDIF unit = 'mum2/cm3' CASE ( 'm_bin1', 'm_bin2', 'm_bin3', 'm_bin4', 'm_bin5', 'm_bin6', & 'm_bin7', 'm_bin8', 'm_bin9', 'm_bin10', 'm_bin11', 'm_bin12', & 'PM2.5', 'PM10', 's_BC', 's_DU', 's_H2O', 's_NH', & 's_NO', 's_OC', 's_SO4', 's_SS' ) IF ( .NOT. salsa ) THEN message_string = 'output of "' // TRIM( var ) // '" requi' // & 'res salsa = .TRUE.' CALL message( 'check_parameters', 'SA0001', 1, 2, 0, 6, 0 ) ENDIF unit = 'kg/m3' CASE DEFAULT unit = 'illegal' END SELECT END SUBROUTINE salsa_check_data_output !------------------------------------------------------------------------------! ! ! Description: ! ------------ !> Subroutine for averaging 3D data !------------------------------------------------------------------------------! SUBROUTINE salsa_3d_data_averaging( mode, variable ) USE control_parameters USE indices USE kinds IMPLICIT NONE CHARACTER (LEN=*) :: mode !< CHARACTER (LEN=*) :: variable !< INTEGER(iwp) :: b !< INTEGER(iwp) :: c !< INTEGER(iwp) :: i !< INTEGER(iwp) :: icc !< INTEGER(iwp) :: j !< INTEGER(iwp) :: k !< REAL(wp) :: df !< For calculating LDSA: fraction of particles !< depositing in the alveolar (or tracheobronchial) !< region of the lung. Depends on the particle size REAL(wp) :: mean_d !< Particle diameter in micrometres REAL(wp) :: nc !< Particle number concentration in units 1/cm**3 REAL(wp) :: temp_bin !< REAL(wp), DIMENSION(:,:,:), POINTER :: to_be_resorted !< points to !< selected output variable temp_bin = 0.0_wp IF ( mode == 'allocate' ) THEN SELECT CASE ( TRIM( variable ) ) CASE ( 'g_H2SO4' ) IF ( .NOT. ALLOCATED( g_H2SO4_av ) ) THEN ALLOCATE( g_H2SO4_av(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ) ENDIF g_H2SO4_av = 0.0_wp CASE ( 'g_HNO3' ) IF ( .NOT. ALLOCATED( g_HNO3_av ) ) THEN ALLOCATE( g_HNO3_av(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ) ENDIF g_HNO3_av = 0.0_wp CASE ( 'g_NH3' ) IF ( .NOT. ALLOCATED( g_NH3_av ) ) THEN ALLOCATE( g_NH3_av(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ) ENDIF g_NH3_av = 0.0_wp CASE ( 'g_OCNV' ) IF ( .NOT. ALLOCATED( g_OCNV_av ) ) THEN ALLOCATE( g_OCNV_av(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ) ENDIF g_OCNV_av = 0.0_wp CASE ( 'g_OCSV' ) IF ( .NOT. ALLOCATED( g_OCSV_av ) ) THEN ALLOCATE( g_OCSV_av(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ) ENDIF g_OCSV_av = 0.0_wp CASE ( 'LDSA' ) IF ( .NOT. ALLOCATED( LDSA_av ) ) THEN ALLOCATE( LDSA_av(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ) ENDIF LDSA_av = 0.0_wp CASE ( 'N_bin1', 'N_bin2', 'N_bin3', 'N_bin4', 'N_bin5', 'N_bin6', & 'N_bin7', 'N_bin8', 'N_bin9', 'N_bin10', 'N_bin11', 'N_bin12' ) IF ( .NOT. ALLOCATED( Nbins_av ) ) THEN ALLOCATE( Nbins_av(nzb:nzt+1,nysg:nyng,nxlg:nxrg,nbins) ) ENDIF Nbins_av = 0.0_wp CASE ( 'Ntot' ) IF ( .NOT. ALLOCATED( Ntot_av ) ) THEN ALLOCATE( Ntot_av(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ) ENDIF Ntot_av = 0.0_wp CASE ( 'm_bin1', 'm_bin2', 'm_bin3', 'm_bin4', 'm_bin5', 'm_bin6', & 'm_bin7', 'm_bin8', 'm_bin9', 'm_bin10', 'm_bin11', 'm_bin12' ) IF ( .NOT. ALLOCATED( mbins_av ) ) THEN ALLOCATE( mbins_av(nzb:nzt+1,nysg:nyng,nxlg:nxrg,nbins) ) ENDIF mbins_av = 0.0_wp CASE ( 'PM2.5' ) IF ( .NOT. ALLOCATED( PM25_av ) ) THEN ALLOCATE( PM25_av(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ) ENDIF PM25_av = 0.0_wp CASE ( 'PM10' ) IF ( .NOT. ALLOCATED( PM10_av ) ) THEN ALLOCATE( PM10_av(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ) ENDIF PM10_av = 0.0_wp CASE ( 's_BC' ) IF ( .NOT. ALLOCATED( s_BC_av ) ) THEN ALLOCATE( s_BC_av(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ) ENDIF s_BC_av = 0.0_wp CASE ( 's_DU' ) IF ( .NOT. ALLOCATED( s_DU_av ) ) THEN ALLOCATE( s_DU_av(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ) ENDIF s_DU_av = 0.0_wp CASE ( 's_H2O' ) IF ( .NOT. ALLOCATED( s_H2O_av ) ) THEN ALLOCATE( s_H2O_av(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ) ENDIF s_H2O_av = 0.0_wp CASE ( 's_NH' ) IF ( .NOT. ALLOCATED( s_NH_av ) ) THEN ALLOCATE( s_NH_av(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ) ENDIF s_NH_av = 0.0_wp CASE ( 's_NO' ) IF ( .NOT. ALLOCATED( s_NO_av ) ) THEN ALLOCATE( s_NO_av(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ) ENDIF s_NO_av = 0.0_wp CASE ( 's_OC' ) IF ( .NOT. ALLOCATED( s_OC_av ) ) THEN ALLOCATE( s_OC_av(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ) ENDIF s_OC_av = 0.0_wp CASE ( 's_SO4' ) IF ( .NOT. ALLOCATED( s_SO4_av ) ) THEN ALLOCATE( s_SO4_av(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ) ENDIF s_SO4_av = 0.0_wp CASE ( 's_SS' ) IF ( .NOT. ALLOCATED( s_SS_av ) ) THEN ALLOCATE( s_SS_av(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ) ENDIF s_SS_av = 0.0_wp CASE DEFAULT CONTINUE END SELECT ELSEIF ( mode == 'sum' ) THEN SELECT CASE ( TRIM( variable ) ) CASE ( 'g_H2SO4', 'g_HNO3', 'g_NH3', 'g_OCNV', 'g_OCSV' ) IF ( TRIM( variable(3:) ) == 'H2SO4' ) THEN icc = 1 to_be_resorted => g_H2SO4_av ELSEIF ( TRIM( variable(3:) ) == 'HNO3' ) THEN icc = 2 to_be_resorted => g_HNO3_av ELSEIF ( TRIM( variable(3:) ) == 'NH3' ) THEN icc = 3 to_be_resorted => g_NH3_av ELSEIF ( TRIM( variable(3:) ) == 'OCNV' ) THEN icc = 4 to_be_resorted => g_OCNV_av ELSEIF ( TRIM( variable(3:) ) == 'OCSV' ) THEN icc = 5 to_be_resorted => g_OCSV_av ENDIF DO i = nxlg, nxrg DO j = nysg, nyng DO k = nzb, nzt+1 to_be_resorted(k,j,i) = to_be_resorted(k,j,i) + & salsa_gas(icc)%conc(k,j,i) ENDDO ENDDO ENDDO CASE ( 'LDSA' ) DO i = nxlg, nxrg DO j = nysg, nyng DO k = nzb, nzt+1 temp_bin = 0.0_wp DO b = 1, nbins ! !-- Diameter in micrometres mean_d = 1.0E+6_wp * Ra_dry(k,j,i,b) * 2.0_wp ! !-- Deposition factor: alveolar (use Ra_dry) df = ( 0.01555_wp / mean_d ) * ( EXP( -0.416_wp * & ( LOG( mean_d ) + 2.84_wp )**2.0_wp ) & + 19.11_wp * EXP( -0.482_wp * & ( LOG( mean_d ) - 1.362_wp )**2.0_wp ) ) ! !-- Number concentration in 1/cm3 nc = 1.0E-6_wp * aerosol_number(b)%conc(k,j,i) ! !-- Lung-deposited surface area LDSA (units mum2/cm3) temp_bin = temp_bin + pi * mean_d**2.0_wp * df * nc ENDDO LDSA_av(k,j,i) = LDSA_av(k,j,i) + temp_bin ENDDO ENDDO ENDDO CASE ( 'N_bin1', 'N_bin2', 'N_bin3', 'N_bin4', 'N_bin5', 'N_bin6', & 'N_bin7', 'N_bin8', 'N_bin9', 'N_bin10', 'N_bin11', 'N_bin12' ) DO i = nxlg, nxrg DO j = nysg, nyng DO k = nzb, nzt+1 DO b = 1, nbins Nbins_av(k,j,i,b) = Nbins_av(k,j,i,b) + & aerosol_number(b)%conc(k,j,i) ENDDO ENDDO ENDDO ENDDO CASE ( 'Ntot' ) DO i = nxlg, nxrg DO j = nysg, nyng DO k = nzb, nzt+1 DO b = 1, nbins Ntot_av(k,j,i) = Ntot_av(k,j,i) + & aerosol_number(b)%conc(k,j,i) ENDDO ENDDO ENDDO ENDDO CASE ( 'm_bin1', 'm_bin2', 'm_bin3', 'm_bin4', 'm_bin5', 'm_bin6', & 'm_bin7', 'm_bin8', 'm_bin9', 'm_bin10', 'm_bin11', 'm_bin12' ) DO i = nxlg, nxrg DO j = nysg, nyng DO k = nzb, nzt+1 DO b = 1, nbins DO c = b, nbins*ncc_tot, nbins mbins_av(k,j,i,b) = mbins_av(k,j,i,b) + & aerosol_mass(c)%conc(k,j,i) ENDDO ENDDO ENDDO ENDDO ENDDO CASE ( 'PM2.5' ) DO i = nxlg, nxrg DO j = nysg, nyng DO k = nzb, nzt+1 temp_bin = 0.0_wp DO b = 1, nbins IF ( 2.0_wp * Ra_dry(k,j,i,b) <= 2.5E-6_wp ) THEN DO c = b, nbins*ncc, nbins temp_bin = temp_bin + aerosol_mass(c)%conc(k,j,i) ENDDO ENDIF ENDDO PM25_av(k,j,i) = PM25_av(k,j,i) + temp_bin ENDDO ENDDO ENDDO CASE ( 'PM10' ) DO i = nxlg, nxrg DO j = nysg, nyng DO k = nzb, nzt+1 temp_bin = 0.0_wp DO b = 1, nbins IF ( 2.0_wp * Ra_dry(k,j,i,b) <= 10.0E-6_wp ) THEN DO c = b, nbins*ncc, nbins temp_bin = temp_bin + aerosol_mass(c)%conc(k,j,i) ENDDO ENDIF ENDDO PM10_av(k,j,i) = PM10_av(k,j,i) + temp_bin ENDDO ENDDO ENDDO CASE ( 's_BC', 's_DU', 's_H2O', 's_NH', 's_NO', 's_OC', 's_SO4', & 's_SS' ) IF ( is_used( prtcl, TRIM( variable(3:) ) ) ) THEN icc = get_index( prtcl, TRIM( variable(3:) ) ) IF ( TRIM( variable(3:) ) == 'BC' ) to_be_resorted => s_BC_av IF ( TRIM( variable(3:) ) == 'DU' ) to_be_resorted => s_DU_av IF ( TRIM( variable(3:) ) == 'NH' ) to_be_resorted => s_NH_av IF ( TRIM( variable(3:) ) == 'NO' ) to_be_resorted => s_NO_av IF ( TRIM( variable(3:) ) == 'OC' ) to_be_resorted => s_OC_av IF ( TRIM( variable(3:) ) == 'SO4' ) to_be_resorted => s_SO4_av IF ( TRIM( variable(3:) ) == 'SS' ) to_be_resorted => s_SS_av DO i = nxlg, nxrg DO j = nysg, nyng DO k = nzb, nzt+1 DO c = ( icc-1 )*nbins+1, icc*nbins to_be_resorted(k,j,i) = to_be_resorted(k,j,i) + & aerosol_mass(c)%conc(k,j,i) ENDDO ENDDO ENDDO ENDDO ENDIF CASE DEFAULT CONTINUE END SELECT ELSEIF ( mode == 'average' ) THEN SELECT CASE ( TRIM( variable ) ) CASE ( 'g_H2SO4', 'g_HNO3', 'g_NH3', 'g_OCNV', 'g_OCSV' ) IF ( TRIM( variable(3:) ) == 'H2SO4' ) THEN icc = 1 to_be_resorted => g_H2SO4_av ELSEIF ( TRIM( variable(3:) ) == 'HNO3' ) THEN icc = 2 to_be_resorted => g_HNO3_av ELSEIF ( TRIM( variable(3:) ) == 'NH3' ) THEN icc = 3 to_be_resorted => g_NH3_av ELSEIF ( TRIM( variable(3:) ) == 'OCNV' ) THEN icc = 4 to_be_resorted => g_OCNV_av ELSEIF ( TRIM( variable(3:) ) == 'OCSV' ) THEN icc = 5 to_be_resorted => g_OCSV_av ENDIF DO i = nxlg, nxrg DO j = nysg, nyng DO k = nzb, nzt+1 to_be_resorted(k,j,i) = to_be_resorted(k,j,i) & / REAL( average_count_3d, KIND=wp ) ENDDO ENDDO ENDDO CASE ( 'LDSA' ) DO i = nxlg, nxrg DO j = nysg, nyng DO k = nzb, nzt+1 LDSA_av(k,j,i) = LDSA_av(k,j,i) & / REAL( average_count_3d, KIND=wp ) ENDDO ENDDO ENDDO CASE ( 'N_bin1', 'N_bin2', 'N_bin3', 'N_bin4', 'N_bin5', 'N_bin6', & 'N_bin7', 'N_bin8', 'N_bin9', 'N_bin10', 'N_bin11', 'N_bin12' ) DO i = nxlg, nxrg DO j = nysg, nyng DO k = nzb, nzt+1 DO b = 1, nbins Nbins_av(k,j,i,b) = Nbins_av(k,j,i,b) & / REAL( average_count_3d, KIND=wp ) ENDDO ENDDO ENDDO ENDDO CASE ( 'Ntot' ) DO i = nxlg, nxrg DO j = nysg, nyng DO k = nzb, nzt+1 Ntot_av(k,j,i) = Ntot_av(k,j,i) & / REAL( average_count_3d, KIND=wp ) ENDDO ENDDO ENDDO CASE ( 'm_bin1', 'm_bin2', 'm_bin3', 'm_bin4', 'm_bin5', 'm_bin6', & 'm_bin7', 'm_bin8', 'm_bin9', 'm_bin10', 'm_bin11', 'm_bin12' ) DO i = nxlg, nxrg DO j = nysg, nyng DO k = nzb, nzt+1 DO b = 1, nbins DO c = b, nbins*ncc, nbins mbins_av(k,j,i,b) = mbins_av(k,j,i,b) & / REAL( average_count_3d, KIND=wp ) ENDDO ENDDO ENDDO ENDDO ENDDO CASE ( 'PM2.5' ) DO i = nxlg, nxrg DO j = nysg, nyng DO k = nzb, nzt+1 PM25_av(k,j,i) = PM25_av(k,j,i) / & REAL( average_count_3d, KIND=wp ) ENDDO ENDDO ENDDO CASE ( 'PM10' ) DO i = nxlg, nxrg DO j = nysg, nyng DO k = nzb, nzt+1 PM10_av(k,j,i) = PM10_av(k,j,i) / & REAL( average_count_3d, KIND=wp ) ENDDO ENDDO ENDDO CASE ( 's_BC', 's_DU', 's_H2O', 's_NH', 's_NO', 's_OC', 's_SO4', & 's_SS' ) IF ( is_used( prtcl, TRIM( variable(3:) ) ) ) THEN icc = get_index( prtcl, TRIM( variable(3:) ) ) IF ( TRIM( variable(3:) ) == 'BC' ) to_be_resorted => s_BC_av IF ( TRIM( variable(3:) ) == 'DU' ) to_be_resorted => s_DU_av IF ( TRIM( variable(3:) ) == 'NH' ) to_be_resorted => s_NH_av IF ( TRIM( variable(3:) ) == 'NO' ) to_be_resorted => s_NO_av IF ( TRIM( variable(3:) ) == 'OC' ) to_be_resorted => s_OC_av IF ( TRIM( variable(3:) ) == 'SO4' ) to_be_resorted => s_SO4_av IF ( TRIM( variable(3:) ) == 'SS' ) to_be_resorted => s_SS_av DO i = nxlg, nxrg DO j = nysg, nyng DO k = nzb, nzt+1 to_be_resorted(k,j,i) = to_be_resorted(k,j,i) / & REAL( average_count_3d, KIND=wp ) ENDDO ENDDO ENDDO ENDIF END SELECT ENDIF END SUBROUTINE salsa_3d_data_averaging !------------------------------------------------------------------------------! ! ! Description: ! ------------ !> Subroutine defining 2D output variables !------------------------------------------------------------------------------! SUBROUTINE salsa_data_output_2d( av, variable, found, grid, mode, local_pf, & two_d, nzb_do, nzt_do ) USE indices USE kinds IMPLICIT NONE CHARACTER (LEN=*) :: grid !< CHARACTER (LEN=*) :: mode !< CHARACTER (LEN=*) :: variable !< CHARACTER (LEN=5) :: vari !< trimmed format of variable INTEGER(iwp) :: av !< INTEGER(iwp) :: b !< running index: size bins INTEGER(iwp) :: c !< running index: mass bins INTEGER(iwp) :: i !< INTEGER(iwp) :: icc !< index of a chemical compound INTEGER(iwp) :: j !< INTEGER(iwp) :: k !< INTEGER(iwp) :: nzb_do !< INTEGER(iwp) :: nzt_do !< LOGICAL :: found !< LOGICAL :: two_d !< flag parameter that indicates 2D variables !< (horizontal cross sections) REAL(wp) :: df !< For calculating LDSA: fraction of particles !< depositing in the alveolar (or tracheobronchial) !< region of the lung. Depends on the particle size REAL(wp) :: fill_value = -9999.0_wp !< value for the _FillValue attribute REAL(wp) :: mean_d !< Particle diameter in micrometres REAL(wp) :: nc !< Particle number concentration in units 1/cm**3 REAL(wp) :: temp_bin !< temporary array for calculating output variables REAL(wp), DIMENSION(nxl:nxr,nys:nyn,nzb_do:nzt_do) :: local_pf !< output REAL(wp), DIMENSION(:,:,:), POINTER :: to_be_resorted !< pointer found = .TRUE. temp_bin = 0.0_wp IF ( TRIM( variable(1:2) ) == 'g_' ) THEN vari = TRIM( variable( 3:LEN( TRIM( variable ) ) - 3 ) ) IF ( av == 0 ) THEN IF ( vari == 'H2SO4') icc = 1 IF ( vari == 'HNO3') icc = 2 IF ( vari == 'NH3') icc = 3 IF ( vari == 'OCNV') icc = 4 IF ( vari == 'OCSV') icc = 5 DO i = nxl, nxr DO j = nys, nyn DO k = nzb_do, nzt_do local_pf(i,j,k) = MERGE( salsa_gas(icc)%conc(k,j,i), & REAL( fill_value, KIND = wp ), & BTEST( wall_flags_0(k,j,i), 0 ) ) ENDDO ENDDO ENDDO ELSE IF ( vari == 'H2SO4' ) to_be_resorted => g_H2SO4_av IF ( vari == 'HNO3' ) to_be_resorted => g_HNO3_av IF ( vari == 'NH3' ) to_be_resorted => g_NH3_av IF ( vari == 'OCNV' ) to_be_resorted => g_OCNV_av IF ( vari == 'OCSV' ) to_be_resorted => g_OCSV_av DO i = nxl, nxr DO j = nys, nyn DO k = nzb_do, nzt_do local_pf(i,j,k) = MERGE( to_be_resorted(k,j,i), & REAL( fill_value, KIND = wp ), & BTEST( wall_flags_0(k,j,i), 0 ) ) ENDDO ENDDO ENDDO ENDIF IF ( mode == 'xy' ) grid = 'zu' ELSEIF ( TRIM( variable(1:4) ) == 'LDSA' ) THEN IF ( av == 0 ) THEN DO i = nxl, nxr DO j = nys, nyn DO k = nzb_do, nzt_do temp_bin = 0.0_wp DO b = 1, nbins ! !-- Diameter in micrometres mean_d = 1.0E+6_wp * Ra_dry(k,j,i,b) * 2.0_wp ! !-- Deposition factor: alveolar df = ( 0.01555_wp / mean_d ) * ( EXP( -0.416_wp * ( LOG( & mean_d ) + 2.84_wp )**2.0_wp ) + 19.11_wp * EXP( & -0.482_wp * ( LOG( mean_d ) - 1.362_wp )**2.0_wp ) ) ! !-- Number concentration in 1/cm3 nc = 1.0E-6_wp * aerosol_number(b)%conc(k,j,i) ! !-- Lung-deposited surface area LDSA (units mum2/cm3) temp_bin = temp_bin + pi * mean_d**2.0_wp * df * nc ENDDO local_pf(i,j,k) = MERGE( temp_bin, REAL( fill_value, KIND = & wp ), BTEST( wall_flags_0(k,j,i), 0 ) ) ENDDO ENDDO ENDDO ELSE DO i = nxl, nxr DO j = nys, nyn DO k = nzb_do, nzt_do local_pf(i,j,k) = MERGE( LDSA_av(k,j,i), REAL( fill_value, & KIND = wp ), BTEST( & wall_flags_0(k,j,i), 0 ) ) ENDDO ENDDO ENDDO ENDIF IF ( mode == 'xy' ) grid = 'zu' ELSEIF ( TRIM( variable(1:5) ) == 'N_bin' ) THEN vari = TRIM( variable( 6:LEN( TRIM( variable ) ) - 3 ) ) IF ( TRIM( vari ) == '1' ) b = 1 IF ( TRIM( vari ) == '2' ) b = 2 IF ( TRIM( vari ) == '3' ) b = 3 IF ( TRIM( vari ) == '4' ) b = 4 IF ( TRIM( vari ) == '5' ) b = 5 IF ( TRIM( vari ) == '6' ) b = 6 IF ( TRIM( vari ) == '7' ) b = 7 IF ( TRIM( vari ) == '8' ) b = 8 IF ( TRIM( vari ) == '9' ) b = 9 IF ( TRIM( vari ) == '10' ) b = 10 IF ( TRIM( vari ) == '11' ) b = 11 IF ( TRIM( vari ) == '12' ) b = 12 IF ( av == 0 ) THEN DO i = nxl, nxr DO j = nys, nyn DO k = nzb_do, nzt_do local_pf(i,j,k) = MERGE( aerosol_number(b)%conc(k,j,i), & REAL( fill_value, KIND = wp ), & BTEST( wall_flags_0(k,j,i), 0 ) ) ENDDO ENDDO ENDDO ELSE DO i = nxl, nxr DO j = nys, nyn DO k = nzb_do, nzt_do local_pf(i,j,k) = MERGE( Nbins_av(k,j,i,b), & REAL( fill_value, KIND = wp ), & BTEST( wall_flags_0(k,j,i), 0 ) ) ENDDO ENDDO ENDDO ENDIF IF ( mode == 'xy' ) grid = 'zu' ELSEIF ( TRIM( variable(1:4) ) == 'Ntot' ) THEN IF ( av == 0 ) THEN DO i = nxl, nxr DO j = nys, nyn DO k = nzb_do, nzt_do temp_bin = 0.0_wp DO b = 1, nbins temp_bin = temp_bin + aerosol_number(b)%conc(k,j,i) ENDDO local_pf(i,j,k) = MERGE( temp_bin, REAL( fill_value, KIND = & wp ), BTEST( wall_flags_0(k,j,i), 0 ) ) ENDDO ENDDO ENDDO ELSE DO i = nxl, nxr DO j = nys, nyn DO k = nzb_do, nzt_do local_pf(i,j,k) = MERGE( Ntot_av(k,j,i), REAL( fill_value, & KIND = wp ), BTEST( & wall_flags_0(k,j,i), 0 ) ) ENDDO ENDDO ENDDO ENDIF IF ( mode == 'xy' ) grid = 'zu' ELSEIF ( TRIM( variable(1:5) ) == 'm_bin' ) THEN vari = TRIM( variable( 6:LEN( TRIM( variable ) ) - 3 ) ) IF ( TRIM( vari ) == '1' ) b = 1 IF ( TRIM( vari ) == '2' ) b = 2 IF ( TRIM( vari ) == '3' ) b = 3 IF ( TRIM( vari ) == '4' ) b = 4 IF ( TRIM( vari ) == '5' ) b = 5 IF ( TRIM( vari ) == '6' ) b = 6 IF ( TRIM( vari ) == '7' ) b = 7 IF ( TRIM( vari ) == '8' ) b = 8 IF ( TRIM( vari ) == '9' ) b = 9 IF ( TRIM( vari ) == '10' ) b = 10 IF ( TRIM( vari ) == '11' ) b = 11 IF ( TRIM( vari ) == '12' ) b = 12 IF ( av == 0 ) THEN DO i = nxl, nxr DO j = nys, nyn DO k = nzb_do, nzt_do temp_bin = 0.0_wp DO c = b, ncc_tot * nbins, nbins temp_bin = temp_bin + aerosol_mass(c)%conc(k,j,i) ENDDO local_pf(i,j,k) = MERGE( temp_bin, REAL( fill_value, & KIND = wp ), BTEST( & wall_flags_0(k,j,i), 0 ) ) ENDDO ENDDO ENDDO ELSE DO i = nxl, nxr DO j = nys, nyn DO k = nzb_do, nzt_do local_pf(i,j,k) = MERGE( mbins_av(k,j,i,b), REAL( fill_value,& KIND = wp ), BTEST( & wall_flags_0(k,j,i), 0 ) ) ENDDO ENDDO ENDDO ENDIF IF ( mode == 'xy' ) grid = 'zu' ELSEIF ( TRIM( variable(1:5) ) == 'PM2.5' ) THEN IF ( av == 0 ) THEN DO i = nxl, nxr DO j = nys, nyn DO k = nzb_do, nzt_do temp_bin = 0.0_wp DO b = 1, nbins IF ( 2.0_wp * Ra_dry(k,j,i,b) <= 2.5E-6_wp ) THEN DO c = b, nbins*ncc, nbins temp_bin = temp_bin + aerosol_mass(c)%conc(k,j,i) ENDDO ENDIF ENDDO local_pf(i,j,k) = MERGE( temp_bin, REAL( fill_value, & KIND = wp ), BTEST( & wall_flags_0(k,j,i), 0 ) ) ENDDO ENDDO ENDDO ELSE DO i = nxl, nxr DO j = nys, nyn DO k = nzb_do, nzt_do local_pf(i,j,k) = MERGE( PM25_av(k,j,i), REAL( fill_value, & KIND = wp ), BTEST( & wall_flags_0(k,j,i), 0 ) ) ENDDO ENDDO ENDDO ENDIF IF ( mode == 'xy' ) grid = 'zu' ELSEIF ( TRIM( variable(1:4) ) == 'PM10' ) THEN IF ( av == 0 ) THEN DO i = nxl, nxr DO j = nys, nyn DO k = nzb_do, nzt_do temp_bin = 0.0_wp DO b = 1, nbins IF ( 2.0_wp * Ra_dry(k,j,i,b) <= 10.0E-6_wp ) THEN DO c = b, nbins*ncc, nbins temp_bin = temp_bin + aerosol_mass(c)%conc(k,j,i) ENDDO ENDIF ENDDO local_pf(i,j,k) = MERGE( temp_bin, REAL( fill_value, & KIND = wp ), BTEST( & wall_flags_0(k,j,i), 0 ) ) ENDDO ENDDO ENDDO ELSE DO i = nxl, nxr DO j = nys, nyn DO k = nzb_do, nzt_do local_pf(i,j,k) = MERGE( PM10_av(k,j,i), REAL( fill_value, & KIND = wp ), BTEST( & wall_flags_0(k,j,i), 0 ) ) ENDDO ENDDO ENDDO ENDIF IF ( mode == 'xy' ) grid = 'zu' ELSEIF ( TRIM( variable(1:2) ) == 's_' ) THEN vari = TRIM( variable( 3:LEN( TRIM( variable ) ) - 3 ) ) IF ( is_used( prtcl, vari ) ) THEN icc = get_index( prtcl, vari ) IF ( av == 0 ) THEN DO i = nxl, nxr DO j = nys, nyn DO k = nzb_do, nzt_do temp_bin = 0.0_wp DO c = ( icc-1 )*nbins+1, icc*nbins, 1 temp_bin = temp_bin + aerosol_mass(c)%conc(k,j,i) ENDDO local_pf(i,j,k) = MERGE( temp_bin, REAL( fill_value, & KIND = wp ), BTEST( & wall_flags_0(k,j,i), 0 ) ) ENDDO ENDDO ENDDO ELSE IF ( vari == 'BC' ) to_be_resorted => s_BC_av IF ( vari == 'DU' ) to_be_resorted => s_DU_av IF ( vari == 'NH' ) to_be_resorted => s_NH_av IF ( vari == 'NO' ) to_be_resorted => s_NO_av IF ( vari == 'OC' ) to_be_resorted => s_OC_av IF ( vari == 'SO4' ) to_be_resorted => s_SO4_av IF ( vari == 'SS' ) to_be_resorted => s_SS_av DO i = nxl, nxr DO j = nys, nyn DO k = nzb_do, nzt_do local_pf(i,j,k) = MERGE( to_be_resorted(k,j,i), & REAL( fill_value, KIND = wp ), & BTEST( wall_flags_0(k,j,i), 0 ) ) ENDDO ENDDO ENDDO ENDIF ELSE local_pf = fill_value ENDIF IF ( mode == 'xy' ) grid = 'zu' ELSE found = .FALSE. grid = 'none' ENDIF END SUBROUTINE salsa_data_output_2d !------------------------------------------------------------------------------! ! ! Description: ! ------------ !> Subroutine defining 3D output variables !------------------------------------------------------------------------------! SUBROUTINE salsa_data_output_3d( av, variable, found, local_pf, nzb_do, & nzt_do ) USE indices USE kinds IMPLICIT NONE CHARACTER (LEN=*), INTENT(in) :: variable !< INTEGER(iwp) :: av !< INTEGER(iwp) :: b !< running index: size bins INTEGER(iwp) :: c !< running index: mass bins INTEGER(iwp) :: i !< INTEGER(iwp) :: icc !< index of a chemical compound INTEGER(iwp) :: j !< INTEGER(iwp) :: k !< INTEGER(iwp) :: nzb_do !< INTEGER(iwp) :: nzt_do !< LOGICAL :: found !< REAL(wp) :: df !< For calculating LDSA: fraction of particles !< depositing in the alveolar (or tracheobronchial) !< region of the lung. Depends on the particle size REAL(wp) :: fill_value = -9999.0_wp !< value for the _FillValue attribute REAL(wp) :: mean_d !< Particle diameter in micrometres REAL(wp) :: nc !< Particle number concentration in units 1/cm**3 REAL(wp) :: temp_bin !< temporary array for calculating output variables REAL(sp), DIMENSION(nxl:nxr,nys:nyn,nzb_do:nzt_do) :: local_pf !< local REAL(wp), DIMENSION(:,:,:), POINTER :: to_be_resorted !< pointer found = .TRUE. temp_bin = 0.0_wp SELECT CASE ( TRIM( variable ) ) CASE ( 'g_H2SO4', 'g_HNO3', 'g_NH3', 'g_OCNV', 'g_OCSV' ) IF ( av == 0 ) THEN IF ( TRIM( variable ) == 'g_H2SO4') icc = 1 IF ( TRIM( variable ) == 'g_HNO3') icc = 2 IF ( TRIM( variable ) == 'g_NH3') icc = 3 IF ( TRIM( variable ) == 'g_OCNV') icc = 4 IF ( TRIM( variable ) == 'g_OCSV') icc = 5 DO i = nxl, nxr DO j = nys, nyn DO k = nzb_do, nzt_do local_pf(i,j,k) = MERGE( salsa_gas(icc)%conc(k,j,i), & REAL( fill_value, KIND = wp ), & BTEST( wall_flags_0(k,j,i), 0 ) ) ENDDO ENDDO ENDDO ELSE IF ( TRIM( variable(3:) ) == 'H2SO4' ) to_be_resorted => g_H2SO4_av IF ( TRIM( variable(3:) ) == 'HNO3' ) to_be_resorted => g_HNO3_av IF ( TRIM( variable(3:) ) == 'NH3' ) to_be_resorted => g_NH3_av IF ( TRIM( variable(3:) ) == 'OCNV' ) to_be_resorted => g_OCNV_av IF ( TRIM( variable(3:) ) == 'OCSV' ) to_be_resorted => g_OCSV_av DO i = nxl, nxr DO j = nys, nyn DO k = nzb_do, nzt_do local_pf(i,j,k) = MERGE( to_be_resorted(k,j,i), & REAL( fill_value, KIND = wp ), & BTEST( wall_flags_0(k,j,i), 0 ) ) ENDDO ENDDO ENDDO ENDIF CASE ( 'LDSA' ) IF ( av == 0 ) THEN DO i = nxl, nxr DO j = nys, nyn DO k = nzb_do, nzt_do temp_bin = 0.0_wp DO b = 1, nbins ! !-- Diameter in micrometres mean_d = 1.0E+6_wp * Ra_dry(k,j,i,b) * 2.0_wp ! !-- Deposition factor: alveolar df = ( 0.01555_wp / mean_d ) * ( EXP( -0.416_wp * & ( LOG( mean_d ) + 2.84_wp )**2.0_wp ) + 19.11_wp & * EXP( -0.482_wp * ( LOG( mean_d ) - 1.362_wp & )**2.0_wp ) ) ! !-- Number concentration in 1/cm3 nc = 1.0E-6_wp * aerosol_number(b)%conc(k,j,i) ! !-- Lung-deposited surface area LDSA (units mum2/cm3) temp_bin = temp_bin + pi * mean_d**2.0_wp * df * nc ENDDO local_pf(i,j,k) = MERGE( temp_bin, & REAL( fill_value, KIND = wp ), & BTEST( wall_flags_0(k,j,i), 0 ) ) ENDDO ENDDO ENDDO ELSE DO i = nxl, nxr DO j = nys, nyn DO k = nzb_do, nzt_do local_pf(i,j,k) = MERGE( LDSA_av(k,j,i), & REAL( fill_value, KIND = wp ), & BTEST( wall_flags_0(k,j,i), 0 ) ) ENDDO ENDDO ENDDO ENDIF CASE ( 'N_bin1', 'N_bin2', 'N_bin3', 'N_bin4', 'N_bin5', 'N_bin6', & 'N_bin7', 'N_bin8', 'N_bin9', 'N_bin10' , 'N_bin11', 'N_bin12' ) IF ( TRIM( variable(6:) ) == '1' ) b = 1 IF ( TRIM( variable(6:) ) == '2' ) b = 2 IF ( TRIM( variable(6:) ) == '3' ) b = 3 IF ( TRIM( variable(6:) ) == '4' ) b = 4 IF ( TRIM( variable(6:) ) == '5' ) b = 5 IF ( TRIM( variable(6:) ) == '6' ) b = 6 IF ( TRIM( variable(6:) ) == '7' ) b = 7 IF ( TRIM( variable(6:) ) == '8' ) b = 8 IF ( TRIM( variable(6:) ) == '9' ) b = 9 IF ( TRIM( variable(6:) ) == '10' ) b = 10 IF ( TRIM( variable(6:) ) == '11' ) b = 11 IF ( TRIM( variable(6:) ) == '12' ) b = 12 IF ( av == 0 ) THEN DO i = nxl, nxr DO j = nys, nyn DO k = nzb_do, nzt_do local_pf(i,j,k) = MERGE( aerosol_number(b)%conc(k,j,i), & REAL( fill_value, KIND = wp ), & BTEST( wall_flags_0(k,j,i), 0 ) ) ENDDO ENDDO ENDDO ELSE DO i = nxl, nxr DO j = nys, nyn DO k = nzb_do, nzt_do local_pf(i,j,k) = MERGE( Nbins_av(k,j,i,b), & REAL( fill_value, KIND = wp ), & BTEST( wall_flags_0(k,j,i), 0 ) ) ENDDO ENDDO ENDDO ENDIF CASE ( 'Ntot' ) IF ( av == 0 ) THEN DO i = nxl, nxr DO j = nys, nyn DO k = nzb_do, nzt_do temp_bin = 0.0_wp DO b = 1, nbins temp_bin = temp_bin + aerosol_number(b)%conc(k,j,i) ENDDO local_pf(i,j,k) = MERGE( temp_bin, & REAL( fill_value, KIND = wp ), & BTEST( wall_flags_0(k,j,i), 0 ) ) ENDDO ENDDO ENDDO ELSE DO i = nxl, nxr DO j = nys, nyn DO k = nzb_do, nzt_do local_pf(i,j,k) = MERGE( Ntot_av(k,j,i), & REAL( fill_value, KIND = wp ), & BTEST( wall_flags_0(k,j,i), 0 ) ) ENDDO ENDDO ENDDO ENDIF CASE ( 'm_bin1', 'm_bin2', 'm_bin3', 'm_bin4', 'm_bin5', 'm_bin6', & 'm_bin7', 'm_bin8', 'm_bin9', 'm_bin10' , 'm_bin11', 'm_bin12' ) IF ( TRIM( variable(6:) ) == '1' ) b = 1 IF ( TRIM( variable(6:) ) == '2' ) b = 2 IF ( TRIM( variable(6:) ) == '3' ) b = 3 IF ( TRIM( variable(6:) ) == '4' ) b = 4 IF ( TRIM( variable(6:) ) == '5' ) b = 5 IF ( TRIM( variable(6:) ) == '6' ) b = 6 IF ( TRIM( variable(6:) ) == '7' ) b = 7 IF ( TRIM( variable(6:) ) == '8' ) b = 8 IF ( TRIM( variable(6:) ) == '9' ) b = 9 IF ( TRIM( variable(6:) ) == '10' ) b = 10 IF ( TRIM( variable(6:) ) == '11' ) b = 11 IF ( TRIM( variable(6:) ) == '12' ) b = 12 IF ( av == 0 ) THEN DO i = nxl, nxr DO j = nys, nyn DO k = nzb_do, nzt_do temp_bin = 0.0_wp DO c = b, ncc_tot * nbins, nbins temp_bin = temp_bin + aerosol_mass(c)%conc(k,j,i) ENDDO local_pf(i,j,k) = MERGE( temp_bin, & REAL( fill_value, KIND = wp ), & BTEST( wall_flags_0(k,j,i), 0 ) ) ENDDO ENDDO ENDDO ELSE DO i = nxl, nxr DO j = nys, nyn DO k = nzb_do, nzt_do local_pf(i,j,k) = MERGE( mbins_av(k,j,i,b), & REAL( fill_value, KIND = wp ), & BTEST( wall_flags_0(k,j,i), 0 ) ) ENDDO ENDDO ENDDO ENDIF CASE ( 'PM2.5' ) IF ( av == 0 ) THEN DO i = nxl, nxr DO j = nys, nyn DO k = nzb_do, nzt_do temp_bin = 0.0_wp DO b = 1, nbins IF ( 2.0_wp * Ra_dry(k,j,i,b) <= 2.5E-6_wp ) THEN DO c = b, nbins * ncc, nbins temp_bin = temp_bin + aerosol_mass(c)%conc(k,j,i) ENDDO ENDIF ENDDO local_pf(i,j,k) = MERGE( temp_bin, & REAL( fill_value, KIND = wp ), & BTEST( wall_flags_0(k,j,i), 0 ) ) ENDDO ENDDO ENDDO ELSE DO i = nxl, nxr DO j = nys, nyn DO k = nzb_do, nzt_do local_pf(i,j,k) = MERGE( PM25_av(k,j,i), & REAL( fill_value, KIND = wp ), & BTEST( wall_flags_0(k,j,i), 0 ) ) ENDDO ENDDO ENDDO ENDIF CASE ( 'PM10' ) IF ( av == 0 ) THEN DO i = nxl, nxr DO j = nys, nyn DO k = nzb_do, nzt_do temp_bin = 0.0_wp DO b = 1, nbins IF ( 2.0_wp * Ra_dry(k,j,i,b) <= 10.0E-6_wp ) THEN DO c = b, nbins * ncc, nbins temp_bin = temp_bin + aerosol_mass(c)%conc(k,j,i) ENDDO ENDIF ENDDO local_pf(i,j,k) = MERGE( temp_bin, & REAL( fill_value, KIND = wp ), & BTEST( wall_flags_0(k,j,i), 0 ) ) ENDDO ENDDO ENDDO ELSE DO i = nxl, nxr DO j = nys, nyn DO k = nzb_do, nzt_do local_pf(i,j,k) = MERGE( PM10_av(k,j,i), & REAL( fill_value, KIND = wp ), & BTEST( wall_flags_0(k,j,i), 0 ) ) ENDDO ENDDO ENDDO ENDIF CASE ( 's_BC', 's_DU', 's_H2O', 's_NH', 's_NO', 's_OC', 's_SO4', 's_SS' ) IF ( is_used( prtcl, TRIM( variable(3:) ) ) ) THEN icc = get_index( prtcl, TRIM( variable(3:) ) ) IF ( av == 0 ) THEN DO i = nxl, nxr DO j = nys, nyn DO k = nzb_do, nzt_do temp_bin = 0.0_wp DO c = ( icc-1 )*nbins+1, icc*nbins temp_bin = temp_bin + aerosol_mass(c)%conc(k,j,i) ENDDO local_pf(i,j,k) = MERGE( temp_bin, & REAL( fill_value, KIND = wp ), & BTEST( wall_flags_0(k,j,i), 0 ) ) ENDDO ENDDO ENDDO ELSE IF ( TRIM( variable(3:) ) == 'BC' ) to_be_resorted => s_BC_av IF ( TRIM( variable(3:) ) == 'DU' ) to_be_resorted => s_DU_av IF ( TRIM( variable(3:) ) == 'NH' ) to_be_resorted => s_NH_av IF ( TRIM( variable(3:) ) == 'NO' ) to_be_resorted => s_NO_av IF ( TRIM( variable(3:) ) == 'OC' ) to_be_resorted => s_OC_av IF ( TRIM( variable(3:) ) == 'SO4' ) to_be_resorted => s_SO4_av IF ( TRIM( variable(3:) ) == 'SS' ) to_be_resorted => s_SS_av DO i = nxl, nxr DO j = nys, nyn DO k = nzb_do, nzt_do local_pf(i,j,k) = MERGE( to_be_resorted(k,j,i), & REAL( fill_value, KIND = wp ), & BTEST( wall_flags_0(k,j,i), 0 ) ) ENDDO ENDDO ENDDO ENDIF ENDIF CASE DEFAULT found = .FALSE. END SELECT END SUBROUTINE salsa_data_output_3d !------------------------------------------------------------------------------! ! ! Description: ! ------------ !> Subroutine defining mask output variables !------------------------------------------------------------------------------! SUBROUTINE salsa_data_output_mask( av, variable, found, local_pf ) USE arrays_3d, & ONLY: tend USE control_parameters, & ONLY: mask_size_l, mask_surface, mid USE surface_mod, & ONLY: get_topography_top_index_ji IMPLICIT NONE CHARACTER(LEN=5) :: grid !< flag to distinquish between staggered grid CHARACTER(LEN=*) :: variable !< CHARACTER(LEN=7) :: vari !< trimmed format of variable INTEGER(iwp) :: av !< INTEGER(iwp) :: b !< loop index for aerosol size number bins INTEGER(iwp) :: c !< loop index for chemical components INTEGER(iwp) :: i !< loop index in x-direction INTEGER(iwp) :: icc !< index of a chemical compound INTEGER(iwp) :: j !< loop index in y-direction INTEGER(iwp) :: k !< loop index in z-direction INTEGER(iwp) :: topo_top_ind !< k index of highest horizontal surface LOGICAL :: found !< LOGICAL :: resorted !< REAL(wp) :: df !< For calculating LDSA: fraction of particles !< depositing in the alveolar (or tracheobronchial) !< region of the lung. Depends on the particle size REAL(wp) :: mean_d !< Particle diameter in micrometres REAL(wp) :: nc !< Particle number concentration in units 1/cm**3 REAL(wp) :: temp_bin !< temporary array for calculating output variables REAL(wp), DIMENSION(mask_size_l(mid,1),mask_size_l(mid,2),mask_size_l(mid,3)) :: local_pf !< REAL(wp), DIMENSION(:,:,:), POINTER :: to_be_resorted !< pointer found = .TRUE. resorted = .FALSE. grid = 's' temp_bin = 0.0_wp SELECT CASE ( TRIM( variable ) ) CASE ( 'g_H2SO4', 'g_HNO3', 'g_NH3', 'g_OCNV', 'g_OCSV' ) vari = TRIM( variable ) IF ( av == 0 ) THEN IF ( vari == 'g_H2SO4') to_be_resorted => salsa_gas(1)%conc IF ( vari == 'g_HNO3') to_be_resorted => salsa_gas(2)%conc IF ( vari == 'g_NH3') to_be_resorted => salsa_gas(3)%conc IF ( vari == 'g_OCNV') to_be_resorted => salsa_gas(4)%conc IF ( vari == 'g_OCSV') to_be_resorted => salsa_gas(5)%conc ELSE IF ( vari == 'g_H2SO4') to_be_resorted => g_H2SO4_av IF ( vari == 'g_HNO3') to_be_resorted => g_HNO3_av IF ( vari == 'g_NH3') to_be_resorted => g_NH3_av IF ( vari == 'g_OCNV') to_be_resorted => g_OCNV_av IF ( vari == 'g_OCSV') to_be_resorted => g_OCSV_av ENDIF CASE ( 'LDSA' ) IF ( av == 0 ) THEN DO i = nxl, nxr DO j = nys, nyn DO k = nzb, nz_do3d temp_bin = 0.0_wp DO b = 1, nbins ! !-- Diameter in micrometres mean_d = 1.0E+6_wp * Ra_dry(k,j,i,b) * 2.0_wp ! !-- Deposition factor: alveolar df = ( 0.01555_wp / mean_d ) * ( EXP( -0.416_wp * & ( LOG( mean_d ) + 2.84_wp )**2.0_wp ) + 19.11_wp & * EXP( -0.482_wp * ( LOG( mean_d ) - 1.362_wp & )**2.0_wp ) ) ! !-- Number concentration in 1/cm3 nc = 1.0E-6_wp * aerosol_number(b)%conc(k,j,i) ! !-- Lung-deposited surface area LDSA (units mum2/cm3) temp_bin = temp_bin + pi * mean_d**2.0_wp * df * nc ENDDO tend(k,j,i) = temp_bin ENDDO ENDDO ENDDO IF ( .NOT. mask_surface(mid) ) THEN DO i = 1, mask_size_l(mid,1) DO j = 1, mask_size_l(mid,2) DO k = 1, mask_size_l(mid,3) local_pf(i,j,k) = tend( mask_k(mid,k), mask_j(mid,j),& mask_i(mid,i) ) ENDDO ENDDO ENDDO ELSE DO i = 1, mask_size_l(mid,1) DO j = 1, mask_size_l(mid,2) topo_top_ind = get_topography_top_index_ji( mask_j(mid,j),& mask_i(mid,i),& grid ) DO k = 1, mask_size_l(mid,3) local_pf(i,j,k) = tend( MIN( topo_top_ind+mask_k(mid,k),& nzt+1 ), & mask_j(mid,j), mask_i(mid,i) ) ENDDO ENDDO ENDDO ENDIF resorted = .TRUE. ELSE to_be_resorted => LDSA_av ENDIF CASE ( 'N_bin1', 'N_bin2', 'N_bin3', 'N_bin4', 'N_bin5', 'N_bin6', & 'N_bin7', 'N_bin8', 'N_bin9', 'N_bin10' , 'N_bin11', 'N_bin12' ) IF ( TRIM( variable(6:) ) == '1' ) b = 1 IF ( TRIM( variable(6:) ) == '2' ) b = 2 IF ( TRIM( variable(6:) ) == '3' ) b = 3 IF ( TRIM( variable(6:) ) == '4' ) b = 4 IF ( TRIM( variable(6:) ) == '5' ) b = 5 IF ( TRIM( variable(6:) ) == '6' ) b = 6 IF ( TRIM( variable(6:) ) == '7' ) b = 7 IF ( TRIM( variable(6:) ) == '8' ) b = 8 IF ( TRIM( variable(6:) ) == '9' ) b = 9 IF ( TRIM( variable(6:) ) == '10' ) b = 10 IF ( TRIM( variable(6:) ) == '11' ) b = 11 IF ( TRIM( variable(6:) ) == '12' ) b = 12 IF ( av == 0 ) THEN IF ( .NOT. mask_surface(mid) ) THEN DO i = 1, mask_size_l(mid,1) DO j = 1, mask_size_l(mid,2) DO k = 1, mask_size_l(mid,3) local_pf(i,j,k) = aerosol_number(b)%conc( mask_k(mid,k),& mask_j(mid,j),& mask_i(mid,i) ) ENDDO ENDDO ENDDO ELSE DO i = 1, mask_size_l(mid,1) DO j = 1, mask_size_l(mid,2) topo_top_ind = get_topography_top_index_ji( mask_j(mid,j),& mask_i(mid,i),& grid ) DO k = 1, mask_size_l(mid,3) local_pf(i,j,k) = aerosol_number(b)%conc( & MIN( topo_top_ind+mask_k(mid,k), & nzt+1 ), & mask_j(mid,j), mask_i(mid,i) ) ENDDO ENDDO ENDDO ENDIF resorted = .TRUE. ELSE to_be_resorted => Nbins_av(:,:,:,b) ENDIF CASE ( 'Ntot' ) IF ( av == 0 ) THEN DO i = nxl, nxr DO j = nys, nyn DO k = nzb, nz_do3d temp_bin = 0.0_wp DO b = 1, nbins temp_bin = temp_bin + aerosol_number(b)%conc(k,j,i) ENDDO tend(k,j,i) = temp_bin ENDDO ENDDO ENDDO IF ( .NOT. mask_surface(mid) ) THEN DO i = 1, mask_size_l(mid,1) DO j = 1, mask_size_l(mid,2) DO k = 1, mask_size_l(mid,3) local_pf(i,j,k) = tend( mask_k(mid,k), mask_j(mid,j),& mask_i(mid,i) ) ENDDO ENDDO ENDDO ELSE DO i = 1, mask_size_l(mid,1) DO j = 1, mask_size_l(mid,2) topo_top_ind = get_topography_top_index_ji( mask_j(mid,j),& mask_i(mid,i),& grid ) DO k = 1, mask_size_l(mid,3) local_pf(i,j,k) = tend( MIN( topo_top_ind+mask_k(mid,k),& nzt+1 ), & mask_j(mid,j), mask_i(mid,i) ) ENDDO ENDDO ENDDO ENDIF resorted = .TRUE. ELSE to_be_resorted => Ntot_av ENDIF CASE ( 'm_bin1', 'm_bin2', 'm_bin3', 'm_bin4', 'm_bin5', 'm_bin6', & 'm_bin7', 'm_bin8', 'm_bin9', 'm_bin10' , 'm_bin11', 'm_bin12' ) IF ( TRIM( variable(6:) ) == '1' ) b = 1 IF ( TRIM( variable(6:) ) == '2' ) b = 2 IF ( TRIM( variable(6:) ) == '3' ) b = 3 IF ( TRIM( variable(6:) ) == '4' ) b = 4 IF ( TRIM( variable(6:) ) == '5' ) b = 5 IF ( TRIM( variable(6:) ) == '6' ) b = 6 IF ( TRIM( variable(6:) ) == '7' ) b = 7 IF ( TRIM( variable(6:) ) == '8' ) b = 8 IF ( TRIM( variable(6:) ) == '9' ) b = 9 IF ( TRIM( variable(6:) ) == '10' ) b = 10 IF ( TRIM( variable(6:) ) == '11' ) b = 11 IF ( TRIM( variable(6:) ) == '12' ) b = 12 IF ( av == 0 ) THEN DO i = nxl, nxr DO j = nys, nyn DO k = nzb, nz_do3d temp_bin = 0.0_wp DO c = b, ncc_tot*nbins, nbins temp_bin = temp_bin + aerosol_mass(c)%conc(k,j,i) ENDDO tend(k,j,i) = temp_bin ENDDO ENDDO ENDDO IF ( .NOT. mask_surface(mid) ) THEN DO i = 1, mask_size_l(mid,1) DO j = 1, mask_size_l(mid,2) DO k = 1, mask_size_l(mid,3) local_pf(i,j,k) = tend( mask_k(mid,k), mask_j(mid,j),& mask_i(mid,i) ) ENDDO ENDDO ENDDO ELSE DO i = 1, mask_size_l(mid,1) DO j = 1, mask_size_l(mid,2) topo_top_ind = get_topography_top_index_ji( mask_j(mid,j),& mask_i(mid,i),& grid ) DO k = 1, mask_size_l(mid,3) local_pf(i,j,k) = tend( MIN( topo_top_ind+mask_k(mid,k),& nzt+1 ), & mask_j(mid,j), mask_i(mid,i) ) ENDDO ENDDO ENDDO ENDIF resorted = .TRUE. ELSE to_be_resorted => mbins_av(:,:,:,b) ENDIF CASE ( 'PM2.5' ) IF ( av == 0 ) THEN DO i = nxl, nxr DO j = nys, nyn DO k = nzb, nz_do3d temp_bin = 0.0_wp DO b = 1, nbins IF ( 2.0_wp * Ra_dry(k,j,i,b) <= 2.5E-6_wp ) THEN DO c = b, nbins * ncc, nbins temp_bin = temp_bin + aerosol_mass(c)%conc(k,j,i) ENDDO ENDIF ENDDO tend(k,j,i) = temp_bin ENDDO ENDDO ENDDO IF ( .NOT. mask_surface(mid) ) THEN DO i = 1, mask_size_l(mid,1) DO j = 1, mask_size_l(mid,2) DO k = 1, mask_size_l(mid,3) local_pf(i,j,k) = tend( mask_k(mid,k), mask_j(mid,j),& mask_i(mid,i) ) ENDDO ENDDO ENDDO ELSE DO i = 1, mask_size_l(mid,1) DO j = 1, mask_size_l(mid,2) topo_top_ind = get_topography_top_index_ji( mask_j(mid,j),& mask_i(mid,i),& grid ) DO k = 1, mask_size_l(mid,3) local_pf(i,j,k) = tend( MIN( topo_top_ind+mask_k(mid,k),& nzt+1 ), & mask_j(mid,j), mask_i(mid,i) ) ENDDO ENDDO ENDDO ENDIF resorted = .TRUE. ELSE to_be_resorted => PM25_av ENDIF CASE ( 'PM10' ) IF ( av == 0 ) THEN DO i = nxl, nxr DO j = nys, nyn DO k = nzb, nz_do3d temp_bin = 0.0_wp DO b = 1, nbins IF ( 2.0_wp * Ra_dry(k,j,i,b) <= 10.0E-6_wp ) THEN DO c = b, nbins * ncc, nbins temp_bin = temp_bin + aerosol_mass(c)%conc(k,j,i) ENDDO ENDIF ENDDO tend(k,j,i) = temp_bin ENDDO ENDDO ENDDO IF ( .NOT. mask_surface(mid) ) THEN DO i = 1, mask_size_l(mid,1) DO j = 1, mask_size_l(mid,2) DO k = 1, mask_size_l(mid,3) local_pf(i,j,k) = tend( mask_k(mid,k), mask_j(mid,j),& mask_i(mid,i) ) ENDDO ENDDO ENDDO ELSE DO i = 1, mask_size_l(mid,1) DO j = 1, mask_size_l(mid,2) topo_top_ind = get_topography_top_index_ji( mask_j(mid,j),& mask_i(mid,i),& grid ) DO k = 1, mask_size_l(mid,3) local_pf(i,j,k) = tend( MIN( topo_top_ind+mask_k(mid,k),& nzt+1 ), & mask_j(mid,j), mask_i(mid,i) ) ENDDO ENDDO ENDDO ENDIF resorted = .TRUE. ELSE to_be_resorted => PM10_av ENDIF CASE ( 's_BC', 's_DU', 's_NH', 's_NO', 's_OC', 's_SO4', 's_SS' ) IF ( av == 0 ) THEN IF ( is_used( prtcl, TRIM( variable(3:) ) ) ) THEN icc = get_index( prtcl, TRIM( variable(3:) ) ) DO i = nxl, nxr DO j = nys, nyn DO k = nzb, nz_do3d temp_bin = 0.0_wp DO c = ( icc-1 )*nbins+1, icc*nbins temp_bin = temp_bin + aerosol_mass(c)%conc(k,j,i) ENDDO tend(k,j,i) = temp_bin ENDDO ENDDO ENDDO ELSE tend = 0.0_wp ENDIF IF ( .NOT. mask_surface(mid) ) THEN DO i = 1, mask_size_l(mid,1) DO j = 1, mask_size_l(mid,2) DO k = 1, mask_size_l(mid,3) local_pf(i,j,k) = tend( mask_k(mid,k), mask_j(mid,j), & mask_i(mid,i) ) ENDDO ENDDO ENDDO ELSE DO i = 1, mask_size_l(mid,1) DO j = 1, mask_size_l(mid,2) topo_top_ind = get_topography_top_index_ji( mask_j(mid,j),& mask_i(mid,i),& grid ) DO k = 1, mask_size_l(mid,3) local_pf(i,j,k) = tend( MIN( topo_top_ind+mask_k(mid,k),& nzt+1 ),& mask_j(mid,j), mask_i(mid,i) ) ENDDO ENDDO ENDDO ENDIF resorted = .TRUE. ELSE IF ( TRIM( variable(3:) ) == 'BC' ) to_be_resorted => s_BC_av IF ( TRIM( variable(3:) ) == 'DU' ) to_be_resorted => s_DU_av IF ( TRIM( variable(3:) ) == 'NH' ) to_be_resorted => s_NH_av IF ( TRIM( variable(3:) ) == 'NO' ) to_be_resorted => s_NO_av IF ( TRIM( variable(3:) ) == 'OC' ) to_be_resorted => s_OC_av IF ( TRIM( variable(3:) ) == 'SO4' ) to_be_resorted => s_SO4_av IF ( TRIM( variable(3:) ) == 'SS' ) to_be_resorted => s_SS_av ENDIF CASE DEFAULT found = .FALSE. END SELECT IF ( .NOT. resorted ) THEN IF ( .NOT. mask_surface(mid) ) THEN ! !-- Default masked output DO i = 1, mask_size_l(mid,1) DO j = 1, mask_size_l(mid,2) DO k = 1, mask_size_l(mid,3) local_pf(i,j,k) = to_be_resorted( mask_k(mid,k), & mask_j(mid,j),mask_i(mid,i) ) ENDDO ENDDO ENDDO ELSE ! !-- Terrain-following masked output DO i = 1, mask_size_l(mid,1) DO j = 1, mask_size_l(mid,2) ! !-- Get k index of highest horizontal surface topo_top_ind = get_topography_top_index_ji( mask_j(mid,j), & mask_i(mid,i), grid ) ! !-- Save output array DO k = 1, mask_size_l(mid,3) local_pf(i,j,k) = to_be_resorted( MIN( topo_top_ind+mask_k(mid,k),& nzt+1 ), & mask_j(mid,j), mask_i(mid,i) ) ENDDO ENDDO ENDDO ENDIF ENDIF END SUBROUTINE salsa_data_output_mask END MODULE salsa_mod