API¶
Simulation¶
Lattice Boltzmann Solver
-
class
lettuce.simulation.
Simulation
(flow, lattice, collision, streaming)[source] Bases:
object
High-level API for simulations.
-
reporters
¶ A list of reporters. Their call functions are invoked after every simulation step (and before the first one).
Type: list
-
initialize
(max_num_steps=500, tol_pressure=0.001)[source] Iterative initialization to get moments consistent with the initial velocity.
Using the initialization does not better TGV convergence. Maybe use a better scheme?
-
initialize_f_neq
()[source] Initialize the distribution function values. The f^(1) contributions are approximated by finite differences. See Krüger et al. (2017).
-
initialize_pressure
(max_num_steps=100000, tol_pressure=1e-06)[source] Reinitialize equilibrium distributions with pressure obtained by a Jacobi solver. Note that this method has to be called before initialize_f_neq.
-
load_checkpoint
(filename)[source] Load f as np.array using pickle module.
-
save_checkpoint
(filename)[source] Write f as np.array using pickle module.
-
step
(num_steps)[source] Take num_steps stream-and-collision steps and return performance in MLUPS.
-
Lattices¶
Stencils and Lattices.
A Stencil, like the D1Q3 class, provides particle velocities (e), weights (w), and speeds of sound. Velocities and weights are stored as numpy arrays.
In contrast, the Lattice lives on the Device (usually a GPU) and its vectors are stored as torch tensors. Its stencil is still accessible trough Lattice.stencil.
-
class
lettuce.lattices.
Lattice
(stencil, device, dtype=torch.float32)[source] Bases:
object
-
D
-
Q
-
classmethod
convert_to_numpy
(tensor)[source]
-
convert_to_tensor
(array)[source]
-
einsum
(equation, fields, **kwargs)[source] Einstein summation on local fields.
-
entropy
(f)[source] entropy according to the H-theorem
-
incompressible_energy
(f)[source] incompressible kinetic energy
-
j
(f)[source] momentum
-
mv
(m, v)[source] matrix-vector multiplication
-
pseudo_entropy_global
(f)[source] pseudo_entropy derived by a Taylor expansion around the weights
-
pseudo_entropy_local
(f)[source] pseudo_entropy derived by a Taylor expansion around the local equilibrium
-
rho
(f)[source] density
-
shear_tensor
(f)[source] computes the shear tensor of a given f in the sense Pi_{lpha eta} = f_i * e_{i lpha} * e_{i eta}
-
u
(f, rho=None, acceleration=None)[source] velocity; the acceleration is used to compute the correct velocity in the presence of a forcing scheme.
-
Stencils¶
-
class
lettuce.stencils.
Stencil
[source] Bases:
object
-
classmethod
D
()[source]
-
classmethod
Q
()[source]
-
classmethod
-
class
lettuce.stencils.
D1Q3
[source] Bases:
lettuce.stencils.Stencil
-
cs
= 0.5773502691896258
-
e
= array([[ 0], [ 1], [-1]])
-
opposite
= [0, 2, 1]
-
w
= array([0.66666667, 0.16666667, 0.16666667])
-
-
class
lettuce.stencils.
D2Q9
[source] Bases:
lettuce.stencils.Stencil
-
cs
= 0.5773502691896258
-
e
= array([[ 0, 0], [ 1, 0], [ 0, 1], [-1, 0], [ 0, -1], [ 1, 1], [-1, 1], [-1, -1], [ 1, -1]])
-
opposite
= [0, 3, 4, 1, 2, 7, 8, 5, 6]
-
w
= array([0.44444444, 0.11111111, 0.11111111, 0.11111111, 0.11111111, 0.02777778, 0.02777778, 0.02777778, 0.02777778])
-
-
class
lettuce.stencils.
D3Q15
[source] Bases:
lettuce.stencils.Stencil
-
cs
= 0.5773502691896258
-
e
= array([[ 0, 0, 0], [ 1, 0, 0], [-1, 0, 0], [ 0, 1, 0], [ 0, -1, 0], [ 0, 0, 1], [ 0, 0, -1], [ 1, 1, 1], [-1, -1, -1], [ 1, 1, -1], [-1, -1, 1], [ 1, -1, 1], [-1, 1, -1], [ 1, -1, -1], [-1, 1, 1]])
-
opposite
= [0, 2, 1, 4, 3, 6, 5, 8, 7, 10, 9, 12, 11, 14, 13]
-
w
= array([0.22222222, 0.11111111, 0.11111111, 0.11111111, 0.11111111, 0.11111111, 0.11111111, 0.01388889, 0.01388889, 0.01388889, 0.01388889, 0.01388889, 0.01388889, 0.01388889, 0.01388889])
-
-
class
lettuce.stencils.
D3Q19
[source] Bases:
lettuce.stencils.Stencil
-
cs
= 0.5773502691896258
-
e
= array([[ 0, 0, 0], [ 1, 0, 0], [-1, 0, 0], [ 0, 1, 0], [ 0, -1, 0], [ 0, 0, 1], [ 0, 0, -1], [ 0, 1, 1], [ 0, -1, -1], [ 0, 1, -1], [ 0, -1, 1], [ 1, 0, 1], [-1, 0, -1], [ 1, 0, -1], [-1, 0, 1], [ 1, 1, 0], [-1, -1, 0], [ 1, -1, 0], [-1, 1, 0]])
-
opposite
= [0, 2, 1, 4, 3, 6, 5, 8, 7, 10, 9, 12, 11, 14, 13, 16, 15, 18, 17]
-
w
= array([0.33333333, 0.05555556, 0.05555556, 0.05555556, 0.05555556, 0.05555556, 0.05555556, 0.02777778, 0.02777778, 0.02777778, 0.02777778, 0.02777778, 0.02777778, 0.02777778, 0.02777778, 0.02777778, 0.02777778, 0.02777778, 0.02777778])
-
-
class
lettuce.stencils.
D3Q27
[source] Bases:
lettuce.stencils.Stencil
-
cs
= 0.5773502691896258
-
e
= array([[ 0, 0, 0], [ 1, 0, 0], [-1, 0, 0], [ 0, 1, 0], [ 0, -1, 0], [ 0, 0, 1], [ 0, 0, -1], [ 0, 1, 1], [ 0, -1, -1], [ 0, 1, -1], [ 0, -1, 1], [ 1, 0, 1], [-1, 0, -1], [ 1, 0, -1], [-1, 0, 1], [ 1, 1, 0], [-1, -1, 0], [ 1, -1, 0], [-1, 1, 0], [ 1, 1, 1], [-1, -1, -1], [ 1, 1, -1], [-1, -1, 1], [ 1, -1, 1], [-1, 1, -1], [ 1, -1, -1], [-1, 1, 1]])
-
opposite
= [0, 2, 1, 4, 3, 6, 5, 8, 7, 10, 9, 12, 11, 14, 13, 16, 15, 18, 17, 20, 19, 22, 21, 24, 23, 26, 25]
-
w
= array([0.2962963 , 0.07407407, 0.07407407, 0.07407407, 0.07407407, 0.07407407, 0.07407407, 0.01851852, 0.01851852, 0.01851852, 0.01851852, 0.01851852, 0.01851852, 0.01851852, 0.01851852, 0.01851852, 0.01851852, 0.01851852, 0.01851852, 0.00462963, 0.00462963, 0.00462963, 0.00462963, 0.00462963, 0.00462963, 0.00462963, 0.00462963])
-
Streaming¶
Streaming Step
Collision¶
Collision models
-
class
lettuce.collision.
BGKCollision
(lattice, tau, force=None)[source] Bases:
object
-
class
lettuce.collision.
KBCCollision2D
(lattice, tau)[source] Bases:
object
Entropic multi-relaxation time model according to Karlin et al. in two dimensions
-
compute_s_seq_from_m
(f, m)[source]
-
kbc_moment_transform
(f)[source] Transforms the f into the KBC moment representation
-
-
class
lettuce.collision.
KBCCollision3D
(lattice, tau)[source] Bases:
object
Entropic multi-relaxation time-relaxation time model according to Karlin et al. in three dimensions
-
compute_s_seq_from_m
(f, m)[source]
-
kbc_moment_transform
(f)[source] Transforms the f into the KBC moment representation
-
-
class
lettuce.collision.
MRTCollision
(lattice, transform, relaxation_parameters)[source] Bases:
object
Multiple relaxation time collision operator
This is an MRT operator in the most general sense of the word. The transform does not have to be linear and can, e.g., be any moment or cumulant transform.
-
class
lettuce.collision.
RegularizedCollision
(lattice, tau)[source] Bases:
object
Regularized LBM according to Jonas Latt and Bastien Chopard (2006)
-
class
lettuce.collision.
SmagorinskyCollision
(lattice, tau, smagorinsky_constant=0.17, force=None)[source] Bases:
object
Smagorinsky large eddy simulation (LES) collision model with BGK operator.
-
class
lettuce.collision.
TRTCollision
(lattice, tau, tau_minus=1.0)[source] Bases:
object
Two relaxation time collision model - standard implementation (cf. Krüger 2017)
-
class
lettuce.collision.
BGKInitialization
(lattice, flow, moment_transformation)[source] Bases:
object
Keep velocity constant.
Observables¶
Observables. Each observable is defined as a callable class. The __call__ function takes f as an argument and returns a torch tensor.
-
class
lettuce.observables.
Observable
(lattice, flow)[source] Bases:
object
-
class
lettuce.observables.
MaximumVelocity
(lattice, flow)[source] Bases:
lettuce.observables.Observable
Maximum velocitiy
-
class
lettuce.observables.
IncompressibleKineticEnergy
(lattice, flow)[source] Bases:
lettuce.observables.Observable
Total kinetic energy of an incompressible flow.
-
class
lettuce.observables.
Enstrophy
(lattice, flow)[source] Bases:
lettuce.observables.Observable
The integral of the vorticity
Notes
The function only works for periodic domains
Reporters¶
Input/output routines. TODO: Logging
-
lettuce.reporters.
write_image
(filename, array2d)[source]
-
lettuce.reporters.
write_vtk
(point_dict, id=0, filename_base='./data/output')[source]
-
class
lettuce.reporters.
VTKReporter
(lattice, flow, interval=50, filename_base='./data/output')[source] Bases:
object
General VTK Reporter for velocity and pressure
-
output_mask
(no_collision_mask)[source] Outputs the no_collision_mask of the simulation object as VTK-file with range [0,1] Usage: vtk_reporter.output_mask(simulation.no_collision_mask)
-
-
class
lettuce.reporters.
ObservableReporter
(observable, interval=1, out=<_io.TextIOWrapper name='<stdout>' mode='w' encoding='UTF-8'>)[source] Bases:
object
A reporter that prints an observable every few iterations.
Examples
Create an Enstrophy reporter.
>>> from lettuce import TaylorGreenVortex3D, Enstrophy, D3Q27, Lattice >>> lattice = Lattice(D3Q27, device="cpu") >>> flow = TaylorGreenVortex(50, 300, 0.1, lattice) >>> enstrophy = Enstrophy(lattice, flow) >>> reporter = ObservableReporter(enstrophy, interval=10) >>> # simulation = ... >>> # simulation.reporters.append(reporter)
-
class
lettuce.reporters.
ErrorReporter
(lattice, flow, interval=1, out=<_io.TextIOWrapper name='<stdout>' mode='w' encoding='UTF-8'>)[source] Bases:
object
Reports numerical errors with respect to analytic solution.
Force¶
Equilibrium¶
-
class
lettuce.equilibrium.
Equilibrium
[source] Bases:
object
-
class
lettuce.equilibrium.
QuadraticEquilibrium
(lattice)[source] Bases:
lettuce.equilibrium.Equilibrium
-
class
lettuce.equilibrium.
IncompressibleQuadraticEquilibrium
(lattice, rho0=1.0)[source] Bases:
lettuce.equilibrium.Equilibrium
-
class
lettuce.equilibrium.
QuadraticEquilibrium_LessMemory
(lattice)[source] Bases:
lettuce.equilibrium.QuadraticEquilibrium
does the same as the normal equilibrium, how ever it uses somewhere around 20% less RAM, but runs about 2% slower on GPU and 11% on CPU
Use this by setting lattice.equilibrium = QuadraticEquilibrium_LessMemory(lattice) before starting your simulation
Boundary¶
Boundary Conditions.
The __call__ function of a boundary defines its application to the distribution functions.
Boundary conditions can define a mask (a boolean numpy array) that specifies the grid points on which the boundary condition operates.
Boundary classes can define two functions make_no_stream_mask and make_no_collision_mask that prevent streaming and collisions on the boundary nodes.
The no-stream mask has the same dimensions as the distribution functions (Q, x, y, (z)) . The no-collision mask has the same dimensions as the grid (x, y, (z)).
-
class
lettuce.boundary.
BounceBackBoundary
(mask, lattice)[source] Bases:
object
Fullway Bounce-Back Boundary
-
make_no_collision_mask
(f_shape)[source]
-
-
class
lettuce.boundary.
AntiBounceBackOutlet
(lattice, direction)[source] Bases:
object
Allows distributions to leave domain unobstructed through this boundary. Based on equations from page 195 of “The lattice Boltzmann method” (2016 by Krüger et al.) Give the side of the domain with the boundary as list [x, y, z] with only one entry nonzero [1, 0, 0] for positive x-direction in 3D; [1, 0] for the same in 2D [0, -1, 0] is negative y-direction in 3D; [0, -1] for the same in 2D
-
make_no_stream_mask
(f_shape)[source]
-
-
class
lettuce.boundary.
EquilibriumBoundaryPU
(mask, lattice, units, velocity, pressure=0)[source] Bases:
object
Sets distributions on this boundary to equilibrium with predefined velocity and pressure. Note that this behavior is generally not compatible with the Navier-Stokes equations. This boundary condition should only be used if no better options are available.
Flows¶
Couette¶
Couette Flow
Poiseuille¶
Poiseuille Flow
Taylor-Green¶
Taylor-Green vortex in 2D and 3D.
Decaying-Turbulence¶
DecayingTurbulence vortex in 2D and 3D. Dimension is set by the stencil. Special Inputs & standard value: wavenumber_energy-peak = 20, initial_energy = 0.5
Additional attributes / properties¶
energy_spectrum: returns a pair [spectrum, wavenumbers]
-
class
lettuce.flows.decayingturbulence.
DecayingTurbulence
(resolution, reynolds_number, mach_number, lattice, k0=20, ic_energy=0.5)[source] Bases:
object
-
analytic_solution
(x, t=0)[source]
-
boundaries
-
energy_spectrum
-
grid
-
initial_solution
(x)[source] Return initial solution. Note: this function sets the characteristic velocity in phyiscal units.
-
Obstacle¶
-
class
lettuce.flows.obstacle.
Obstacle
(shape, reynolds_number, mach_number, lattice, domain_length_x, char_length=1, char_velocity=1)[source] Bases:
object
Flow class to simulate the flow around an object (mask). It consists of one inflow (equilibrium boundary) and one outflow (anti-bounce-back-boundary), leading to a flow in positive x direction.
Parameters: - shape (Tuple[int]:) – Grid resolution.
- domain_length_x (float) – Length of the domain in physical units.
-
mask
¶ Boolean mask to define the obstacle. The shape of this object is the shape of the grid. Initially set to zero (no obstacle).
Type: np.array with dtype = np.bool
Examples
Initialization of flow around a cylinder:
>>> from lettuce import Lattice, D2Q9 >>> flow = Obstacle( >>> shape=(101, 51), >>> reynolds_number=100, >>> mach_number=0.1, >>> lattice=lattice, >>> domain_length_x=10.1 >>> ) >>> x, y = flow.grid >>> condition = np.sqrt((x-2.5)**2+(y-2.5)**2) < 1. >>> flow.mask[np.where(condition)] = 1
-
boundaries
-
grid
-
initial_solution
(x)[source]
-
mask
-
lettuce.flows.obstacle.
Obstacle2D
(resolution_x, resolution_y, reynolds_number, mach_number, lattice, char_length_lu)[source]
-
lettuce.flows.obstacle.
Obstacle3D
(resolution_x, resolution_y, resolution_z, reynolds_number, mach_number, lattice, char_length_lu)[source]
Utility¶
Utility functions.
-
exception
lettuce.util.
LettuceException
[source] Bases:
Exception
-
exception
lettuce.util.
LettuceWarning
[source] Bases:
UserWarning
-
exception
lettuce.util.
InefficientCodeWarning
[source] Bases:
lettuce.util.LettuceWarning
-
exception
lettuce.util.
ExperimentalWarning
[source] Bases:
lettuce.util.LettuceWarning
-
lettuce.util.
get_subclasses
(classname, module)[source]
-
lettuce.util.
torch_gradient
(f, dx=1, order=2)[source] Function to calculate the first derivative of tensors. Orders O(h²); O(h⁴); O(h⁶) are implemented.
Notes
See [1]. The function only works for periodic domains
References
[1] Fornberg B. (1988) Generation of Finite Difference Formulas on Arbitrarily Spaced Grids, Mathematics of Computation 51, no. 184 : 699-706. PDF.
-
lettuce.util.
torch_jacobi
(f, p, dx, device, dim, tol_abs=1e-10, max_num_steps=100000)[source] Jacobi solver for the Poisson pressure equation
-
lettuce.util.
grid_fine_to_coarse
(lattice, f_fine, tau_fine, tau_coarse)[source]
-
lettuce.util.
pressure_poisson
(units, u, rho0, tol_abs=1e-10, max_num_steps=100000)[source] Solve the pressure poisson equation using a jacobi scheme.
Parameters: - units (lettuce.UnitConversion) – The flow instance.
- u (torch.Tensor) – The velocity tensor.
- rho0 (torch.Tensor) – Initial guess for the density (i.e., pressure).
- tol_abs (float) – The tolerance for pressure convergence.
Returns: rho – The converged density (i.e., pressure).
Return type: torch.Tensor
-
lettuce.util.
append_axes
(array, n)[source]
Command-Line Interface¶
Console script for lettuce. To get help for terminal commands, open a console and type:
>>> lettuce --help