.

Non-equilibrium cloud formation

CloudMicrophysics.MicrophysicsNonEqModule
Non-equilibrium bulk microphysics scheme, which includes:
  • condensation and evaporation of cloud liquid water and deposition and sublimation of cloud ice (relaxation to equilibrium)
source
CloudMicrophysics.MicrophysicsNonEq.τ_relaxFunction
τ_relax(liquid)
τ_relax(ice)
  • liquid or ice - a type for cloud liquid water or ice

Returns the relaxation timescale for condensation and evaporation of cloud liquid water or the relaxation timescale for sublimation and deposition of cloud ice.

source
CloudMicrophysics.MicrophysicsNonEq.conv_q_vap_to_q_liq_iceFunction
conv_q_vap_to_q_liq_ice(liquid, q_sat, q)
conv_q_vap_to_q_liq_ice(ice, q_sat, q)
  • liquid or ice - a struct with cloud water or ice free parameters
  • q_sat - PhasePartition at equilibrium
  • q - current PhasePartition

Returns the cloud water tendency due to condensation and evaporation or cloud ice tendency due to sublimation and vapor deposition. The tendency is obtained assuming a relaxation to equilibrium with a constant timescale and is based on the difference between specific humidity in equilibrium at the current temperature and the current cloud condensate.

source
CloudMicrophysics.MicrophysicsNonEq.conv_q_vap_to_q_liq_ice_MM2015Function
conv_q_vap_to_q_liq_ice_MM2015(liquid, tps, q, ρ, T)
conv_q_vap_to_q_liq_ice_MM2015(ice, tps, q, ρ, T)
  • liquid or ice - a struct with cloud water or ice free parameters
  • tps - thermodynamics parameters struct
  • q - current PhasePartition
  • ρ - air density [kg/m3]
  • T - air temperature [K]

Returns the cloud water tendency due to condensation and evaporation or cloud ice tendency due to sublimation and vapor deposition. The formulation is based on Morrison and Grabowski 2008 and Morrison and Milbrandt 2015

source
CloudMicrophysics.MicrophysicsNonEq.terminal_velocityFunction
terminal_velocity(sediment, vel, ρ, q)
  • sediment - a struct with sedimentation type (cloud liquid or ice)
  • vel - a struct with terminal velocity parameters
  • ρₐ - air density
  • q - cloud liquid or ice specific content

Returns the mass weighted average terminal velocity assuming a monodisperse size distribution with prescribed number concentration. The fall velocity of individual particles is parameterized following Chen et. al 2022, DOI: 10.1016/j.atmosres.2022.106171

source

0-moment precipitation microphysics

CloudMicrophysics.Microphysics0MModule
Microphysics0M

Zero-moment bulk microphysics scheme that instantly removes moisture above certain threshold. This is equivalent to instanteneous conversion of cloud condensate into precipitation and precipitation fallout with infinite terminal velocity.

source
CloudMicrophysics.Microphysics0M.remove_precipitationFunction
remove_precipitation(params_0M, q; q_vap_sat)
  • params_0M - a struct with 0-moment parameters
  • q - current PhasePartition
  • q_vap_sat - specific humidity at saturation

Returns the q_tot tendency due to the removal of precipitation. The tendency is obtained assuming a relaxation with a constant timescale to a state with precipitable water removed. The threshold for when to remove q_tot is defined either by the condensate specific content or supersaturation. The thresholds and the relaxation timescale are defined in ClimaParams.

source

1-moment precipitation microphysics

CloudMicrophysics.Microphysics1MModule
One-moment bulk microphysics scheme, which includes:
  • terminal velocity of precipitation
  • autoconversion of cloud liquid water into rain and of cloud ice into snow
  • accretion due to collisions between categories of condensed species
  • evaporation and sublimation of hydrometeors
  • melting of snow into rain
source
CloudMicrophysics.Microphysics1M.get_v0Function
get_v0(v, ρ)
  • v - a struct with bulk 1-moment terminal velocity parameters
  • ρ - air density (only for Rain)

Returns the proportionality coefficient in terminal velocity(r/r0).

source
CloudMicrophysics.Microphysics1M.get_n0Function
get_n0(pdf, q_sno, ρ)
  • pdf - a struct with parameters for snow, ice, and rain size distribution
  • q_sno - snow specific content (only for Snow)
  • ρ - air density (only for Snow)

Returns the intercept parameter of the assumed Marshall-Palmer distribution

source
CloudMicrophysics.Microphysics1M.lambda_inverseFunction
lambda_inverse(pdf, mass, q, ρ)
  • pdf, mass - structs with particle size distribution and mass parameters
  • q - specific content of rain, ice or snow
  • ρ - air density

Returns the inverse of rate parameter of the assumed size distribution of particles (rain drops, ice crystals, snow crystals). The value is clipped at 1e-8.

source
CloudMicrophysics.Microphysics1M.terminal_velocityFunction
terminal_velocity(precip, vel, ρ, q)
  • precip - a struct with precipitation type (rain or snow)
  • vel - a struct with terminal velocity parameters
  • ρ - air density
  • q - rain or snow specific content

Returns the mass weighted average terminal velocity assuming a Marshall-Palmer (1948) distribution of particles. Fall velocity of individual rain drops is parameterized:

  • assuming an empirical power-law relations for velocity == Blk1MVelType
  • following Chen et. al 2022, DOI: 10.1016/j.atmosres.2022.106171, for velocity == Chen2022VelType
source
CloudMicrophysics.Microphysics1M.conv_q_liq_to_q_raiFunction
conv_q_liq_to_q_rai(acnv, q_liq, smooth_transition)
  • acnv - 1M autoconversion parameters
  • q_liq - liquid water specific content
  • smooth_transition - a flag to switch on smoothing

Returns the q_rai tendency due to collisions between cloud droplets (autoconversion), parametrized following Kessler (1995).

source
CloudMicrophysics.Microphysics1M.conv_q_ice_to_q_sno_no_supersatFunction
conv_q_ice_to_q_sno_no_supersat(acnv, q_ice, smooth_transition)
  • acnv - 1M autoconversion parameters
  • q_ice - cloud ice specific content
  • smooth_transition - a flag to switch on smoothing

Returns the q_sno tendency due to autoconversion from ice. This is a simplified version of a snow autoconversion rate that can be used in simulations where there is no supersaturation (for example in TC.jl when using saturation adjustment).

source
CloudMicrophysics.Microphysics1M.conv_q_ice_to_q_snoFunction
conv_q_ice_to_q_sno(ice, aps, tps, q, ρ, T)
  • ice - a struct with ice parameters
  • aps - a struct with air properties
  • tps - a struct with thermodynamics parameters
  • q - phase partition
  • ρ - air density
  • T - air temperature

Returns the q_sno tendency due to autoconversion from ice. Parameterized following Harrington et al. (1996) and Kaul et al. (2015).

source
CloudMicrophysics.Microphysics1M.accretionFunction
accretion(cloud, precip, vel, ce, q_clo, q_pre, ρ)
  • cloud - type for cloud water or cloud ice
  • precip - type for rain or snow
  • vel - a struct with terminal velocity parameters
  • ce - collision efficiency parameters
  • q_clo - cloud water or cloud ice specific content
  • q_pre - rain water or snow specific content
  • ρ - air density

Returns the source of precipitating water (rain or snow) due to collisions with cloud water (liquid or ice).

source
CloudMicrophysics.Microphysics1M.accretion_rain_sinkFunction
accretion_rain_sink(rain, ice, vel, ce, q_ice, q_rai, ρ)
  • rain - rain type parameters
  • ice - ice type parameters
  • vel - terminal velocity parameters for rain
  • ce - collision efficiency parameters
  • q_ice - cloud ice specific content
  • q_rai - rain water specific content
  • ρ - air density

Returns the sink of rain water (partial source of snow) due to collisions with cloud ice.

source
CloudMicrophysics.Microphysics1M.accretion_snow_rainFunction
accretion_snow_rain(ce, type_i, type_j, blk1m_type_i, blk1m_type_j, q_i, q_j, ρ)
  • ce - collision efficiency parameters
  • i - snow for temperatures below freezing or rain for temperatures above freezing
  • j - rain for temperatures below freezing or snow for temperatures above freezing
  • type_i, type_j - a type for snow or rain
  • blk1mveltype_ti, blk1mveltype_tj - 1M terminal velocity parameters
  • q_ - specific content of snow or rain
  • ρ - air density

Returns the accretion rate between rain and snow. Collisions between rain and snow result in snow at temperatures below freezing and in rain at temperatures above freezing.

source
CloudMicrophysics.Microphysics1M.evaporation_sublimationFunction
evaporation_sublimation(rain, vel, aps, tps, q, q_rai, ρ, T)
evaporation_sublimation(snow, vel, aps, tps, q, q_sno, ρ, T)
  • rain - a struct with rain parameters
  • snow - a struct with snow parameters
  • vel - a struct with terminal velocity parameters
  • aps - a struct with air parameters
  • tps - a struct with thermodynamics parameters
  • q - phase partition
  • q_rai - rain specific content
  • q_sno - snow specific content
  • ρ - air density
  • T - air temperature

Returns the tendency due to rain evaporation or snow sublimation.

source
CloudMicrophysics.Microphysics1M.snow_meltFunction
snow_melt(snow, vel, aps, tps, q_sno, ρ, T)
  • snow - snow parameters
  • vel - terminal velocity parameters
  • aps - air properties
  • tps - thermodynamics parameters
  • q_sno - snow water specific content
  • ρ - air density
  • T - air temperature

Returns the tendency due to snow melt.

source

2-moment precipitation microphysics

CloudMicrophysics.Microphysics2MModule

Double-moment bulk microphysics parametrizations including:

  • autoconversion, accretion, self-collection, breakup, mean terminal velocity of raindrops and rain evaporation rates from Seifert and Beheng 2006,
  • additional double-moment bulk microphysics autoconversion and accretion rates from: Khairoutdinov and Kogan 2000, Beheng 1994, Tripoli and Cotton 1980, and Liu and Daum 2004.
source
CloudMicrophysics.Microphysics2M.pdf_cloud_parametersFunction
pdf_cloud_parameters(pdf_c, qₗ, ρₐ, Nₗ)
  • pdf_c - a struct with SB2006 cloud droplets size distribution parameters

  • qₗ - cloud water specific content

  • ρₐ - air density

  • Nₗ cloud droplet number concentration

    Returns the mean mass of cloud droplets xc [μg], a multiplier needed to convert xc to base SI units χ [-], and the two cloud droplet mass distribution parameters Ac [m^-3, μg^-3] and Bc [μg^-1], It also returns four cloud droplet size distribution parameters Cc [(1/mm3)^4], Ec [1/mm3], ϕc [-], ψc [-] that are consistent with the assumed cloud droplet mass distribution.

source
CloudMicrophysics.Microphysics2M.pdf_rain_parametersFunction
pdf_rain_parameters(pdf_r, qᵣ, ρₐ, Nᵣ)
  • pdf_r - a struct with SB2006 raindrops size distribution parameters

  • qᵣ - rain water specific content

  • ρₐ - air density

  • Nᵣ rain drop number concentration

    Returns the parameters of the rain drop diameter distribution λr [1/m], N₀r [1/m4], optionally limited within prescribed ranges. Also returns the mean mass of rain drops xr [kg], and the two rain drop mass distribution parameters Ar [1/m3 * (1/kg)^(νr+1)], Br [(1/kg)^μr].

source
CloudMicrophysics.Microphysics2M.size_distributionFunction
size_distribution(pdf, D, q, ρ, N)
  • pdf - struct containing size distribution parameters of cloud or rain

  • D - particle size (i.e. maximum dimension of particle)

  • q - cloud or rain water specific content

  • ρₐ - density of air

  • N - cloud or rain water number concentration

    Returns the size distribution value for rain or cloud particles for a given particle size

source
CloudMicrophysics.Microphysics2M.get_size_distribution_boundFunction
get_size_distribution_bound(pdf, q, N, ρₐ, tolerance)
  • pdf_r - struct containing size distribution parameters for cloud or rain

  • q - mass mixing ratio of cloud or rain water

  • N - number mixing ratio of cloud or rain

  • ρₐ - density of air

  • tolerance - tolerance for integration error

    Returns Dmax value such that (1 - tolerance) = 1/N * ∫ N'(D) dD from 0 to Dmax. All inputs and output Dmax are in base SI units. For rain size distribution Dmax is obtained analytically. For cloud size distribution D_max is calculated through a linear approximation of the bounds from numerical solutions.

source
CloudMicrophysics.Microphysics2M.autoconversionFunction
autoconversion(scheme, q_liq, q_rai, ρ, N_liq)
  • acnv, pdf_c - structs with autoconversion and cloud size distribution parameters
  • q_liq - cloud water specific content
  • q_rai - rain water specific content
  • ρ - air density
  • N_liq - cloud droplet number density

Returns a LiqRaiRates object containing q_liq, N_liq, q_rai, N_rai tendencies due to collisions between cloud droplets (autoconversion) for scheme == SB2006Type

source
CloudMicrophysics.Microphysics2M.accretionFunction
accretion(scheme, q_liq, q_rai, ρ, N_liq)
  • scheme - type for 2-moment accretion parameterization
  • q_liq - cloud water specific content
  • q_rai - rain water specific content
  • ρ - air density
  • N_liq - cloud droplet number density

Returns a LiqRaiRates object containing q_liq, N_liq, q_rai, N_rai tendencies due to collisions between raindrops and cloud droplets (accretion) for scheme == SB2006Type

source
accretion(accretion_scheme, q_liq, q_rai, ρ)
  • accretion_scheme - type for 2-moment rain accretion parameterization
  • q_liq - cloud water specific content
  • q_rai - rain water specific content
  • ρ - air density (for KK2000Type and Beheng1994Type)

Returns the accretion rate of rain, parametrized following

  • Khairoutdinov and Kogan (2000) for scheme == KK2000Type
  • Beheng (1994) for scheme == B1994Type
  • Tripoli and Cotton (1980) for scheme == TC1980Type
source
CloudMicrophysics.Microphysics2M.liquid_self_collectionFunction
liquid_self_collection(scheme, q_liq, ρ, dN_liq_dt_au)
  • scheme - type for 2-moment liquid self-collection parameterization
  • q_liq - cloud water specific content
  • ρ - air density
  • dN_liq_dt_au - rate of change of cloud droplets number density due to autoconversion

Returns the cloud droplets number density tendency due to collisions of cloud droplets that produce larger cloud droplets (self-collection) for scheme == SB2006Type

source
CloudMicrophysics.Microphysics2M.autoconversion_and_liquid_self_collectionFunction
autoconversion_and_liquid_self_collection(scheme, q_liq, q_rai, ρ, N_liq)
  • scheme - type for 2-moment rain autoconversion parameterization
  • q_liq - cloud water specific content
  • q_rai - rain water specific content
  • ρ - air density
  • N_liq - cloud droplet number density

Returns a named tupple containing a LiqRaiRates object for the autoconversion rate and the liquid self-collection rate for scheme == SB2006Type

source
CloudMicrophysics.Microphysics2M.rain_self_collectionFunction
rain_self_collection(scheme, q_rai, ρ, N_rai)
  • scheme - type for 2-moment rain self-collection parameterization
  • q_rai - rain water specific content
  • ρ - air density
  • N_rai - raindrops number density

Returns the raindrops number density tendency due to collisions of raindrops that produce larger raindrops (self-collection) for scheme == SB2006Type

source
CloudMicrophysics.Microphysics2M.rain_breakupFunction
rain_breakup(scheme, q_rai, ρ, dN_rai_dt_sc)
  • scheme - type for 2-moment liquid self-collection parameterization
  • q_rai - rain water specific content
  • ρ - air density
  • N_rai - raindrops number density
  • dN_rai_dt_sc - rate of change of raindrops number density due to self-collection

Returns the raindrops number density tendency due to breakup of raindrops that produce smaller raindrops for scheme == SB2006Type

source
CloudMicrophysics.Microphysics2M.rain_self_collection_and_breakupFunction
rain_self_collection_and_breakup(SB2006, q_rai, ρ, N_rai)
  • SB2006 - a struct with SB2006 parameters for raindrops size distribution, self collection, and breakup
  • q_rai - rain water specific content
  • ρ - air density
  • N_rai - raindrops number density

Returns a named tupple containing the raindrops self-collection and breakup rates for scheme == SB2006Type

source
CloudMicrophysics.Microphysics2M.rain_terminal_velocityFunction
rain_terminal_velocity(SB2006, vel, q_rai, ρ, N_rai)
  • SB2006 - a struct with SB2006 rain size distribution parameters
  • vel - a struct with terminal velocity parameters
  • q_rai - rain water specific content [kg/kg]
  • ρ - air density [kg/m^3]
  • N_rai - raindrops number density [1/m^3]

Returns a tuple containing the number and mass weigthed mean fall velocities of raindrops in [m/s]. Assuming an exponential size distribution from Seifert and Beheng 2006 for scheme == SB2006Type Fall velocity of individual rain drops is parameterized:

  • assuming an empirical relation similar to Rogers (1993) for velo_scheme == SB2006VelType
  • following Chen et. al 2022, DOI: 10.1016/j.atmosres.2022.106171 for velo_scheme == Chen2022Type
source
CloudMicrophysics.Microphysics2M.rain_evaporationFunction
rain_evaporation(evap, aps, tps, q, q_rai, ρ, N_rai, T)
  • evap - evaporation parameterization scheme
  • aps - air properties
  • tps - thermodynamics parameters
  • q - phase partition
  • q_rai - rain specific content
  • ρ - air density
  • N_rai - raindrops number density
  • T - air temperature

Returns a named tuple containing the tendency of raindrops number density and rain water specific content due to rain rain_evaporation, assuming a power law velocity relation for fall velocity of individual drops and an exponential size distribution, for scheme == SB2006Type

source
CloudMicrophysics.Microphysics2M.conv_q_liq_to_q_raiFunction
conv_q_liq_to_q_rai(acnv, q_liq, ρ, N_d; smooth_transition)
  • acnv - 2-moment rain autoconversion parameterization
  • q_liq - cloud water specific content
  • ρ - air density
  • N_d - prescribed cloud droplet number concentration

Returns the q_rai tendency due to collisions between cloud droplets (autoconversion), parametrized following:

  • Khairoutdinov and Kogan (2000) for scheme == KK2000Type
  • Beheng (1994) for scheme == B1994Type
  • Tripoli and Cotton (1980) for scheme == TC1980Type
  • Liu and Daum (2004) for scheme ==LD2004Type

The Beheng1994Type, TC1980Type and LD2004Type of schemes additionally accept smooth_transition flag that smoothes their thershold behaviour if set to true. The default value is false.

source

P3 scheme

Construct parameterization set

CloudMicrophysics.Parameters.ParametersP3Type
ParametersP3{FT}

Parameters for P3 bulk microphysics scheme.

From Morrison and Milbrandt 2015 [20]

Fields

  • mass: Mass-size relation

  • area: Area-size relation

  • slope: Slope relation

  • vent: Ventilation relation

  • ρ_i: Cloud ice density [kg m⁻³]

  • ρ_l: Cloud liquid water density [kg m⁻³]

  • T_freeze: Water freeze temperature [K]

source
CloudMicrophysics.Parameters.ParametersP3Method
ParametersP3(FT)

Create a ParametersP3 object from a floating point type FT.

Examples

julia> import CloudMicrophysics.Parameters as CMP

julia> CMP.ParametersP3(Float64)
ParametersP3{Float64}
├── mass: MassPowerLaw
│   ├── α_va = 0.018537721864540644 [kg μm^(-β_va)]
│   └── β_va = 1.9 [-]
├── area: AreaPowerLaw
│   ├── γ = 0.2285 [μm^(2-σ)]
│   └── σ = 1.88 [-]
├── slope: SlopePowerLaw
│   ├── a = 0.00191 [m^b]
│   ├── b = 0.8 [-]
│   ├── c = 2.0 [-]
│   └── μ_max = 6.0 [-]
├── vent: VentilationSB2005
│   ├── vent_a = 0.78 [-]
│   └── vent_b = 0.308 [-]
├── ρ_i = 916.7 [kg m⁻³]
├── ρ_l = 1000.0 [kg m⁻³]
└── T_freeze = 273.15 [K]
source
CloudMicrophysics.Parameters.ParametersP3Method
ParametersP3(toml_dict::CP.AbstractTOMLDict; [slope_law = :powerlaw])

Create a ParametersP3 object from a ClimaParams TOML dictionary.

Arguments

  • toml_dict::CP.AbstractTOMLDict: A ClimaParams TOML dictionary
  • slope_law: Slope law to use (:constant or, by default, :powerlaw)
source

Sub-parameterizations

CloudMicrophysics.Parameters.MassPowerLawType
MassPowerLaw{FT}

Parameters for mass(size) relation.

From measurements of mass grown by vapor diffusion and aggregation in midlatitude cirrus by Brown and Francis (1995) doi: 10.1175/1520-0426(1995)012<0410:IMOTIW>2.0.CO;2

A part of the ParametersP3 parameter set.

Note

The BF1995_mass_coeff_alpha parameter is provided in units of [g μm^(-β_va)] but the α_va field is stored in SI-like units of [kg m^(-β_va)] for consistency with the rest of the code.

Fields

  • α_va: Coefficient in mass(size) relation [kg m^(-β_va)]

  • β_va: Coefficient in mass(size) relation [-]

source
CloudMicrophysics.Parameters.AreaPowerLawType
AreaPowerLaw{FT}

Parameters for area(size) relation.

A part of the ParametersP3 parameter set.

Fields

  • γ: Coefficient in area(size) for ice side plane, column, bullet, and planar polycrystal aggregates from Mitchell 1996 [μm^(2-σ)]

  • σ: Coefficient in area(size) for ice side plane, column, bullet, and planar polycrystal aggregates from Mitchell 1996 [-]

source
CloudMicrophysics.Parameters.SlopePowerLawType
SlopePowerLaw{FT}

Slope parameter μ as a power law in shape parameter λ:

\[μ(λ) = a λ^b - c\]

and is limited to:

\[0 ≤ μ ≤ μ_{max}\]

See also Eq. 3 in Morrison and Milbrandt 2015.

A part of the ParametersP3 parameter set.

Fields

  • a: Scale [m^b]

  • b: Power [-]

  • c: Offset [-]

  • μ_max: Upper limiter [-]

source
CloudMicrophysics.Parameters.VentilationSB2005Type
VentilationSB2005{FT}

Parameters for ventilation factor:

\[F(r) = a_{vent} + b_{vent} N_{Sc}^{1/3} N_{Re}(r)^{1/2}\]

From Seifert and Beheng 2005, doi: 10.1007/s00703-005-0112-4

A part of the ParametersP3 parameter set.

Fields

  • vent_a: Ventilation factor a [-]

  • vent_b: Ventilation factor b [-]

source

Obtain particle state

CloudMicrophysics.P3Scheme.P3StateType
P3State{FT}

State of the P3 scheme.

This struct bundles the P3 parameterizations params, the provided rime state (F_rim, ρ_r), and the derived threshold variables (D_th, D_gr, D_cr, ρ_g).

To obtain a P3State object, use the get_state function.

Fields

  • params: CMP.ParametersP3 object

  • F_rim: Rime mass fraction

  • ρ_r: Rime density

  • ρ_g: Graupel density

  • D_th: Spherical ice threshold

  • D_gr: Graupel threshold

  • D_cr: Partially rimed ice threshold

source
CloudMicrophysics.P3Scheme.get_stateFunction
get_state(params; F_rim, ρ_r)

Create a P3State object from a CMP.ParametersP3 object and rime state parameters.

Arguments

Examples

julia> import CloudMicrophysics.Parameters as CMP, CloudMicrophysics.P3Scheme as P3

julia> FT = Float32;

julia> params = CMP.ParametersP3(FT);

julia> state = P3.get_state(params; F_rim = FT(0.5), ρ_r = FT(916.7))
P3State{Float32}
├── params = {MassPowerLaw, AreaPowerLaw, SlopePowerLaw, VentilationSB2005}
├── F_rim = 0.5 [-]
├── ρ_r = 916.7 [kg/m^3]
├── ρ_g = 702.8062 [kg/m^3]
├── D_th = 9.728088e-5 [m]
├── D_gr = 0.00012385941 [m]
└── D_cr = 0.00023259086 [m]
source
get_state(dist::P3Distribution)

Return the particle state from a P3Distribution object.

source

State relationships

CloudMicrophysics.P3Scheme.get_ρ_dFunction
get_ρ_d(mass::MassPowerLaw, F_rim, ρ_r)

Exact solution for the density of the unrimed portion of the particle as function of the rime mass fraction F_rim, mass power law parameters mass, and rime density ρ_r.

Arguments

Returns

  • ρ_d: density of the unrimed portion of the particle [kg/m³]

Examples

julia> import CloudMicrophysics.Parameters as CMP, 
              ClimaParams as CP, 
              CloudMicrophysics.P3Scheme as P3

julia> FT = Float64;

julia> mass = CMP.MassPowerLaw(CP.create_toml_dict(FT));

julia> F_rim, ρ_r = FT(0.5), FT(916.7);

julia> ρ_d = P3.get_ρ_d(mass, F_rim, ρ_r)
488.9120789986412
source
CloudMicrophysics.P3Scheme.get_ρ_gFunction
get_ρ_g(ρ_r, F_rim, ρ_d)

Return the density of total (deposition + rime) ice mass for graupel [kg/m³]

Arguments

  • ρ_r: rime density (L_rim/B_rim) [kg/m³]
  • F_rim: rime mass fraction (L_rim / L_ice) [-]
  • ρ_d: density of the unrimed portion of the particle [kg/m³], see get_ρ_d

Returns

  • ρ_g: density of total (deposition + rime) ice mass for graupel [kg/m³]

Examples

julia> import CloudMicrophysics.P3Scheme as P3

julia> FT = Float64; 

julia> ρ_r, F_rim, ρ_d = FT(916.7), FT(0.5), FT(916.7);

julia> ρ_g = P3.get_ρ_g(ρ_r, F_rim, ρ_d)
916.7

Notes:

See Eq. 16 in [20].

source

Fetch from state

CloudMicrophysics.P3Scheme.threshold_tupleFunction
threshold_tuple(state)

Return a tuple of the thresholds for the current state.

This function is useful for providing thresholds to quadgk, for example.

Note

If the state is unrimed, there is only one threshold, D_th. Otherwise (rimed state), there are three thresholds, D_th < D_gr < D_cr.

source

Derived quantities

Main particle methods

These methods are a function of the particle diameter, D.

CloudMicrophysics.P3Scheme.ice_densityFunction
ice_density(state::P3State, D)

Return the density of a particle based on where it falls in the particle-size-based properties regime.

Arguments

  • state: P3State object
  • D: maximum particle dimension [m]

The density of nonspherical particles is assumed to be the particle mass divided by the volume of a sphere with the same D [20]. Needed for aspect ratio calculation, so we assume zero liquid fraction.

source
CloudMicrophysics.P3Scheme.∂ice_mass_∂DFunction
∂ice_mass_∂D(state, D)

Calculate the derivative of the ice mass with respect to the maximum particle dimension.

Used in snow melting calculations.

Arguments

  • state: P3State object
  • D: maximum dimension of the particle [m]
source
CloudMicrophysics.P3Scheme.ice_areaFunction
ice_area(state, D)

Return the cross-sectional area of a particle based on where it falls in the particle-size-based properties regime.

Arguments

  • state: P3State object
  • D: maximum particle dimension [m]
source
CloudMicrophysics.P3Scheme.ϕᵢFunction
ϕᵢ(state, D)

Returns the aspect ratio (ϕ) for an ice particle

Arguments

  • state: P3State object
  • D: maximum dimension of ice particle [m]

The density of nonspherical particles is assumed to be equal to the particle mass divided by the volume of a spherical particle with the same D_max [20]. Assuming zero liquid fraction and oblate shape.

source

Supporting methods

CloudMicrophysics.P3Scheme.mass_sphericalFunction
mass_spherical(ρ, D)

Calculate the mass as a function of size for spherical particles, used for small ice, liquid, and graupel.

Arguments

  • ρ: bulk density [kg m⁻³]
  • D: maximum particle dimension [m]
source

Distribution parameters

CloudMicrophysics.P3Scheme.P3DistributionType
P3Distribution{FT}

Distribution of ice particles in size.

Fields

  • state: Particle state, see P3State

  • L: The mass concentration [kg/m³]

  • N: The number concentration [1/m³]

  • log_λ: Logarithm of the slope parameter

  • log_N₀: Logarithm of the intercept parameter

source
CloudMicrophysics.P3Scheme.get_distribution_parametersFunction
get_distribution_parameters(state; L, N, [log_λ_min, log_λ_max])

Solve for the distribution parameters given the state, and the mass (L) and number (N) concentrations.

The assumed distribution is of the form

\[N(D) = N₀ D^μ e^{-λD}\]

where N(D) is the number concentration at diameter D and μ is the slope parameter. The slope parameter is parameterized, e.g. CMP.SlopePowerLaw or CMP.SlopeConstant.

This algorithm solves for log_λ = log(λ) and log_N₀ = log(N₀) given L and N by solving the equations:

\[\begin{align*} \log(L) &= \log ∫_0^∞ m(D) N(D)\ \mathrm{d}D, \\ \log(N) &= \log ∫_0^∞ N(D)\ \mathrm{d}D, \\ \end{align*}\]

where m(D) is the mass of a particle at diameter D (see ice_mass). The procedure is decribed in detail in the P3 docs.

Arguments

Keyword arguments

  • L: The mass concentration [kg/m³]
  • N: The number concentration [1/m³]
  • search_bound_min: The minimum value of the search bounds [log(1/m)], default is log(1e1)
  • search_bound_max: The maximum value of the search bounds [log(1/m)], default is log(1e7)

Returns

  • P3Distribution object with the distribution parameters log_λ and log_N₀

Examples

julia> import CloudMicrophysics.Parameters as CMP,
              CloudMicrophysics.P3Scheme   as P3

julia> params = CMP.ParametersP3(Float64);

julia> state = P3.get_state(params; F_rim = 0.0, ρ_r = 400.0);

julia> dist = P3.get_distribution_parameters(state; L = 1e-3, N = 1e3)
P3Distribution{Float64}
├── state: is unrimed
├── log_λ = 5.4897008376530385 [log(1/m)]
└── log_N₀ = 12.397456116635176 [log(1/m^3)]
source

Distribution relationships

CloudMicrophysics.P3Scheme.log_N_div_N₀Function
log_N_div_N₀(state, log_λ)
log_N_div_N₀(slope, log_λ)

Compute log(N/N₀) given either the state or the slope parameterization, and log_λ

Note: This function is equivalent to log_integrate_moment_psd(0, Inf, 1, 0, μ, log_λ)

source
CloudMicrophysics.P3Scheme.get_log_N₀Function
get_log_N₀(state; N, log_λ)
get_log_N₀(params; N, log_λ)

Compute log(N₀) given the state, N, and log_λ

Arguments

Keyword arguments

  • N: The number concentration [1/m³]
  • log_λ: The log of the slope parameter [log(1/m)]
source

Derived integral quantities

These methods integrate over the particle size distribution.

CloudMicrophysics.P3Scheme.ice_particle_terminal_velocityFunction
ice_particle_terminal_velocity(state, D, Chen2022, ρₐ, use_aspect_ratio)

Returns the terminal velocity of a single ice particle as a function of its size (maximum dimension, D) using the Chen 2022 parametrization.

Arguments

  • state: The P3State object
  • D: Maximum particle dimension
  • Chen2022: The CMP.Chen2022VelType object with terminal velocity parameters
  • ρₐ: Air density
  • use_aspect_ratio: Bool flag set to true if we want to consider the effects of particle aspect ratio on its terminal velocity (default: true)
source
CloudMicrophysics.P3Scheme.ice_terminal_velocityFunction
ice_terminal_velocity(dist, Chen2022, ρₐ, use_aspect_ratio)

Compute the mass and number weighted fall speeds for ice

See Eq. C10 of [20] and use the Chen 2022 terminal velocity scheme.

Arguments

  • dist: The P3Distribution object
  • Chen2022: The CMP.Chen2022VelType object
  • ρₐ: The density of air
  • use_aspect_ratio: Bool flag set to true if we want to consider the effects of particle aspect ratio on its terminal velocity (default: true)
  • accurate: Set to true to perform a more accurate numerical integration see ∫fdD for details. Default is false.

Returns

  • v_n: The number weighted fall speed
  • v_m: The mass weighted fall speed
source
CloudMicrophysics.P3Scheme.het_ice_nucleationFunction
het_ice_nucleation(pdf_c, p3, tps, q, N, T, ρₐ, p, aerosol)
  • aerosol - aerosol parameters (supported types: desert dust, illite, kaolinite)
  • tps - thermodynamics parameters
  • qₚ - phase partition
  • N_liq - cloud water number concentration
  • RH - relative humidity
  • T - temperature
  • ρₐ - air density
  • dt - model time step

Returns a named tuple with ice number concentration and ice content hetergoeneous freezing rates from cloud droplets.

source
CloudMicrophysics.P3Scheme.ice_meltFunction
ice_melt(p3, Chen2022, aps, tps, Tₐ, ρₐ, dt; ∫kwargs...)

Arguments

  • dist: a P3Distribution object
  • Chen2022: struct containing Chen 2022 velocity parameters
  • aps: air properties
  • tps: thermodynamics parameters
  • Tₐ: temperature (K)
  • ρₐ: air density
  • dt: model time step (for limiting the tendnecy)

Keyword arguments

  • ∫kwargs: Named tuple of keyword arguments passed to ∫fdD

Returns the melting rate of ice (QIMLT in Morrison and Mildbrandt (2015)).

source

Supporting integral methods

CloudMicrophysics.P3Scheme.∫fdDFunction
∫fdD(f, state::P3State; [D_max = 1], kwargs...)

Integrate the function f over the size distribution of the ice particles.

Usage

This function is useful for integrating functions over the size distribution of the ice particles. It is a light wrapper around QGK.quadgk that automatically inserts the calculated size distribution thresholds (i.e. D_th, D_gr, D_cr) as integration limits. The upper integration limit is by default set to D_max = 1 m, with the presumption that all integrals will decay to zero much quicker than this. Formally, integration is to infinity.

This method calls ∫fdD_error, which returns both the value of the integral and the estimated error. Since the error is typically not of interest, this method only returns the value of the integral.

Arguments

  • f: The function to integrate
  • state: The P3State object
  • D_max: The maximum diameter to integrate to [m]. Default: D_max = 1 m.
  • kwargs: Additional optional keyword arguments to pass to QGK.quadgk
    • rtol: The relative tolerance for the integration, default: rtol = sqrt(eps(FT))
    • atol: The absolute tolerance for the integration, default: atol = 0
    • maxevals: The maximum number of function evaluations, default: maxevals = 10^7
    • order: The order of the quadrature rule, default: order = 7

Returns

  • value: The value of the integral
Integral accuracy

To achieve highest accuracy, which can be challenging when integrating the N′ice function, it is recommended to increase the order of the quadrature rule and set rtol = 0. Experimentally, order = 44 has been found to be sufficient for most cases.

For convenience, passing accurate = true will set rtol = 0 and order = 44.

Examples

julia> import CloudMicrophysics.Parameters as CMP,
              CloudMicrophysics.P3Scheme   as P3

julia> params = CMP.ParametersP3(Float64);

julia> state = P3.get_state(params; F_rim = 0.0, ρ_r = 400.0);

julia> f(D) = D^3;  # Define a function to integrate

julia> P3.∫fdD(f, state)    # Integrate the function
0.25

julia> P3.∫fdD(state) do D  # Integrate with a `do`-block
           P3.ice_mass(state, D)
       end
0.006392317884295288
source
CloudMicrophysics.P3Scheme.∫fdD_errorFunction
∫fdD_error(f, state::P3State; [D_max = 1], kwargs...)

Integrate the function f over the size distribution of the ice particles.

Returns

  • value: The value of the integral
  • error: The estimated error of the integral

Notes

See ∫fdD, which only returns the value of the integral and not the error, for details.

source

Aerosol model

CloudMicrophysics.AerosolModelModule
A container for information on aerosol size distribution
and chemical properties.

The size distribution is a sum of lognormal internally mixed modes.
The chemical composition can be expressed using kappa parameter
or hygroscopicity parameter B.
source
CloudMicrophysics.AerosolModel.Mode_BType
Mode_B

Represents the sizes and chemical composition of aerosol particles in one size distribution mode. The mode is assumed to be made up of internally mixed components and follow a lognormal size distribution. The chemical composition of aerosol particles in this mode is described using the parameters from Abdul-Razzak and Ghan 2000.

source
CloudMicrophysics.AerosolModel.Mode_κType
Mode_κ

Represents the sizes and chemical composition of aerosol particles in one size distribution mode. The mode is assumed to be made up of internally mixed components and follow a lognormal size distribution. The chemical composition of aerosol particles in this mode is described using the parameters from Petters and Kreidenweis 2007.

source
CloudMicrophysics.AerosolModel.AerosolDistributionType
AerosolDistribution

Represents the aerosol size distribution as a tuple with different modes. All modes have to either be of type ModeB (Abdul-Razzak and Ghan 2000) or of type Modeκ (Petters and Kreidenweis 2007).

Constructors

AerosolDistribution(modes::T)
source

Aerosol activation

CloudMicrophysics.AerosolActivationModule
Aerosol activation scheme, which includes:
  • mean hygroscopicity for each mode of the aerosol size distribution
  • critical supersaturation for each mode of the aerosol size distribution
  • maximum supersaturation
  • total number of particles actived
  • total mass of particles actived
source
CloudMicrophysics.AerosolActivation.mean_hygroscopicity_parameterFunction
mean_hygroscopicity_parameter(ap, ad)
  • ap - a struct with aerosol activation parameters
  • ad - a struct with aerosol distribution (B or κ based)

Returns a tuple of hygroscopicity parameters (one tuple element for each aerosol size distribution mode). The tuple is computed either as mass-weighted B parameters (Abdul-Razzak and Ghan 2000) or volume weighted kappa parameters (Petters and Kreidenweis 2007). Implemented via a dispatch based on aerosol distribution mode type.

source
CloudMicrophysics.AerosolActivation.max_supersaturationFunction
max_supersaturation(ap, ad, aip, tps, T, p, w, q)
  • ap - a struct with aerosol activation parameters
  • ad - a struct with aerosol distribution
  • aip - a struct with air parameters
  • tps - a struct with thermodynamics parameters
  • T - air temperature
  • p - air pressure
  • w - vertical velocity
  • q - phase partition

Returns the maximum supersaturation.

source
CloudMicrophysics.AerosolActivation.N_activated_per_modeFunction
N_activated_per_mode(ap, ad, aip, tps, T, p, w, q)
  • ap - a struct with aerosol activation parameters
  • ad - aerosol distribution struct
  • aip - a struct with air parameters
  • tps - a struct with thermodynamics parameters
  • T - air temperature
  • p - air pressure
  • w - vertical velocity
  • q - phase partition

Returns the number of activated aerosol particles in each aerosol size distribution mode.

source
CloudMicrophysics.AerosolActivation.M_activated_per_modeFunction
M_activated_per_mode(ap, ad, aip, tps, T, p, w, q)
  • ap - a struct with aerosol activation parameters
  • ad - a struct with aerosol distribution parameters
  • aip - a struct with air parameters
  • tps - a struct with thermodynamics parameters
  • T - air temperature
  • p - air pressure
  • w - vertical velocity
  • q - phase partition

Returns the mass of activated aerosol particles per mode of the aerosol size distribution.

source
CloudMicrophysics.AerosolActivation.total_N_activatedFunction
total_N_activated(ap, ad, aip, tps, T, p, w, q)
  • ap - a struct with aerosol activation parameters
  • ad - aerosol distribution struct
  • aip - a struct with air properties
  • tps - a struct with thermodynamics parameters
  • T - air temperature
  • p - air pressure
  • w - vertical velocity
  • q - phase partition

Returns the total number of activated aerosol particles.

source
CloudMicrophysics.AerosolActivation.total_M_activatedFunction
total_M_activated(ap, ad, aip, tps, T, p, w, q)
  • ap - a struct with aerosol activation parameters
  • ad - aerosol distribution struct
  • aip - a struct with air properties
  • tps - a struct with thermodynamics parameters
  • T - air temperature
  • p - air pressure
  • w - vertical velocity
  • q - phase partition

Returns the total mass of activated aerosol particles.

source

Artifact calling

Heterogeneous ice nucleation

CloudMicrophysics.HetIceNucleation.dust_activated_number_fractionFunction
dust_activated_number_fraction(dust, ip, Si, T)
  • dust - a struct with dust parameters
  • ip - a struct with ice nucleation parameters
  • Si - ice saturation ratio
  • T - air temperature [K]

Returns the number fraction of mineral dust particles acting as deposition nuclei (n ice nuclei / n dust particles). From Mohler et al 2006 Table 2 (averages from different measurements excluding those where a was not measured)

source
CloudMicrophysics.HetIceNucleation.MohlerDepositionRateFunction
MohlerDepositionRate(dust, ip, S_i, T, dSi_dt, N_aer)
  • dust - a struct with dust parameters
  • ip - a struct with ice nucleation parameters
  • Si - ice saturation
  • T - ambient temperature
  • dSi_dt - change in ice saturation over time
  • N_aer - number of unactivated aerosols

Returns the ice nucleation rate from deposition. From Mohler et al 2006 equation 5.

source
CloudMicrophysics.HetIceNucleation.deposition_JFunction
deposition_J(dust, Δa_w)
  • dust - a struct with dust parameters
  • Δa_w - change in water activity [unitless].

Returns the deposition nucleation rate coefficient, J, in m^-2 s^-1 for different minerals in liquid droplets. The free parameters m and c are derived from China et al (2017) see DOI: 10.1002/2016JD025817 Returns zero for unsupported aerosol types.

source
CloudMicrophysics.HetIceNucleation.ABIFM_JFunction
ABIFM_J(dust, Δa_w)
  • dust - a struct with dust parameters
  • Δa_w - change in water activity [unitless].

Returns the immersion freezing nucleation rate coefficient, J, in m^-2 s^-1 for different minerals in liquid droplets. The free parameters m and c are taken from Knopf & Alpert 2013 see DOI: 10.1039/C3FD00035D Returns zero for unsupported aerosol types.

source
CloudMicrophysics.HetIceNucleation.P3_deposition_N_iFunction
P3_deposition_N_i(ip, T)
  • ip - a struct with ice nucleation parameters,
  • T - air temperature [K].

Returns the number of ice nucleated via deposition nucleation with units of m^-3. From Thompson et al 2004 eqn 2 as used in Morrison & Milbrandt 2015.

source
CloudMicrophysics.HetIceNucleation.P3_het_N_iFunction
P3_het_N_i(ip, T, N_l, B, V_l, a, Δt)
  • ip - a struct with ice nucleation parameters,
  • T - air temperature [K],
  • N_l - number of droplets [m^-3],
  • B - water-type dependent parameter [cm^-3 s^-1],
  • V_l - volume of droplets to be frozen [m^3],
  • a - empirical parameter [C^-1],
  • Δt - timestep.

Returns the number of ice nucleated within Δt via heterogeneous freezing with units of m^-3. From Pruppacher & Klett 1997 eqn (9-51) as used in Morrison & Milbrandt 2015.

source
CloudMicrophysics.HetIceNucleation.INP_concentration_frequencyFunction
INP_concentration_frequency(params,INPC,T)
  • params - a struct with INPC(T) distribution parameters
  • INPC - concentration of ice nucleating particles [m^-3]
  • T - air temperature [K]

Returns the relative frequency of a given INP concentration, depending on the temperature. Based on Frostenberg et al., 2023. See DOI: 10.5194/acp-23-10883-2023

source

Homogeneous ice nucleation

CloudMicrophysics.HomIceNucleation.homogeneous_J_cubicFunction
homogeneous_J_cubic(ip, Δa_w)
  • ip - a struct with ice nucleation parameters
  • Δa_w - change in water activity [-].

Returns the homogeneous freezing nucleation rate coefficient, J, in m^-3 s^-1 for sulphuric acid solutions. Parameterization based on Koop 2000, DOI: 10.1038/35020537.

source
CloudMicrophysics.HomIceNucleation.homogeneous_J_linearFunction
homogeneous_J_linear(ip, Δa_w)
  • ip - a struct with ice nucleation parameters
  • Δa_w - change in water activity [-].

Returns the homogeneous freezing nucleation rate coefficient, J, in m^-3 s^-1 for sulphuric acid solutions. Parameterization derived from a linear fit of the Koop 2000 parameterization, DOI: 10.1038/35020537.

source

Cloud diagnostics

CloudMicrophysics.CloudDiagnostics.radar_reflectivity_1MFunction
radar_reflectivity_1M(precip, q, ρ)
  • precip - struct with 1-moment microphysics rain free parameters
  • q - specific content of rain
  • ρ - air density

Returns logarithmic radar reflectivity for the 1-moment microphysics based on the assumed rain particle size distribution. Normalized by the reflectivty of 1 millimiter drop in a volume of 1m3. The values are clipped at -150 dBZ.

source
CloudMicrophysics.CloudDiagnostics.radar_reflectivity_2MFunction
radar_reflectivity_2M(structs, q_liq, q_rai, N_liq, N_rai, ρ_air)

- `structs` - structs microphysics 2-moment with SB2006 cloud droplets
              and raindrops size distributions parameters
- `q_liq` - cloud water specific content
- `q_rai` - rain water specific content
- `N_liq` - cloud droplet number density
- `N_rai` - rain droplet number density
- `ρ_air` - air density

Returns logarithmic radar reflectivity for the 2-moment microphysics SB2006 based on the assumed cloud and rain particle size distribuions. Normalized by the reflectivty of 1 millimiter drop in a volume of 1m3. The values are clipped at -150 dBZ.

source
CloudMicrophysics.CloudDiagnostics.effective_radius_Liu_Hallet_97Function
effective_radius_Liu_Hallet_97(q_liq, q_rai, N_liq, N_rai, ρ_air, ρ_w)

- `q_liq` - cloud water specific content
- `q_rai` - rain water specific content
- `N_liq` - cloud droplet number density
- `N_rai` - rain droplet number density
- `ρ_air` - air density
- `ρ_w` - water density

Returns effective radius using the "1/3" power law from Liu and Hallett (1997). If not provided by the user, it is assumed that there is no rain present and that the cloud droplet number concentration is 100 1/cm3.

source
CloudMicrophysics.CloudDiagnostics.effective_radius_2MFunction
effective_radius_2M(structs, q_liq, q_rai, N_liq, N_rai, ρ_air)

- `structs` - structs with SB2006 cloud droplets and raindrops
            size distribution parameters
- `q_liq` - cloud water specific content
- `q_rai` - rain water specific content
- `N_liq` - cloud droplet number density
- `N_rai` - rain droplet number density
- `ρ_air` - air density

Returns effective radius for the 2-moment microphysics scheme. Computed based on the assumed cloud and rain particle size distributions.

source

Common utility functions

CloudMicrophysics.Common.G_funcFunction
G_func(air_props, tps, T, Liquid())
G_func(air_props, tps, T, Ice())
  • air_props - struct with air parameters
  • tps - struct with thermodynamics parameters
  • T - air temperature
  • Liquid(), Ice() - liquid or ice phase to dispatch over.

Utility function combining thermal conductivity and vapor diffusivity effects.

source
CloudMicrophysics.Common.logistic_functionFunction
logistic_function(x, x_0, k)
  • x - independent variable
  • x_0 - threshold value for x
  • k - growth rate of the curve, characterizing steepness of the transition

Returns the value of the logistic function for smooth transitioning at thresholds. This is a normalized curve changing from 0 to 1 while x varies from 0 to Inf (for positive k). For x < 0 the value at x = 0 (zero) is returned. For x_0 = 0 H(x) is returned.

source
CloudMicrophysics.Common.logistic_function_integralFunction
logistic_function_integral(x, x_0, k)
  • x - independent variable
  • x_0 - threshold value for x
  • k - growth rate of the logistic function, characterizing steepness of the transition

Returns the value of the indefinite integral of the logistic function, for smooth transitioning of piecewise linear profiles at thresholds. This curve smoothly transition from y = 0 for 0 < x < x0 to y = x - x0 for x_0 < x.

source
CloudMicrophysics.Common.H2SO4_soln_saturation_vapor_pressureFunction
H2SO4_soln_saturation_vapor_pressure(prs, x, T)
  • prs - a struct with H2SO4 solution free parameters
  • x - wt percent sulphuric acid [unitless]
  • T - air temperature [K].

Returns the saturation vapor pressure above a sulphuric acid solution droplet in Pa. x is, for example, 0.1 if droplets are 10 percent sulphuric acid by weight

source
CloudMicrophysics.Common.a_w_xTFunction
a_w_xT(H2SO4_prs, tps, x, T)
  • H2SO4_prs - a struct with H2SO4 solution free parameters
  • tps - a struct with thermodynamics parameters
  • x - wt percent sulphuric acid [unitless]
  • T - air temperature [K].

Returns water activity of H2SO4 containing droplet. x is, for example, 0.1 if droplets are 10 percent sulphuric acid by weight.

source
CloudMicrophysics.Common.a_w_eTFunction
a_w_eT(tps, e, T)
  • tps - struct with thermodynamics parameters
  • e - partial pressure of water [Pa]
  • T - air temperature [K].

Returns water activity of pure water droplet. Valid when droplet is in equilibrium with surroundings.

source
CloudMicrophysics.Common.Chen2022_monodisperse_pdfFunction
Chen2022_monodisperse_pdf(a, b, c, D)
  • a, b, c, - free parameters defined in Chen etl al 2022
  • D - droplet diameter

Returns the addends of the bulk fall speed of rain or ice particles following Chen et al 2022 DOI: 10.1016/j.atmosres.2022.106171 in [m/s]. Assuming monodisperse droplet distribution.

source
CloudMicrophysics.Common.Chen2022_exponential_pdfFunction
Chen2022_exponential_pdf(a, b, c, λ_inv, k)
  • a, b, c, - free parameters defined in Chen etl al 2022
  • λ_inv - inverse of the size distribution parameter
  • k - size distribution moment for which we compute the bulk fall speed

Returns the addends of the bulk fall speed of rain or ice particles following Chen et al 2022 DOI: 10.1016/j.atmosres.2022.106171 in [m/s]. Assuming exponential size distribution and hence μ=0.

source
CloudMicrophysics.Common.Chen2022_vel_coeffs_B1Function
Chen2022_vel_coeffs_B1(coeffs, ρₐ)
  • coeffs - a struct with terminal velocity free parameters
  • ρₐ - air density

Returns the coefficients from Table B1 Appendix B in Chen et al 2022 DOI: 10.1016/j.atmosres.2022.106171

source
CloudMicrophysics.Common.Chen2022_vel_coeffs_B2Function
Chen2022_vel_coeffs_B2(coeffs, ρₐ, ρᵢ)
  • coeffs - a struct with terminal velocity free parameters
  • ρₐ - air density
  • ρᵢ - apparent density of ice particles

Returns the coefficients from Table B2 Appendix B in Chen et al 2022 DOI: 10.1016/j.atmosres.2022.106171

source
CloudMicrophysics.Common.Chen2022_vel_coeffs_B4Function
Chen2022_vel_coeffs_B4(coeffs, ρₐ, ρᵢ)
  • coeffs - a struct with terminal velocity free parameters
  • ρₐ - air density
  • ρᵢ - apparent density of ice particles

Returns the coefficients from Table B4 Appendix B in Chen et al 2022 DOI: 10.1016/j.atmosres.2022.106171

source

Parameters

CloudMicrophysics.Parameters.ArizonaTestDustType
ArizonaTestDust{FT}

Parameters for Arizona Test Dust from Mohler et al, 2006. DOI: 10.5194/acp-6-3007-2006

Fields

  • S₀_warm: S₀ for T > T_thr [-]

  • S₀_cold: S₀ for T < T_thr [-]

  • a_warm: a for T > T_thr [-]

  • a_cold: a for T < T_thr [-]

  • deposition_m: m coefficient for deposition nucleation J [-]

  • deposition_c: c coefficient for deposition nucleation J [-]

  • ABIFM_m: m coefficient for immersion freezing J [-]

  • ABIFM_c: c coefficient for immersion freezing J [-]

source
CloudMicrophysics.Parameters.DesertDustType
DesertDust{FT}

Parameters for desert dust from Knopf and Alpert 2013 DOI: 10.1039/C3FD00035D and from Mohler et al, 2006 DOI: 10.5194/acp-6-3007-2006

Fields

  • S₀_warm: S₀ for T > T_thr [-]

  • S₀_cold: S₀ for T < T_thr [-]

  • a_warm: a for T > T_thr [-]

  • a_cold: a for T < T_thr [-]

  • ABIFM_m: m coefficient for immersion freezing J [-]

  • ABIFM_c: c coefficient for immersion freezing J [-]

source
CloudMicrophysics.Parameters.AsianDustType
AsianDust{FT}

Parameters for Asian Dust

Fields

  • deposition_m: m coefficient for deposition nucleation J [-]

  • deposition_c: c coefficient for deposition nucleation J [-]

  • ABIFM_m: m coefficient for immersion freezing J [-]

  • ABIFM_c: c coefficient for immersion freezing J [-]

source
CloudMicrophysics.Parameters.IlliteType
Illite{FT}

Parameters for illite from Knopf and Alpert 2013 DOI: 10.1039/C3FD00035D

Fields

  • deposition_m: m coefficient for deposition nucleation J [-]

  • deposition_c: c coefficient for deposition nucleation J [-]

  • ABIFM_m: m coefficient for immersion freezing J [-]

  • ABIFM_c: c coefficient for immersion freezing J [-]

source
CloudMicrophysics.Parameters.KaoliniteType
Kaolinite{FT}

Parameters for kaolinite from Knopf and Alpert 2013 DOI: 10.1039/C3FD00035D and China et al 2017 DOI: 10.1002/2016JD025817

Fields

  • deposition_m: m coefficient for deposition nucleation J [-]

  • deposition_c: c coefficient for deposition nucleation J [-]

  • ABIFM_m: m coefficient for immersion freezing J [-]

  • ABIFM_c: c coefficient for immersion freezing J [-]

source
CloudMicrophysics.Parameters.FeldsparType
Feldspar{FT}

Parameters for Feldspar from Alpert et al 2022 DOI: 10.1039/D1EA00077B

Fields

  • deposition_m: m coefficient for deposition nucleation J [-]

  • deposition_c: c coefficient for deposition nucleation J [-]

source
CloudMicrophysics.Parameters.FerrihydriteType
Ferrihydrite{FT}

Parameters for Ferrihydrite from Alpert et al 2022 DOI: 10.1039/D1EA00077B

Fields

  • deposition_m: m coefficient for deposition nucleation J [-]

  • deposition_c: c coefficient for deposition nucleation J [-]

source
CloudMicrophysics.Parameters.DustType
Dust{FT}

Parameters for generic dust.

Fields

  • deposition_m: m coefficient for deposition nucleation J [-]

  • deposition_c: c coefficient for deposition nucleation J [-]

  • ABIFM_m: m coefficient for immersion freezing J [-]

  • ABIFM_c: c coefficient for immersion freezing J [-]

source
CloudMicrophysics.Parameters.SeasaltType
Seasalt{FT}

Parameters for seasalt

Fields

  • M: molar mass [kg/mol]

  • ρ: density [kg/m3]

  • ϕ: osmotic coefficient [-]

  • ν: ion number [-]

  • ϵ: water soluble mass fraction [-]

  • κ: hygroscopicity parameter [-]

source
CloudMicrophysics.Parameters.SulfateType
Sulfate{FT}

Parameters for sulfate aerosol

Fields

  • M: molar mass [kg/mol]

  • ρ: density [kg/m3]

  • ϕ: osmotic coefficient [-]

  • ν: ion number [-]

  • ϵ: water soluble mass fraction [-]

  • κ: hygroscopicity parameter [-]

source
CloudMicrophysics.Parameters.AerosolActivationParametersType
AerosolActivationParameters{FT}

Parameters for Abdul-Razzak and Ghan 2000 aerosol activation scheme DOI: 10.1029/1999JD901161

Fields

  • M_w: molar mass of water [kg/mol]

  • R: gas constant [J/mol/K]

  • ρ_w: cloud water density [kg/m3]

  • σ: surface tension of water [N/m]

  • g: gravitational_acceleration [m/s2]

  • f1: scaling coefficient in Abdul-Razzak and Ghan 2000 [-]

  • f2: scaling coefficient in Abdul-Razzak and Ghan 2000 [-]

  • g1: scaling coefficient in Abdul-Razzak and Ghan 2000 [-]

  • g2: scaling coefficient in Abdul-Razzak and Ghan 2000 [-]

  • p1: power of (zeta / eta) in Abdul-Razzak and Ghan 2000 [-]

  • p2: power of (S_m^2 / (zeta + 3 * eta)) in Abdul-Razzak and Ghan 2000 [-]

source
CloudMicrophysics.Parameters.H2SO4SolutionParametersType
H2SO4SolutionParameters{FT}

Parameters for water activity of H2SO4 solutions from Luo et al 1995. DOI: 10.1029/94GL02988

Fields

  • T_max: max temperature for which the parameterization is valid [K]

  • T_min: min temperature for which the parameterization is valid [K]

  • w_2: coefficient [-]

  • c1: coefficient [-]

  • c2: coefficient [-]

  • c3: coefficient [-]

  • c4: coefficient [-]

  • c5: coefficient [-]

  • c6: coefficient [-]

  • c7: coefficient [-]

source
CloudMicrophysics.Parameters.Mohler2006Type
Mohler2006{FT}

Parameters for ice nucleation from Mohler et al 2006 DOI: 10.5194/acp-6-3007-2006

Fields

  • Sᵢ_max: max allowed supersaturation [-]

  • T_thr: threshold temperature [K]

source
CloudMicrophysics.Parameters.Koop2000Type
Koop2000{FT}

Parameters for ice nucleation from Koop et al 2000 DOI: 10.1038/35020537

Fields

  • Δa_w_min: min Δaw [-]

  • Δa_w_max: max Δaw [-]

  • c₁: coefficient [-]

  • c₂: coefficient [-]

  • c₃: coefficient [-]

  • c₄: coefficient [-]

  • linear_c₁: coefficient [-]

  • linear_c₂: coefficient [-]

source
CloudMicrophysics.Parameters.Parameters0MType
Parameters0M{FT}

Parameters for zero-moment bulk microphysics scheme

Fields

  • τ_precip: precipitation timescale [s]

  • qc_0: condensate specific content precipitation threshold [-]

  • S_0: supersaturation precipitation threshold [-]

source
CloudMicrophysics.Parameters.ParticleMassType
ParticleMass{FT}

A struct with coefficients of the assumed mass(size) relationship for particles

m(r) = m0 χm (r/r0)^(me + Δm)

Fields

  • r0: particle length scale [m]

  • m0: mass size relation coefficient [kg]

  • me: mass size relation coefficient [-]

  • Δm: mass size relation coefficient [-]

  • χm: mass size relation coefficient [-]

source
CloudMicrophysics.Parameters.ParticleAreaType
ParticleArea{FT}

A struct with coefficients of the assumed crosssectionarea(size) relationship for particles

a(r) = a0 χa (r/r0)^(ae + Δa)

Fields

  • a0: cross section size relation coefficient [m2]

  • ae: cross section size relation coefficient [-]

  • Δa: cross section size relation coefficient [-]

  • χa: cross section size relation coefficient [-]

source
CloudMicrophysics.Parameters.Acnv1MType
Acnv1M{FT}

A struct with autoconversion parameters

Fields

  • τ: autoconversion timescale [s]

  • q_threshold: condensate specific content autoconversion threshold [-]

  • k: threshold smooth transition steepness [-]

source
CloudMicrophysics.Parameters.CloudLiquidType
CloudLiquid{FT}

The parameters and type for cloud liquid water condensate

Fields

  • τ_relax: condensation evaporation non_equil microphysics relaxation timescale [s]

  • ρw: water density [kg/m3]

  • r_eff: effective radius [m]

source
CloudMicrophysics.Parameters.CloudIceType
CloudIce{FT, MS}

The parameters and type for cloud ice condensate

Fields

  • pdf: a struct with size distribution parameters

  • mass: a struct with mass size relation parameters

  • r0: particle length scale [m]

  • r_ice_snow: ice snow threshold radius [m]

  • τ_relax: deposition sublimation non_equil microphysics relaxation timescale [s]

  • ρᵢ: cloud ice apparent density [kg/m3]

  • r_eff: effective radius [m]

source
CloudMicrophysics.Parameters.RainType
Rain{FT, MS, AR, VT, AC}

The parameters and type for rain

Fields

  • pdf: a struct with size distribution parameters

  • mass: a struct with mass size relation parameters

  • area: a struct with cross section size relation parameers

  • vent: a struct with ventilation coefficients

  • acnv1M: a struct with cloud water to rain autoconversion parameters

  • r0: particle length scale [m]

source
CloudMicrophysics.Parameters.SnowType
Snow{FT, PD, MS, AR, VT, AP, AC}

The parameters and type for snow

Fields

  • pdf: a struct with size distribution parameters

  • mass: a struct with mass size relation parameters

  • area: a struct with cross section size relation parameers

  • vent: a struct with ventilation coefficients

  • aspr: a struct with aspect ratio parameters

  • acnv1M: a struct with ice to snow autoconversion parameters

  • r0: particle length scale [m]

  • T_freeze: freezing temperature of water [K]

  • ρᵢ: snow apparent density [kg/m3]

source
CloudMicrophysics.Parameters.CollisionEffType
CollisionEff{FT}

Collision efficiency parameters for the 1-moment scheme

Fields

  • e_liq_rai: cloud liquid-rain collision efficiency [-]

  • e_liq_sno: cloud liquid-snow collision efficiency [-]

  • e_ice_rai: cloud ice-rain collision efficiency [-]

  • e_ice_sno: cloud ice-snow collision efficiency [-]

  • e_rai_sno: rain-snow collision efficiency [-]

source
CloudMicrophysics.Parameters.KK2000Type
KK2000

The type and parameters for 2-moment precipitation formation by Khairoutdinov and Kogan (2000)

DOI:10.1175/1520-0493(2000)128<0229:ANCPPI>2.0.CO;2

Fields

  • acnv: Autoconversion parameters

  • accr: Accretion parameters

source
CloudMicrophysics.Parameters.AcnvKK2000Type
AcnvKK2000

Khairoutdinov and Kogan (2000) autoconversion parameters

Fields

  • A: Autoconversion coefficient A

  • a: Autoconversion coefficient a

  • b: Autoconversion coefficient b

  • c: Autoconversion coefficient c

source
CloudMicrophysics.Parameters.B1994Type
B1994

The type and parameter for 2-moment precipitation formation by Beheng (1994) DOI: 10.1016/0169-8095(94)90020-5

Fields

  • acnv: Autoconversion coeff C

  • accr: Autoconversion coeff a

source
CloudMicrophysics.Parameters.AcnvB1994Type
AcnvB1994

Beheng (1994) autoconversion parameters

Fields

  • C: Autoconversion coeff C

  • a: Autoconversion coeff a

  • b: Autoconversion coeff b

  • c: Autoconversion coeff c

  • N_0: Autoconversion coeff N_0

  • d_low: Autoconversion coeff d_low

  • d_high: Autoconversion coeff d_high

  • k: Threshold for smooth tranistion steepness

source
CloudMicrophysics.Parameters.TC1980Type
TC1980

The type and parameters for 2-moment precipitation formation by Tripoli and Cotton (1980)

DOI: 10.1175/1520-0450(1980)019<1037:ANIOSF>2.0.CO;2

Fields

  • acnv: Autoconversion parameters

  • accr: Accretion parameters

source
CloudMicrophysics.Parameters.AcnvTC1980Type
AcnvTC1980

Tripoli and Cotton (1980) autoconversion parameters

Fields

  • a: Autoconversion coefficient a

  • b: Autoconversion coefficient b

  • D: Autoconversion coefficient D

  • r_0: Autoconversion coefficient r_0

  • me_liq: Autoconversion coefficient me_liq

  • m0_liq_coeff: Autoconversion coefficient m0liqcoeff

  • k: Threshold for smooth tranistion steepness

source
CloudMicrophysics.Parameters.LD2004Type
LD2004

The type and parameters for 2-moment precipitation formation by Liu and Daum (2004)

DOI: 10.1175/1520-0469(2004)061<1539:POTAPI>2.0.CO;2

Fields

  • R_6C_0: Autoconversion coefficient R6C0

  • E_0: Autoconversion coefficient E_0

  • ρ_w: liquid water density [kg/m3]

  • k: Threshold for smooth tranistion steepness

source
CloudMicrophysics.Parameters.VarTimescaleAcnvType
VarTimescaleAcnv

The type for 2-moment precipitation formation based on the 1-moment parameterization with variable time scale Azimi et al (2023)

Fields

  • τ: Timescale [s]

  • α: Powerlaw coefficient [-]

source
CloudMicrophysics.Parameters.SB2006Type
SB2006

The type and parameters for 2-moment precipitation formation by Seifert and Beheng (2006). The pdf_r type choses between running with or without limiters on raindrop size distribution parameters

DOI: 10.1007/s00703-005-0112-4

Fields

  • pdf_c: Cloud particle size distribution parameters

  • pdf_r: Rain particle size distribution parameters

  • acnv: Autoconversion parameters

  • accr: Accretion parameters

  • self: Rain selfcollection parameters

  • brek: Rain breakup parameters

  • evap: Rain evaporation parameters

source
CloudMicrophysics.Parameters.RainParticlePDF_SB2006Type
RainParticlePDF_SB2006

Rain size distribution parameters from SB2006 but without the limiters

Fields

  • νr: Raindrop size distribution coefficient νr

  • μr: Raindrop size distribution coefficient μr

  • ρw: Cloud liquid water density [kg/m3]

  • ρ0: Reference air density [kg/m3]

  • xr_min: Raindrop minimal mass

source
CloudMicrophysics.Parameters.CloudParticlePDF_SB2006Type
CloudParticlePDF_SB2006

Cloud droplets size distribution parameters from SB2006

Fields

  • νc: Cloud droplet size distribution coefficient νc

  • μc: Cloud droplet size distribution coefficient μc

  • ρw: Cloud liquid water density [kg/m3]

source
CloudMicrophysics.Parameters.AcnvSB2006Type
AcnvSB2006

Autoconversion parameters from SB2006

Fields

  • kcc: Collection kernel coefficient

  • x_star: Minimum mass of rain droplets

  • ρ0: Reference air density [kg/m3]

  • A: Autoconversion correcting function coeff A

  • a: Autoconversion correcting function coeff a

  • b: Autoconversion correcting function coeff b

source
CloudMicrophysics.Parameters.AccrSB2006Type
AccrSB2006

Accretion parameters from SB2006

Fields

  • kcr: Collection kernel coefficient Kcr

  • τ0: Accretion correcting function coefficient τ_0

  • ρ0: Reference air density [kg/m3]

  • c: Accretion correcting function coefficient c

source
CloudMicrophysics.Parameters.SelfColSB2006Type
SelfColSB2006

Rain selfcollection parameters from SB2006

Fields

  • krr: Collection kernel coefficient krr

  • κrr: Collection kernel coefficient kappa rr

  • d: Raindrop self collection coefficient d

source
CloudMicrophysics.Parameters.BreakupSB2006Type
BreakupSB2006

Rain breakup parameters from SB2006

Fields

  • Deq: Raindrop equilibrium mean diamater

  • Dr_th: Raindrop breakup mean diamater threshold

  • kbr: Raindrops breakup coefficient kbr

  • κbr: Raindrops breakup coefficient kappa br

source
CloudMicrophysics.Parameters.EvaporationSB2006Type
EvaporationSB2006

Rain evaporation parameters from SB2006

Fields

  • av: Ventillation coefficient a

  • bv: Ventillation coefficient b

  • α: Rain evapoartion coefficient α

  • β: Rain evapoartion coefficient β

  • ρ0: Reference air density [kg/m3]

source
CloudMicrophysics.Parameters.Blk1MVelTypeRainType
Blk1MVelTypeRain

The type for precipitation terminal velocity from the simple 1-moment scheme for rain

Fields

  • r0: particle length scale [m]

  • ve: rain terminal velocity size relation coefficient [-]

  • Δv: rain terminal velocity size relation coefficient [-]

  • χv: rain terminal velocity size relation coefficient [-]

  • ρw: cloud water density [kg/m3]

  • C_drag: rain drop drag coefficient [-]

  • grav: gravitational acceleration [m/s2]

source
CloudMicrophysics.Parameters.Blk1MVelTypeSnowType
Blk1MVelTypeSnow

The type for precipitation terminal velocity from the simple 1-moment scheme for snow

Fields

  • r0: particle length scale [m]

  • ve: snow terminal velocity size relation coefficient [-]

  • Δv: snow terminal velocity size relation coefficient [-]

  • χv: snow terminal velocity size relation coefficient [-]

  • v0: snow terminal velocity size relation coefficient [m/s]

source
CloudMicrophysics.Parameters.Chen2022VelTypeSmallIceType
Chen2022VelTypeSmallIce

The type for precipitation terminal velocity from Chen et al 2022 for small ice. See Table B3 for parameter definitions. DOI: 10.1016/j.atmosres.2022.106171

Fields

  • A

  • B

  • C

  • E

  • F

  • G

  • cutoff: cutoff for small vs large ice particle dimension [m]

source
CloudMicrophysics.Parameters.Chen2022VelTypeLargeIceType
Chen2022VelTypeLargeIce

The type for precipitation terminal velocity from Chen et al 2022 for large ice. See Table B4 for parameter definitions. DOI: 10.1016/j.atmosres.2022.106171

Fields

  • A

  • B

  • C

  • E

  • F

  • G

  • H

  • cutoff: cutoff for small vs large ice particle dimension [m]

source

Precipitation susceptibility

CloudMicrophysics.PrecipitationSusceptibility.precipitation_susceptibility_autoconversionFunction
precipitation_susceptibility_autoconversion(param_set, scheme, q_liq, q_rai, ρ, N_liq)
  • scheme - type for 2-moment rain autoconversion parameterization
  • q_liq - cloud water specific content
  • q_rai - rain water specific content
  • ρ - air density
  • N_liq - cloud droplet number density

Returns the precipitation susceptibility rates due to autoconversion as a precip_susceptibility_rates object, using automatic differentiation. Works for any 2-moment scheme, as long as autoconversion is defined for it.

source
CloudMicrophysics.PrecipitationSusceptibility.precipitation_susceptibility_accretionFunction
precipitation_susceptibility_accretion(param_set, scheme, q_liq, q_rai, ρ, N_liq)
  • scheme - type for 2-moment rain autoconversion parameterization
  • q_liq - cloud water specific content
  • q_rai - rain water specific content
  • ρ - air density
  • N_liq - cloud droplet number density

Returns the precipitation susceptibility rates due to accretion as a precip_susceptibility_rates object, using automatic differentiation. Works for any 2-moment scheme, as long as accretion is defined for it.

source