Collision models (lbmpy.methods)

This module defines the classes defining all types of lattice Boltzmann methods available in lbmpy, together with a set of factory functions used to create their instances. The factory functions are organized in three levels of abstraction. Objects of the method classes should be created only through these factory functions, never manually.

Methods in lbmpy can be distinguished into three categories:
  • Raw Moment-based methods, including the classical single relaxation-time (SRT, BGK), two relaxation-time (TRT) and multiple relaxation-time (MRT) methods, as well as entropic variants of the TRT method.

  • Central Moment-based methods, which are multiple relaxation-time methods using relaxation in central moment space.

  • Cumulant-based methods, multiple relaxation-time methods using relaxation in cumulant space.

Abstract LB Method Base Class

class LbmCollisionRule(lb_method, *args, **kwargs)

A pystencils AssignmentCollection that additionally holds an AbstractLbMethod

class AbstractLbMethod(stencil)

Abstract base class for all LBM methods.

property stencil

Discrete set of velocities, represented by lbmpy.stencils.LBStencil

property dim

The method’s spatial dimensionality

property pre_collision_pdf_symbols

Tuple of symbols representing the pdf values before collision

property post_collision_pdf_symbols

Tuple of symbols representing the pdf values after collision

abstract property relaxation_rates

Tuple containing the relaxation rates of each moment

property relaxation_matrix

Returns a qxq diagonal matrix which contains the relaxation rate for each moment on the diagonal

property symbolic_relaxation_matrix

Returns a qxq diagonal matrix which contains the relaxation rate for each moment on the diagonal. In contrast to the normal relaxation matrix all numeric values are replaced by sympy symbols

property subs_dict_relaxation_rate

returns a dictionary which maps the replaced numerical relaxation rates to its original numerical value

abstract property conserved_quantity_computation

Returns an instance of class lbmpy.methods.AbstractConservedQuantityComputation

abstract property weights

Returns a sequence of weights, one for each lattice direction

abstract get_equilibrium()

Returns equation collection, to compute equilibrium values. The equations have the post collision symbols as left hand sides and are functions of the conserved quantities

abstract get_collision_rule()

Returns an LbmCollisionRule i.e. an equation collection with a reference to the method. This collision rule defines the collision operator.

Raw Moment-based methods

These methods are represented by instances of lbmpy.methods.momentbased.MomentBasedLbMethod and will derive collision equations either in population space (SRT, TRT, entropic TRT), or in raw moment space (MRT variants).

Creation Functions

The following factory functions create raw moment-based methods using variants of the regular hydrodynamic maxwellian equilibrium.

create_srt(stencil, relaxation_rate, continuous_equilibrium=True, **kwargs)

Creates a single relaxation time (SRT) lattice Boltzmann model also known as BGK model.

Internally calls either create_with_discrete_maxwellian_equilibrium() or create_with_continuous_maxwellian_equilibrium(), depending on the value of continuous_equilibrium.

If not specified otherwise, collision equations will be derived in population space.

Parameters:
  • stencil – instance of lbmpy.stencils.LBStencil

  • relaxation_rate – relaxation rate (inverse of the relaxation time) usually called \(\omega\) in LBM literature

  • continuous_equilibrium – determines if the discrete or continuous maxwellian equilibrium is used to compute the equilibrium moments

Returns:

lbmpy.methods.momentbased.MomentBasedLbMethod instance

create_trt(stencil, relaxation_rate_even_moments, relaxation_rate_odd_moments, continuous_equilibrium=True, **kwargs)

Creates a two relaxation time (TRT) lattice Boltzmann model, where even and odd moments are relaxed differently. In the SRT model the exact wall position of no-slip boundaries depends on the viscosity, the TRT method does not have this problem.

Parameters are similar to lbmpy.methods.create_srt(), but instead of one relaxation rate there are two relaxation rates: one for even moments (determines viscosity) and one for odd moments. If unsure how to choose the odd relaxation rate, use the function lbmpy.methods.create_trt_with_magic_number()

Internally calls either create_with_discrete_maxwellian_equilibrium() or create_with_continuous_maxwellian_equilibrium(), depending on the value of continuous_equilibrium.

If not specified otherwise, collision equations will be derived in population space.

Returns:

lbmpy.methods.momentbased.MomentBasedLbMethod instance

create_trt_with_magic_number(stencil, relaxation_rate, magic_number=3 / 16, **kwargs)

Creates a two relaxation time (TRT) lattice Boltzmann method, where the relaxation time for odd moments is determines from the even moment relaxation time and a “magic number”. For possible parameters see lbmpy.methods.create_trt()

Internally calls either create_with_discrete_maxwellian_equilibrium() or create_with_continuous_maxwellian_equilibrium(), depending on the value of continuous_equilibrium.

If not specified otherwise, collision equations will be derived in population space.

Parameters:
  • stencil – instance of lbmpy.stencils.LBStencil

  • relaxation_rate – relaxation rate (inverse of the relaxation time) usually called \(\omega\) in LBM literature

  • magic_number – magic number which is used to calculate the relxation rate for the odd moments.

Returns:

lbmpy.methods.momentbased.MomentBasedLbMethod instance

create_mrt_orthogonal(stencil, relaxation_rates, continuous_equilibrium=True, weighted=None, nested_moments=None, conserved_moments=True, **kwargs)

Returns an orthogonal multi-relaxation time model for the stencils D2Q9, D3Q15, D3Q19 and D3Q27. These MRT methods are just one specific version - there are many MRT methods possible for all these stencils which differ by the linear combination of moments and the grouping into equal relaxation times. To create a generic MRT method use create_with_discrete_maxwellian_equilibrium

Internally calls either create_with_discrete_maxwellian_equilibrium() or create_with_continuous_maxwellian_equilibrium(), depending on the value of continuous_equilibrium.

If not specified otherwise, collision equations will be derived in raw moment space.

Parameters:
  • stencil – instance of lbmpy.stencils.LBStencil

  • relaxation_rates – relaxation rates for the moments

  • continuous_equilibrium – determines if the discrete or continuous maxwellian equilibrium is used to compute the equilibrium moments

  • weighted – whether to use weighted or unweighted orthogonality

  • nested_moments – a list of lists of modes, grouped by common relaxation times. If this argument is not provided, Gram-Schmidt orthogonalization of the default modes is performed. The default modes equal the raw moments except for the separation of the shear and bulk viscosity.

  • conserved_moments – If lower order moments are conserved or not.

create_trt_kbc(dim, shear_relaxation_rate, higher_order_relaxation_rate, method_name='KBC-N4', continuous_equilibrium=True, **kwargs)

Creates a method with two relaxation rates, one for lower order moments which determines the viscosity and one for higher order moments. In entropic models this second relaxation rate is chosen subject to an entropy condition. Which moments are relaxed by which rate is determined by the method_name

Internally calls either create_with_discrete_maxwellian_equilibrium() or create_with_continuous_maxwellian_equilibrium(), depending on the value of continuous_equilibrium.

If not specified otherwise, collision equations will be derived in population space.

Parameters:
  • dim – 2 or 3, leads to stencil D2Q9 or D3Q27

  • shear_relaxation_rate – relaxation rate that determines viscosity

  • higher_order_relaxation_rate – relaxation rate for higher order moments

  • method_name – string ‘KBC-Nx’ where x can be an number from 1 to 4, for details see “Karlin 2015: Entropic multi relaxation lattice Boltzmann models for turbulent flows”

  • continuous_equilibrium – determines if the discrete or continuous maxwellian equilibrium is used to compute the equilibrium moments.

Class

class MomentBasedLbMethod(stencil, equilibrium, relaxation_dict, conserved_quantity_computation=None, force_model=None, zero_centered=False, fraction_field=None, moment_transform_class=<class 'lbmpy.moment_transforms.rawmomenttransforms.PdfsToMomentsByChimeraTransform'>)

Moment based LBM is a class to represent the single (SRT), two (TRT) and multi relaxation time (MRT) methods. These methods work by transforming the pdfs into moment space using a linear transformation. In the moment-space each component (moment) is relaxed to an equilibrium moment by a certain relaxation rate. These equilibrium moments can e.g. be determined by taking the equilibrium moments of the continuous Maxwellian.

Parameters:
  • stencil – see lbmpy.stencils.LBStencil

  • equilibrium – Instance of lbmpy.equilibrium.AbstractEquilibrium, defining the equilibrium distribution used by this method.

  • relaxation_dict – a dictionary mapping moments in either tuple or polynomial formulation to their relaxation rate.

  • conserved_quantity_computation – instance of lbmpy.methods.AbstractConservedQuantityComputation. This determines how conserved quantities are computed, and defines the symbols used in the equilibrium moments like e.g. density and velocity.

  • force_model – Instance of lbmpy.forcemodels.AbstractForceModel, or None if no forcing terms are required

  • zero_centered – Determines the PDF storage format, regular or centered around the equilibrium’s background distribution.

  • moment_transform_class – transformation class to transform PDFs to the moment space (subclass of lbmpy.moment_transforms.AbstractRawMomentTransform), or None if equations are to be derived in population space.

property force_model

Force model employed by this method.

property relaxation_info_dict

Dictionary mapping this method’s moments to their relaxation rates and equilibrium values. Beware: Changes to this dictionary are not reflected in the method. For changing relaxation rates, use relaxation_rate_dict instead.

property conserved_quantity_computation

Returns an instance of class lbmpy.methods.AbstractConservedQuantityComputation

property equilibrium_distribution

Returns this method’s equilibrium distribution (see lbmpy.equilibrium.AbstractEquilibrium

property moment_transform_class

The transform class (subclass of lbmpy.moment_transforms.AbstractRawMomentTransform defining the transformation of populations to moment space.

property moment_space_collision

Returns whether collision is derived in terms of moments or in terms of populations only.

property moments

Moments relaxed by this method.

property moment_equilibrium_values

Equilibrium values of this method’s moments.

property relaxation_rates

Relaxation rates for this method’s moments.

property relaxation_rate_dict

Dictionary mapping moments to relaxation rates. Changes are reflected by the method.

property zeroth_order_equilibrium_moment_symbol

Returns a symbol referring to the zeroth-order moment of this method’s equilibrium distribution, which is the area under it’s curve (see lbmpy.equilibrium.AbstractEquilibrium.zeroth_order_moment_symbol).

property first_order_equilibrium_moment_symbols

Returns a vector of symbols referring to the first-order moment of this method’s equilibrium distribution, which is its mean value. (see lbmpy.equilibrium.AbstractEquilibrium.first_order_moment_symbols).

property weights

Returns a sequence of weights, one for each lattice direction

override_weights(weights)

Manually set this method’s lattice weights.

get_equilibrium(conserved_quantity_equations=None, include_force_terms=False, pre_simplification=False, subexpressions=False, keep_cqc_subexpressions=True)

Returns equation collection, to compute equilibrium values. The equations have the post collision symbols as left-hand sides and are functions of the conserved quantities

Parameters:
  • conserved_quantity_equations (Optional[AssignmentCollection]) – equations to compute conserved quantities.

  • include_force_terms (bool) – if set to True the equilibrium is shifted by forcing terms coming from the force model of the method.

  • pre_simplification (bool) – with or without pre-simplifications for the calculation of the collision

  • subexpressions (bool) – if set to false all subexpressions of the equilibrium assignments are plugged into the main assignments

  • keep_cqc_subexpressions (bool) – if equilibrium is returned without subexpressions keep_cqc_subexpressions determines if also subexpressions to calculate conserved quantities should be plugged into the main assignments

Return type:

LbmCollisionRule

get_equilibrium_terms()

Returns this method’s equilibrium populations as a vector.

get_collision_rule(conserved_quantity_equations=None, pre_simplification=True)

Returns an LbmCollisionRule i.e. an equation collection with a reference to the method. This collision rule defines the collision operator.

Return type:

LbmCollisionRule

set_zeroth_moment_relaxation_rate(relaxation_rate)

Alters the relaxation rate of the zeroth-order moment.

set_first_moment_relaxation_rate(relaxation_rate)

Alters the relaxation rates of the first-order moments.

set_conserved_moments_relaxation_rate(relaxation_rate)

Alters the relaxation rates of the zeroth- and first-order moments.

set_force_model(force_model)

Updates this method’s force model.

Central Moment-based methods

These methods are represented by instances of lbmpy.methods.momentbased.CentralMomentBasedLbMethod and will derive collision equations in central moment space.

Creation Functions

The following factory functions create central moment-based methods using variants of the regular hydrodynamic maxwellian equilibrium.

create_central_moment(stencil, relaxation_rates, nested_moments=None, continuous_equilibrium=True, conserved_moments=True, fraction_field=None, **kwargs)

Creates moment based LB method where the collision takes place in the central moment space.

Internally calls either create_with_discrete_maxwellian_equilibrium() or create_with_continuous_maxwellian_equilibrium(), depending on the value of continuous_equilibrium.

Parameters:
  • stencil – instance of lbmpy.stencils.LBStencil

  • relaxation_rates – relaxation rates (inverse of the relaxation times) for each moment

  • nested_moments – a list of lists of modes, grouped by common relaxation times.

  • continuous_equilibrium – determines if the discrete or continuous maxwellian equilibrium is used to compute the equilibrium moments.

  • conserved_moments – If lower order moments are conserved or not.

  • fraction_field – fraction field for the PSM method

Returns:

lbmpy.methods.momentbased.CentralMomentBasedLbMethod instance

Class

class CentralMomentBasedLbMethod(stencil, equilibrium, relaxation_dict, conserved_quantity_computation=None, force_model=None, zero_centered=False, central_moment_transform_class=<class 'lbmpy.moment_transforms.centralmomenttransforms.BinomialChimeraTransform'>)

Central Moment based LBM is a class to represent the single (SRT), two (TRT) and multi relaxation time (MRT) methods, where the collision is performed in the central moment space. These methods work by transforming the pdfs into moment space using a linear transformation and then shiftig them into the central moment space. In the central moment space each component (moment) is relaxed to an equilibrium moment by a certain relaxation rate. These equilibrium moments can e.g. be determined by taking the equilibrium moments of the continuous Maxwellian.

Parameters:
  • stencil – see lbmpy.stencils.LBStencil

  • equilibrium – Instance of lbmpy.equilibrium.AbstractEquilibrium, defining the equilibrium distribution used by this method.

  • relaxation_dict – a dictionary mapping moments in either tuple or polynomial formulation to their relaxation rate.

  • conserved_quantity_computation – instance of lbmpy.methods.AbstractConservedQuantityComputation. This determines how conserved quantities are computed, and defines the symbols used in the equilibrium moments like e.g. density and velocity.

  • force_model – Instance of lbmpy.forcemodels.AbstractForceModel, or None if no forcing terms are required

  • zero_centered – Determines the PDF storage format, regular or centered around the equilibrium’s background distribution.

  • central_moment_transform_class – transformation class to transform PDFs to central moment space (subclass of lbmpy.moment_transforms.AbstractCentralMomentTransform)

property force_model

Force model employed by this method.

property relaxation_info_dict

Dictionary mapping this method’s moments to their relaxation rates and equilibrium values. Beware: Changes to this dictionary are not reflected in the method. For changing relaxation rates, use relaxation_rate_dict instead.

property conserved_quantity_computation

Returns an instance of class lbmpy.methods.AbstractConservedQuantityComputation

property equilibrium_distribution

Returns this method’s equilibrium distribution (see lbmpy.equilibrium.AbstractEquilibrium

property central_moment_transform_class

The transform class (subclass of lbmpy.moment_transforms.AbstractCentralMomentTransform defining the transformation of populations to central moment space.

property moments

Central moments relaxed by this method.

property moment_equilibrium_values

Equilibrium values of this method’s moments.

property relaxation_rates

Relaxation rates for this method’s moments.

property relaxation_rate_dict

Dictionary mapping moments to relaxation rates. Changes are reflected by the method.

property zeroth_order_equilibrium_moment_symbol

Returns a symbol referring to the zeroth-order moment of this method’s equilibrium distribution, which is the area under it’s curve (see lbmpy.equilibrium.AbstractEquilibrium.zeroth_order_moment_symbol).

property first_order_equilibrium_moment_symbols

Returns a vector of symbols referring to the first-order moment of this method’s equilibrium distribution, which is its mean value. (see lbmpy.equilibrium.AbstractEquilibrium.first_order_moment_symbols).

property weights

Returns a sequence of weights, one for each lattice direction

property relaxation_matrix: MutableDenseMatrix

Returns a qxq diagonal matrix which contains the relaxation rate for each moment on the diagonal

get_equilibrium(conserved_quantity_equations=None, subexpressions=False, pre_simplification=False, keep_cqc_subexpressions=True, include_force_terms=False)

Returns equation collection, to compute equilibrium values. The equations have the post collision symbols as left-hand sides and are functions of the conserved quantities

Parameters:
  • conserved_quantity_equations (Optional[AssignmentCollection]) – equations to compute conserved quantities.

  • subexpressions (bool) – if set to false all subexpressions of the equilibrium assignments are plugged into the main assignments

  • pre_simplification (bool) – with or without pre_simplifications for the calculation of the collision

  • keep_cqc_subexpressions (bool) – if equilibrium is returned without subexpressions keep_cqc_subexpressions determines if also subexpressions to calculate conserved quantities should be plugged into the main assignments

  • include_force_terms (bool) – if set to True the equilibrium is shifted by forcing terms coming from the force model of the method.

Return type:

LbmCollisionRule

get_collision_rule(conserved_quantity_equations=None, pre_simplification=False)

Returns an LbmCollisionRule i.e. an equation collection with a reference to the method. This collision rule defines the collision operator.

Return type:

LbmCollisionRule

Cumulant-based methods

These methods are represented by instances of lbmpy.methods.cumulantbased.CumulantBasedLbMethod and will derive collision equations in cumulant space.

Creation Functions

The following factory functions create cumulant-based methods using the regular continuous hydrodynamic maxwellian equilibrium.

create_cumulant(stencil, relaxation_rates, cumulant_groups, conserved_moments=True, fraction_field=None, **kwargs)

Creates a cumulant-based lattice Boltzmann method.

Parameters:
  • stencil – instance of lbmpy.stencils.LBStencil

  • relaxation_rates – relaxation rates for each cumulant group. If None are provided a list of symbolic relaxation rates is created and used. If only a list with one entry is provided this relaxation rate is used for determine the viscosity of the simulation. All other cumulants are relaxed with one. If a cumulant force model is provided the first order cumulants are relaxed with two to ensure that the force is applied correctly to the momentum groups

  • cumulant_groups – Nested sequence of polynomial expressions defining the cumulants to be relaxed. All cumulants within one group are relaxed with the same relaxation rate.

  • conserved_moments – If lower order moments are conserved or not.

  • kwargs – See create_with_continuous_maxwellian_equilibrium()

Returns:

lbmpy.methods.cumulantbased.CumulantBasedLbMethod instance

create_with_monomial_cumulants(stencil, relaxation_rates, conserved_moments=True, **kwargs)

Creates a cumulant lattice Boltzmann model using the given stencil’s canonical monomial cumulants.

Parameters:
  • stencil – instance of lbmpy.stencils.LBStencil

  • relaxation_rates – relaxation rates for each cumulant group. If None are provided a list of symbolic relaxation rates is created and used. If only a list with one entry is provided this relaxation rate is used for determine the viscosity of the simulation. All other cumulants are relaxed with one. If a cumulant force model is provided the first order cumulants are relaxed with two to ensure that the force is applied correctly to the momentum groups

  • conserved_moments – If lower order moments are conserved or not.

  • kwargs – See create_cumulant()

Returns:

lbmpy.methods.cumulantbased.CumulantBasedLbMethod instance

create_with_default_polynomial_cumulants(stencil, relaxation_rates, conserved_moments=True, **kwargs)

Creates a cumulant lattice Boltzmann model based on the default polynomial set of [GSchonherrPK15].

Args: See create_cumulant().

Returns:

lbmpy.methods.cumulantbased.CumulantBasedLbMethod instance

Class

class CumulantBasedLbMethod(stencil, equilibrium, relaxation_dict, conserved_quantity_computation=None, force_model=None, zero_centered=False, central_moment_transform_class=<class 'lbmpy.moment_transforms.centralmomenttransforms.BinomialChimeraTransform'>, cumulant_transform_class=<class 'lbmpy.moment_transforms.cumulanttransforms.CentralMomentsToCumulantsByGeneratingFunc'>)

This class implements cumulant-based lattice boltzmann methods which relax all the non-conserved quantities as either monomial or polynomial cumulants. It is mostly inspired by the work presented in [GSchonherrPK15].

This method is implemented modularly as the transformation from populations to central moments to cumulants is governed by subclasses of lbmpy.moment_transforms.AbstractMomentTransform which can be specified by constructor argument. This allows the selection of the most efficient transformation for a given setup.

Parameters:
  • stencil – see lbmpy.stencils.LBStencil

  • equilibrium – Instance of lbmpy.equilibrium.AbstractEquilibrium, defining the equilibrium distribution used by this method.

  • relaxation_dict – a dictionary mapping cumulants in either tuple or polynomial formulation to their relaxation rate.

  • conserved_quantity_computation – instance of lbmpy.methods.AbstractConservedQuantityComputation. This determines how conserved quantities are computed, and defines the symbols used in the equilibrium moments like e.g. density and velocity.

  • force_model – Instance of lbmpy.forcemodels.AbstractForceModel, or None if no forcing terms are required

  • zero_centered – Determines the PDF storage format, regular or centered around the equilibrium’s background distribution.

  • central_moment_transform_class – transformation class to transform PDFs to central moment space (subclass of lbmpy.moment_transforms.AbstractCentralMomentTransform)

  • cumulant_transform_class – transform class to get from the central moment space to the cumulant space

property force_model

Force model employed by this method.

property relaxation_info_dict

Dictionary mapping this method’s cumulants to their relaxation rates and equilibrium values. Beware: Changes to this dictionary are not reflected in the method. For changing relaxation rates, use relaxation_rate_dict instead.

property conserved_quantity_computation

Returns an instance of class lbmpy.methods.AbstractConservedQuantityComputation

property equilibrium_distribution

Returns this method’s equilibrium distribution (see lbmpy.equilibrium.AbstractEquilibrium

property central_moment_transform_class

The transform class (subclass of lbmpy.moment_transforms.AbstractCentralMomentTransform defining the transformation of populations to central moment space.

property cumulant_transform_class

The transform class defining the transform from central moment to cumulant space.

property cumulants

Cumulants relaxed by this method.

property cumulant_equilibrium_values

Equilibrium values of this method’s cumulants.

property relaxation_rates

Relaxation rates for this method’s cumulants.

property relaxation_rate_dict

Dictionary mapping cumulants to relaxation rates. Changes are reflected by the method.

property zeroth_order_equilibrium_moment_symbol

Returns a symbol referring to the zeroth-order moment of this method’s equilibrium distribution, which is the area under its curve (see lbmpy.equilibrium.AbstractEquilibrium.zeroth_order_moment_symbol).

property first_order_equilibrium_moment_symbols

Returns a vector of symbols referring to the first-order moment of this method’s equilibrium distribution, which is its mean value. (see lbmpy.equilibrium.AbstractEquilibrium.first_order_moment_symbols).

property weights

Returns a sequence of weights, one for each lattice direction

get_equilibrium(conserved_quantity_equations=None, subexpressions=False, pre_simplification=False, keep_cqc_subexpressions=True, include_force_terms=False)

Returns equation collection, to compute equilibrium values. The equations have the post collision symbols as left-hand sides and are functions of the conserved quantities

Parameters:
  • conserved_quantity_equations (Optional[AssignmentCollection]) – equations to compute conserved quantities.

  • subexpressions (bool) – if set to false all subexpressions of the equilibrium assignments are plugged into the main assignments

  • pre_simplification (bool) – with or without pre_simplifications for the calculation of the collision

  • keep_cqc_subexpressions (bool) – if equilibrium is returned without subexpressions keep_cqc_subexpressions determines if also subexpressions to calculate conserved quantities should be plugged into the main assignments

  • include_force_terms (bool) – if set to True the equilibrium is shifted by forcing terms coming from the force model of the method.

Return type:

LbmCollisionRule

get_collision_rule(conserved_quantity_equations=None, pre_simplification=False)

Returns an LbmCollisionRule i.e. an equation collection with a reference to the method. This collision rule defines the collision operator.

Return type:

AssignmentCollection

Default Moment sets

The following functions provide default sets of polynomial raw moments, central moments and cumulants to be used in MRT-type methods.

cascaded_moment_sets_literature(stencil)

Returns default groups of central moments or cumulants to be relaxed with common relaxation rates as stated in literature. Groups are ordered like this:

  • First group is density

  • Second group are the momentum modes

  • Third group are the shear modes

  • Fourth group is the bulk mode

  • Remaining groups do not govern hydrodynamic properties

Parameters:

stencil – instance of lbmpy.stencils.LBStencil. Can be D2Q9, D3Q7, D3Q15, D3Q19 or D3Q27

mrt_orthogonal_modes_literature(stencil, is_weighted)

Returns a list of lists of modes, grouped by common relaxation times. This is for commonly used MRT models found in literature.

Parameters:
  • stencil – instance of lbmpy.stencils.LBStencil. Can be D2Q9, D3Q15, D3Q19 or D3Q27

  • is_weighted – whether to use weighted or unweighted orthogonality

MRT schemes as described in the following references are used

Low-Level Method Creation Interface

The following classes and factory functions constitute the lower levels of abstraction in method creation. They are called from the higher-level functions described above, or, in special cases, by the user directly.

Custom method variants in population space, raw and central moment space based on the hydrodynamic Maxwellian equilibrium may be created using either lbmpy.methods.creationfunctions.create_with_discrete_maxwellian_equilibrium() or create_with_continuous_maxwellian_equilibrium().

Methods may also be created using custom equilibrium distributions using lbmpy.methods.creationfunctions.create_from_equilibrium().

The desired collision space, and the transform classes defining the manner of transforming populations to that space and back, are defined using lbmpy.enums.CollisionSpace and lbmpy.methods.CollisionSpaceInfo.

Collision Space Info

Low-Level Creation Functions

create_with_discrete_maxwellian_equilibrium(stencil, moment_to_relaxation_rate_dict, compressible=False, zero_centered=False, delta_equilibrium=False, force_model=None, equilibrium_order=2, c_s_sq=1 / 3, **kwargs)

Creates a moment-based LBM by taking a dictionary of moments with corresponding relaxation rates.

These moments are relaxed against the moments of the discrete Maxwellian distribution (see lbmpy.equilibrium.DiscreteHydrodynamicMaxwellian).

Internally, this method calls lbmpy.methods.create_from_equilibrium().

Parameters:
  • stencil – instance of lbmpy.stencils.LBStenil

  • moment_to_relaxation_rate_dict – dict that has as many entries as the stencil. Each moment, which can be represented by an exponent tuple or in polynomial form (see lbmpy.moments), is mapped to a relaxation rate.

  • compressible – If False, the incompressible equilibrium formulation is used to better approximate the incompressible Navier-Stokes equations. Otherwise, the default textbook equilibrium is used.

  • zero_centered – If True, the zero-centered storage format for PDFs is used, storing only their deviation from the background distribution (given by the lattice weights).

  • delta_equilibrium – Takes effect only if zero-centered storage is used. If True, the equilibrium distribution is also given only as its deviation from the background distribution.

  • force_model – instance of lbmpy.forcemodels.AbstractForceModel, or None if no external forces are present.

  • equilibrium_order – approximation order of macroscopic velocity \(\mathbf{u}\) in the equilibrium

  • c_s_sq – Speed of sound squared

  • kwargs – See lbmpy.methods.create_from_equilibrium()

Returns:

Instance of a subclass of lbmpy.methods.AbstractLbMethod.

create_with_continuous_maxwellian_equilibrium(stencil, moment_to_relaxation_rate_dict, compressible=False, zero_centered=False, delta_equilibrium=False, force_model=None, equilibrium_order=2, c_s_sq=1 / 3, **kwargs)

Creates a moment-based LBM by taking a dictionary of moments with corresponding relaxation rates. These moments are relaxed against the moments of the continuous Maxwellian distribution (see lbmpy.equilibrium.ContinuousHydrodynamicMaxwellian).

Internally, this method calls lbmpy.methods.create_from_equilibrium().

Parameters:
  • stencil – instance of lbmpy.stencils.LBStenil

  • moment_to_relaxation_rate_dict – dict that has as many entries as the stencil. Each moment, which can be represented by an exponent tuple or in polynomial form (see lbmpy.moments), is mapped to a relaxation rate.

  • compressible – If False, the incompressible equilibrium formulation is used to better approximate the incompressible Navier-Stokes equations. Otherwise, the default textbook equilibrium is used.

  • zero_centered – If True, the zero-centered storage format for PDFs is used, storing only their deviation from the background distribution (given by the lattice weights).

  • delta_equilibrium – Takes effect only if zero-centered storage is used. If True, the equilibrium distribution is also given only as its deviation from the background distribution.

  • force_model – Instance of lbmpy.forcemodels.AbstractForceModel, or None if no external forces are present.

  • equilibrium_order – approximation order of macroscopic velocity \(\mathbf{u}\) in the equilibrium

  • c_s_sq – Speed of sound squared

  • kwargs – See lbmpy.methods.create_from_equilibrium()

Returns:

Instance of a subclass of lbmpy.methods.AbstractLbMethod.

create_from_equilibrium(stencil, equilibrium, conserved_quantity_computation, moment_to_relaxation_rate_dict, collision_space_info=CollisionSpaceInfo(collision_space=<CollisionSpace.POPULATIONS: 1>, raw_moment_transform_class=None, central_moment_transform_class=None, cumulant_transform_class=None), zero_centered=False, force_model=None, fraction_field=None)

Creates a lattice Boltzmann method in either population, moment, or central moment space, from a given discrete velocity set and equilibrium distribution.

This function takes a stencil, an equilibrium distribution, an appropriate conserved quantity computation instance, a dictionary mapping moment polynomials to their relaxation rates, and a collision space info determining the desired collision space. It returns a method instance relaxing the given moments against their equilibrium values computed from the given distribution, in the given collision space.

Parameters:
  • stencil – Instance of lbmpy.stencils.LBStencil

  • equilibrium – Instance of a subclass of lbmpy.equilibrium.AbstractEquilibrium.

  • conserved_quantity_computation – Instance of a subclass of lbmpy.methods.AbstractConservedQuantityComputation, which must provide equations for the conserved quantities used in the equilibrium.

  • moment_to_relaxation_rate_dict – Dictionary mapping moment polynomials to relaxation rates.

  • collision_space_info – Instance of CollisionSpaceInfo, defining the method’s desired collision space and the manner of transformation to this space. Cumulant-based methods are not supported yet.

  • zero_centered – Whether or not the zero-centered storage format should be used. If True, the given equilibrium must either be a deviation-only formulation, or must provide a background distribution for PDFs to be centered around.

  • force_model – Instance of lbmpy.forcemodels.AbstractForceModel, or None if no external forces are present.

Conserved Quantity Computation

The classes of the conserved quantity computation (CQC) submodule define an LB Method’s conserved quantities and the equations to compute them. For hydrodynamic methods, lbmpy.methods.DensityVelocityComputation is the typical choice. For custom methods, however, a custom CQC class might have to be created.

class AbstractConservedQuantityComputation

This class defines how conserved quantities are computed as functions of the pdfs. Conserved quantities are used for output and as input to the equilibrium in the collision step

Depending on the method they might also be computed slightly different, e.g. due to a force model.

An additional method describes how to get the conserved quantities for the equilibrium for initialization. In most cases the inputs can be used directly, but for some methods they have to be altered slightly. For example in zero centered hydrodynamic schemes with force model, the density has to be decreased by one, and the given velocity has to be shifted dependent on the force.

../_images/moment_shift.svg
abstract property conserved_quantities

Dict, mapping names (symbol) to dimensionality (int) For example: {‘density’ : 1, ‘velocity’ : 3} The naming strings can be used in output_equations_from_pdfs() and equilibrium_input_equations_from_init_values()

defined_symbols(order='all')

Returns a dict, mapping names of conserved quantities to their symbols

abstract property default_values

Returns a dict of symbol to default value, where “default” means that the equilibrium simplifies to the weights if these values are inserted. Hydrodynamic example: rho=1, u_i = 0

abstract equilibrium_input_equations_from_pdfs(pdfs, **kwargs)

Returns an equation collection that defines all necessary quantities to compute the equilibrium as functions of the pdfs. For hydrodynamic LBM schemes this is usually the density and velocity.

Parameters:

pdfs – values or symbols for the pdf values

abstract output_equations_from_pdfs(pdfs, output_quantity_names_to_symbols, **kwargs)

Returns an equation collection that defines conserved quantities for output. These conserved quantities might be slightly different that the ones used as input for the equilibrium e.g. due to a force model.

Parameters:
  • pdfs – values for the pdf entries

  • output_quantity_names_to_symbols – dict mapping of conserved quantity names (See conserved_quantities()) to symbols or field accesses where they should be written to

abstract equilibrium_input_equations_from_init_values(**kwargs)

Returns an equation collection that defines all necessary quantities to compute the equilibrium as function of given conserved quantities. Parameters can be names that are given by symbol names of conserved_quantities(). For all parameters not specified each implementation should use sensible defaults. For example hydrodynamic schemes use density=1 and velocity=0.

class DensityVelocityComputation(stencil, compressible, zero_centered, force_model=None, background_density=1, density_symbol=rho, density_deviation_symbol=delta_rho, velocity_symbols=(u_0, u_1, u_2), c_s_sq=c_s**2, second_order_moment_symbols=(p_0, p_1, p_2, p_3, p_4, p_5, p_6, p_7, p_8), zeroth_order_moment_symbol=None, first_order_moment_symbols=None)

This class emits equations for in- and output of the conserved quantities of regular hydrodynamic lattice Boltzmann methods, which are density and velocity. The following symbolic quantities are manages by this class:

  • Density \(\rho\), background density \(\rho_0\) (typically set to \(1\)) and the density deviation \(\delta \rho\). They are connected through \(\rho = \rho_0 + \delta\rho\).

  • Velocity \(\mathbf{u} = (u_0, u_1, u_2)\).

Furthermore, this class provides output functionality for the second-order moment tensor \(p\).

Parameters:
  • stencil – see lbmpy.stencils.LBStencil

  • compressibleTrue indicates the usage of a compressible equilibrium (see lbmpy.equilibrium.ContinuousHydrodynamicMaxwellian), and sets the reference density to \(\rho\). False indicates an incompressible equilibrium, using \(\rho\) as reference density.

  • zero_centered – Indicates whether or not PDFs are stored in regular or zero-centered format.

  • force_model – Employed force model. See lbmpy.forcemodels.

property conserved_quantities

Dict, mapping names (symbol) to dimensionality (int) For example: {‘density’ : 1, ‘velocity’ : 3} The naming strings can be used in output_equations_from_pdfs() and equilibrium_input_equations_from_init_values()

property compressible

Indicates whether a compressible or incompressible equilibrium is used.

defined_symbols(order='all')

Returns a dict, mapping names of conserved quantities to their symbols

property zero_centered_pdfs

Whether regular or zero-centered storage is employed.

property zeroth_order_moment_symbol

Symbol corresponding to the zeroth-order moment of the stored PDFs, i.e. density_deviation_symbol if zero-centered, else density_symbol.

property first_order_moment_symbols

Symbol corresponding to the first-order moment vector of the stored PDFs, divided by the reference density.

property density_symbol

Symbol for the density.

property density_deviation_symbol

Symbol for the density deviation.

property background_density

Symbol or value of the background density.

property velocity_symbols

Symbols for the velocity.

property default_values

Returns a dict of symbol to default value, where “default” means that the equilibrium simplifies to the weights if these values are inserted. Hydrodynamic example: rho=1, u_i = 0

equilibrium_input_equations_from_pdfs(pdfs, force_substitution=True)

Returns an equation collection that defines all necessary quantities to compute the equilibrium as functions of the pdfs. For hydrodynamic LBM schemes this is usually the density and velocity.

Parameters:

pdfs – values or symbols for the pdf values

equilibrium_input_equations_from_init_values(density=1, velocity=(0, 0, 0), force_substitution=True)

Returns an equation collection that defines all necessary quantities to compute the equilibrium as function of given conserved quantities. Parameters can be names that are given by symbol names of conserved_quantities(). For all parameters not specified each implementation should use sensible defaults. For example hydrodynamic schemes use density=1 and velocity=0.

output_equations_from_pdfs(pdfs, output_quantity_names_to_symbols, force_substitution=True)

Returns an equation collection that defines conserved quantities for output. These conserved quantities might be slightly different that the ones used as input for the equilibrium e.g. due to a force model.

Parameters:
  • pdfs – values for the pdf entries

  • output_quantity_names_to_symbols – dict mapping of conserved quantity names (See conserved_quantities()) to symbols or field accesses where they should be written to