Welcome to lettuce’s documentation!

https://raw.githubusercontent.com/lettucecfd/lettuce/master/.source/img/logo_lettuce_typo.png CI Status Documentation Status https://img.shields.io/lgtm/grade/python/g/lettucecfd/lettuce.svg?logo=lgtm&logoWidth=18

GPU-accelerated Lattice Boltzmann Simulations in Python

Lettuce is a Computational Fluid Dynamics framework based on the lattice Boltzmann method (LBM).

It provides

  • GPU-accelerated computation based on PyTorch
  • Rapid Prototyping in 2D and 3D
  • Usage of neural networks and automatic differentiation within LBM

Resources

When using lettuce please cite:

@inproceedings{bedrunka2021lettuce,
  title={Lettuce: PyTorch-Based Lattice Boltzmann Framework},
  author={Bedrunka, Mario Christopher and Wilde, Dominik and Kliemank, Martin and Reith, Dirk and Foysi, Holger and Kr{\"a}mer, Andreas},
  booktitle={High Performance Computing: ISC High Performance Digital 2021 International Workshops, Frankfurt am Main, Germany, June 24--July 2, 2021, Revised Selected Papers},
  pages={40},
  organization={Springer Nature}
}

Getting Started

The following Python code will run a two-dimensional Taylor-Green vortex on a GPU:

import torch
from lettuce import BGKCollision, StandardStreaming, Lattice, D2Q9, TaylorGreenVortex2D, Simulation

device = "cuda:0"   # for running on cpu: device = "cpu"
dtype = torch.float32

lattice = Lattice(D2Q9, device, dtype)
flow = TaylorGreenVortex2D(resolution=256, reynolds_number=10, mach_number=0.05, lattice=lattice)
collision = BGKCollision(lattice, tau=flow.units.relaxation_parameter_lu)
streaming = StandardStreaming(lattice)
simulation = Simulation(flow=flow, lattice=lattice,  collision=collision, streaming=streaming)
mlups = simulation.step(num_steps=1000)

print("Performance in MLUPS:", mlups)

More advanced examples are available as jupyter notebooks:

Installation

  • Install the anaconda package manager from www.anaconda.org

  • Create a new conda environment and install all dependencies:

    conda create -n lettuce -c pytorch -c conda-forge\
         "pytorch>=1.2" matplotlib pytest click cudatoolkit "pyevtk>=1.2"
    
  • Activate the conda environment:

    conda activate lettuce
    
  • Clone this repository from github

  • Change into the cloned directory

  • Run the install script:

    python setup.py install
    
  • Run the test cases:

    python setup.py test
    
  • Check out the convergence order, running on CPU:

    lettuce --no-cuda convergence
    
  • For running a CUDA-driven LBM simulation on one GPU omit the –no-cuda. If CUDA is not found, make sure that cuda drivers are installed and compatible with the installed cudatoolkit (see conda install command above).

  • Check out the performance, running on GPU:

    lettuce benchmark
    

Credits

We use the following third-party packages:

This package was created with Cookiecutter and the audreyr/cookiecutter-pypackage project template.

License

  • Free software: MIT license, as found in the LICENSE file.

Installation

From sources

The sources for lettuce can be downloaded from the Github repo.

You can either clone the public repository:

$ git clone git://github.com/lettucecfd/lettuce

Or download the tarball:

$ curl  -OL https://github.com/lettucecfd/lettuce/tarball/master

Once you have a copy of the source, you can install it with:

$ python setup.py install

Usage

To use lettuce in a project:

import lettuce

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]
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

class lettuce.streaming.StandardStreaming(lattice)[source]

Bases: object

Standard Streaming step on a regular grid.

no_stream_mask

Boolean mask with the same shape as the distribution function f. If None, stream all (also around all boundaries).

Type:torch.Tensor
no_stream_mask

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

class lettuce.observables.EnergySpectrum(lattice, flow)[source]

Bases: lettuce.observables.Observable

The kinetic energy spectrum

spectrum_from_u(u)[source]

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

class lettuce.force.Guo(lattice, tau, acceleration)[source]

Bases: object

source_term(u)[source]
u_eq(f)[source]
ueq_scaling_factor
class lettuce.force.ShanChen(lattice, tau, acceleration)[source]

Bases: object

source_term(u)[source]
u_eq(f)[source]
ueq_scaling_factor

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.

class lettuce.boundary.EquilibriumOutletP(lattice, direction, rho0=1.0)[source]

Bases: lettuce.boundary.AntiBounceBackOutlet

Equilibrium outlet with constant pressure.

make_no_collision_mask(f_shape)[source]
make_no_stream_mask(f_shape)[source]

Flows

Couette

Couette Flow

class lettuce.flows.couette.CouetteFlow2D(resolution, reynolds_number, mach_number, lattice)[source]

Bases: object

analytic_solution(x, t=0)[source]
boundaries
grid
initial_solution(x)[source]

Poiseuille

Poiseuille Flow

class lettuce.flows.poiseuille.PoiseuilleFlow2D(resolution, reynolds_number, mach_number, lattice, initialize_with_zeros=True)[source]

Bases: object

acceleration
analytic_solution(grid)[source]
boundaries
grid
initial_solution(grid)[source]

Taylor-Green

Taylor-Green vortex in 2D and 3D.

class lettuce.flows.taylorgreen.TaylorGreenVortex2D(resolution, reynolds_number, mach_number, lattice)[source]

Bases: object

analytic_solution(x, t=0)[source]
boundaries
grid
initial_solution(x)[source]
class lettuce.flows.taylorgreen.TaylorGreenVortex3D(resolution, reynolds_number, mach_number, lattice)[source]

Bases: object

boundaries
grid
initial_solution(x)[source]

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

Contributing

Contributions are welcome, and they are greatly appreciated! Every little bit helps, and credit will always be given.

You can contribute in many ways:

Types of Contributions

Report Bugs

Report bugs at https://github.com/lettucecfd/lettuce/issues.

If you are reporting a bug, please include:

  • Your operating system name and version.
  • Any details about your local setup that might be helpful in troubleshooting.
  • Detailed steps to reproduce the bug.

Fix Bugs

Look through the GitHub issues for bugs. Anything tagged with “bug” and “help wanted” is open to whoever wants to implement it.

Implement Features

Look through the GitHub issues for features. Anything tagged with “enhancement” and “help wanted” is open to whoever wants to implement it.

Write Documentation

lettuce could always use more documentation, whether as part of the official lettuce docs, in docstrings, or even on the web in blog posts, articles, and such.

Submit Feedback

The best way to send feedback is to file an issue at https://github.com/lettucecfd/lettuce/issues.

If you are proposing a feature:

  • Explain in detail how it would work.
  • Keep the scope as narrow as possible, to make it easier to implement.
  • Remember that this is a volunteer-driven project, and that contributions are welcome :)

Get Started!

Ready to contribute? Here’s how to set up lettuce for local development.

  1. Fork the lettuce repo on GitHub.

  2. Clone your fork locally:

    $ git clone git@github.com:your_name_here/lettuce.git
    
  3. Install your local copy into a virtualenv. Assuming you have virtualenvwrapper installed, this is how you set up your fork for local development:

    $ mkvirtualenv lettuce
    $ cd lettuce/
    $ python setup.py develop
    
  4. Create a branch for local development:

    $ git checkout -b name-of-your-bugfix-or-feature
    

    Now you can make your changes locally.

  5. When you’re done making changes, check that your changes pass flake8 and the tests, including testing other Python versions with tox:

    $ flake8 lettuce tests
    $ python setup.py test or py.test
    $ tox
    

    To get flake8 and tox, just pip install them into your virtualenv.

  6. Commit your changes and push your branch to GitHub:

    $ git add .
    $ git commit -m "Your detailed description of your changes."
    $ git push origin name-of-your-bugfix-or-feature
    
  7. Submit a pull request through the GitHub website.

Pull Request Guidelines

Before you submit a pull request, check that it meets these guidelines:

  1. The pull request should include tests.
  2. If the pull request adds functionality, the docs should be updated. Put your new functionality into a function with a docstring, and add the feature to the list in README.rst.
  3. The pull request should work for Python 2.7, 3.4, 3.5 and 3.6, and for PyPy. Check https://travis-ci.org/lettucecfd/lettuce/pull_requests and make sure that the tests pass for all supported Python versions.

Tips

To run a subset of tests:

$ py.test tests.test_lettuce

Deploying

A reminder for the maintainers on how to deploy. Make sure all your changes are committed (including an entry in HISTORY.rst). Then run:

$ bumpversion patch # possible: major / minor / patch
$ git push
$ git push --tags

Travis will then deploy to PyPI if tests pass.

Credits

Development Lead

Contributors

  • Dominik Wilde
  • Mario Bedrunka
  • You?

History

0.2.2 (2022-01-12)

Move to lettucecfd organization

Fixes:

  • move CI to github actions
  • readme with links to presentation and manuscript

Enhancements:

  • D3Q15 stencil
  • D3Q27 moment transform
  • Porous media example
  • Code reformatting

0.2.1 (2021-03-30)

Manuscript submission

0.2.0 (2020-04-20)

0.1.0 (2019-05-02)

  • First release

Indices and tables