Equilibrium Distributions (lbmpy.equilibrium)

This module contains various classes encapsulating equilibrium distributions used in the lattice Boltzmann method. These include both the continuous and the discretized variants of the Maxwellian equilibrium of hydrodynamics. Furthermore, a lightweight wrapper class for custom discrete equilibria is provided. Custom equilibria may also be implemented by manually overriding the abstract base class lbmpy.equilibrium.AbstractEquilibrium.

Abstract Base Class

class AbstractEquilibrium(dim=3)

Abstract Base Class for description of equilibrium distribution functions used in lattice Boltzmann methods.

Equilibrium Representation:

This class provides the common interface for describing equilibrium distribution functions, which is then used by the various method classes in the derivation of collision equations. An equilibrium distribution is defined by either its continuous equation (see continuous_equation) or a set of discrete populations (see discrete_populations and lbmpy.equilibrium.GenericDiscreteEquilibrium). The distribution function may be given either in its regular, absolute form; or only as its deviation from the rest state, represented by the background distribution (see background_distribution).

Computation of Statistical Modes:

The major computational task of an equilbrium class is the computation of the distribution’s statistical modes. For discrete distributions, the subclass lbmpy.equilibrium.GenericDiscreteEquilibrium provides a generic implementation for their computation. For continuous distributions, computation of raw moments, central moments, and cumulants is more complicated, but may also be simplified using special tricks.

As the computation of statistical modes is a time-consuming process, the abstract base class provides caching functionality to avoid recomputing quantities that are already known.

Instructions to Override:

If you wish to model a simple custom discrete distribution, just using the class lbmpy.equilibrium.GenericDiscreteEquilibrium might already be sufficient. If, however, you need to implement more specific functionality, custom properties, a background distribution, etc., or if you wish to model a continuous distribution, you will have to set up a custom subclass of AbstractEquilibrium.

A subclass must implement all abstract properties according to their docstrings. For computation of statistical modes, a large part of the infrastructure is already given in the abstract base class. The public interface for computing e.g. raw moments reduces the computation of polynomial moments to their contained monomials (for details on how moments are represented in lbmpy, see lbmpy.moments). The values of both polynomial and monomial moments, once computed, will be cached per instance of the equilibrium class. To take full advantage of the caching functionality, you will have to override only _monomial_raw_moment() and its central moment and cumulant counterparts. These methods will be called only once for each monomial quantity when it is required for the first time. Afterward, the cached value will be used.

property dim

This distribution’s spatial dimensionality.

abstract property deviation_only

Whether or not this equilibrium distribution is represented only by its deviation from the background distribution.

abstract property continuous_equation

Returns the continuous equation defining this equilibrium distribution, or None if no such equation is available.

abstract property discrete_populations

Returns the discrete populations of this equilibrium distribution as a tuple, or None if none are available.

abstract property background_distribution

Returns this equilibrium distribution’s background distribution, which is the distribution the discrete populations are centered around in the case of zero-centered storage. If no background distribution is available, None must be returned.

abstract property zeroth_order_moment_symbol

Returns a symbol referring to the zeroth-order moment of this distribution, which is the area under it’s curve.

abstract property first_order_moment_symbols

Returns a vector of symbols referring to the first-order moment of this distribution, which is its mean value.

moment(exponent_tuple_or_polynomial)

Returns this equilibrium distribution’s moment specified by exponent_tuple_or_polynomial.

Parameters:

exponent_tuple_or_polynomial – Moment specification, see lbmpy.moments.

moments(exponent_tuples_or_polynomials)

Returns a tuple of this equilibrium distribution’s moments specified by ‘exponent_tuple_or_polynomial’.

Parameters:

exponent_tuples_or_polynomials – Sequence of moment specifications, see lbmpy.moments.

central_moment(exponent_tuple_or_polynomial, frame_of_reference)

Returns this equilibrium distribution’s central moment specified by exponent_tuple_or_polynomial, computed according to the given frame_of_reference.

Parameters:
  • exponent_tuple_or_polynomial – Moment specification, see lbmpy.moments.

  • frame_of_reference – The frame of reference with respect to which the central moment should be computed.

central_moments(exponent_tuples_or_polynomials, frame_of_reference)

Returns a list this equilibrium distribution’s central moments specified by exponent_tuples_or_polynomials, computed according to the given frame_of_reference.

Parameters:
  • exponent_tuples_or_polynomials – Sequence of moment specifications, see lbmpy.moments.

  • frame_of_reference – The frame of reference with respect to which the central moment should be computed.

cumulant(exponent_tuple_or_polynomial, rescale=True)

Returns this equilibrium distribution’s cumulant specified by exponent_tuple_or_polynomial.

Parameters:
  • exponent_tuple_or_polynomial – Moment specification, see lbmpy.moments.

  • rescale – If True, the cumulant value should be multiplied by the zeroth-order moment.

cumulants(exponent_tuples_or_polynomials, rescale=True)

Returns a list of this equilibrium distribution’s cumulants specified by exponent_tuples_or_polynomial.

Parameters:
  • exponent_tuples_or_polynomials – Sequence of moment specifications, see lbmpy.moments.

  • rescale – If True, the cumulant value should be multiplied by the zeroth-order moment.

abstract _monomial_raw_moment(exponents)

See lbmpy.equilibrium.AbstractEquilibrium.moment().

abstract _monomial_central_moment(exponents, frame_of_reference)

See lbmpy.equilibrium.AbstractEquilibrium.central_moment().

abstract _monomial_cumulant(exponents, rescale)

See lbmpy.equilibrium.AbstractEquilibrium.cumulant().

Generic Discrete Equilibria

Use the following class for custom discrete equilibria.

class GenericDiscreteEquilibrium(stencil, equilibrium_pdfs, zeroth_order_moment_symbol, first_order_moment_symbols, deviation_only=False)

Class for encapsulating arbitrary discrete equilibria, given by their equilibrium populations.

This class takes both a stencil and a sequence of populations modelling a discrete distribution function and provides basic functionality for computing and caching that distribution’s statistical modes.

Parameters:
  • stencil – Discrete velocity set, see lbmpy.stencils.LBStencil.

  • equilibrium_pdfs – List of q populations, describing the particle distribution on the discrete velocity set given by the stencil.

  • zeroth_order_moment_symbol – Symbol corresponding to the distribution’s zeroth-order moment, the area under it’s curve (see zeroth_order_moment_symbol).

  • first_order_moment_symbols – Sequence of symbols corresponding to the distribution’s first-order moment, the vector of its mean values (see first_order_moment_symbols).

  • deviation_only – Set to True if the given populations model only the deviation from a rest state, to be used in junction with the zero-centered storage format.

property deviation_only

Whether or not this equilibrium distribution is represented only by its deviation from the background distribution.

property continuous_equation

Always returns None.

property discrete_populations

Returns the discrete populations of this equilibrium distribution as a tuple, or None if none are available.

property background_distribution

Always returns None. To specify a background distribution, override this class.

property zeroth_order_moment_symbol

Returns a symbol referring to the zeroth-order moment of this distribution, which is the area under it’s curve.

property first_order_moment_symbols

Returns a vector of symbols referring to the first-order moment of this distribution, which is its mean value.

Maxwellian Equilibria for Hydrodynamics

The following classes represent the continuous and the discrete variant of the Maxwellian equilibrium for hydrodynamics.

class ContinuousHydrodynamicMaxwellian(dim=3, compressible=True, deviation_only=False, order=None, rho=rho, rho_background=1, delta_rho=delta_rho, u=(u_0, u_1, u_2), v=(v_0, v_1, v_2), c_s_sq=c_s**2)

The standard continuous Maxwellian equilibrium distribution for hydrodynamics.

This class represents the Maxwellian equilibrium distribution of hydrodynamics in its continuous form in \(d\) dimensions [KrugerKK+17]:

\[\Psi \left( \rho, \mathbf{u}, \mathbf{\xi} \right) = \rho \left( \frac{1}{2 \pi c_s^2} \right)^{d/2} \exp \left( \frac{- (\mathbf{\xi} - \mathbf{u})^2 }{2 c_s^2} \right)\]

Beyond this classic, ‘compressible’ form of the equilibrium, an alternative form known as the incompressible equilibrium of the LBM can be obtained by setting the flag compressible=False. The continuous incompressible equilibrium can be expressed as ([GruszczynskiML+20, HL97]):

\[\Psi^{\mathrm{incomp}} \left( \rho, \mathbf{u}, \mathbf{\xi} \right) = \Psi \left( \rho_0, \mathbf{u}, \mathbf{\xi} \right) + \Psi \left( \delta\rho, \mathbf{0}, \mathbf{\xi} \right)\]

Here, \(\rho_0\) (typically \(\rho_0 = 1\)) denotes the background density, and \(\delta\rho\) is the density deviation, such that the total fluid density amounts to \(\rho = \rho_0 + \delta\rho\).

To simplify computations when the zero-centered storage format is used for PDFs, both equilibrium variants can also be expressed in a deviation-only or delta-equilibrium form, which is obtained by subtracting the constant background distribution \(\Psi (\rho_0, \mathbf{0})\). The delta form expresses the equilibrium distribution only by its deviation from the rest state:

\[\begin{split}\delta\Psi \left( \rho, \mathbf{u}, \mathbf{\xi} \right) &= \Psi \left( \rho, \mathbf{u}, \mathbf{\xi} \right) - \Psi \left( \rho_0, \mathbf{0}, \mathbf{\xi} \right) \\ \delta\Psi^{\mathrm{incomp}} \left( \rho, \mathbf{u}, \mathbf{\xi} \right) &= \Psi^{\mathrm{incomp}} \left( \rho, \mathbf{u}, \mathbf{\xi} \right) - \Psi \left( \rho_0, \mathbf{0}, \mathbf{\xi} \right)\end{split}\]
Parameters:
  • dim – Spatial dimensionality

  • compressible – If False, the incompressible equilibrium is created

  • deviation_only – If True, the delta-equilibrium is created

  • order – The discretization order in velocity to which computed statistical modes should be truncated

  • rho – Symbol or value for the density

  • rho_background – Symbol or value for the background density

  • delta_rho – Symbol or value for the density deviation

  • u – Sequence of symbols for the macroscopic velocity

  • v – Sequence of symbols for the particle velocity \(\xi\)

  • c_s_sq – Symbol or value for the squared speed of sound

property deviation_only

Whether or not this equilibrium distribution is represented only by its deviation from the background distribution.

property continuous_equation

Returns the continuous equation defining this equilibrium distribution, or None if no such equation is available.

property zeroth_order_moment_symbol

Returns a symbol referring to the zeroth-order moment of this distribution, which is the area under it’s curve.

property first_order_moment_symbols

Returns a vector of symbols referring to the first-order moment of this distribution, which is its mean value.

property background_distribution

Returns this equilibrium distribution’s background distribution, which is the distribution the discrete populations are centered around in the case of zero-centered storage. If no background distribution is available, None must be returned.

property discrete_populations

Returns the discrete populations of this equilibrium distribution as a tuple, or None if none are available.

central_moment(exponent_tuple_or_polynomial, velocity=None)

Returns this equilibrium distribution’s central moment specified by exponent_tuple_or_polynomial, computed according to the given frame_of_reference.

Parameters:
  • exponent_tuple_or_polynomial – Moment specification, see lbmpy.moments.

  • frame_of_reference – The frame of reference with respect to which the central moment should be computed.

central_moments(exponent_tuples_or_polynomials, velocity=None)

Returns a list this equilibrium distribution’s central moments specified by exponent_tuples_or_polynomials, computed according to the given frame_of_reference.

Parameters:
  • exponent_tuples_or_polynomials – Sequence of moment specifications, see lbmpy.moments.

  • frame_of_reference – The frame of reference with respect to which the central moment should be computed.

cumulant(exponent_tuple_or_polynomial, rescale=True)

Returns this equilibrium distribution’s cumulant specified by exponent_tuple_or_polynomial.

Parameters:
  • exponent_tuple_or_polynomial – Moment specification, see lbmpy.moments.

  • rescale – If True, the cumulant value should be multiplied by the zeroth-order moment.

class DiscreteHydrodynamicMaxwellian(stencil, compressible=True, deviation_only=False, order=2, rho=rho, delta_rho=delta_rho, u=(u_0, u_1, u_2), c_s_sq=c_s**2)

The textbook discretization of the Maxwellian equilibrium distribution of hydrodynamics.

This class represents the default discretization of the Maxwellian in velocity space, computed from the distribution’s expansion in Hermite polynomials (cf. [KrugerKK+17]). In \(d\) dimensions, its populations \(f_i\) on a given stencil \((\mathbf{c}_i)_{i=0,\dots,q-1}\) are given by

\[f_i (\rho, \mathbf{u}) = w_i \rho \left( 1 + \frac{\mathbf{c}_i \cdot \mathbf{u}}{c_s^2} + \frac{(\mathbf{c}_i \cdot \mathbf{u})^2}{2 c_s^4} - \frac{\mathbf{u} \cdot \mathbf{u}}{2 c_s^2} \right).\]

Here \(w_i\) denote the Hermite integration weights, also called lattice weights. The incompressible variant of this distribution [HL97] can be written as

\[f_i^{\mathrm{incomp}} (\rho, \mathbf{u}) = w_i \rho + w_i \rho_0 \left( \frac{\mathbf{c}_i \cdot \mathbf{u}}{c_s^2} + \frac{(\mathbf{c}_i \cdot \mathbf{u})^2}{2 c_s^4} - \frac{\mathbf{u} \cdot \mathbf{u}}{2 c_s^2} \right).\]

Again, for usage with zero-centered PDF storage, both distributions may be expressed in a delta-form by subtracting their values at the background rest state at \(\rho = \rho_0\), \(\mathbf{u} = \mathbf{0}\), which are exactly the lattice weights:

\[\begin{split}\delta f_i &= f_i - w_i \\ \delta f_i^{\mathrm{incomp}} &= f_i^{\mathrm{incomp}} - w_i \\\end{split}\]
Parameters:
  • stencil – Discrete velocity set for the discretization, see lbmpy.stencils.LBStencil

  • compressible – If False, the incompressible equilibrium is created

  • deviation_only – If True, the delta-equilibrium is created

  • order – The discretization order in velocity to which computed statistical modes should be truncated

  • rho – Symbol or value for the density

  • delta_rho – Symbol or value for the density deviation

  • u – Sequence of symbols for the macroscopic velocity

  • c_s_sq – Symbol or value for the squared speed of sound

property deviation_only

Whether or not this equilibrium distribution is represented only by its deviation from the background distribution.

property background_distribution

Returns the discrete Maxwellian background distribution, which amounts exactly to the lattice weights.

central_moment(exponent_tuple_or_polynomial, velocity=None)

Returns this equilibrium distribution’s central moment specified by exponent_tuple_or_polynomial, computed according to the given frame_of_reference.

Parameters:
  • exponent_tuple_or_polynomial – Moment specification, see lbmpy.moments.

  • frame_of_reference – The frame of reference with respect to which the central moment should be computed.

central_moments(exponent_tuples_or_polynomials, velocity=None)

Returns a list this equilibrium distribution’s central moments specified by exponent_tuples_or_polynomials, computed according to the given frame_of_reference.

Parameters:
  • exponent_tuples_or_polynomials – Sequence of moment specifications, see lbmpy.moments.

  • frame_of_reference – The frame of reference with respect to which the central moment should be computed.

cumulant(exponent_tuple_or_polynomial, rescale=True)

Returns this equilibrium distribution’s cumulant specified by exponent_tuple_or_polynomial.

Parameters:
  • exponent_tuple_or_polynomial – Moment specification, see lbmpy.moments.

  • rescale – If True, the cumulant value should be multiplied by the zeroth-order moment.