Aerosol condensation process and its verification test: error of using Taylor expansion to the solutions of $\rm H_2SO_4/SO_4$

$\color{red}{\text{Last update: 07-25-2019, test strategy, result and presentation need revision.}}$

1 Introduction

Condensation represents the process that the vapor of a compound condenses onto the surface of existing particle population. MAM currently considers two trace gases: sulfur acid gas ($\rm H_2SO_4$, and its condensed particles $\rm SO_4$) and a single lumped gas species of secondary organic aerosol (SOAG, and its condensed particles SOA \cite{liu-2012-gmd}).

By assuming that the particles are spherical and the size distribution of the particle population is described by a few log-normal modes, the governing equations for the condensation of $\rm H_2SO_4$ read:

\begin{align} C_{_{H_2SO_4,i}} &= \int_{-\infty}^{\infty} 2 \pi D_{_p} \mathbb{D}_{_{g,H_2SO_4}} f(K_{_{n,H_2SO_4}},\alpha_{_{H_2SO_4}}) n_{_{i}}(\ln D_{_{p}}) d\ln D_{_{p}} \,, \label{kmt-h2so4} \\ \frac{dq_{_{v,H_2SO_4}}}{dt} &= P_{_{chem}} - q_{_{H_2SO_4}} \sum_{i=1}^{I} C_{_{H_2SO_4,i}} \,, \label{qgh2so4} \\ \frac{dq_{_{m,SO_4,i}}}{dt} &= C_{_{H_2SO_4,i}} \times q_{_{v,H_2SO_4}} \,, i = 1,2,...,I \,, \label{qaso4} \end{align}

where the definitions of variables in Eq. \eqref{kmt-h2so4} can be found in a separate model documentation; $q_{_{v,H_2SO_4}}$ is the mass mixing ratio of $\rm H_2SO_4$ (unit: kmol of $\rm H_2SO_4$ per kmol of dry air); $P_{_{chem}}$ is the net chemical production rate of $\rm H_2SO_4$ due to the gas- and aqueous-phase chemistry (unit: kmol of $\rm H_2SO_4$ per kmol of dry air per second); $q_{_{m,SO_4,i}}$ is the mass mixing ratio of $\rm SO_4$ in mode $i$.

To solve for $q_{_{v,H_2SO_4}}$ and $q_{_{m,SO_4,i}}$, MAM assumes that $P_{_{chem}}$ and $C_{_{H_2SO_4,i}}$ are fixed within one time step $\Delta t_{_{cond}}$ ($\Delta t_{_{cond}} = t_{_1} - t_{_0}$). Therefore, semi-analytical solutions can be derived for $q_{_{v,H_2SO_4}}$ and $q_{_{m,SO_4,i}}$. More details about the basics notations, model parameterization, assumptions and numerical implementation can be found in a separate ''MAM_Basics'' and model documentation.

2 Goals

Though there are semi-analytical solutions for $q_{_{v,H_2SO_4}}$ and $q_{_{m,SO_4,i}}$, MAM applies the Taylor expansion to the semi-analytical solutions with respect to time and truncates at the second-order term, when $\Delta t_{_{cond}} \sum\limits_{i=1}^{I}C_{_{H_2SO_4,i}} \le 10^{-3}$. The main purpose to use Taylor expansion here is to avoid the division of a small number in the formula of semi-analytical solution (see Eq. \eqref{semi-ana-h2so4} and \eqref{semi-ana-so4} in Section 3.1). However, the magnitude of error introduced by the Taylor expansion remains unknown. Thus in this test, we aim at examing the difference between semi-analytical solution and Taylor expansion of semi-analytical solutions of $\rm H_2SO_4/SO_4$. Note that in this test, $P_{_{chem}}$ and $C_{_{H_2SO_4,i}}$ are fixed over the whole simulation time and we do not consider the existence of SOAG/SOA.

3 Steps to establish the verification test

This verification test uses a box model that inherits the data structure and condensation subroutines from MAM in E3SM. The box model serves as a driver program to specify the necessary parameters in order to calculate the MAM and reference solutions.

3.1 Reference solution

The reference solutions of this verification test take the semi-analytical formulae for $\rm H_2SO_4/SO_4$, based on the assumptions mentioned in Section 1. The semi-analytical solutions of Eq. \eqref{qgh2so4} and \eqref{qaso4} at $t = t_{_1}$ read:

\begin{align} q_{_{v,H_2SO_4}}(t_{_1}) &= \left( q_{_{v,H_2SO_4}}(t_{_0}) - \frac{P_{_{chem}} (t_{_0})} {C_{_{sum,H_2SO_4}} (t_{_0})} \right) \exp^{-C_{_{sum,H_2SO_4}} (t_{_0}) \Delta t_{_{cond}}} + \frac{P_{_{chem}} (t_{_0})} {C_{_{sum,H_2SO_4}} (t_{_0})} \,, \label{semi-ana-h2so4} \\ q_{_{m,SO_4,i}}(t_{_1}) &= q_{_{m,SO_4,i}}(t_{_0}) + \frac{C_{_{H_2SO_4,i}} (t_{_0})} {C_{_{sum,H_2SO_4}} (t_{_0})} \left[ \left( q_{_{v,H_2SO_4}}(t_{_0}) - \frac{P_{_{chem}} (t_{_0})} {C_{_{sum,H_2SO_4}} (t_{_0})} \right) \left( 1-\exp^{-C_{_{sum,H_2SO_4}} (t_{_0}) \Delta t_{_{cond}}} \right) + P_{_{chem}} (t_{_0}) \Delta t_{_{cond}} \right] \,. \label{semi-ana-so4} \end{align}

where $C_{_{sum,H_2SO_4}} (t_{_0}) = \sum\limits_{i=1}^{I}C_{_{H_2SO_4,i}} (t_{_0})$.

3.2 MAM solution

When $10^{-20} \le C_{_{sum,H_2SO_4}} (t_{_0}) \Delta t_{_{cond}} \le 10^{-3}$, MAM applies the Taylor expansion to the semi-analytical solution of $q_{_{v,H_2SO_4}}(t)$ with respect to time and truncates at the second-order term. The approximated solution of $q_{_{v,H_2SO_4}}(t)$ at $t = t_{_1}$ reads:

\begin{align} q_{_{v,H_2SO_4}}(t_{_1}) &= q_{_{v,H_2SO_4}}(t_{_0}) \left[ 1 - C_{_{sum,H_2SO_4}} (t_{_0}) \Delta t_{_{cond}} + \frac{\left( C_{_{sum,H_2SO_4}} (t_{_0}) \Delta t_{_{cond}} \right)^2}{2} \right] \nonumber \\ &+ P_{_{chem}} (t_{_0}) \Delta t_{_{cond}} \left[ 1 - \frac{C_{_{sum,H_2SO_4}} (t_{_0}) \Delta t_{_{cond}}} {2} + \frac{\left( C_{_{sum,H_2SO_4}} (t_{_0}) \Delta t_{_{cond}} \right)^2}{6} \right] \,. \label{mam-h2so4} \end{align}

The approximated solution of $q_{_{m,SO_4,i}}(t)$ at $t = t_{_1}$ reads:

\begin{equation} \label{mam-so4} q_{_{m,SO_4,i}}(t_{_1}) = q_{_{m,SO_4,i}}(t_{_0}) + \frac{C_{_{H_2SO_4,i}} (t_{_0})}{C_{_{sum,H_2SO_4}} (t_{_0})} \left( q_{_{v,H_2SO_4}}(t_{_0}) + P_{_{chem}} (t_{_0}) \Delta t_{_{cond}} - q_{_{v,H_2SO_4}}(t_{_1}) \right) \,. \end{equation}

3.3 Parameters required for this verification test

To perform this verification test, the following quantities are specified: %

\begin{enumerate} \item All the parameters used to calculate the $C_{_{H_2SO_4,i}}$. These include: % \begin{itemize} \item Accommodation coefficient $\alpha_{_{H_2SO_4}} = 0.65$ (unitless, \cite{poschl-1998-jpca}). \item Temperature $T$, for which we use 273 K as an example. \item Pressure $P$, for which we use 1000 hPa as an exmaple. \item Molecular diffusion volume (unitless) $\mathbb{V}_{_{D,H_2SO_4}} = 42.88$ and $\mathbb{V}_{_{D,air}} = 20.1$. \item Molecular weight $M_{_{w,H_2SO_4}} = 98.0783997$ and $M_{_{w,air}} = 28.966$. \item Number concentration of particle population in mode $i$ ($N_{_{t,i}}$, unit: $\# \cdot m^{-3}$). We set two groups of $N_{_{t,i}}$ so that $10^{-20} \le C_{_{sum,H_2SO_4}} (t_{_0}) \Delta t_{_{cond}} \le 10^{-3}$ ($\Delta t_{_{cond}}$ is set correspondingly below): \begin{align} \rm Group1&: N_{_{t,1}} = 10^8,& N_{_{t,2}} = 10^9,& &N_{_{t,3}} &= 10^5,& &N_{_{t,4}} = 2 \times 10^8 \,, \nonumber \\ \rm Group2&: N_{_{t,1}} = 10^4,& N_{_{t,2}} = 10^5,& &N_{_{t,3}} &= 10,& &N_{_{t,4}} = 2 \times 10^4 \,. \nonumber \end{align} Group1 is referred to as a situation with typical aerosol number concentrations in reality while Group2 is referred to as a situation with low aerosol number concentrations. \item Geometric mean dry diameter ($D_{_{gn,d,i}}$) and geometric standard deviation ($\sigma_{_{g,i}}$) for each mode. The values of $D_{_{gn,d,i}}$ and $\sigma_{_{g,i}}$ are chosen to represent the midrange values for tropospheric aerosol observations \cite{dick-2004-jgr, liu-2012-gmd}, and are listed in the ''MAM_Basics'' documentation. The geometric mean wet diameter ($D_{_{gn,w,i}}$) is the actual value used to calculate the $C_{_{l,i}}$ and is assumed to be 1.2 times $D_{_{gn,d,i}}$ for simplicity. \item Gauss-Hermite quadrature points and their associated weights. For this test, both MAM and reference solutions use 2 Gauss-Hermite quadrature points originally in MAM (less accurate \cite{zwillinger2002crc}). \end{itemize} % \item Simulation length $\Delta t_{_{tot}}$, which is set to 1800 seconds. \item Time step $\Delta t_{_{cond}}$ (unit: second). We set two groups of $\Delta t_{_{cond}}$ corresponding to the values of $N_{_{t,i}}$ above so that $10^{-20} \le C_{_{H_2SO_4,i}} (t_{_0}) \Delta t_{_{cond}} \le 10^{-3}$: \begin{align} \rm Group1&: \Delta t_{_{cond}} = 1 \,, \nonumber \\ \rm Group2&: \Delta t_{_{cond}} = 1800 \,. \nonumber \end{align} \item Density (unit: $\rm kg \cdot m^{-3}$) and molecular weight (unit: $\rm g \cdot mol^{-1}$) of aerosol species ($\rho_{_{l}}$ and $M_{_{w,l}}$), for which the following values are used: $$[\rho_{_{SO_4}}:1770 \ ; \rho_{_{POA}}:1000 \ ; \rho_{_{SOA}}:1000 \ ; \rho_{_{BC}}:1700 \ ; \rho_{_{SS}}:1900 \ ; \rho_{_{DST}}:2600 \ ; \rho_{_{MOA}}:1601] \,, \nonumber$$ $$[M_{_{w,SO_4}}:115 \ ; M_{_{w,POA}}:12 \ ; M_{_{w,SOA}}:12 \ ; M_{_{w,BC}}:12 \ ; M_{_{w,SS}}:58.5 \ ; M_{_{w,DST}}:135 \ ; M_{_{w,MOA}}:250092] \,. \nonumber$$ \item Initial mass fractions of aerosol species in mode $i$ ($f_{_{m,l,i}}$), whose values are shown below (numbers in the parentheses) as the example setting. \begin{align} \text{ Accumulation}:& &\text{SS(50%)},& &\text{$\rm SO_4$(50%)}, \nonumber \\ \text{ Aitken}:& &\text{SS(50%)},& &\text{$\rm SO_4$(50%)}, \nonumber \\ \text{ Coarse}:& &\text{DST(50%)},& &\text{$\rm SO_4$(50%)}, \nonumber \\ \text{Primary carbon}:& &\text{POA(50%)},& &\text{BC(50%)}. \nonumber \end{align} No initial SOA and $\rm SO_4$ are specified for the primary carbon mode because that mode only contains BC, POA and MOA. By specifying the $T$, $p$, $R$, $N_{_{t,i}}$, $D_{_{gn,d,i}}$, $\sigma_{_{g,i,}}$, $\rho_{_{l}}$, $M_{_{w,l}}$ and $f_{_{m,l,i}}$, we could then calculate the initial $q_{_{m,l,i}}$ as: \begin{equation} q_{_{m,l,i}} = \frac{\frac{\pi}{6} N_{_{t,i}} D_{_{gn,d,i}}^3 \exp \left( \frac{9}{2} \ln^2 \sigma_{_{g,i}} \right)} {\sum\limits_{l}\frac{f_{_{m,l,i}}}{\rho_{_{l}}}} \cdot f_{_{m,l,i}} \cdot \frac{\rho_{_{l}}RT}{1000 M_{_{w,l}} p} \,. \end{equation} \item Initial mass mixing ratio of $q_{_{v,H_2SO_4}}$, which is set to $10^{-13}$, close to the global mean value at the surface layer. \item The net chemical production rate $P_{_{chem}}$ for $\rm H_2SO_4$, which is set to zero for simplicity. \end{enumerate}

With these parameters given, totally 10 (= 5 (4 modes for aerosols + 1 gas species) $\times$ 2 (groups)) verification tests are performed.

3.4 PASS/FAIL criterion

The PASS/FAIL criterion is determined based on two situations:

\begin{itemize} \item If $\Delta t_{_{cond}} = \Delta t_{_{tot}}$, we calculate the MAM solutions at $t = t_{_1}$ based on Eq. \eqref{mam-h2so4} and \eqref{mam-so4} (denoted as $q_{_{v,MAM}}$ and $q_{_{m,MAM,i}}$), and the reference solutions at $t = t_{_1}$ based on Eq. \eqref{semi-ana-h2so4} and \eqref{semi-ana-so4} (denoted as $q_{_{v,REF}}$ and $q_{_{m,REF,i}}$). Then we calculate the relative error $\phi$ as (note $P_{_{chem}}$ is assumed as zero): \begin{align} \phi_{_v} &= \frac{\left| q_{_{v,MAM}} - q_{_{v,REF}} \right|}{q_{_{v,REF}}} = \frac{\left| \left[ 1 - C_{_{sum,H_2SO_4}} (t_{_0}) \Delta t_{_{cond}} + \frac{\left( C_{_{sum,H_2SO_4}} (t_{_0}) \Delta t_{_{cond}} \right)^2}{2} \right] - \exp^{-C_{_{sum,H_2SO_4}} (t_{_0}) \Delta t_{_{cond}}} \right|} {\exp^{-C_{_{sum,H_2SO_4}} (t_{_0}) \Delta t_{_{cond}}}} \nonumber \\ &\lt \frac{\left( C_{_{sum,H_2SO_4}} (t_{_0}) \Delta t_{_{cond}} \right)^3}{6} \exp^{C_{_{sum,H_2SO_4}} (t_{_0}) \Delta t_{_{cond}}} \,, \label{phi_gas} \\ \phi_{_{m,i}} &= \frac{\left| q_{_{m,MAM,i}} - q_{_{m,REF,i}} \right|}{q_{_{m,REF,i}}} = \frac{\frac{C_{_{H_2SO_4,i}} (t_{_0})}{C_{_{sum,H_2SO_4}} (t_{_0})} \left| q_{_{v,MAM}} - q_{_{v,REF}} \right|}{q_{_{m,SO_4,i}}(t_{_0}) + \frac{C_{_{H_2SO_4,i}} (t_{_0})}{C_{_{sum,H_2SO_4}} (t_{_0})} \left( q_{_{v,H_2SO_4}}(t_{_0}) - q_{_{v,REF}} \right)} \lt \frac{\left| q_{_{v,MAM}} - q_{_{v,REF}} \right|}{q_{_{v,H_2SO_4}}(t_{_0}) - q_{_{v,REF}}} = \frac{\phi_{_v}}{\exp^{C_{_{sum,H_2SO_4}} (t_{_0}) \Delta t_{_{cond}}} - 1} \,. \label{phi_aer} \end{align} \item If $\Delta t_{_{cond}} < \Delta t_{_{tot}}$, the $\phi$ estimated for one $\Delta t_{_{cond}}$ is the same as those in Eq. \eqref{phi_gas} and \eqref{phi_aer}. The $\phi$ estimated for one $\Delta t$ is then bounded as: \begin{align} \phi_{_v} &= \frac{\left| q_{_{v,MAM}} - q_{_{v,REF}} \right|}{q_{_{v,REF}}} = \frac{\left| \left[ 1 - C_{_{sum,H_2SO_4}} (t_{_0}) \Delta t_{_{cond}} + \frac{\left( C_{_{sum,H_2SO_4}} (t_{_0}) \Delta t_{_{cond}} \right)^2}{2} \right]^\frac{\Delta t_{_{tot}}}{\Delta t_{_{cond}}} - \exp^{-C_{_{sum,H_2SO_4}} (t_{_0}) \Delta t_{_{tot}}} \right|}{\exp^{-C_{_{sum,H_2SO_4}} (t_{_0}) \Delta t_{_{tot}}}} \nonumber \\ &\lt \frac{\left| \left[ \exp^{-C_{_{sum,H_2SO_4}} (t_{_0}) \Delta t_{_{cond}}} + \frac{\left( C_{_{sum,H_2SO_4}} (t_{_0}) \Delta t_{_{cond}} \right)^3}{6} \right]^\frac{\Delta t_{_{tot}}}{\Delta t_{_{cond}}} - \exp^{-C_{_{sum,H_2SO_4}} (t_{_0}) \Delta t_{_{tot}}} \right|}{\exp^{-C_{_{sum,H_2SO_4}} (t_{_0}) \Delta t_{_{tot}}}} \,, \label{phi_gas2} \\ \phi_{_{m,i}} &= \frac{\left| q_{_{m,MAM,i}} - q_{_{m,REF,i}} \right|}{q_{_{m,REF,i}}} = \frac{\frac{C_{_{H_2SO_4,i}} (t_{_0})}{C_{_{sum,H_2SO_4}} (t_{_0})} \left| q_{_{v,MAM}} - q_{_{v,REF}} \right|}{q_{_{m,SO_4,i}}(t_{_0}) + \frac{C_{_{H_2SO_4,i}} (t_{_0})}{C_{_{sum,H_2SO_4}} (t_{_0})} \left( q_{_{v,H_2SO_4}}(t_{_0}) - q_{_{v,REF}} \right)} \nonumber \\ &\lt \frac{\left| q_{_{v,MAM}} - q_{_{v,REF}} \right|}{q_{_{v,H_2SO_4}}(t_{_0}) - q_{_{v,REF}}} = \frac{\phi_{_v}}{\exp^{C_{_{sum,H_2SO_4}} (t_{_0}) \Delta t_{_{tot}}}-1} \,. \label{phi_aer2} \end{align} \end{itemize}

Therefore, an overall PASS is issued if the following criteria are met for all the verification tests:

\begin{itemize} \item $\phi_{_v}$ satisfies the Eq. \eqref{phi_gas} and $\phi_{_{m,i}}$ satisfies the Eq. \eqref{phi_aer} when $\Delta t_{_{cond}} = \Delta t_{_{tot}}$; \item $\phi_{_v}$ satisfies the Eq. \eqref{phi_gas2} and $\phi_{_{m,i}}$ satisfies the Eq. \eqref{phi_aer2} when $\Delta t_{_{cond}} < \Delta t_{_{tot}}$; \end{itemize}

4 Source codes and configuration files of verification test

The verification test requires a test driver to perform the verification test and the configuration files for compilation and execution of cmdv-test-runner (e.g., the yaml file). If we go the root directory of CMDV-Verification folder that is cloned from GitHub Repo, the test driver code (driver.F90) is under $\rm \color{blue}{tests/mam/restructed\_tests/F90\_SRC/test\_drivers/condensation/h2so4\_taylor}$ and the yaml file is under $\rm \color{blue}{tests/mam/restructed\_tests/yaml\_files\_for\_test\_runner/condensation/h2so4\_taylor.test.yaml}$.

Note that this test is not expected to run nightly and thus not included in the latest CMDV-Verification GitHub Repo. In order to run this test, please checkout the version with hash 5e0c49ae4b6866742414547a995e085ba92d60d8 first.

5 Verification Results

If the verification test runs successfully, we should obtain the netcdf file named "condensation_verification_test.nc" with the MAM and reference solutions for $q_{_{v,H_2SO_4}}$ and $q_{_{m,SO_4,i}}$. Assuming that the netcdf file is put under the "postprocess_input" folder, the following cells would calculate the relative errors and check whether a PASS/FAIL should be issued.

In [1]:
# load the library
import numpy as np
import xarray as xr
In [2]:
# define the file path and name
nc_file   = "./postprocess_input/condensation_verification_test.nc"

# read in the MAM and reference solutions
ds        = xr.open_dataset(nc_file)
pver      = ds['pver'].size
mode      = ds['mode'].size
mam_so4   = ds['mam_so4']
ref_so4   = ds['ref_so4']
mam_h2so4 = ds['mam_h2so4']
ref_h2so4 = ds['ref_h2so4']
uptkrate  = ds['uptkrate']
dt        = ds['dt']
In [3]:
# calculate the relative errors of H2SO4/SO4
phi_gas                  = np.empty(pver)
phi_aer                  = np.empty(shape=[mode,pver])
for k in np.arange(pver):
    if  np.abs(dt[k] - 1800.)/1800. < 1.e-6:
        print("dt = %.f" % dt[k])
        csum             = uptkrate[k]*dt[k]
        phi_gas[k]       = np.abs(mam_h2so4[k] - ref_h2so4[k]) / ref_h2so4[k]
        for n in np.arange(mode):
            phi_aer[n,k] = np.abs(mam_so4[n,k] - ref_so4[n,k]) / ref_so4[n,k]
            print("Relative error of "+u'SO\u2084'+":   %.4e" % phi_aer[n,k])
            phi_aer[n,k] = phi_aer[n,k] < phi_gas[k] / (np.exp(csum) - 1)
        print("Relative error of "+u'H\u2082SO\u2084'+": %.4e" % phi_gas[k])
        phi_gas[k]       = phi_gas[k] < (csum**3/6.)*np.exp(csum)
        print("csum*dt = %.4e" % csum.values)
    else:
        print("dt = %.f" % dt[k])
        csum             = uptkrate[k]*dt[k]
        csumt            = uptkrate[k]*1800.
        phi_gas[k]       = np.abs(mam_h2so4[k] - ref_h2so4[k]) / ref_h2so4[k]
        for n in np.arange(mode):
            phi_aer[n,k] = np.abs(mam_so4[n,k] - ref_so4[n,k]) / ref_so4[n,k]
            print("Relative error of "+u'SO\u2084'+":   %.4e" % phi_aer[n,k])
            phi_aer[n,k] = phi_aer[n,k] < phi_gas[k] / (np.exp(csumt) - 1)
        print("Relative error of "+u'H\u2082SO\u2084'+": %.4e" % phi_gas[k])
        phi_gas[k]       = phi_gas[k] < np.abs(np.power((np.exp(-csum) + csum**3/6.),1800./dt[k])*np.exp(csumt)-1)
        print("csum*dt = %.4e" % csum.values)
dt = 1
Relative error of SO₄:   1.9820e-11
Relative error of SO₄:   1.5113e-10
Relative error of SO₄:   1.4852e-13
Relative error of SO₄:   3.5504e-08
Relative error of H₂SO₄: 7.4587e-08
csum*dt = 6.2870e-04
dt = 1800
Relative error of SO₄:   1.9885e-12
Relative error of SO₄:   1.5128e-11
Relative error of SO₄:   1.5159e-14
Relative error of SO₄:   2.1337e-09
Relative error of H₂SO₄: 2.4148e-13
csum*dt = 1.1317e-04
In [4]:
# run this cell to issue PASS/FAIL
num = np.count_nonzero(phi_gas)
if  num == pver:
    num = np.count_nonzero(phi_aer)
    if  num == pver * mode:
        print("PASS")
    else:
        print("FAIL")
else:
    print("FAIL")
PASS

References

[1] Liu X., Easter R. C., Ghan S. J. et al., ``Toward a minimal representation of aerosols in climate models: description and evaluation in the Community Atmosphere Model CAM5'', Geoscientific Model Development, vol. 5, number 3, pp. 709--739, 2012. online

[2] P{\"o}schl Ulrich, Canagaratna Manjula, Jayne John T et al., ``Mass accommodation coefficient of H2SO4 vapor on aqueous sulfuric acid surfaces and gaseous diffusion coefficient of H2SO4 in N2/H2O'', The Journal of Physical Chemistry A, vol. 102, number 49, pp. 10082--10089, 1998.

[3] Easter Richard C., Ghan Steven J., Zhang Yang et al., ``MIRAGE: Model description and evaluation of aerosols and trace gases'', Journal of Geophysical Research: Atmospheres, vol. 109, number D20, pp. , 2004.

[4] Daniel Zwillinger, ``CRC standard mathematical tables and formulae'', 2002.