PySDM is a package for simulating the dynamics of population of particles. It is intended to serve as a building block for simulation systems modelling fluid flows involving a dispersed phase, with PySDM being responsible for representation of the dispersed phase. Currently, the development is focused on atmospheric cloud physics applications, in particular on modelling the dynamics of particles immersed in moist air using the particle-based (a.k.a. super-droplet) approach to represent aerosol/cloud/rain microphysics. The package core is a Pythonic high-performance implementation of the Super-Droplet Method (SDM) Monte-Carlo algorithm for representing collisional growth (Shima et al. 2009), hence the name.
PySDM has two alternative parallel number-crunching backends
available: multi-threaded CPU backend based on Numba
and GPU-resident backend built on top of ThrustRTC.
The Numba backend named CPU
is the default, and features multi-threaded parallelism for
multi-core CPUs.
It uses the just-in-time compilation technique based on the LLVM infrastructure.
The ThrustRTC backend named GPU
offers GPU-resident operation of PySDM
leveraging the SIMT
parallelisation model.
Using the GPU
backend requires nVidia hardware and CUDA driver.
For an overview paper on PySDM v1 (and the preferred item to cite if using PySDM), see Bartman et al. 2021 arXiv e-print (submitted to JOSS). For a list of talks and other materials on PySDM, see the project wiki.
PySDM dependencies are: Numpy, Numba, SciPy, Pint, chempy, ThrustRTC and CURandRTC.
To install PySDM using pip
, use: pip install git+https://github.com/atmos-cloud-sim-uj/PySDM.git
.
For development purposes, we suggest cloning the repository and installing it using pip -e
.
Test-time dependencies are listed in the requirements.txt
file.
PySDM examples listed below are hosted in a separate repository and constitute
a separate PySDM_examples
Python package.
The examples have additional dependencies listed in PySDM_examples
package setup.py
file.
Running the examples requires the the PySDM_examples
package to be installed (or within PYTHONPATH
).
Since the examples package includes Jupyter notebooks, the suggested install and launch steps are:
git clone https://github.com/atmos-cloud-sim-uj/PySDM-examples.git
cd PySDM-examples
pip install -e .
jupyter-notebook
- Shima et al. 2009 Fig. 2
(Box model, coalescence only, test case employing Golovin analytical solution) - Berry 1967 Figs. 5, 8 & 10
(Box model, coalescence only, test cases for realistic kernels)
- Arabas & Shima 2017 Fig. 5
(Adiabatic parcel, monodisperse size spectrum activation/deactivation test case) - Yang et al. 2018 Fig. 2:
(Adiabatic parcel, polydisperse size spectrum activation/deactivation test case)
- Kreidenweis et al. 2003 Fig 1:
(Adiabatic parcel, polydisperse size spectrum, aqueous‐phase SO2 oxidation test case) - Jaruga and Pawlowska 2018:
(same test case as above, different numerical settings)
- Shipway & Hill 2012
(Fig. 1 with thermodynamics & condensation only - no particle displacement)
2D kinematic (prescribed-flow) Sc-mimicking aerosol collisional processing (warm-rain) examples (CPU backend only, stay tuned...)
- Arabas et al. 2015 Figs. 8 & 9:
(interactive web-GUI with product selection, parameter sliders and netCDF/plot export buttons) - Bartman et al. 2021 (in preparation):
(default-settings based script generating a netCDF file and loading it subsequently to create the animation below)
In order to depict the PySDM API with a practical example, the following listings provide sample code roughly reproducing the Figure 2 from Shima et al. 2009 paper using PySDM from Python, Julia and Matlab. It is a coalescence-only set-up in which the initial particle size spectrum is exponential and is deterministically sampled to match the condition of each super-droplet having equal initial multiplicity:
Julia (click to expand)
using Pkg
Pkg.add("PyCall")
Pkg.add("Plots")
using PyCall
si = pyimport("PySDM.physics").si
ConstantMultiplicity = pyimport("PySDM.initialisation.spectral_sampling").ConstantMultiplicity
Exponential = pyimport("PySDM.initialisation.spectra").Exponential
n_sd = 2^15
initial_spectrum = Exponential(norm_factor=8.39e12, scale=1.19e5 * si.um^3)
attributes = Dict()
attributes["volume"], attributes["n"] = ConstantMultiplicity(spectrum=initial_spectrum).sample(n_sd)
Matlab (click to expand)
si = py.importlib.import_module('PySDM.physics').si;
ConstantMultiplicity = py.importlib.import_module('PySDM.initialisation.spectral_sampling').ConstantMultiplicity;
Exponential = py.importlib.import_module('PySDM.initialisation.spectra').Exponential;
n_sd = 2^15;
initial_spectrum = Exponential(pyargs(...
'norm_factor', 8.39e12, ...
'scale', 1.19e5 * si.um ^ 3 ...
));
tmp = ConstantMultiplicity(initial_spectrum).sample(int32(n_sd));
attributes = py.dict(pyargs('volume', tmp{1}, 'n', tmp{2}));
Python
from PySDM.physics import si
from PySDM.initialisation.spectral_sampling import ConstantMultiplicity
from PySDM.initialisation.spectra import Exponential
n_sd = 2**15
initial_spectrum = Exponential(norm_factor=8.39e12, scale=1.19e5 * si.um**3)
attributes = {}
attributes['volume'], attributes['n'] = ConstantMultiplicity(initial_spectrum).sample(n_sd)
The key element of the PySDM interface is the Core
class which instances are used to manage the system state and control the simulation.
Instantiation of the Core
class is handled by the Builder
as exemplified below:
Julia (click to expand)
Builder = pyimport("PySDM").Builder
Box = pyimport("PySDM.environments").Box
Coalescence = pyimport("PySDM.dynamics").Coalescence
Golovin = pyimport("PySDM.physics.coalescence_kernels").Golovin
CPU = pyimport("PySDM.backends").CPU
ParticlesVolumeSpectrum = pyimport("PySDM.products.state").ParticlesVolumeSpectrum
builder = Builder(n_sd=n_sd, backend=CPU)
builder.set_environment(Box(dt=1 * si.s, dv=1e6 * si.m^3))
builder.add_dynamic(Coalescence(kernel=Golovin(b=1.5e3 / si.s)))
products = [ParticlesVolumeSpectrum()]
particles = builder.build(attributes, products)
Matlab (click to expand)
Builder = py.importlib.import_module('PySDM').Builder;
Box = py.importlib.import_module('PySDM.environments').Box;
Coalescence = py.importlib.import_module('PySDM.dynamics').Coalescence;
Golovin = py.importlib.import_module('PySDM.physics.coalescence_kernels').Golovin;
CPU = py.importlib.import_module('PySDM.backends').CPU;
ParticlesVolumeSpectrum = py.importlib.import_module('PySDM.products.state').ParticlesVolumeSpectrum;
builder = Builder(pyargs('n_sd', int32(n_sd), 'backend', CPU));
builder.set_environment(Box(pyargs('dt', 1 * si.s, 'dv', 1e6 * si.m ^ 3)));
builder.add_dynamic(Coalescence(pyargs('kernel', Golovin(1.5e3 / si.s))));
products = py.list({ ParticlesVolumeSpectrum() });
particles = builder.build(attributes, products);
Python
from PySDM import Builder
from PySDM.environments import Box
from PySDM.dynamics import Coalescence
from PySDM.physics.coalescence_kernels import Golovin
from PySDM.backends import CPU
from PySDM.products.state import ParticlesVolumeSpectrum
builder = Builder(n_sd=n_sd, backend=CPU)
builder.set_environment(Box(dt=1 * si.s, dv=1e6 * si.m**3))
builder.add_dynamic(Coalescence(kernel=Golovin(b=1.5e3 / si.s)))
products = [ParticlesVolumeSpectrum()]
particles = builder.build(attributes, products)
The backend
argument may be set to CPU
or GPU
what translates to choosing the multi-threaded backend or the
GPU-resident computation mode, respectively.
The employed Box
environment corresponds to a zero-dimensional framework
(particle positions are not considered).
The vectors of particle multiplicities n
and particle volumes v
are
used to initialise super-droplet attributes.
The Coalescence
Monte-Carlo algorithm (Super Droplet Method) is registered as the only
dynamic in the system (other available dynamics representing
condensational growth and particle displacement).
Finally, the build()
method is used to obtain an instance
of Core
which can then be used to control time-stepping and
access simulation state.
The run(nt)
method advances the simulation by nt
timesteps.
In the listing below, its usage is interleaved with plotting logic
which displays a histogram of particle mass distribution
at selected timesteps:
Julia (click to expand)
rho_w = pyimport("PySDM.physics.constants").rho_w
using Plots
radius_bins_edges = 10 .^ range(log10(10*si.um), log10(5e3*si.um), length=32)
for step = 0:1200:3600
particles.run(step - particles.n_steps)
plot!(
radius_bins_edges[1:end-1] / si.um,
particles.products["dv/dlnr"].get(radius_bins_edges) * rho_w / si.g,
linetype=:steppost,
xaxis=:log,
xlabel="particle radius [µm]",
ylabel="dm/dlnr [g/m^3/(unit dr/r)]",
label="t = $step s"
)
end
savefig("plot.svg")
Matlab (click to expand)
rho_w = py.importlib.import_module('PySDM.physics.constants').rho_w;
radius_bins_edges = logspace(log10(10 * si.um), log10(5e3 * si.um), 32);
for step = 0:1200:3600
particles.run(int32(step - particles.n_steps))
x = radius_bins_edges / si.um;
y = particles.products{"dv/dlnr"}.get(py.numpy.array(radius_bins_edges)) * rho_w / si.g;
stairs(...
x(1:end-1), ...
double(py.array.array('d',py.numpy.nditer(y))), ...
'DisplayName', sprintf("t = %d s", step) ...
);
hold on
end
hold off
set(gca,'XScale','log');
xlabel('particle radius [µm]')
ylabel("dm/dlnr [g/m^3/(unit dr/r)]")
legend()
Python
from PySDM.physics.constants import rho_w
from matplotlib import pyplot
import numpy as np
radius_bins_edges = np.logspace(np.log10(10 * si.um), np.log10(5e3 * si.um), num=32)
for step in [0, 1200, 2400, 3600]:
particles.run(step - particles.n_steps)
pyplot.step(x=radius_bins_edges[:-1] / si.um,
y=particles.products['dv/dlnr'].get(radius_bins_edges) * rho_w / si.g,
where='post', label=f"t = {step}s")
pyplot.xscale('log')
pyplot.xlabel('particle radius [µm]')
pyplot.ylabel("dm/dlnr [g/m$^3$/(unit dr/r)]")
pyplot.legend()
pyplot.savefig('readme.svg')
The resultant plot looks as follows:
Julia (click to expand)
using PyCall
using Plots
si = pyimport("PySDM.physics").si
spectral_sampling = pyimport("PySDM.initialisation").spectral_sampling
multiplicities = pyimport("PySDM.initialisation").multiplicities
spectra = pyimport("PySDM.initialisation").spectra
r_wet_init = pyimport("PySDM.initialisation").r_wet_init
CPU = pyimport("PySDM.backends").CPU
AmbientThermodynamics = pyimport("PySDM.dynamics").AmbientThermodynamics
Condensation = pyimport("PySDM.dynamics").Condensation
Parcel = pyimport("PySDM.environments").Parcel
Builder = pyimport("PySDM").Builder
products = pyimport("PySDM.products")
env = Parcel(
dt=.25 * si.s,
mass_of_dry_air=1e3 * si.kg,
p0=1122 * si.hPa,
q0=20 * si.g / si.kg,
T0=300 * si.K,
w= 2.5 * si.m / si.s
)
spectrum=spectra.Lognormal(norm_factor=1e4/si.mg, m_mode=50*si.nm, s_geom=1.4)
kappa = .5 * si.dimensionless
cloud_range = (.5 * si.um, 25 * si.um)
output_interval = 4
output_points = 40
n_sd = 256
builder = Builder(backend=CPU, n_sd=n_sd)
builder.set_environment(env)
builder.add_dynamic(AmbientThermodynamics())
builder.add_dynamic(Condensation(kappa=kappa))
r_dry, specific_concentration = spectral_sampling.Logarithmic(spectrum).sample(n_sd)
r_wet = r_wet_init(r_dry, env, kappa)
attributes = Dict()
attributes["n"] = multiplicities.discretise_n(specific_concentration * env.mass_of_dry_air)
attributes["dry volume"] = builder.formulae.trivia.volume(radius=r_dry)
attributes["volume"] = builder.formulae.trivia.volume(radius=r_wet)
particles = builder.build(attributes, products=[
products.PeakSupersaturation(),
products.CloudDropletEffectiveRadius(radius_range=cloud_range),
products.CloudDropletConcentration(radius_range=cloud_range),
products.WaterMixingRatio(radius_range=cloud_range)
])
cell_id=1
output = Dict("z" => Array{Float32}(undef, output_points+1))
for (_, product) in particles.products
output[product.name] = Array{Float32}(undef, output_points+1)
output[product.name][1] = product.get()[cell_id]
end
output["z"][1] = env.__getitem__("z")[cell_id]
for step = 2:output_points+1
particles.run(steps=output_interval)
for (_, product) in particles.products
output[product.name][step] = product.get()[cell_id]
end
output["z"][step]=env.__getitem__("z")[cell_id]
end
plots = []
for (_, product) in particles.products
append!(plots, [plot(output[product.name], output["z"], ylabel="z [m]", xlabel=product.unit, title=product.name)])
end
plot(plots..., layout=(1,4))
savefig("parcel.svg")
Matlab (click to expand)
si = py.importlib.import_module('PySDM.physics').si;
spectral_sampling = py.importlib.import_module('PySDM.initialisation').spectral_sampling;
multiplicities = py.importlib.import_module('PySDM.initialisation').multiplicities;
spectra = py.importlib.import_module('PySDM.initialisation').spectra;
r_wet_init = py.importlib.import_module('PySDM.initialisation').r_wet_init;
CPU = py.importlib.import_module('PySDM.backends').CPU;
AmbientThermodynamics = py.importlib.import_module('PySDM.dynamics').AmbientThermodynamics;
Condensation = py.importlib.import_module('PySDM.dynamics').Condensation;
Parcel = py.importlib.import_module('PySDM.environments').Parcel;
Builder = py.importlib.import_module('PySDM').Builder;
products = py.importlib.import_module('PySDM.products');
env = Parcel(pyargs( ...
'dt', .25 * si.s, ...
'mass_of_dry_air', 1e3 * si.kg, ...
'p0', 1122 * si.hPa, ...
'q0', 20 * si.g / si.kg, ...
'T0', 300 * si.K, ...
'w', 2.5 * si.m / si.s ...
));
spectrum = spectra.Lognormal(pyargs('norm_factor', 1e4/si.mg, 'm_mode', 50 * si.nm, 's_geom', 1.4));
kappa = .5;
cloud_range = py.tuple({.5 * si.um, 25 * si.um});
output_interval = 1;
output_points = 40;
n_sd = 256;
builder = Builder(pyargs('backend', CPU, 'n_sd', int32(n_sd)));
builder.set_environment(env);
builder.add_dynamic(AmbientThermodynamics())
builder.add_dynamic(Condensation(pyargs('kappa', kappa)))
tmp = spectral_sampling.Logarithmic(spectrum).sample(int32(n_sd));
r_dry = tmp{1};
specific_concentration = tmp{2};
r_wet = r_wet_init(r_dry, env, kappa);
attributes = py.dict(pyargs( ...
'n', multiplicities.discretise_n(specific_concentration * env.mass_of_dry_air), ...
'dry volume', builder.formulae.trivia.volume(r_dry), ...
'volume', builder.formulae.trivia.volume(r_wet) ...
));
particles = builder.build(attributes, py.list({ ...
products.PeakSupersaturation(), ...
products.CloudDropletEffectiveRadius(pyargs('radius_range', cloud_range)), ...
products.CloudDropletConcentration(pyargs('radius_range', cloud_range)), ...
products.WaterMixingRatio(pyargs('radius_range', cloud_range)) ...
}));
cell_id = int32(0);
output_size = [output_points+1, 1 + length(py.list(particles.products.keys()))];
output_types = repelem({'double'}, output_size(2));
output_names = ['z', cellfun(@string, cell(py.list(particles.products.keys())))];
output = table(...
'Size', output_size, ...
'VariableTypes', output_types, ...
'VariableNames', output_names ...
);
for pykey = py.list(keys(particles.products))
get = py.getattr(particles.products{pykey{1}}.get(), '__getitem__');
key = string(pykey{1});
output{1, key} = get(cell_id);
end
get = py.getattr(env, '__getitem__');
zget = py.getattr(get('z'), '__getitem__');
output{1, 'z'} = zget(cell_id);
for i=2:output_points+1
particles.run(pyargs('steps', int32(output_interval)));
for pykey = py.list(keys(particles.products))
get = py.getattr(particles.products{pykey{1}}.get(), '__getitem__');
key = string(pykey{1});
output{i, key} = get(cell_id);
end
output{i, 'z'} = zget(cell_id);
end
i=1;
for pykey = py.list(keys(particles.products))
product = particles.products{pykey{1}};
subplot(1, width(output)-1, i);
plot(output{:, string(pykey{1})}, output.z);
title(string(product.name));
xlabel(string(product.unit));
ylabel('z [m]');
i=i+1;
end
Python
from matplotlib import pyplot
from PySDM.physics import si
from PySDM.initialisation import spectral_sampling, multiplicities, spectra, r_wet_init
from PySDM.backends import CPU
from PySDM.dynamics import AmbientThermodynamics, Condensation
from PySDM.environments import Parcel
from PySDM import Builder, products
env = Parcel(
dt=.25 * si.s,
mass_of_dry_air=1e3 * si.kg,
p0=1122 * si.hPa,
q0=20 * si.g / si.kg,
T0=300 * si.K,
w=2.5 * si.m / si.s
)
spectrum = spectra.Lognormal(norm_factor=1e4/si.mg, m_mode=50*si.nm, s_geom=1.5)
kappa = .5 * si.dimensionless
cloud_range = (.5 * si.um, 25 * si.um)
output_interval = 4
output_points = 40
n_sd = 256
builder = Builder(backend=CPU, n_sd=n_sd)
builder.set_environment(env)
builder.add_dynamic(AmbientThermodynamics())
builder.add_dynamic(Condensation(kappa=kappa))
r_dry, specific_concentration = spectral_sampling.Logarithmic(spectrum).sample(n_sd)
r_wet = r_wet_init(r_dry, env, kappa)
attributes = {
'n': multiplicities.discretise_n(specific_concentration * env.mass_of_dry_air),
'dry volume': builder.formulae.trivia.volume(radius=r_dry),
'volume': builder.formulae.trivia.volume(radius=r_wet)
}
particles = builder.build(attributes, products=[
products.PeakSupersaturation(),
products.CloudDropletEffectiveRadius(radius_range=cloud_range),
products.CloudDropletConcentration(radius_range=cloud_range),
products.WaterMixingRatio(radius_range=cloud_range)
])
cell_id = 0
output = {product.name: [product.get()[cell_id]] for product in particles.products.values()}
output['z'] = [env['z'][cell_id]]
for step in range(output_points):
particles.run(steps=output_interval)
for product in particles.products.values():
output[product.name].append(product.get()[cell_id])
output['z'].append(env['z'][cell_id])
fig, axs = pyplot.subplots(1, len(particles.products), sharey="all")
for i, (key, product) in enumerate(particles.products.items()):
axs[i].plot(output[key], output['z'], marker='.')
axs[i].set_title(product.name)
axs[i].set_xlabel(product.unit)
axs[i].grid()
pyplot.savefig('parcel.svg')
- backends:
- CPU=Numba: multi-threaded CPU backend using LLVM-powered just-in-time compilation
- GPU=ThrustRTC: GPU-resident backend using NVRTC runtime compilation library for CUDA
- initialisation:
- multiplicities: integer-valued discretisation with sanity checks for errors due to type casting
- r_wet_init: kappa-Keohler-based equilibrium in unsaturated conditions (RH=1 used in root-finding above saturation)
- spatial_sampling: pseudorandom sampling using NumPy's default RNG
- spectra: Exponential and Lognormal classes
- spectral_sampling: linear, logarithmic and constant_multiplicity classes
- physics:
- coalescence_kernels (selected)
- constants: physical constants partly imported from SciPy and chempy packages
- dimensional_analysis: tool for enabling dimensional analysis of the code for unit tests (based on pint)
- formulae: physical formulae boosted with Numba (e.g., for initialisation)
- environments:
- Box: bare zero-dimensional framework
- Parcel: zero-dimensional adiabatic parcel framework
- Kinematic1D: single-column time-varying-updraft framework with moisture advection handled by PyMPDATA
- Kinematic2D: two-dimensional single-eddy prescribed-flow framework with moisture and heat advection handled by PyMPDATA
- dynamics:
- Coalescence: SDM implementation with adaptive timestepping
- Condensation: bespoke solver with implicit-in-particle-size integration and adaptive timestepping (Numba only as of now, soon on all backends)
- AqueousChemistry: aqueous-phase chemistry (incl. SO2 oxidation)
- Displacement: includes advection with the flow & sedimentation
- EulerianAdvection: wrapper class for triggering integration in the Eulerian advection solver underlying used by the selected environment
- Attributes (selected):
- Products (selected):
Development of PySDM is supported by the EU through a grant of the Foundation for Polish Science (POIR.04.04.00-00-5E1C/18).
copyright: Jagiellonian University
licence: GPL v3
- https://patents.google.com/patent/US7756693B2
- https://patents.google.com/patent/EP1847939A3
- https://patents.google.com/patent/JP4742387B2
- https://patents.google.com/patent/CN101059821B
- SCALE-SDM (Fortran):
https://github.com/Shima-Lab/SCALE-SDM_BOMEX_Sato2018/blob/master/contrib/SDM/sdm_coalescence.f90 - Pencil Code (Fortran):
https://github.com/pencil-code/pencil-code/blob/master/src/particles_coagulation.f90 - PALM LES (Fortran):
https://palm.muk.uni-hannover.de/trac/browser/palm/trunk/SOURCE/lagrangian_particle_model_mod.f90 - libcloudph++ (C++):
https://github.com/igfuw/libcloudphxx/blob/master/src/impl/particles_impl_coal.ipp - LCM1D (Python)
https://github.com/SimonUnterstrasser/ColumnModel - superdroplet (Cython/Numba/C++11/Fortran 2008/Julia)
https://github.com/darothen/superdroplet - NTLP (FORTRAN)
https://github.com/Folca/NTLP/blob/SuperDroplet/les.F
- PartMC (Fortran):
https://github.com/compdyn/partmc
- pyrcel: https://github.com/darothen/pyrcel
- PyBox: https://github.com/loftytopping/PyBox
- py-cloud-parcel-model: https://github.com/emmasimp/py-cloud-parcel-model