Package PyPartMC

Expand source code
# pylint: disable=invalid-name,wrong-import-position
import os
from collections import namedtuple
from contextlib import contextmanager
from pathlib import Path


# https://github.com/diegoferigo/cmake-build-extension/blob/master/src/cmake_build_extension/__init__.py
@contextmanager
def __build_extension_env():
    cookies = []
    # https://docs.python.org/3/whatsnew/3.8.html#bpo-36085-whatsnew
    if hasattr(os, "add_dll_directory"):
        basepath = os.path.dirname(os.path.abspath(__file__))
        dllspath = os.path.join(basepath, "..")
        os.environ["PATH"] = dllspath + os.pathsep + os.environ["PATH"]
        for path in os.environ.get("PATH", "").split(os.pathsep):
            if path and Path(path).is_absolute() and Path(path).is_dir():
                cookies.append(os.add_dll_directory(path))
    try:
        yield
    finally:
        for cookie in cookies:
            cookie.close()


def __generate_si():
    prefixes = {
        "T": 1e12,
        "G": 1e9,
        "M": 1e6,
        "k": 1e3,
        "h": 1e2,
        "da": 1e1,
        "": 1e0,
        "d": 1e-1,
        "c": 1e-2,
        "m": 1e-3,
        "u": 1e-6,
        "n": 1e-9,
        "p": 1e-12,
    }
    units = {
        "m": 1e0,
        "g": 1e-3,
        "s": 1e0,
        "K": 1e0,
        "Pa": 1e0,
        "mol": 1e0,
        "W": 1e0,
        "J": 1e0,
        "N": 1e0,
    }
    return namedtuple("SI", [prefix + unit for prefix in prefixes for unit in units])(
        **{
            prefix_k + unit_k: prefix_v * unit_v
            for prefix_k, prefix_v in prefixes.items()
            for unit_k, unit_v in units.items()
        }
    )


si = __generate_si()
""" 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! """

with __build_extension_env():
    import _PyPartMC
    from _PyPartMC import *
    from _PyPartMC import __all__ as _PyPartMC_all  # pylint: disable=no-name-in-module
    from _PyPartMC import __version__, __versions_of_build_time_dependencies__

    __all__ = tuple([*_PyPartMC_all, "si"])

Global variables

var si

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!

Functions

def condense_equilib_particle(arg0: _PyPartMC.EnvState, arg1: _PyPartMC.AeroData, arg2: _PyPartMC.AeroParticle)

Determine the water equilibrium state of a single particle.

def condense_equilib_particles(arg0: _PyPartMC.EnvState, arg1: _PyPartMC.AeroData, arg2: _PyPartMC.AeroState)

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

def diam2rad(arg0: float) ‑> float

Convert diameter (m) to radius (m).

def histogram_1d(arg0: _PyPartMC.BinGrid, arg1: List[float], arg2: List[float]) ‑> List[float]

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

def histogram_2d(arg0: _PyPartMC.BinGrid, arg1: List[float], arg2: _PyPartMC.BinGrid, arg3: List[float], arg4: List[float]) ‑> List[List[float]]

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

def input_state(arg0: str, arg1: _PyPartMC.AeroData, arg2: _PyPartMC.AeroState, arg3: _PyPartMC.GasData, arg4: _PyPartMC.GasState, arg5: _PyPartMC.EnvState)

Read current state from netCDF output file.

def loss_rate(arg0: _PyPartMC.Scenario, arg1: float, arg2: float, arg3: _PyPartMC.AeroData, arg4: _PyPartMC.EnvState) ‑> float

Evaluate a loss rate function.

def loss_rate_dry_dep(arg0: float, arg1: float, arg2: _PyPartMC.AeroData, arg3: _PyPartMC.EnvState) ‑> float

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

def output_state(arg0: str, arg1: _PyPartMC.AeroData, arg2: _PyPartMC.AeroState, arg3: _PyPartMC.GasData, arg4: _PyPartMC.GasState, arg5: _PyPartMC.EnvState)

Output current state to netCDF file.

def pow2_above(arg0: int) ‑> int

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

def rad2diam(arg0: float) ‑> float

Convert radius (m) to diameter (m).

def rand_init(arg0: int)

Initializes 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(arg0: float, arg1: float) ‑> float

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

def run_part(arg0: _PyPartMC.Scenario, arg1: _PyPartMC.EnvState, arg2: _PyPartMC.AeroData, arg3: _PyPartMC.AeroState, arg4: _PyPartMC.GasData, arg5: _PyPartMC.GasState, arg6: _PyPartMC.RunPartOpt, arg7: _PyPartMC.CampCore, arg8: _PyPartMC.Photolysis)

Do a particle-resolved Monte Carlo simulation.

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

Do a time block

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

Do a single time step

def sphere_rad2vol(arg0: float) ‑> float

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

def sphere_vol2rad(arg0: float) ‑> float

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

Classes

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.

init(self: _PyPartMC.AeroData, arg0: json) -> None

Ancestors

  • pybind11_builtins.pybind11_object

Instance variables

var densities

Return array of aerosol species densities

var frac_dim

Volume fractal dimension (1)

var kappa

Returns array of aerosol species hygroscopicity parameter kappa

var molecular_weights

Returns array of aerosol species molecular weights

var n_source

Number of aerosol sources

var prime_radius

Radius of primary particles (m)

var sources

return list of source names

var species

returns list of aerosol species names

var vol_fill_factor

Volume filling factor (1)

Methods

def density(self: _PyPartMC.AeroData, arg0: str) ‑> float

Return density of an aerosol species

def diam2vol(self: _PyPartMC.AeroData, arg0: float) ‑> float

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

def rad2vol(self: _PyPartMC.AeroData, arg0: float) ‑> float

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

def spec_by_name(self: _PyPartMC.AeroData, arg0: str) ‑> int

Returns the number of the species in AeroData with the given name

def vol2diam(self: _PyPartMC.AeroData, arg0: float) ‑> float

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

def vol2rad(self: _PyPartMC.AeroData, arg0: float) ‑> float

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

class AeroDist (...)

init(self: _PyPartMC.AeroDist, arg0: _PyPartMC.AeroData, arg1: json) -> None

Ancestors

  • pybind11_builtins.pybind11_object

Instance variables

var n_mode

Number of aerosol modes

var num_conc

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

Methods

def mode(self: _PyPartMC.AeroDist, arg0: int) ‑> _PyPartMC.AeroMode

returns the mode of a given index

class AeroMode (...)

init(self: _PyPartMC.AeroMode, arg0: _PyPartMC.AeroData, arg1: json) -> None

Ancestors

  • pybind11_builtins.pybind11_object

Instance variables

var char_radius

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

var gsd

Geometric standard deviation

var name

Mode name, used to track particle sources

var num_conc

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

var type

Mode type (given by module constants)

var vol_frac

Species fractions by volume

var vol_frac_std

Species fraction standard deviation

Methods

def num_dist(self: _PyPartMC.AeroMode, arg0: _PyPartMC.BinGrid, arg1: _PyPartMC.AeroData) ‑> List[float]

returns the binned number concenration of a mode

def set_sample(self: _PyPartMC.AeroMode, arg0: List[float], arg1: List[float])
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.

init(self: _PyPartMC.AeroParticle, arg0: _PyPartMC.AeroData, arg1: List[float]) -> None

Ancestors

  • pybind11_builtins.pybind11_object

Instance variables

var absorb_cross_sect

Absorption cross-section (m^-2).

var asymmetry

Asymmetry parameter (1).

var density

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

var diameter

Total diameter of the particle (m).

var dry_diameter

Total dry diameter of the particle (m).

var dry_radius

Total dry radius of the particle (m).

var dry_volume

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

var greatest_create_time

Last time a constituent was created (s).

var id

Unique ID number.

var least_create_time

First time a constituent was created (s).

var mass

Total mass of the particle (kg).

var moles

Total moles in the particle (1).

var radius

Total radius of the particle (m).

var refract_core

Refractive index of the core (1).

var refract_shell

Refractive index of the shell (1).

var scatter_cross_sect

Scattering cross-section (m^-2).

var solute_kappa

Returns the average of the solute kappas (1).

var sources

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

var species_masses

Mass of all species in the particle (kg).

var volume

Total volume of the particle (m^3).

var volumes

Constituent species volumes (m^3)

Methods

def approx_crit_rel_humid(self: _PyPartMC.AeroParticle, arg0: _PyPartMC.EnvState) ‑> float

Returns the approximate critical relative humidity (1).

def coagulate(self: _PyPartMC.AeroParticle, arg0: _PyPartMC.AeroParticle) ‑> _PyPartMC.AeroParticle

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

def crit_diameter(self: _PyPartMC.AeroParticle, arg0: _PyPartMC.EnvState) ‑> float

Returns the critical diameter (m).

def crit_rel_humid(self: _PyPartMC.AeroParticle, arg0: _PyPartMC.EnvState) ‑> float

Returns the critical relative humidity (1).

def mobility_diameter(self: _PyPartMC.AeroParticle, arg0: _PyPartMC.EnvState) ‑> float

Mobility diameter of the particle (m).

def set_vols(self: _PyPartMC.AeroParticle, arg0: List[float])

Sets the aerosol particle volumes.

def species_mass(*args, **kwargs)

Overloaded function.

  1. species_mass(self: _PyPartMC.AeroParticle, arg0: int) -> float

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

  1. species_mass(self: _PyPartMC.AeroParticle, arg0: str) -> float

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

def species_volume(*args, **kwargs)

Overloaded function.

  1. species_volume(self: _PyPartMC.AeroParticle, arg0: int) -> float

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

  1. species_volume(self: _PyPartMC.AeroParticle, arg0: str) -> float

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

def zero(self: _PyPartMC.AeroParticle)

Resets an aero_particle to be zero.

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.

init(self: _PyPartMC.AeroState, arg0: _PyPartMC.AeroData, arg1: float, arg2: str) -> None

Ancestors

  • pybind11_builtins.pybind11_object

Instance variables

var dry_diameters

returns the dry diameter of each particle in the population

var ids

returns the IDs of all particles.

var num_concs

returns the number concentration of each particle in the population

var total_mass_conc

returns the total mass concentration of the population

var total_num_conc

returns the total number concentration of the population

Methods

def add(self: _PyPartMC.AeroState, arg0: _PyPartMC.AeroState)

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

def add_particle(self: _PyPartMC.AeroState, arg0: _PyPartMC.AeroParticle)

add a particle to an AeroState

def add_particles(self: _PyPartMC.AeroState, arg0: _PyPartMC.AeroState)

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

def bin_average_comp(self: _PyPartMC.AeroState, arg0: _PyPartMC.BinGrid)

composition-averages population using BinGrid

def copy_weight(self: _PyPartMC.AeroState, arg0: _PyPartMC.AeroState)

copy weighting from another AeroState

def crit_rel_humids(self: _PyPartMC.AeroState, arg0: _PyPartMC.EnvState) ‑> List[float]

returns the critical relative humidity of each particle in the population

def diameters(self: _PyPartMC.AeroState, include: Optional[List[str]] = None, exclude: Optional[List[str]] = None) ‑> List[float]

returns the diameter of each particle in the population

def dist_sample(self: _PyPartMC.AeroState, AeroDist: _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 make_dry(self: _PyPartMC.AeroState)

Make all particles dry (water set to zero).

def masses(self: _PyPartMC.AeroState, include: Optional[List[str]] = None, exclude: Optional[List[str]] = None) ‑> List[float]

returns the total mass of each particle in the population

def mixing_state(self: _PyPartMC.AeroState, include: Optional[List[str]] = None, exclude: Optional[List[str]] = None, group: Optional[List[str]] = None) ‑> Tuple[float, float, float]

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

def mobility_diameters(self: _PyPartMC.AeroState, arg0: _PyPartMC.EnvState) ‑> List[float]

returns the mobility diameter of each particle in the population

def particle(self: _PyPartMC.AeroState, arg0: int) ‑> _PyPartMC.AeroParticle

returns the particle of a given index

def rand_particle(self: _PyPartMC.AeroState) ‑> _PyPartMC.AeroParticle

returns a random particle from the population

def remove_particle(self: _PyPartMC.AeroState, arg0: int)

remove particle of a given index

def sample(self: _PyPartMC.AeroState, arg0: _PyPartMC.AeroState, arg1: float)

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(self: _PyPartMC.AeroState, arg0: _PyPartMC.AeroState, arg1: float)

!> 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 volumes(self: _PyPartMC.AeroState, include: Optional[List[str]] = None, exclude: Optional[List[str]] = None) ‑> List[float]

returns the volume of each particle in the population

def zero(self: _PyPartMC.AeroState)

remove all particles from an AeroState

class BinGrid (...)

init(self: _PyPartMC.BinGrid, arg0: float, arg1: str, arg2: float, arg3: float) -> None

Ancestors

  • pybind11_builtins.pybind11_object

Instance variables

var centers

Bin centers

var edges

Bin edges

class CampCore (...)

An interface between PartMC and the CAMP

init(self: _PyPartMC.CampCore) -> None

Ancestors

  • pybind11_builtins.pybind11_object
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.

init(self: _PyPartMC.EnvState, arg0: json) -> None

Ancestors

  • pybind11_builtins.pybind11_object

Instance variables

var air_density

Air density (kg m^{-3})

var elapsed_time

returns time since start_time (s).

var height

Box height (m)

var pressure

Ambient pressure (Pa)

var rh

returns the current relative humidity of the environment state

var start_time

returns start time (s since 00:00 UTC on start_day)

var temp

returns the current temperature of the environment state

Methods

def set_temperature(self: _PyPartMC.EnvState, arg0: float)

sets the temperature of the environment state

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).

init(self: _PyPartMC.GasData, arg0: tuple) -> None

Ancestors

  • pybind11_builtins.pybind11_object

Instance variables

var n_spec

(arg0: _PyPartMC.GasData) -> int

var species

returns list of gas species names

Methods

def spec_by_name(self: _PyPartMC.GasData, arg0: str) ‑> int

returns the number of the species in gas with the given name

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.

init(self: _PyPartMC.GasState, arg0: _PyPartMC.GasData) -> None

instantiates and initializes based on GasData

Ancestors

  • pybind11_builtins.pybind11_object

Instance variables

var mix_rats

provides access (read of write) to the array of mixing ratios

var n_spec

returns number of gas species

Methods

def mix_rat(self: _PyPartMC.GasState, arg0: str) ‑> float

returns the mixing ratio of a gas species

def set_size(self: _PyPartMC.GasState)

sets the GasState to the size of GasData

class Photolysis (...)

PartMC interface to a photolysis module

init(self: _PyPartMC.Photolysis) -> None

Ancestors

  • pybind11_builtins.pybind11_object
class RunPartOpt (...)

Options controlling the execution of run_part().

init(self: _PyPartMC.RunPartOpt, arg0: json) -> None

Ancestors

  • pybind11_builtins.pybind11_object

Instance variables

var del_t

time step

var t_max

total simulation time

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.

init(self: _PyPartMC.Scenario, arg0: _PyPartMC.GasData, arg1: _PyPartMC.AeroData, arg2: json) -> None

instantiates and initializes from a JSON object

Ancestors

  • pybind11_builtins.pybind11_object

Instance variables

var aero_dilution_n_times

returns the number of times specified for dilution

var aero_dilution_rate

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

var aero_dilution_time

Aerosol-background dilution set-point times (s)

var aero_emissions_n_times

returns the number of times specified for emissions

var aero_emissions_rate_scale

Aerosol emission rate scales at set-points (1)

var aero_emissions_time

(arg0: _PyPartMC.Scenario) -> List[float]

Methods

def aero_background(self: _PyPartMC.Scenario, arg0: _PyPartMC.AeroData, arg1: int) ‑> _PyPartMC.AeroDist

returns aero_background AeroDists at a given index

def aero_emissions(self: _PyPartMC.Scenario, arg0: _PyPartMC.AeroData, arg1: int) ‑> _PyPartMC.AeroDist

returns aero_emissions AeroDists at a given index

def init_env_state(self: _PyPartMC.Scenario, arg0: _PyPartMC.EnvState, arg1: float)

initializes the EnvState