Boundary Conditions

class LbBoundary(name=None, calculate_force_on_boundary=False)

Base class that all boundaries should derive from.

Parameters:

name – optional name of the boundary.

property additional_data

Return a list of (name, type) tuples for additional data items required in this boundary These data items can either be initialized in separate kernel see additional_data_kernel_init or by Python callbacks - see additional_data_callback

property additional_data_init_callback

Return a callback function called with a boundary data setter object and returning a dict of data-name to data for each element that should be initialized

get_additional_code_nodes(lb_method)

Return a list of code nodes that will be added in the generated code before the index field loop.

Parameters:

lb_method – lattice Boltzmann method. See lbmpy.creationfunctions.create_lb_method()

class NoSlip(name=None, calculate_force_on_boundary=False)

No-Slip, (half-way) simple bounce back boundary condition, enforcing zero velocity at obstacle. Populations leaving the boundary node \(\mathbf{x}_b\) at time \(t\) are reflected back with \(\mathbf{c}_{\overline{i}} = -\mathbf{c}_{i}\)

\[f_{\overline{i}}(\mathbf{x}_b, t + \Delta t) = f^{\star}_{i}(\mathbf{x}_b, t)\]
Parameters:
  • name – optional name of the boundary.

  • calculate_force_on_boundary – stores the force for each PDF at the boundary in a force vector

get_additional_code_nodes(lb_method)

Return a list of code nodes that will be added in the generated code before the index field loop.

Parameters:

lb_method – lattice Boltzmann method. See lbmpy.creationfunctions.create_lb_method()

class NoSlipLinearBouzidi(name=None, init_wall_distance=None, data_type='double', calculate_force_on_boundary=False)

No-Slip, (half-way) simple bounce back boundary condition with interpolation to increase accuracy: [BFL01]. In order to make the boundary condition work properly a Python callback function needs to be provided to calculate the distance from the wall for each cell near to the boundary. If this is not done the boundary condition will fall back to a simple NoSlip boundary. Furthermore, for this boundary condition a second fluid cell away from the wall is needed. If the second fluid cell is not available (e.g. because it is marked as boundary as well), the boundary condition should fall back to a NoSlip boundary as well.

Parameters:
  • name – optional name of the boundary.

  • init_wall_distance – Python callback function to calculate the wall distance for each cell near to the boundary

  • data_type – data type of the wall distance q

property additional_data

Used internally only. For the NoSlipLinearBouzidi boundary the distance to the obstacle of every direction is needed. This information is stored in the index vector.

get_additional_code_nodes(lb_method)

Return a list of code nodes that will be added in the generated code before the index field loop.

Parameters:

lb_method – lattice Boltzmann method. See lbmpy.creationfunctions.create_lb_method()

property additional_data_init_callback

Return a callback function called with a boundary data setter object and returning a dict of data-name to data for each element that should be initialized

class QuadraticBounceBack(relaxation_rate, name=None, init_wall_distance=None, data_type='double', calculate_force_on_boundary=False)

Second order accurate bounce back boundary condition. Implementation details are provided in a demo notebook here: https://pycodegen.pages.i10git.cs.fau.de/lbmpy/notebooks/demo_interpolation_boundary_conditions.html

Parameters:
  • relaxation_rate – relaxation rate to realise a BGK scheme for recovering the pre collision PDF value.

  • name – optional name of the boundary.

  • init_wall_distance – Python callback function to calculate the wall distance for each cell near to the boundary

  • data_type – data type of the wall distance q

property additional_data

Used internally only. For the NoSlipLinearBouzidi boundary the distance to the obstacle of every direction is needed. This information is stored in the index vector.

property additional_data_init_callback

Return a callback function called with a boundary data setter object and returning a dict of data-name to data for each element that should be initialized

get_additional_code_nodes(lb_method)

Return a list of code nodes that will be added in the generated code before the index field loop.

Parameters:

lb_method – Lattice Boltzmann method. See lbmpy.creationfunctions.create_lb_method()

Returns:

list containing LbmWeightInfo

class FreeSlip(stencil, normal_direction=None, name=None)

Free-Slip boundary condition, which enforces a zero normal fluid velocity \(u_n = 0\) but places no restrictions on the tangential fluid velocity \(u_t\).

Parameters:
  • stencil – LBM stencil which is used for the simulation

  • normal_direction – optional normal direction pointing from wall to fluid. If the Free slip boundary is applied to a certain side in the domain it is not necessary to calculate the normal direction since it can be stated for all boundary cells. This reduces the memory space for the index array significantly.

  • name – optional name of the boundary.

property additional_data

Used internally only. For the FreeSlip boundary the information of the normal direction for each pdf direction is needed. This information is stored in the index vector.

property additional_data_init_callback

Return a callback function called with a boundary data setter object and returning a dict of data-name to data for each element that should be initialized

get_additional_code_nodes(lb_method)

Return a list of code nodes that will be added in the generated code before the index field loop.

Parameters:

lb_method – lattice Boltzmann method. See lbmpy.creationfunctions.create_lb_method()

class WallFunctionBounce(lb_method, pdfs, normal_direction, wall_function_model, mean_velocity=None, sampling_shift=1, maronga_sampling_shift=None, dt=1, dy=1, y=0.5, target_friction_velocity=None, weight_method=WeightMethod.LATTICE_WEIGHT, name=None, data_type='double')

Wall function based on the bounce back idea, cf. [HOK21]. Its implementation is extended to the D3Q27 stencil, whereas different weights of the drag distribution are proposed.

Parameters:
  • lb_method – LB method which is used for the simulation

  • pdfs – Symbolic representation of the particle distribution functions.

  • normal_direction – Normal direction of the wall. Currently, only straight and axis-aligned walls are supported.

  • wall_function_model – Wall function that is used to retrieve the wall stress \(tau_w\) during the simulation. See lbmpy.boundaries.wall_treatment.WallFunctionModel for more details

  • mean_velocity – Optional field or field access for the mean velocity. As wall functions are typically defined in terms of the mean velocity, it is recommended to provide this variable. Per default, the instantaneous velocity obtained from pdfs is used for the wall function.

  • sampling_shift – Optional sampling shift for the velocity sampling. Can be provided as symbolic variable or integer. In both cases, the user must assure that the sampling shift is at least 1, as sampling in boundary cells is not physical. Per default, a sampling shift of 1 is employed which corresponds to a sampling in the first fluid cell normal to the wall. For lower friction Reynolds numbers, choosing a sampling shift >1 has shown to improve the results for higher resolutions. Mutually exclusive with the Maronga sampling shift.

  • maronga_sampling_shift – Optionally, apply a correction factor to the wall shear stress proposed by Maronga et al. [MKR20]. Has only been tested and validated for the MOST wall function. No guarantee is given that it also works with other wall functions. Mutually exclusive with the standard sampling shift.

  • dt – time discretisation. Usually one in LB units

  • dy – space discretisation. Usually one in LB units

  • y – distance from the wall

  • target_friction_velocity – A target friction velocity can be given if an estimate is known a priori. This target friction velocity will be used as initial guess for implicit wall functions to ensure convergence of the Newton algorithm.

  • weight_method – The extension of the WFB to a D3Q27 stencil is non-unique. Different weights can be chosen to define the drag distribution onto the pdfs. Per default, weights corresponding to the weights in the D3Q27 stencil are chosen.

  • name – Optional name of the boundary.

  • data_type – Floating-point precision. Per default, double.

class WeightMethod(value)

An enumeration.

get_additional_code_nodes(lb_method)

Return a list of code nodes that will be added in the generated code before the index field loop.

Parameters:

lb_method – lattice Boltzmann method. See lbmpy.creationfunctions.create_lb_method()

class UBB(velocity, density=None, adapt_velocity_to_force=False, dim=None, name=None, data_type='double')
Velocity bounce back boundary condition, enforcing specified velocity at obstacle. Furthermore, a density

at the wall can be implied. The boundary condition is implemented with the following formula:

\[f_{\overline{i}}(\mathbf{x}_b, t + \Delta t) = f^{\star}_{i}(\mathbf{x}_b, t) - 2 w_{i} \rho_{w} \frac{\mathbf{c}_i \cdot \mathbf{u}_w}{c_s^2}\]
Parameters:
  • velocity – Prescribe the fluid velocity \(\mathbf{u}_w\) at the wall. Can either be a constant, an access into a field, or a callback function. The callback functions gets a numpy record array with members, x, y, z, dir (direction) and velocity which has to be set to the desired velocity of the corresponding link

  • density – Prescribe the fluid density \(\rho_{w}\) at the wall. If not prescribed the density is calculated from the PDFs at the wall. The density can only be set constant.

  • adapt_velocity_to_force – adapts the velocity to the correct equilibrium when the lattice Boltzmann method holds a forcing term. If no forcing term is set and adapt_velocity_to_force is set to True it has no effect.

  • dim – number of spatial dimensions

  • name – optional name of the boundary.

property additional_data

In case of the UBB boundary additional data is a velocity vector. This vector is added to each cell to realize velocity profiles for the inlet.

property additional_data_init_callback

Initialise additional data of the boundary. For an example see tutorial 02 or lbmpy.geometry.add_pipe_inflow_boundary

get_additional_code_nodes(lb_method)

Return a list of code nodes that will be added in the generated code before the index field loop.

Parameters:

lb_method – Lattice Boltzmann method. See lbmpy.creationfunctions.create_lb_method()

Returns:

list containing LbmWeightInfo and NeighbourOffsetArrays

property velocity_is_callable

Returns True is velocity is callable. This means the velocity should be initialised via a callback function. This is useful if the inflow velocity should have a certain profile for instance

class SimpleExtrapolationOutflow(normal_direction, stencil, name=None)

Simple Outflow boundary condition [GSchonherrPK15], equation F.1 (listed below). This boundary condition extrapolates missing populations from the last layer of fluid cells onto the boundary by copying them in the normal direction.

\[f_{\overline{1}jkxyzt} = f_{\overline{1}jk(x - \Delta x)yzt}\]
Parameters:
  • normal_direction – direction vector normal to the outflow

  • stencil – stencil used by the lattice boltzmann method

  • name – optional name of the boundary.

get_additional_code_nodes(lb_method)

Return a list of code nodes that will be added in the generated code before the index field loop.

Parameters:

lb_method – Lattice Boltzmann method. See lbmpy.creationfunctions.create_lb_method()

Returns:

list containing NeighbourOffsetArrays

class ExtrapolationOutflow(normal_direction, lb_method, dt=1, dx=1, name=None, streaming_pattern='pull', zeroth_timestep=Timestep.BOTH, initial_density=None, initial_velocity=None, data_type='double')

Outflow boundary condition [GSchonherrPK15], equation F.2, with u neglected (listed below). This boundary condition interpolates populations missing on the boundary in normal direction. For this interpolation, the PDF values of the last time step are used. They are interpolated between fluid cell and boundary cell. To get the PDF values from the last time step an index array is used which stores them.

\[f_{\overline{1}jkxyzt} = f_{\overline{1}jk(x - \Delta x)yz(t - \Delta t)} c \theta^{\frac{1}{2}} \frac{\Delta t}{\Delta x} + \left(1 - c \theta^{\frac{1}{2}} \frac{\Delta t}{\Delta x} \right) f_{\overline{1}jk(x - \Delta x)yzt}\]
Parameters:
  • normal_direction – direction vector normal to the outflow

  • lb_method – the lattice Boltzmann method to be used in the simulation

  • dt – lattice time step size

  • dx – lattice spacing distance

  • name – optional name of the boundary.

  • streaming_pattern – Streaming pattern to be used in the simulation

  • zeroth_timestep – for in-place patterns, whether the initial setup corresponds to an even or odd time step

  • initial_density – floating point constant or callback taking spatial coordinates (x, y [,z]) as positional arguments, specifying the initial density on boundary nodes

  • initial_velocity – tuple of floating point constants or callback taking spatial coordinates (x, y [,z]) as positional arguments, specifying the initial velocity on boundary nodes

property additional_data

Used internally only. For the ExtrapolationOutflow information of the previous PDF values is needed. This information is stored in the index vector.

property additional_data_init_callback

Return a callback function called with a boundary data setter object and returning a dict of data-name to data for each element that should be initialized

get_additional_code_nodes(lb_method)

Return a list of code nodes that will be added in the generated code before the index field loop.

Parameters:

lb_method – Lattice Boltzmann method. See lbmpy.creationfunctions.create_lb_method()

Returns:

list containing NeighbourOffsetArrays

class FixedDensity(density, name=None)
Boundary condition for prescribing a density at the wall. Through \(p = c_s^2 \rho\) this boundary condition

can also function as a pressure boundary condition.

\[f_{\overline{i}}(\mathbf{x}_b, t + \Delta t) = - f^{\star}_{i}(\mathbf{x}_b, t) + 2 w_{i} \rho_{w} (1 + \frac{(\mathbf{c}_i \cdot \mathbf{u}_w)^2}{2c_s^4} + \frac{\mathbf{u}_w^2}{2c_s^2})\]
Parameters:
  • density – value of the density which should be set.

  • name – optional name of the boundary.

class DiffusionDirichlet(concentration, velocity_field=None, name=None, data_type='double')

Concentration boundary which is used for concentration or thermal boundary conditions of convection-diffusion equation Base on https://doi.org/10.1103/PhysRevE.85.016701.

Parameters:
  • concentration – can either be a constant, an access into a field, or a callback function. The callback functions gets a numpy record array with members, x, y, z, dir (direction) and concentration which has to be set to the desired velocity of the corresponding link

  • velocity_field – if velocity field is given the boundary value is approximated by using the discrete equilibrium.

  • name – optional name of the boundary.

  • data_type – data type of the concentration value. default is double

property additional_data

In case of the UBB boundary additional data is a velocity vector. This vector is added to each cell to realize velocity profiles for the inlet.

property additional_data_init_callback

Initialise additional data of the boundary. For an example see tutorial 02 or lbmpy.geometry.add_pipe_inflow_boundary

get_additional_code_nodes(lb_method)

Return a list of code nodes that will be added in the generated code before the index field loop.

Parameters:

lb_method – Lattice Boltzmann method. See lbmpy.creationfunctions.create_lb_method()

Returns:

list containing LbmWeightInfo

class NeumannByCopy(name=None, calculate_force_on_boundary=False)

Neumann boundary condition which is implemented by coping the PDF values to achieve similar values at the fluid and the boundary node

get_additional_code_nodes(lb_method)

Return a list of code nodes that will be added in the generated code before the index field loop.

Parameters:

lb_method – Lattice Boltzmann method. See lbmpy.creationfunctions.create_lb_method()

Returns:

list containing NeighbourOffsetArrays

class StreamInConstant(constant, name=None)

Boundary condition that takes a constant and overrides the boundary PDFs with this value. This is used for debugging mainly.

Parameters:
  • constant – value which should be set for the PDFs at the boundary cell.

  • name – optional name of the boundary.

get_additional_code_nodes(lb_method)

Return a list of code nodes that will be added in the generated code before the index field loop.

Parameters:

lb_method – Lattice Boltzmann method. See lbmpy.creationfunctions.create_lb_method()

Returns:

list containing NeighbourOffsetArrays