PyPartMC

logo

PyPartMC

PyPartMC is a Python interface to PartMC, a particle-resolved Monte-Carlo code for atmospheric aerosol simulation. Development of PyPartMC has been intended to remove limitations to the use of Fortran-implemented PartMC. PyPartMC facilitates the dissemination of computational research results by streamlining independent execution of PartMC simulations (also during peer-review processes). Additionally, the ability to easily package examples, simple simulations, and results in a web-based notebook allows PyPartMC to support the efforts of many members of the scientific community, including researchers, instructors, and students, with nominal software and hardware requirements.

Documentation of PyPartMC is hosted at https://open-atmos.github.io/PyPartMC. PyPartMC is implemented in C++ and it also constitutes a C++ API to the PartMC Fortran internals. The Python API can facilitate using PartMC from other environments - see, e.g., Julia and Matlab examples below.

For an outline of the project, rationale, architecture, and features, refer to: D'Aquino et al., 2024 (SoftwareX) (please cite if PyPartMC is used in your research). For a list of talks and other relevant resources, please see project Wiki. If interested in contributing to PyPartMC, please have a look a the notes for developers.

US Funding PL Funding

License: GPL v3 Copyright tests+pypi API docs codecov DOI PyPI version Project Status: Active – The project has reached a stable, usable state and is being actively developed. pyOpenSci Peer-Reviewed

Python 3 Linux OK macOS OK Windows OK Jupyter

Installation

Using the command-line pip tool (also applies to conda environments)

pip install PyPartMC

Note that, depending on the environment (OS, hardware, Python version), the pip-install invocation may either trigger a download of a pre-compiled binary, or trigger compilation of PyPartMC. In the latter case, a Fortran compiler and some development tools includiong CMake, m4 and perl are required (while all non-Python dependencies are included in the PyPartMC source archive). In both cases, all Python dependencies will be resolved by pip.

In a Jupyter notebook cell (also on Colab or jupyter-hub instances)

! pip install PyPartMC
import PyPartMC

Jupyter notebooks with examples

Note: clicking the badges below redirects to cloud-computing platforms. The mybinder.org links allow anonymous execution, Google Colab requires logging in with a Google account, ARM JupyerHub requires logging in with an ARM account (and directing Jupyter to a particular notebook within the examples folder).

The example notebooks feature additional dependencies that can be installed with:

pip install PyPartMC[examples]
  • Urban plume scenario demo (as in PartMC):
    View notebook Open In Colab Binder ARM JupyterHub
  • Dry-Wet Particle Size Equilibration with PartMC and PySDM:
    View notebook Open In Colab Binder ARM JupyterHub Voila
  • Simulation output processing example (loading from netCDF files using PyPartMC):
    View notebook Open In Colab Binder ARM JupyterHub
  • Optical properties calculation using external Python package (PyMieScatt):
    View notebook Open In Colab Binder ARM JupyterHub
  • Cloud parcel example featuring supersaturation-evolution-coupled CCN activation and drop growth:
    View notebook Open In Colab Binder ARM JupyterHub
  • Immersion freezing example:
    View notebook Open In Colab Binder ARM JupyterHub
  • Coagulation model intercomparison for additive (Golovin) kernel with: PyPartMC, PySDM, Droplets.jl and dustpy:
    View notebook Open In Colab Binder ARM JupyterHub
  • Particle simulation with multiphase chemistry handled using CAMP (without coagulation):
    View notebook Open In Colab Binder ARM JupyterHub

Features

  • works on Linux, macOS and Windows (compatibility assured with CI builds)
  • hassle-free installation using pip (prior PartMC installation not needed)
  • works out of the box on mybinder.org, Google Colab and alike
  • ships with a set of examples maintained in a form of Jupyter notebooks
  • Pythonic API (but retaining PartMC jargon) incl. Python GC deallocation of Fortran objects
  • specification of parameters using native Python datatypes (lists, dicts) in place of PartMC spec files
  • code snippets in README depicting how to use PyPartMC from Julia and Matlab (also executed on CI)
  • auto-generated API docs on the web
  • support for [de]serialization of selected wrapped structures using JSON
  • based on unmodified PartMC code
  • does not use or require shell or any pre-installed libraries
  • aiming at 100% unit test coverage

Usage examples

The listings below depict how the identical task of randomly sampling particles from an aerosol size distribution in PartMC can be done in different programming languages.

For a Fortran equivalent of the Python, Julia and Matlab programs below, see the readme_fortran folder.

Python

import numpy as np

import PyPartMC as ppmc
from PyPartMC import si

aero_data = ppmc.AeroData((
    #      [density, ions in solution, molecular weight, kappa, abifm_m, abifm_c]
    {"OC": [1000 *si.kg/si.m**3, 0, 1e-3 *si.kg/si.mol, 0.001, 0, 0]},
    {"BC": [1800 *si.kg/si.m**3, 0, 1e-3 *si.kg/si.mol, 0, 0 , 0]},
))

aero_dist = ppmc.AeroDist(
    aero_data,
    [{
        "cooking": {
            "mass_frac": [{"OC": [1]}],
            "diam_type": "geometric",
            "mode_type": "log_normal",
            "num_conc": 3200 / si.cm**3,
            "geom_mean_diam": 8.64 * si.nm,
            "log10_geom_std_dev": 0.28,
        },
        "diesel": {
            "mass_frac": [{"OC": [0.3]}, {"BC": [0.7]}],
            "diam_type": "geometric",
            "mode_type": "log_normal",
            "num_conc": 2900 / si.cm**3,
            "geom_mean_diam": 50 * si.nm,
            "log10_geom_std_dev": 0.24,
        }
    }],
)

n_part = 100
aero_state = ppmc.AeroState(aero_data, n_part, "nummass_source")
aero_state.dist_sample(aero_dist)
print(np.dot(aero_state.masses(), aero_state.num_concs), "# kg/m3")

Julia (using PyCall.jl)

using Pkg
Pkg.add("PyCall")

using PyCall
ppmc = pyimport("PyPartMC")
si = ppmc["si"]

aero_data = ppmc.AeroData((
  #       (density, ions in solution, molecular weight, kappa, abifm_m, abifm_c)
  Dict("OC"=>(1000 * si.kg/si.m^3, 0, 1e-3 * si.kg/si.mol, 0.001, 0, 0)),
  Dict("BC"=>(1800 * si.kg/si.m^3, 0, 1e-3 * si.kg/si.mol, 0, 0, 0))
))

aero_dist = ppmc.AeroDist(aero_data, (
  Dict( 
    "cooking" => Dict(
      "mass_frac" => (Dict("OC" => (1,)),),
      "diam_type" => "geometric",
      "mode_type" => "log_normal",
      "num_conc" => 3200 / si.cm^3,
      "geom_mean_diam" => 8.64 * si.nm,
      "log10_geom_std_dev" => .28,
    ),
    "diesel" => Dict(
      "mass_frac" => (Dict("OC" => (.3,)), Dict("BC" => (.7,))),
      "diam_type" => "geometric",
      "mode_type" => "log_normal",
      "num_conc" => 2900 / si.cm^3,
      "geom_mean_diam" => 50 * si.nm,
      "log10_geom_std_dev" => .24,
    )
  ),
))

n_part = 100
aero_state = ppmc.AeroState(aero_data, n_part, "nummass_source")
aero_state.dist_sample(aero_dist)
print(aero_state.masses()'aero_state.num_concs, "# kg/m3")

Matlab (using Matlab's built-in Python interface)

notes (see the PyPartMC Matlab CI workflow for an example on how to achieve it on Ubuntu 20):

  • Matlab ships with convenience copies of C, C++ and Fortran runtime libraries which are dlopened() by default; one way to make PyPartMC OK with it is to [pip-]install by compiling from source using the very same version of GCC that Matlab borrowed these libraries from (e.g., GCC 9 for Matlab R2022a, etc);
  • Matlab needs to use the same Python interpretter/venv as the pip invocation used to install PyPartMC;
ppmc = py.importlib.import_module('PyPartMC');
si = py.importlib.import_module('PyPartMC').si;

aero_data = ppmc.AeroData(py.tuple({ ...
  py.dict(pyargs("OC", py.tuple({1000 * si.kg/si.m^3, 0, 1e-3 * si.kg/si.mol, 0.001, 0, 0}))), ...
  py.dict(pyargs("BC", py.tuple({1800 * si.kg/si.m^3, 0, 1e-3 * si.kg/si.mol, 0, 0, 0}))) ...
}));

aero_dist = ppmc.AeroDist(aero_data, py.tuple({ ...
  py.dict(pyargs( ...
    "cooking", py.dict(pyargs( ...
      "mass_frac", py.tuple({py.dict(pyargs("OC", py.tuple({1})))}), ...
      "diam_type", "geometric", ...
      "mode_type", "log_normal", ...
      "num_conc", 3200 / si.cm^3, ...
      "geom_mean_diam", 8.64 * si.nm, ...
      "log10_geom_std_dev", .28 ...
    )), ...
    "diesel", py.dict(pyargs( ...
      "mass_frac", py.tuple({ ...
        py.dict(pyargs("OC", py.tuple({.3}))), ...
        py.dict(pyargs("BC", py.tuple({.7}))), ...
      }), ...
      "diam_type", "geometric", ...
      "mode_type", "log_normal", ...
      "num_conc", 2900 / si.cm^3, ...
      "geom_mean_diam", 50 * si.nm, ...
      "log10_geom_std_dev", .24 ...
    )) ...
  )) ...
}));

n_part = 100;
aero_state = ppmc.AeroState(aero_data, n_part, "nummass_source");
aero_state.dist_sample(aero_dist);
masses = cell(aero_state.masses());
num_concs = cell(aero_state.num_concs);
fprintf('%g # kg/m3\n', dot([masses{:}], [num_concs{:}]))

usage in other projects

PyPartMC is used within the test workflow of the PySDM project.

Other packages with relevant feature scope

  • aerosolGDEFoam: OpenFOAM CFD-coupled aerosol dynamics including nucleation, coagulation, and surface growth
  • AIOMFAC and AIOMFAC-web: Fortran-implemented aerosol thermodynamic model for calculation of activity coefficients in organic-inorganic mixtures – from simple binary solutions to complex multicomponent systems
  • DustPy: Python package for modelling dust evolution in protoplanetary disks (differences: focus on astrophysical applications vs. atmospheric aerosol)
  • multilayerpy: kinetic multi-layer model for aerosol particles and films
  • PyBox: aerosol simulation model featuring gas and particle chamistry (differences: PyBox focuses on chemical mechanisms; PyPartMC is an interface to PartMC which focuses on physics - e.g., collisions of aerosol particles - while chemical processes are handled with external software, e.g., CAMP or MOSAIC)
  • PyCHAM: CHemistry with Aerosol Microphysics in Python Box Model for modelling of indoor environments, including aerosol chambers
  • PySDM: particle-based Monte-Carlo aerosol-cloud simulation package (differences: PySDM focuses on growth and breakup processes relevant to cloud droplets; PyPartMC focuses on processes relevant to air pollutants and their chemical and physical transformations)
  • SSH-aerosol: C++/Fortran package for simulating evolution of primary and secondary atmospheric aerosols

FAQ

  • Q: How to install PyPartMC with MOSAIC enabled?
    A: Installation can be done using pip, however, pip needs to be instructed not to use binary packages available at pypi.org but rather to compile from source (pip will download the source from pip.org), and the path to compiled MOSAIC library needs to be provided at compile-time; the following command should convey it:
MOSAIC_HOME=<<PATH_TO_MOSAIC_LIB>> pip install --force-reinstall --no-binary=PyPartMC PyPartMC
  • Q: Why pip install PyPartMC triggers compilation on my brand new Apple machine, while it quickly downloads and installs binary packages when executed on older Macs, Windows or Linux?
    A: We are providing binary wheels on PyPI for Apple-silicon (arm64) machines for selected macOS version made available by Github. In case the macOS version you are using is newer, compilation from source is triggered.

  • Q: Why some of the constructors expect data to be passed as lists of single-entry dictionaries instead of multi-element dictionaries?
    A: This is intentional and related with PartMC relying on the order of elements within spec-file input; while Python dictionaries preserve ordering (insertion order), JSON format does not, and we intend to make these data structures safe to be [de]serialized using JSON.

  • Q: How to check the version of PartMC that PyPartMC was compiled against?
    A: Version numbers of compile-time dependencies of PyPartMC, including PartMC, can be accessed as follows:

import PyPartMC
PyPartMC.__versions_of_build_time_dependencies__['PartMC']
  • Q: Why m4 and perl are required at compile time?
    A: PyPartMC includes parts of netCDF and HDF5 codebases which depend on m4 and perl, respectively, for generating source files before compilation.

Troubleshooting

Common installation issues

error: [Errno 2] No such file or directory: 'cmake'

Try rerunning after installing CMake, e.g., using apt-get install cmake (Ubuntu/Debian), brew install cmake (homebrew on macOS) or using MSYS2 on Windows.

No CMAKE_Fortran_COMPILER could be found.

Try installing a Fortran compiler (e.g., brew reinstall gcc with Homebrew on macOS or using MSYS2 on Windows).

Could not find NC_M4 using the following names: m4, m4.exe

Try installing m4 (e.g., using MSYS2 on Windows).

Acknowledgement and citations

We would greatly appreciate citation of the PartMC model description paper (Riemer et al., 2009) and the PyPartMC description paper (D’Aquino et al., 2024) if PyPartMC was used in your study. The citations are:

  • Riemer, N., M. West, R. A. Zaveri, R. C. Easter: Simulating the evolution of soot mixing-state with a particle-resolved aerosol model
    J. Geophys. Res., 114, D09202, 2009, DOI: 10.1029/2008JD011073
  • D’Aquino, Z., S. Arabas, J. H. Curtis, A. Vaishnav, N. Riemer, M. West: PyPartMC: A pythonic interfact to a particle-resolved, Monte Carlo aerosol simulation framework
    SoftwareX, 25, 101613, 2024, DOI: 10.1016/j.softx.2023.101613

The following paragraph provides a more substantial description of PartMC (text released into the public domain and can be freely copied by anyone for any purpose):

PartMC is a stochastic, particle-resolved aerosol box model. It tracks the composition of many computational particles (104 to 106) within a well-mixed air volume, each represented by a composition vector that evolves based on physical and chemical processes. The physical processes—including Brownian coagulation, new particle formation, emissions, dilution, and deposition—are simulated using a stochastic Monte Carlo approach via a Poisson process while chemical processes are simulated deterministically for each computational particle. The weighted flow algorithm (DeVille, Riemer, and West, 2011, 2019) enhances efficiency and reduces ensemble variance. Detailed numerical methods are described in Riemer et al. (2009), DeVille et al. (2011, 2019), and Curtis et al. (2016). PartMC is open-source under the GNU GPL v2 and available at github.com/compdyn/partmc.

References:

  • Curtis, J. H., M. D. Michelotti, N. Riemer, M. T. Heath, M. West: Accelerated simulation of stochastic particle removal processes in particle-resolved aerosol models, J. Computational Phys., 322, 21-32, 2016, DOI: 10.1016/j.jcp.2016.06.029
  • DeVille, L., N. Riemer, M. West, Convergence of a generalized weighted flow algorithm for stochastic particle coagulation, J. Computational Dynamics, 6, 69-94, 2019, DOI: 10.3934/jcd.2019003
  • DeVille, R. E. L., N. Riemer, M. West, The Weighted Flow Algorithm (WFA) for stochastic particle coagulation, J. Computational Phys., 230, 8427-8451, 2011, DOI: 10.1016/j.jcp.2011.07.027
  • Riemer, N., M. West, R. A. Zaveri, R. C. Easter, Simulating the evolution of soot mixing-state with a particle-resolved aerosol model, J. Geophys. Res., 114, D09202, 2009., DOI: 10.1029/2008JD011073

Credits

PyPartMC:

authors: PyPartMC developers
funding: US Department of Energy Atmospheric System Research programme, Polish National Science Centre
copyright: University of Illinois at Urbana-Champaign
licence: GPL v3

PartMC:

authors: Nicole Riemer, Matthew West, Jeff Curtis et al.
licence: GPL v2 or later

  1"""
  2.. include::../../README.md
  3"""
  4
  5#  pylint: disable=invalid-name
  6import importlib.metadata
  7import inspect
  8
  9# pylint: disable=invalid-name,wrong-import-position
 10import os
 11from collections import namedtuple
 12from contextlib import contextmanager
 13from pathlib import Path
 14
 15import nanobind
 16
 17
 18# https://github.com/diegoferigo/cmake-build-extension/blob/master/src/cmake_build_extension/__init__.py
 19@contextmanager
 20def __build_extension_env():
 21    cookies = []
 22    # https://docs.python.org/3/whatsnew/3.8.html#bpo-36085-whatsnew
 23    if hasattr(os, "add_dll_directory"):
 24        basepath = os.path.dirname(os.path.abspath(__file__))
 25        dllspath = os.path.join(basepath, "..")
 26        os.environ["PATH"] = dllspath + os.pathsep + os.environ["PATH"]
 27        for path in os.environ.get("PATH", "").split(os.pathsep):
 28            if path and Path(path).is_absolute() and Path(path).is_dir():
 29                cookies.append(os.add_dll_directory(path))
 30    try:
 31        yield
 32    finally:
 33        for cookie in cookies:
 34            cookie.close()
 35
 36
 37def __generate_si():
 38    prefixes = {
 39        "T": 1e12,
 40        "G": 1e9,
 41        "M": 1e6,
 42        "k": 1e3,
 43        "h": 1e2,
 44        "da": 1e1,
 45        "": 1e0,
 46        "d": 1e-1,
 47        "c": 1e-2,
 48        "m": 1e-3,
 49        "u": 1e-6,
 50        "n": 1e-9,
 51        "p": 1e-12,
 52    }
 53    units = {
 54        "m": 1e0,
 55        "g": 1e-3,
 56        "s": 1e0,
 57        "K": 1e0,
 58        "Pa": 1e0,
 59        "mol": 1e0,
 60        "W": 1e0,
 61        "J": 1e0,
 62        "N": 1e0,
 63    }
 64    return namedtuple("SI", [prefix + unit for prefix in prefixes for unit in units])(
 65        **{
 66            prefix_k + unit_k: prefix_v * unit_v
 67            for prefix_k, prefix_v in prefixes.items()
 68            for unit_k, unit_v in units.items()
 69        }
 70    )
 71
 72
 73si = __generate_si()
 74""" a utility namedtuple aimed at clrifying physics-related code by providing
 75    SI-prefix-aware unit multipliers, resulting in e.g.: `p = 1000 * si.hPa`
 76    notation. Note: no dimensional analysis is done! """
 77
 78from ._PyPartMC import *  # pylint: disable=import-error
 79from ._PyPartMC import (  # pylint: disable=import-error
 80    __versions_of_build_time_dependencies__,
 81)
 82
 83# Hacky workaround for missing docs in pdoc auto-generated documentation.
 84# After the switch to nanobind, the docs became empty despite "__doc__" being
 85# accessible in all of PyPartMC's objects. The code below manually populates
 86# the "__all__" atrribute of the package. Additionally, functions in the generated
 87# docs would be listed as nanobind objects with no additional documentation.
 88# To solve that, dummy functions of the same name are created, and their "__doc__"
 89# attribute is manually set to the "original" objects' "__doc__"
 90if os.getenv("PDOC_GENERATE_PYPARTMC_DOCS") == "1":
 91    all_items = []
 92    for name, obj in inspect.getmembers(
 93        _PyPartMC  # pylint: disable=undefined-variable
 94    ):
 95        if callable(obj):
 96            if not inspect.isclass(obj):
 97                exec(name + " = lambda : 0")  # pylint: disable=exec-used
 98                temp = "_PyPartMC." + name + ".__doc__"
 99                setattr(eval(name), "__doc__", eval(temp))  # pylint: disable=eval-used
100            all_items.append(name)
101
102    __all__ = tuple([*all_items, "si"])
103
104__version__ = importlib.metadata.version(__package__)
105
106# workaround for MATLAB bindings
107# pylint: disable=undefined-variable
108setattr(nanobind, "nb_type_0", type(_PyPartMC.AeroData))
class AeroBinned:

Aerosol number and volume distributions stored per size bin. These quantities are densities both in volume (per m^3) and in radius (per log_width).

AeroBinned()

__init__(self, arg: PyPartMC._PyPartMC.AeroData, /) -> None __init__(self, arg0: PyPartMC._PyPartMC.AeroData, arg1: PyPartMC._PyPartMC.BinGrid, /) -> None

num_conc

(self) -> [std::valarray]

Returns the number concentration of each bin (#/m^3/log_width)

vol_conc

(self) -> [std::valarray]

Returns the volume concentration per bin per species (m^3/m^3/log_width)

def add_aero_dist(unknown):

add_aero_dist(self, arg0: PyPartMC._PyPartMC.BinGrid, arg1: PyPartMC._PyPartMC.AeroDist, /) -> None

Add an AeroDist to an AeroBinned.

class AeroData:

Aerosol material properties and associated data.

The data in this structure is constant, as it represents physical quantities that cannot change over time.

Each aerosol species is identified by an index i = 1,...,aero_data_n_spec(aero_data). Then \c name(i) is the name of that species, \c density(i) is its density, etc. The ordering of the species is arbitrary and should not be relied upon (currently it is the order in the species data file). The only exception is that it is possible to find out which species is water from the \c i_water variable.

The names of the aerosol species are not important to PartMC, as only the material properties are used. The names are used for input and output, and also for communication with MOSAIC. For the MOSAIC interface to work correctly the species must be named the same, but without the \c _a suffix.

AeroData()

__init__(self, arg: PyPartMC._PyPartMC.CampCore, /) -> None __init__(self, arg: [JSON], /) -> None

def spec_by_name(unknown):

spec_by_name(self, arg: str, /) -> int

Return the number of the species in AeroData with the given name.

n_source

(self) -> int

Number of aerosol sources.

sources

(self) -> tuple

Return list of source names.

frac_dim

(self) -> float

Volume fractal dimension (1).

vol_fill_factor

(self) -> float

Volume filling factor (1).

prime_radius

(self) -> float

Radius of primary particles (m).

densities

(self) -> [std::valarray]

Return array of aerosol species densities.

kappa

(self) -> [std::valarray]

Returns array of aerosol species hygroscopicity parameter kappa.

molecular_weights

(self) -> [std::valarray]

Return array of aerosol species molecular weights.

def density(unknown):

density(self, arg: str, /) -> float

Return density of an aerosol species.

def rad2vol(unknown):

rad2vol(self, arg: float, /) -> float

Convert geometric radius (m) to mass-equivalent volume (m^3).

def vol2rad(unknown):

vol2rad(self, arg: float, /) -> float

Convert mass-equivalent volume (m^3) to geometric radius (m).

def diam2vol(unknown):

diam2vol(self, arg: float, /) -> float

Convert geometric diameter (m) to mass-equivalent volume (m^3).

def vol2diam(unknown):

vol2diam(self, arg: float, /) -> float

Convert mass-equivalent volume (m^3) to geometric diameter (m).

species

(self) -> tuple

Return list of aerosol species names.

class AeroDist:
AeroDist()

__init__(self, arg0: PyPartMC._PyPartMC.AeroData, arg1: [JSON], /) -> None

n_mode

(self) -> int

Number of aerosol modes in the distribution.

num_conc

(self) -> float

Total number concentration of a distribution (#/m^3).

def mode(unknown):

mode(self, arg: int, /) -> PyPartMC._PyPartMC.AeroMode

Return the mode of a given index.

class AeroMode:

An aerosol size distribution mode.

AeroMode()

__init__(self, arg0: PyPartMC._PyPartMC.AeroData, arg1: [JSON], /) -> None

num_conc

(self) -> float

Provide access (read or write) to the total number concentration of a mode.

def num_dist(unknown):

num_dist(self, arg0: PyPartMC._PyPartMC.BinGrid, arg1: PyPartMC._PyPartMC.AeroData, /) -> [std::valarray]

Return the binned number concenration of a mode.

vol_frac

(self) -> [std::valarray]

Species fractions by volume.

vol_frac_std

(self) -> [std::valarray]

Species fraction standard deviation.

char_radius

(self) -> float

Characteristic radius, with meaning dependent on mode type (m).

gsd

(self) -> float

Geometric standard deviation.

def set_sample(unknown):

set_sample(self, arg0: [std::valarray], arg1: [std::valarray], /) -> None

sample_num_conc

(self) -> [std::valarray]

Sample bin number concentrations (m^{-3}).

sample_radius

(self) -> [std::valarray]

Sample bin radii (m).

type

(self) -> str

Mode type (given by module constants).

name

(self) -> str

Mode name, used to track particle sources.

class AeroParticle:

Single aerosol particle data structure.

The \c vol array stores the total volumes of the different species that make up the particle. This array must have length equal to aero_data%%n_spec, so that \c vol(i) is the volume (in m^3) of the i'th aerosol species.

AeroParticle()

__init__(self, arg0: PyPartMC._PyPartMC.AeroData, arg1: [std::valarray], /) -> None

volumes

(self) -> [std::valarray]

Constituent species volumes (m^3)

volume

(self) -> float

Total volume of the particle (m^3).

def species_volume(unknown):

species_volume(self, arg: int, /) -> float species_volume(self, arg: str, /) -> float

Volume of a single species in the particle (m^3).

dry_volume

(self) -> float

Total dry volume of the particle (m^3).

radius

(self) -> float

Total radius of the particle (m).

dry_radius

(self) -> float

Total dry radius of the particle (m).

diameter

(self) -> float

Total diameter of the particle (m).

dry_diameter

(self) -> float

Total dry diameter of the particle (m).

mass

(self) -> float

Total mass of the particle (kg).

def species_mass(unknown):

species_mass(self, arg: int, /) -> float species_mass(self, arg: str, /) -> float

Mass of a single species in the particle (kg).

species_masses

(self) -> [std::valarray]

Mass of all species in the particle (kg).

solute_kappa

(self) -> float

Return the average of the solute kappas (1).

moles

(self) -> float

Total moles in the particle (1).

absorb_cross_sect

(self) -> float

Absorption cross-section (m^-2).

scatter_cross_sect

(self) -> float

Scattering cross-section (m^-2).

asymmetry

(self) -> float

Asymmetry parameter (1).

refract_shell

(self) -> complex

Refractive index of the shell (1).

refract_core

(self) -> complex

Refractive index of the core (1).

sources

(self) -> [std::valarray]

Number of original particles from each source that coagulated to form particle.

least_create_time

(self) -> float

First time a constituent was created (s).

greatest_create_time

(self) -> float

Last time a constituent was created (s).

id

(self) -> int

Unique ID number.

is_frozen

(self) -> int

Frozen status - particle is ice if 1.

def mobility_diameter(unknown):

mobility_diameter(self, arg: PyPartMC._PyPartMC.EnvState, /) -> float

Mobility diameter of the particle (m).

density

(self) -> float

Average density of the particle (kg/m^3)

def approx_crit_rel_humid(unknown):

approx_crit_rel_humid(self, arg: PyPartMC._PyPartMC.EnvState, /) -> float

Return the approximate critical relative humidity (1).

def crit_rel_humid(unknown):

crit_rel_humid(self, arg: PyPartMC._PyPartMC.EnvState, /) -> float

Return the critical relative humidity (1).

def crit_diameter(unknown):

crit_diameter(self, arg: PyPartMC._PyPartMC.EnvState, /) -> float

Return the critical diameter (m).

def coagulate(unknown):

coagulate(self, arg: PyPartMC._PyPartMC.AeroParticle, /) -> PyPartMC._PyPartMC.AeroParticle

Coagulate two particles together to make a new one. The new particle will not have its ID set.

def zero(unknown):

zero(self) -> None

Reset an aero_particle to be zero.

def set_vols(unknown):

set_vols(self, arg: [std::valarray], /) -> None

Set the aerosol particle volumes.

class AeroState:

The current collection of aerosol particles.

The particles in \c aero_state_t are stored in a single flat array (the \c apa data member), with a sorting into size bins and weight groups/classes possibly stored in the \c aero_sorted data member (if \c valid_sort is true).

Every time we remove particles we keep track of the particle ID and the action performed in the aero_info_array_t structure. This is typically cleared each time we output data to disk.

AeroState()

__init__(self, arg0: PyPartMC._PyPartMC.AeroData, arg1: float, arg2: str, /) -> None __init__(self, arg0: PyPartMC._PyPartMC.AeroData, arg1: float, arg2: str, arg3: PyPartMC._PyPartMC.CampCore, /) -> None

total_num_conc

(self) -> float

Return the total number concentration of the population.

total_mass_conc

(self) -> float

Return the total mass concentration of the population.

frozen_fraction

(self) -> float

Return the fraction of ice particles in the population.

num_concs

(self) -> [std::valarray]

Return the number concentration of each particle in the population.

def masses(unknown):

masses(self, include: collections.abc.Sequence[str] | None | None = None, exclude: collections.abc.Sequence[str] | None | None = None) -> [std::valarray]

Return the total mass of each particle in the population.

def volumes(unknown):

volumes(self, include: collections.abc.Sequence[str] | None | None = None, exclude: collections.abc.Sequence[str] | None | None = None) -> [std::valarray]

Return the volume of each particle in the population.

dry_diameters

(self) -> [std::valarray]

Return the dry diameter of each particle in the population.

def mobility_diameters(unknown):

mobility_diameters(self, arg: PyPartMC._PyPartMC.EnvState, /) -> [std::valarray]

Return the mobility diameter of each particle in the population.

def diameters(unknown):

diameters(self, include: collections.abc.Sequence[str] | None | None = None, exclude: collections.abc.Sequence[str] | None | None = None) -> [std::valarray]

Return the diameter of each particle in the population.

def crit_rel_humids(unknown):

crit_rel_humids(self, arg: PyPartMC._PyPartMC.EnvState, /) -> [std::valarray]

Return the critical relative humidity of each particle in the population.

def make_dry(unknown):

make_dry(self) -> None

Make all particles dry (water set to zero).

ids

(self) -> [std::valarray]

Return the IDs of all particles.

def mixing_state(unknown):

mixing_state(self, include: collections.abc.Sequence[str] | None | None = None, exclude: collections.abc.Sequence[str] | None | None = None, group: collections.abc.Sequence[str] | None | None = None) -> tuple[float, float, float]

returns the mixing state parameters (d_alpha, d_gamma, chi) of the population

def bin_average_comp(unknown):

bin_average_comp(self, arg: PyPartMC._PyPartMC.BinGrid, /) -> None

Composition-averages population using BinGrid.

def particle(unknown):

particle(self, arg: int, /) -> PyPartMC._PyPartMC.AeroParticle

Return the particle of a given index.

def rand_particle(unknown):

rand_particle(self) -> PyPartMC._PyPartMC.AeroParticle

Return a random particle from the population.

def dist_sample(unknown):

dist_sample(self, AeroDist: PyPartMC._PyPartMC.AeroDist, sample_prop: float = 1.0, create_time: float = 0.0, allow_doubling: bool = True, allow_halving: bool = True) -> int

sample particles for AeroState from an AeroDist

def add_particle(unknown):

add_particle(self, arg: PyPartMC._PyPartMC.AeroParticle, /) -> None

Add a particle to an AeroState.

def add(unknown):

add(self, arg: PyPartMC._PyPartMC.AeroState, /) -> None

aero_state += aero_state_delta, including combining the weights, so the new concentration is the weighted average of the two concentrations.

def add_particles(unknown):

add_particles(self, arg: PyPartMC._PyPartMC.AeroState, /) -> None

aero_state += aero_state_delta, with the weight left unchanged so the new concentration is the sum of the two concentrations.

def sample(unknown):

sample(self, arg0: PyPartMC._PyPartMC.AeroState, arg1: float, /) -> None

Generates a random sample by removing particles from aero_state_from and adding them to aero_state_to, transfering weight as well. This is the equivalent of aero_state_add().

def sample_particles(unknown):

sample_particles(self, arg0: PyPartMC._PyPartMC.AeroState, arg1: float, /) -> None

!> Generates a random sample by removing particles from aero_state_from and adding them to aero_state_to, which must be already allocated (and should have its weight set).

       None of the weights are altered by this sampling, making this the
       equivalent of aero_state_add_particles().
def copy_weight(unknown):

copy_weight(self, arg: PyPartMC._PyPartMC.AeroState, /) -> None

Copy weighting from another AeroState.

def remove_particle(unknown):

remove_particle(self, arg: int, /) -> None

Remove particle of a given index.

def zero(unknown):

zero(self) -> None

Remove all particles from an AeroState.

class BinGrid:

1D grid, either logarithmic or linear.

BinGrid()

__init__(self, arg0: int, arg1: str, arg2: float, arg3: float, /) -> None

edges

(self) -> [std::valarray]

Return bin edges.

centers

(self) -> [std::valarray]

Return bin centers.

widths

(self) -> [std::valarray]

Return bin widths.

class CampCore:

An interface between PartMC and the CAMP

CampCore()

__init__(self) -> None __init__(self, arg: str, /) -> None

class EnvState:

Current environment state.

All quantities are instantaneous, describing the state at a particular instant of time. Constant data and other data not associated with the current environment state is stored in scenario_t.

EnvState()

__init__(self, arg: [JSON], /) -> None

def set_temperature(unknown):

set_temperature(self, arg: float, /) -> None

Set the temperature of the environment state.

temp

(self) -> float

Return the current temperature of the environment state.

rh

(self) -> float

Return the current relative humidity of the environment state.

elapsed_time

(self) -> float

Return time since start_time (s).

start_time

(self) -> float

Return the simulation start time (s since 00:00 UTC on start_day).

height

(self) -> float

Box height (m).

pressure

(self) -> float

Ambient pressure (Pa).

air_density

(self) -> float

Air density (kg m^{-3}).

additive_kernel_coefficient

(self) -> float

Scaling coefficient for additive coagulation kernel.

class GasData:

Constant gas data.

Each gas species is identified by an integer \c i between 1 and \c gas_data_n_spec(gas_data). Species \c i has name \c gas_data%%name(i). The variable gas data describing the current mixing ratios is stored in the gas_state_t structure, so the mixing ratio of species \c i is gas_state%%mix_rat(i).

GasData()

__init__(self, arg: PyPartMC._PyPartMC.CampCore, /) -> None __init__(self, arg: tuple, /) -> None

n_spec

(self) -> int

def spec_by_name(unknown):

spec_by_name(self, arg: str, /) -> int

Return the number of the species in gas with the given name.

species

(self) -> tuple

Return list of gas species names.

class GasState:

Current state of the gas mixing ratios in the system.

The gas species are defined by the gas_data_t structure, so that \c gas_state%%mix_rat(i) is the current mixing ratio of the gas with name \c gas_data%%name(i), etc.

By convention, if gas_state_is_allocated() return \c .false., then the gas_state is treated as zero for all operations on it. This will be the case for new \c gas_state_t structures.

GasState()

__init__(self, arg: PyPartMC._PyPartMC.GasData, /) -> None

Instantiate and initialize based on GasData.

n_spec

(self) -> int

Return number of gas species.

def set_size(unknown):

set_size(self) -> None

Set the GasState to the size of GasData.

def mix_rat(unknown):

mix_rat(self, arg: str, /) -> float

Return the mixing ratio of a gas species.

mix_rats

(self) -> [std::valarray]

Provide access (read or write) to the array of mixing ratios.

class Photolysis:

PartMC interface to a photolysis module

Photolysis()

__init__(self) -> None __init__(self, arg: PyPartMC._PyPartMC.CampCore, /) -> None

class RunExactOpt:

Options controlling the execution of run_exact().

RunExactOpt()

__init__(self, arg0: [JSON], arg1: PyPartMC._PyPartMC.EnvState, /) -> None

t_max

(self) -> float

Total simulation time.

class RunPartOpt:

Options controlling the execution of run_part().

RunPartOpt()

__init__(self, arg: [JSON], /) -> None

t_max

(self) -> float

Total simulation time.

del_t

(self) -> float

Time step.

class RunSectOpt:

Options controlling the execution of run_sect().

RunSectOpt()

__init__(self, arg0: [JSON], arg1: PyPartMC._PyPartMC.EnvState, /) -> None

t_max

(self) -> float

Total simulation time.

del_t

(self) -> float

Time step.

class Scenario:

This is everything needed to drive the scenario being simulated.

The temperature, pressure, emissions and background states are profiles prescribed as functions of time by giving a number of times and the corresponding data. Simple data such as temperature and pressure is linearly interpolated between times, with constant interpolation outside of the range of times. Gases and aerosols are interpolated with gas_state_interp_1d() and aero_dist_interp_1d(), respectively.

Scenario()

__init__(self, arg0: PyPartMC._PyPartMC.GasData, arg1: PyPartMC._PyPartMC.AeroData, arg2: [JSON], /) -> None

instantiates and initializes from a JSON object

def init_env_state(unknown):

init_env_state(self, arg0: PyPartMC._PyPartMC.EnvState, arg1: float, /) -> None

Initialize the EnvState.

def aero_emissions(unknown):

aero_emissions(self, arg0: PyPartMC._PyPartMC.AeroData, arg1: int, /) -> PyPartMC._PyPartMC.AeroDist

Return aero_emissions AeroDists at a given index.

aero_emissions_n_times

(self) -> int

Return the number of times specified for emissions.

aero_emissions_rate_scale

(self) -> [std::valarray]

Aerosol emission rate scales at set-points (1).

aero_emissions_time

(self) -> [std::valarray]

def aero_background(unknown):

aero_background(self, arg0: PyPartMC._PyPartMC.AeroData, arg1: int, /) -> PyPartMC._PyPartMC.AeroDist

Return aero_background AeroDists at a given index.

aero_dilution_n_times

(self) -> int

Return the number of times specified for dilution.

aero_dilution_rate

(self) -> [std::valarray]

Aerosol-background dilution rates at set-points (s^{-1}).

aero_dilution_time

(self) -> [std::valarray]

Aerosol-background dilution set-point times (s).

def condense_equilib_particle():

condense_equilib_particle(arg0: PyPartMC._PyPartMC.EnvState, arg1: PyPartMC._PyPartMC.AeroData, arg2: PyPartMC._PyPartMC.AeroParticle, /) -> None

Determine the water equilibrium state of a single particle.

def condense_equilib_particles():

condense_equilib_particles(arg0: PyPartMC._PyPartMC.EnvState, arg1: PyPartMC._PyPartMC.AeroData, arg2: PyPartMC._PyPartMC.AeroState, /) -> None

Call condense_equilib_particle() on each particle in the aerosol to ensure that every particle has its water content in equilibrium.

def diam2rad():

diam2rad(arg: float, /) -> float

Convert diameter (m) to radius (m).

def histogram_1d():

histogram_1d(arg0: PyPartMC._PyPartMC.BinGrid, arg1: [std::valarray], arg2: [std::valarray], /) -> [std::valarray]

Return a 1D histogram with of the given weighted data, scaled by the bin sizes.

def histogram_2d():

histogram_2d(arg0: PyPartMC._PyPartMC.BinGrid, arg1: [std::valarray], arg2: PyPartMC._PyPartMC.BinGrid, arg3: [std::valarray], arg4: [std::valarray], /) -> list[list[float]]

Return a 2D histogram with of the given weighted data, scaled by the bin sizes.

def input_exact():
def input_sectional():
def input_state():
def loss_rate():

loss_rate(arg0: PyPartMC._PyPartMC.Scenario, arg1: float, arg2: float, arg3: PyPartMC._PyPartMC.AeroData, arg4: PyPartMC._PyPartMC.EnvState, /) -> float

Evaluate a loss rate function.

def loss_rate_dry_dep():

loss_rate_dry_dep(arg0: float, arg1: float, arg2: PyPartMC._PyPartMC.AeroData, arg3: PyPartMC._PyPartMC.EnvState, /) -> float

Compute and return the dry deposition rate for a given particle.

def output_state():

output_state(arg0: str, arg1: PyPartMC._PyPartMC.AeroData, arg2: PyPartMC._PyPartMC.AeroState, arg3: PyPartMC._PyPartMC.GasData, arg4: PyPartMC._PyPartMC.GasState, arg5: PyPartMC._PyPartMC.EnvState, /) -> None

Output current state to netCDF file.

def pow2_above():

pow2_above(arg: int, /) -> int

Return the least power-of-2 that is at least equal to n.

def rad2diam():

rad2diam(arg: float, /) -> float

Convert radius (m) to diameter (m).

def rand_init():

rand_init(arg: int, /) -> None

Initialize the random number generator to the state defined by the given seed. If the seed is 0 then a seed is auto-generated from the current time.

def rand_normal():

rand_normal(arg0: float, arg1: float, /) -> float

Generate a normally distributed random number with the given mean and standard deviation.

def run_exact():
def run_part():
def run_part_timeblock():

run_part_timeblock(arg0: PyPartMC._PyPartMC.Scenario, arg1: PyPartMC._PyPartMC.EnvState, arg2: PyPartMC._PyPartMC.AeroData, arg3: PyPartMC._PyPartMC.AeroState, arg4: PyPartMC._PyPartMC.GasData, arg5: PyPartMC._PyPartMC.GasState, arg6: PyPartMC._PyPartMC.RunPartOpt, arg7: PyPartMC._PyPartMC.CampCore, arg8: PyPartMC._PyPartMC.Photolysis, arg9: int, arg10: int, arg11: float, arg12: float, arg13: float, arg14: int, /) -> tuple[float, float, int]

Do multiple time steps over a block of time.

def run_part_timestep():

run_part_timestep(arg0: PyPartMC._PyPartMC.Scenario, arg1: PyPartMC._PyPartMC.EnvState, arg2: PyPartMC._PyPartMC.AeroData, arg3: PyPartMC._PyPartMC.AeroState, arg4: PyPartMC._PyPartMC.GasData, arg5: PyPartMC._PyPartMC.GasState, arg6: PyPartMC._PyPartMC.RunPartOpt, arg7: PyPartMC._PyPartMC.CampCore, arg8: PyPartMC._PyPartMC.Photolysis, arg9: int, arg10: float, arg11: float, arg12: float, arg13: int, /) -> tuple[float, float, int]

Do a single time step.

def run_sect():
def sphere_rad2vol():

sphere_rad2vol(arg: float, /) -> float

Convert geometric radius (m) to mass-equivalent volume for spherical particles.

def sphere_vol2rad():

sphere_vol2rad(arg: float, /) -> float

Convert mass-equivalent volume (m^3) to geometric radius (m) for spherical particles.

si = SI(Tm=1000000000000.0, Tg=1000000000.0, Ts=1000000000000.0, TK=1000000000000.0, TPa=1000000000000.0, Tmol=1000000000000.0, TW=1000000000000.0, TJ=1000000000000.0, TN=1000000000000.0, Gm=1000000000.0, Gg=1000000.0, Gs=1000000000.0, GK=1000000000.0, GPa=1000000000.0, Gmol=1000000000.0, GW=1000000000.0, GJ=1000000000.0, GN=1000000000.0, Mm=1000000.0, Mg=1000.0, Ms=1000000.0, MK=1000000.0, MPa=1000000.0, Mmol=1000000.0, MW=1000000.0, MJ=1000000.0, MN=1000000.0, km=1000.0, kg=1.0, ks=1000.0, kK=1000.0, kPa=1000.0, kmol=1000.0, kW=1000.0, kJ=1000.0, kN=1000.0, hm=100.0, hg=0.1, hs=100.0, hK=100.0, hPa=100.0, hmol=100.0, hW=100.0, hJ=100.0, hN=100.0, dam=10.0, dag=0.01, das=10.0, daK=10.0, daPa=10.0, damol=10.0, daW=10.0, daJ=10.0, daN=10.0, m=1.0, g=0.001, s=1.0, K=1.0, Pa=1.0, mol=1.0, W=1.0, J=1.0, N=1.0, dm=0.1, dg=0.0001, ds=0.1, dK=0.1, dPa=0.1, dmol=0.1, dW=0.1, dJ=0.1, dN=0.1, cm=0.01, cg=1e-05, cs=0.01, cK=0.01, cPa=0.01, cmol=0.01, cW=0.01, cJ=0.01, cN=0.01, mm=0.001, mg=1e-06, ms=0.001, mK=0.001, mPa=0.001, mmol=0.001, mW=0.001, mJ=0.001, mN=0.001, um=1e-06, ug=1e-09, us=1e-06, uK=1e-06, uPa=1e-06, umol=1e-06, uW=1e-06, uJ=1e-06, uN=1e-06, nm=1e-09, ng=1.0000000000000002e-12, ns=1e-09, nK=1e-09, nPa=1e-09, nmol=1e-09, nW=1e-09, nJ=1e-09, nN=1e-09, pm=1e-12, pg=1e-15, ps=1e-12, pK=1e-12, pPa=1e-12, pmol=1e-12, pW=1e-12, pJ=1e-12, pN=1e-12)

a utility namedtuple aimed at clrifying physics-related code by providing SI-prefix-aware unit multipliers, resulting in e.g.: p = 1000 * si.hPa notation. Note: no dimensional analysis is done!