Moment Transforms (lbmpy.moment_transforms)

The moment_transforms module provides an ecosystem for transformation of quantities between lattice Boltzmann population space and other spaces of observable quantities. Currently, transforms to raw and central moment space are available.

The common base class lbmpy.moment_transforms.AbstractMomentTransform defines the interface all transform classes share. This interface, together with the fundamental principles all transforms adhere to, is explained in the following.

Principles of Collision Space Transforms

Each class of this module implements a bijective map \(\mathcal{M}\) from population space to raw moment or central moment space, capable of transforming the particle distribution function \(\mathbf{f}\) to (almost) arbitrary user-defined moment sets.

Monomial and Polynomial Moments

We discriminate monomial and polynomial moments. The monomial moments are the canonical basis of their corresponding space. Polynomial moments are defined as linear combinations of this basis. Monomial moments are typically expressed by exponent tuples \((\alpha, \beta, \gamma)\). The monomial raw moments are defined as

\[m_{\alpha \beta \gamma} = \sum_{i = 0}^{q - 1} f_i c_{i,x}^{\alpha} c_{i,y}^{\beta} c_{i,z}^{\gamma}\]

and the monomial central moments are given by

\[\kappa_{\alpha \beta \gamma} = \sum_{i = 0}^{q - 1} f_i (c_{i,x} - u_x)^{\alpha} (c_{i,y} - u_y)^{\beta} (c_{i,z} - u_z)^{\gamma}.\]

Polynomial moments are, in turn, expressed by polynomial expressions \(p\) in the coordinate symbols \(x\), \(y\) and \(z\). An exponent tuple \((\alpha, \beta, \gamma)\) corresponds directly to the mixed monomial expression \(x^{\alpha} y^{\beta} z^{\gamma}\). Polynomial moments are then constructed as linear combinations of these monomials. For example, the polynomial

\[p(x,y,z) = x^2 + y^2 + z^2 + 1.\]

defines both the polynomial raw moment

\[M = m_{200} + m_{020} + m_{002}\]

and central moment

\[K = \kappa_{200} + \kappa_{020} + \kappa_{002}.\]

Transformation Matrices

The collision space basis for an MRT LB method on a \(DdQq\) lattice is spanned by a set of \(q\) polynomial quantities, given by polynomials \(\left( p_i \right)_{i=0, \dots, q-1}\). Through the polynomials, a polynomialization matrix \(P\) is defined such that

\[\begin{split}\mathbf{M} &= P \mathbf{m} \\ \mathbf{K} &= P \mathbf{\kappa}\end{split}\]

Both rules can also be written using matrix multiplication, such that

\[\mathbf{m} = M \mathbf{f} \qquad \mathbf{\kappa} = K \mathbf{f}.\]

Further, there exists a mapping from raw to central moment space using (monomial or polynomial) shift matrices (see set_up_shift_matrix) such that

\[\mathbf{\kappa} = N \mathbf{m} \qquad \mathbf{K} = N^P \mathbf{M}.\]

Working with the transformation matrices, those mappings can easily be inverted. This module provides classes that derive efficient implementations of these transformations.

Moment Aliasing

For a collision space transform to be invertible, exactly \(q\) independent polynomial quantities of the collision space must be chosen. In that case, since all transforms discussed here are linear, the space of populations is isomorphic to the chosen moment space. The important word here is ‘independent’. Even if \(q\) distinct moment polynomials are chosen, due to discrete lattice artifacts, they might not span the entire collision space if chosen wrongly. The discrete lattice structure gives rise to moment aliasing effects. The most simple such effect occurs in the monomial raw moments, where are all nonzero even and all odd exponents are essentially the same. For example, we have \(m_{400} = m_{200}\) or \(m_{305} = m_{101}\). This rule holds on every direct-neighborhood stencil. Sparse stencils, like \(D3Q15\), may introduce additional aliases. On the \(D3Q15\) stencil, due to its structure, we have also \(m_{112} = m_{110}\) and even \(m_{202} = m_{220} = m_{022} = m_{222}\).

Including aliases in a set of monomial moments will lead to a non-invertible transform and is hence not possible. In polynomials, however, aliases may occur without breaking the transform. Several established sets of polynomial moments, like in orthogonal raw moment space MRT methods, will comprise \(q\) polynomials \(\mathbf{M}\) that in turn are built out of \(r > q\) monomial moments \(\mathbf{m}\). In that set of monomials, non-independent moments have to be included by definition. Since the full transformation matrix \(M^P = PM\) is still invertible, this is not a problem in general. If, however, the two steps of the transform are computed separately, it becomes problematic, as the matrices \(M \in \mathbb{R}^{r \times q}\) and \(P \in \mathbb{R}^{q \times r}\) are not invertible on their own.

But there is a remedy. By eliminating aliases from the moment polynomials, they can be reduced to a new set of polynomials based on a non-aliased reduced vector of monomial moments \(\tilde{\mathbf{m}} \in \mathbb{R}^{q}\), expressed through the reduced matrix \(\tilde{P}\).

Interfaces and Usage

Construction

All moment transform classes expect either a sequence of exponent tuples or a sequence of polynomials that define the set of quantities spanning the destination space. If polynomials are given, the exponent tuples of the underlying set of monomials are extracted automatically. Depending on the implementation, moment aliases may be eliminated in the process, altering both the polynomials and the set of monomials.

Forward and Backward Transform

The core functionality of the transform classes is given by the methods forward_transform and backward_transform. They return assignment collections containing the equations to compute moments from populations, and vice versa.

Symbols Used

Unless otherwise specified by the user, monomial and polynomial quantities occur in the transformation equations according to the naming conventions listed in the following table:

Monomial

Polynomial

Pre-Collision Raw Moment

\(m_{\alpha \beta \gamma}\)

\(M_{i}\)

Post-Collision Raw Moment

\(m_{post, \alpha \beta \gamma}\)

\(M_{post, i}\)

Pre-Collision Central Moment

\(\kappa_{\alpha \beta \gamma}\)

\(K_{i}\)

Post-Collision Central Moment

\(\kappa_{post, \alpha \beta \gamma}\)

\(K_{post, i}\)

These symbols are also exposed by the member properties pre_collision_symbols, post_collision_symbols, pre_collision_monomial_symbols and post_collision_monomial_symbols.

Forward Transform

Implementations of the lbmpy.moment_transforms.AbstractMomentTransform.forward_transform method derive equations for the transform from population to moment space, to compute pre-collision moments from pre-collision populations. The returned AssignmentCollection has the pre_collision_symbols as left-hand sides of its main assignments, computed from the given pdf_symbols and, potentially, the macroscopic density and velocity. Depending on the implementation, the pre_collision_monomial_symbols may be part of the assignment collection in the form of subexpressions.

Backward Transform

Implementations of lbmpy.moment_transforms.AbstractMomentTransform.backward_transform contain the post-collision polynomial quantities as free symbols of their equation right-hand sides, and compute the post-collision populations from those. The PDF symbols are the left-hand sides of the main assignments.

Absorption of Conserved Quantity Equations

Transformations from the population space to any space of observable quantities may absorb the equations defining the macroscopic quantities entering the equilibrium (typically the density \(\rho\) and the velocity \(\mathbf{u}\)). This means that the forward_transform will possibly rewrite the assignments given in the constructor argument conserved_quantity_equations to reduce the total operation count. For example, in the transformation step from populations to raw moments (see lbmpy.moment_transforms.PdfsToMomentsByChimeraTransform), \(\rho\) can be aliased as the zeroth-order moment \(m_{000}\). Assignments to the conserved quantities will then be part of the AssignmentCollection returned by forward_transform and need not be added to the collision rule separately.

Simplification

Both forward_transform and backward_transform expect a keyword argument simplification which can be used to direct simplification steps applied during the derivation of the transformation equations. Possible values are:

  • False or 'none': No simplification is to be applied

  • True or 'default': A default simplification strategy specific to the implementation is applied.

    The actual simplification steps depend strongly on the nature of the equations. They are defined by the implementation. It is the responsibility of the implementation to select the most effective simplification strategy.

  • 'default_with_cse': Same as 'default', but with an additional pass of common subexpression elimination.

Working With Monomials

In certain situations, we want the forward_transform to yield equations for the monomial symbols \(m_{\alpha \beta \gamma}\) and \(\kappa_{\alpha \beta \gamma}\) only, and the backward_transform to return equations that also expect these symbols as input. In this case, it is not sufficient to pass a set of monomials or exponent tuples to the constructor, as those are still treated as polynomials internally. Instead, both transform methods expose keyword arguments return_monomials and start_from_monomials, respectively. If set to true, equations in the monomial moments are returned. They are best used only together with the exponent_tuples constructor argument to have full control over the monomials. If polynomials are passed to the constructor, the behaviour of these flags is generally not well-defined, especially in the presence of aliases.

The Transform Classes

Abstract Base Class

class AbstractMomentTransform(stencil, equilibrium_density, equilibrium_velocity, moment_exponents=None, moment_polynomials=None, conserved_quantity_equations=None, background_distribution=None, pre_collision_symbol_base=None, post_collision_symbol_base=None, pre_collision_monomial_symbol_base=None, post_collision_monomial_symbol_base=None)

Abstract Base Class for classes providing transformations between moment spaces.

property pre_collision_symbols

List of symbols corresponding to the pre-collision quantities that will be the left-hand sides of assignments returned by forward_transform().

property post_collision_symbols

List of symbols corresponding to the post-collision quantities that are input to the right-hand sides of assignments returned by:func:backward_transform.

property pre_collision_monomial_symbols

List of symbols corresponding to the pre-collision monomial quantities that might exist as left-hand sides of subexpressions in the assignment collection returned by forward_transform().

property post_collision_monomial_symbols

List of symbols corresponding to the post-collision monomial quantities that might exist as left-hand sides of subexpressions in the assignment collection returned by backward_transform().

abstract forward_transform(*args, **kwargs)

Implemented in a subclass, will return the forward transform equations.

abstract backward_transform(*args, **kwargs)

Implemented in a subclass, will return the backward transform equations.

property absorbs_conserved_quantity_equations

Whether or not the given conserved quantity equations will be included in the assignment collection returned by forward_transform(), possibly in simplified form.

Moment Space Transforms

class PdfsToMomentsByMatrixTransform(stencil, moment_polynomials, equilibrium_density, equilibrium_velocity, conserved_quantity_equations=None, **kwargs)

Transform between populations and moment space spanned by a polynomial basis, using matrix-vector multiplication.

property absorbs_conserved_quantity_equations

Whether or not the given conserved quantity equations will be included in the assignment collection returned by forward_transform(), possibly in simplified form.

forward_transform(pdf_symbols, simplification=True, subexpression_base='sub_f_to_M', return_monomials=False)

Returns an assignment collection containing equations for pre-collision polynomial moments, expressed in terms of the pre-collision populations by matrix-multiplication.

The moment transformation matrix \(M\) provided by lbmpy.moments.moment_matrix() is used to compute the pre-collision moments as \(\mathbf{M} = M \cdot \mathbf{f}\), which is returned element-wise.

Parameters:
  • pdf_symbols – List of symbols that represent the pre-collision populations

  • simplification – Simplification specification. See AbstractMomentTransform

  • subexpression_base – The base name used for any subexpressions of the transformation.

  • return_monomials – Return equations for monomial moments. Use only when specifying moment_exponents in constructor!

backward_transform(pdf_symbols, simplification=True, subexpression_base='sub_k_to_f', start_from_monomials=False)

Returns an assignment collection containing equations for post-collision populations, expressed in terms of the post-collision polynomial moments by matrix-multiplication.

The moment transformation matrix \(M\) provided by lbmpy.moments.moment_matrix() is inverted and used to compute the pre-collision moments as \(\mathbf{f}^{\ast} = M^{-1} \cdot \mathbf{M}_{\mathrm{post}}\), which is returned element-wise.

Simplifications

If simplification is enabled, the equations for populations \(f_i\) and \(f_{\bar{i}}\) of opposite stencil directions \(\mathbf{c}_i\) and \(\mathbf{c}_{\bar{i}} = - \mathbf{c}_i\) are split into their symmetric and antisymmetric parts \(f_i^{\mathrm{sym}}, f_i^{\mathrm{anti}}\), such that

\[ \begin{align}\begin{aligned}f_i = f_i^{\mathrm{sym}} + f_i^{\mathrm{anti}}\\f_{\bar{i}} = f_i^{\mathrm{sym}} - f_i^{\mathrm{anti}}\end{aligned}\end{align} \]
Parameters:
  • pdf_symbols – List of symbols that represent the post-collision populations

  • simplification – Simplification specification. See AbstractMomentTransform

  • subexpression_base – The base name used for any subexpressions of the transformation.

  • start_from_monomials – Return equations for monomial moments. Use only when specifying moment_exponents in constructor!

class PdfsToMomentsByChimeraTransform(stencil, moment_polynomials, equilibrium_density, equilibrium_velocity, conserved_quantity_equations=None, **kwargs)

Transform between populations and moment space spanned by a polynomial basis, using the raw-moment chimera transform in the forward direction and matrix-vector multiplication in the backward direction.

property absorbs_conserved_quantity_equations

Whether or not the given conserved quantity equations will be included in the assignment collection returned by forward_transform(), possibly in simplified form.

get_cq_to_moment_symbols_dict(moment_symbol_base)

Returns a dictionary mapping the density and velocity symbols to the correspondig zeroth- and first-order raw moment symbols

forward_transform(pdf_symbols, simplification=True, subexpression_base='sub_f_to_m', return_monomials=False)

Returns an assignment collection containing equations for pre-collision polynomial moments, expressed in terms of the pre-collision populations, using the raw moment chimera transform.

The chimera transform for raw moments is given by [GSchonherrPK15] :

\[\begin{split}f_{xyz} &:= f_i \text{ such that } c_i = (x,y,z)^T \\ m_{xy|\gamma} &:= \sum_{z \in \{-1, 0, 1\} } f_{xyz} \cdot z^{\gamma} \\ m_{x|\beta \gamma} &:= \sum_{y \in \{-1, 0, 1\}} m_{xy|\gamma} \cdot y^{\beta} \\ m_{\alpha \beta \gamma} &:= \sum_{x \in \{-1, 0, 1\}} m_{x|\beta \gamma} \cdot x^{\alpha}\end{split}\]

The obtained raw moments are afterward combined to the desired polynomial moments.

Conserved Quantity Equations

If given, this transform absorbs the conserved quantity equations and simplifies them using the monomial raw moment equations, if simplification is enabled.

De-Aliasing

If more than \(q\) monomial moments are extracted from the polynomial set, the polynomials are de-aliased by eliminating aliases according to the stencil using lbmpy.moments.non_aliased_polynomial_raw_moments.

Simplification

If simplification is enabled, the absorbed conserved quantity equations are - if possible - rewritten using the monomial symbols. If the conserved quantities originate somewhere else than in the lower-order moments (like from an external field), they are not affected by this simplification. Furthermore, aliases and constants are propagated in the chimera equations.

Parameters:
  • pdf_symbols – List of symbols that represent the pre-collision populations

  • simplification – Simplification specification. See AbstractMomentTransform

  • subexpression_base – The base name used for any subexpressions of the transformation.

  • return_monomials – Return equations for monomial moments. Use only when specifying moment_exponents in constructor!

backward_transform(pdf_symbols, simplification=True, subexpression_base='sub_k_to_f', start_from_monomials=False)

Returns an assignment collection containing equations for post-collision populations, expressed in terms of the post-collision polynomial moments by matrix-multiplication.

The post-collision monomial moments \(\mathbf{m}_{\mathrm{post}}\) are first obtained from the polynomials. Then, the monomial transformation matrix \(M_r\) provided by lbmpy.moments.moment_matrix() is inverted and used to compute the post-collision populations as \(\mathbf{f}_{\mathrm{post}} = M_r^{-1} \cdot \mathbf{m}_{\mathrm{post}}\).

De-Aliasing

See PdfsToMomentsByChimeraTransform.forward_transform.

Simplifications

If simplification is enabled, the equations for populations \(f_i\) and \(f_{\bar{i}}\) of opposite stencil directions \(\mathbf{c}_i\) and \(\mathbf{c}_{\bar{i}} = - \mathbf{c}_i\) are split into their symmetric and antisymmetric parts \(f_i^{\mathrm{sym}}, f_i^{\mathrm{anti}}\), such that

\[ \begin{align}\begin{aligned}f_i = f_i^{\mathrm{sym}} + f_i^{\mathrm{anti}}\\f_{\bar{i}} = f_i^{\mathrm{sym}} - f_i^{\mathrm{anti}}\end{aligned}\end{align} \]
Parameters:
  • pdf_symbols – List of symbols that represent the post-collision populations

  • simplification – Simplification specification. See AbstractMomentTransform

  • subexpression_base – The base name used for any subexpressions of the transformation.

  • start_from_monomials – Return equations for monomial moments. Use only when specifying moment_exponents in constructor!

Central Moment Space Transforms

class PdfsToCentralMomentsByMatrix(stencil, moment_polynomials, equilibrium_density, equilibrium_velocity, **kwargs)

Transform from populations to central moment space by matrix-vector multiplication.

forward_transform(pdf_symbols, simplification=True, subexpression_base='sub_f_to_k', return_monomials=False)

Returns an assignment collection containing equations for pre-collision polynomial central moments, expressed in terms of the pre-collision populations by matrix-multiplication.

The central moment transformation matrix \(K\) provided by lbmpy.moments.moment_matrix() is used to compute the pre-collision moments as \(\mathbf{K} = K \cdot \mathbf{f}\), which are returned element-wise.

Parameters:
  • pdf_symbols – List of symbols that represent the pre-collision populations

  • simplification – Simplification specification. See AbstractMomentTransform

  • subexpression_base – The base name used for any subexpressions of the transformation.

  • return_monomials – Return equations for monomial moments. Use only when specifying moment_exponents in constructor!

backward_transform(pdf_symbols, simplification=True, subexpression_base='sub_k_to_f', start_from_monomials=False)

Returns an assignment collection containing equations for post-collision populations, expressed in terms of the post-collision polynomial central moments by matrix-multiplication.

The moment transformation matrix \(K\) provided by lbmpy.moments.moment_matrix() is inverted and used to compute the pre-collision moments as \(\mathbf{f}^{\ast} = K^{-1} \cdot \mathbf{K}_{\mathrm{post}}\), which is returned element-wise.

Parameters:
  • pdf_symbols – List of symbols that represent the post-collision populations

  • simplification – Simplification specification. See AbstractMomentTransform

  • subexpression_base – The base name used for any subexpressions of the transformation.

  • start_from_monomials – Return equations for monomial moments. Use only when specifying moment_exponents in constructor!

class BinomialChimeraTransform(stencil, moment_polynomials, equilibrium_density, equilibrium_velocity, conserved_quantity_equations=None, **kwargs)

Transform from populations to central moments using a chimera transform implementing the binomial expansion.

property absorbs_conserved_quantity_equations

Whether or not the given conserved quantity equations will be included in the assignment collection returned by forward_transform(), possibly in simplified form.

forward_transform(pdf_symbols, simplification=True, subexpression_base='sub_f_to_k', return_monomials=False)

Returns equations for polynomial central moments, computed from pre-collision populations through a cascade of three steps.

First, the monomial raw moment vector \(\mathbf{m}\) is computed using the raw-moment chimera transform (see lbmpy.moment_transforms.PdfsToMomentsByChimeraTransform).

Second, we obtain monomial central moments from monomial raw moments using the binomial chimera transform:

\[\begin{split}\kappa_{ab|\gamma} &:= \sum_{c = 0}^{\gamma} \binom{\gamma}{c} v_z^{\gamma - c} m_{abc} \\ \kappa_{a|\beta\gamma} &:= \sum_{b = 0}^{\beta} \binom{\beta}{b} v_z^{\beta - b} \kappa_{ab|\gamma} \\ \kappa_{\alpha\beta\gamma} &:= \sum_{a = 0}^{\alpha} \binom{\alpha}{a} v_z^{\alpha - a} \kappa_{a|\beta\gamma} \\\end{split}\]

Lastly, the polynomial central moments are computed using the polynomialization matrix as \(\mathbf{K} = P \mathbf{\kappa}\).

Conserved Quantity Equations

If given, this transform absorbs the conserved quantity equations and simplifies them using the raw moment equations, if simplification is enabled.

Simplification

If simplification is enabled, the absorbed conserved quantity equations are - if possible - rewritten using the monomial symbols. If the conserved quantities originate somewhere else than in the lower-order moments (like from an external field), they are not affected by this simplification.

The raw moment chimera transform is simplified by propagation of aliases.

The equations of the binomial chimera transform are simplified by expressing conserved raw moments in terms of the conserved quantities, and subsequent propagation of aliases, constants, and any expressions that are purely products of conserved quantities.

De-Aliasing

If more than \(q\) monomial moments are extracted from the polynomial set, they are de-aliased and reduced to a set of only \(q\) moments using the same rules as for raw moments. For polynomialization, a special reduced matrix \(\tilde{P}\) is used, which is computed using lbmpy.moments.central_moment_reduced_monomial_to_polynomial_matrix.

Parameters:
  • pdf_symbols – List of symbols that represent the pre-collision populations

  • simplification – Simplification specification. See AbstractMomentTransform

  • subexpression_base – The base name used for any subexpressions of the transformation.

  • return_monomials – Return equations for monomial moments. Use only when specifying moment_exponents in constructor!

backward_transform(pdf_symbols, simplification=True, subexpression_base='sub_k_to_f', start_from_monomials=False)

Returns an assignment collection containing equations for post-collision populations, expressed in terms of the post-collision polynomial central moments by three steps.

The post-collision monomial central moments \(\mathbf{\kappa}_{\mathrm{post}}\) are first obtained from the polynomials through multiplication with \(P^{-1}\).

Afterward, monomial post-collision raw moments are obtained from monomial central moments using the binomial chimera transform:

\[\begin{split}m^{\ast}_{ab|\gamma} &:= \sum_{c = 0}^{\gamma} \binom{\gamma}{c} v_z^{\gamma - c} \kappa^{\ast}_{abc} \\ m^{\ast}_{a|\beta\gamma} &:= \sum_{b = 0}^{\beta} \binom{\beta}{b} v_z^{\beta - b} m^{\ast}_{ab|\gamma} \\ m^{\ast}_{\alpha\beta\gamma} &:= \sum_{a = 0}^{\alpha} \binom{\alpha}{a} v_z^{\alpha - a} m^{\ast}_{a|\beta\gamma} \\\end{split}\]

Finally, the monomial raw moment transformation matrix \(M_r\) provided by lbmpy.moments.moment_matrix() is inverted and used to compute the pre-collision moments as \(\mathbf{f}_{\mathrm{post}} = M_r^{-1} \cdot \mathbf{m}_{\mathrm{post}}\).

De-Aliasing:

See PdfsToCentralMomentsByShiftMatrix.forward_transform.

Simplifications

If simplification is enabled, the inverse shift matrix equations are simplified by recursively inserting lower-order moments into equations for higher-order moments. To this end, these equations are factored recursively by the velocity symbols.

The equations of the binomial chimera transform are simplified by propagation of aliases.

Further, the equations for populations \(f_i\) and \(f_{\bar{i}}\) of opposite stencil directions \(\mathbf{c}_i\) and \(\mathbf{c}_{\bar{i}} = - \mathbf{c}_i\) are split into their symmetric and antisymmetric parts \(f_i^{\mathrm{sym}}, f_i^{\mathrm{anti}}\), such that

\[ \begin{align}\begin{aligned}f_i = f_i^{\mathrm{sym}} + f_i^{\mathrm{anti}}\\f_{\bar{i}} = f_i^{\mathrm{sym}} - f_i^{\mathrm{anti}}\end{aligned}\end{align} \]
Parameters:
  • pdf_symbols – List of symbols that represent the post-collision populations

  • simplification – Simplification specification. See AbstractMomentTransform

  • subexpression_base – The base name used for any subexpressions of the transformation.

  • start_from_monomials – Return equations for monomial moments. Use only when specifying moment_exponents in constructor!

class FastCentralMomentTransform(stencil, moment_polynomials, equilibrium_density, equilibrium_velocity, conserved_quantity_equations=None, **kwargs)

Transform from populations to central moments, using the fast central-moment transform equations introduced by [GSchonherrPK15].

Attention: The fast central moment transform has originally been designed for the D3Q27 stencil, and is also tested and safely usable with D2Q9 and D3Q19. While the forward- transform does not pose any problems, the backward equations may be inefficient, or even not cleanly derivable for other stencils. Use with care!

forward_transform(pdf_symbols, simplification=True, subexpression_base='sub_f_to_k', return_monomials=False)

Returns an assignment collection containing equations for pre-collision polynomial central moments, expressed in terms of the pre-collision populations.

The monomial central moments are computed from populations through the central-moment chimera transform:

\[\begin{split}f_{xyz} &:= f_i \text{ such that } c_i = (x,y,z)^T \\ \kappa_{xy|\gamma} &:= \sum_{z \in \{-1, 0, 1\} } f_{xyz} \cdot (z - u_z)^{\gamma} \\ \kappa_{x|\beta \gamma} &:= \sum_{y \in \{-1, 0, 1\}} \kappa_{xy|\gamma} \cdot (y - u_y)^{\beta} \\ \kappa_{\alpha \beta \gamma} &:= \sum_{x \in \{-1, 0, 1\}} \kappa_{x|\beta \gamma} \cdot (x - u_x)^{\alpha}\end{split}\]

The polynomial moments are afterward computed from the monomials by matrix-multiplication using the polynomialization matrix \(P\).

De-Aliasing

If more than \(q\) monomial moments are extracted from the polynomial set, they are de-aliased and reduced to a set of only \(q\) moments using the same rules as for raw moments. For polynomialization, a special reduced matrix \(\tilde{P}\) is used, which is computed using lbmpy.moments.central_moment_reduced_monomial_to_polynomial_matrix.

Parameters:
  • pdf_symbols – List of symbols that represent the pre-collision populations

  • simplification – Simplification specification. See AbstractMomentTransform

  • subexpression_base – The base name used for any subexpressions of the transformation.

  • return_monomials – Return equations for monomial moments. Use only when specifying moment_exponents in constructor!

backward_transform(pdf_symbols, simplification=True, subexpression_base='sub_k_to_f', start_from_monomials=False)

Returns an assignment collection containing equations for post-collision populations, expressed in terms of the post-collision polynomial central moments using the backward fast central moment transform.

First, monomial central moments are obtained from the polynomial moments by multiplication with \(P^{-1}\). Then, the elementwise equations of the matrix multiplication \(K^{-1} \cdot \mathbf{K}\) with the monomial central moment matrix (see PdfsToCentralMomentsByMatrix) are recursively simplified by extracting certain linear combinations of velocities, to obtain equations similar to the ones given in [GSchonherrPK15].

The backward transform is designed for D3Q27, inherently generalizes to D2Q9, and is tested for D3Q19. It also returns correct equations for D3Q15, whose efficiency is however questionable.

De-Aliasing:

See FastCentralMomentTransform.forward_transform.

Parameters:
  • pdf_symbols – List of symbols that represent the post-collision populations

  • simplification – Simplification specification. See AbstractMomentTransform

  • subexpression_base – The base name used for any subexpressions of the transformation.

  • start_from_monomials – Return equations for monomial moments. Use only when specifying moment_exponents in constructor!

class PdfsToCentralMomentsByShiftMatrix(stencil, moment_polynomials, equilibrium_density, equilibrium_velocity, conserved_quantity_equations=None, **kwargs)

Transform from populations to central moments using a shift matrix.

property absorbs_conserved_quantity_equations

Whether or not the given conserved quantity equations will be included in the assignment collection returned by forward_transform(), possibly in simplified form.

forward_transform(pdf_symbols, simplification=True, subexpression_base='sub_f_to_k', return_monomials=False)

Returns equations for polynomial central moments, computed from pre-collision populations through a cascade of three steps.

First, the monomial raw moment vector \(\mathbf{m}\) is computed using the raw-moment chimera transform (see lbmpy.moment_transforms.PdfsToMomentsByChimeraTransform). Then, the monomial shift matrix \(N\) provided by lbmpy.moments.set_up_shift_matrix is used to compute the monomial central moment vector as \(\mathbf{\kappa} = N \mathbf{m}\). Lastly, the polynomial central moments are computed using the polynomialization matrix as \(\mathbf{K} = P \mathbf{\kappa}\).

Conserved Quantity Equations

If given, this transform absorbs the conserved quantity equations and simplifies them using the raw moment equations, if simplification is enabled.

Simplification

If simplification is enabled, the absorbed conserved quantity equations are - if possible - rewritten using the monomial symbols. If the conserved quantities originate somewhere else than in the lower-order moments (like from an external field), they are not affected by this simplification.

The relations between conserved quantities and raw moments are used to simplify the equations obtained from the shift matrix. Further, these equations are simplified by recursively inserting lower-order moments into equations for higher-order moments.

De-Aliasing

If more than \(q\) monomial moments are extracted from the polynomial set, they are de-aliased and reduced to a set of only \(q\) moments using the same rules as for raw moments. For polynomialization, a special reduced matrix \(\tilde{P}\) is used, which is computed using lbmpy.moments.central_moment_reduced_monomial_to_polynomial_matrix.

Parameters:
  • pdf_symbols – List of symbols that represent the pre-collision populations

  • simplification – Simplification specification. See AbstractMomentTransform

  • subexpression_base – The base name used for any subexpressions of the transformation.

  • return_monomials – Return equations for monomial moments. Use only when specifying moment_exponents in constructor!

backward_transform(pdf_symbols, simplification=True, subexpression_base='sub_k_to_f', start_from_monomials=False)

Returns an assignment collection containing equations for post-collision populations, expressed in terms of the post-collision polynomial central moments by matrix-multiplication including the shift matrix.

The post-collision monomial central moments \(\mathbf{\kappa}_{\mathrm{post}}\) are first obtained from the polynomials through multiplication with \(P^{-1}\). The shift-matrix is inverted as well, to obtain the monomial raw moments as \(\mathbf{m}_{post} = N^{-1} \mathbf{\kappa}_{post}\). Finally, the monomial raw moment transformation matrix \(M_r\) provided by lbmpy.moments.moment_matrix() is inverted and used to compute the pre-collision moments as \(\mathbf{f}_{\mathrm{post}} = M_r^{-1} \cdot \mathbf{m}_{\mathrm{post}}\).

De-Aliasing:

See PdfsToCentralMomentsByShiftMatrix.forward_transform.

Simplifications

If simplification is enabled, the inverse shift matrix equations are simplified by recursively inserting lower-order moments into equations for higher-order moments. To this end, these equations are factored recursively by the velocity symbols.

Further, the equations for populations \(f_i\) and \(f_{\bar{i}}\) of opposite stencil directions \(\mathbf{c}_i\) and \(\mathbf{c}_{\bar{i}} = - \mathbf{c}_i\) are split into their symmetric and antisymmetric parts \(f_i^{\mathrm{sym}}, f_i^{\mathrm{anti}}\), such that

\[ \begin{align}\begin{aligned}f_i = f_i^{\mathrm{sym}} + f_i^{\mathrm{anti}}\\f_{\bar{i}} = f_i^{\mathrm{sym}} - f_i^{\mathrm{anti}}\end{aligned}\end{align} \]
Parameters:
  • pdf_symbols – List of symbols that represent the post-collision populations

  • simplification – Simplification specification. See AbstractMomentTransform

  • subexpression_base – The base name used for any subexpressions of the transformation.

  • start_from_monomials – Return equations for monomial moments. Use only when specifying moment_exponents in constructor!

Cumulant Space Transforms

class CentralMomentsToCumulantsByGeneratingFunc(stencil, cumulant_polynomials, equilibrium_density, equilibrium_velocity, cumulant_exponents=None, pre_collision_symbol_base='C', post_collision_symbol_base='C_post', pre_collision_monomial_symbol_base='c', post_collision_monomial_symbol_base='c_post', **kwargs)
property required_central_moments

The required central moments as a sorted list of exponent tuples

forward_transform(central_moment_base='kappa', simplification=True, subexpression_base='sub_k_to_C', return_monomials=False)

Implemented in a subclass, will return the forward transform equations.

backward_transform(central_moment_base='kappa_post', simplification=True, subexpression_base='sub_C_to_k', start_from_monomials=False)

Implemented in a subclass, will return the backward transform equations.