[1]:
import sympy as sp
from lbmpy.chapman_enskog import ChapmanEnskogAnalysis
from lbmpy.creationfunctions import create_lb_method
sp.init_printing()

Demo: Automatic Chapman Enskog Analysis

First, we create a SRT lattice Boltzmann method. It is defined as the set of moments, together with one relaxation rate per moment.

[2]:
method = create_lb_method(method='trt', stencil='D3Q19', compressible=False)
method
[2]:
Moment Eq. Value Relaxation Rate
$1$ $\rho$ $\omega_{0}$
$x$ $u_{0}$ $\omega_{1}$
$y$ $u_{1}$ $\omega_{1}$
$z$ $u_{2}$ $\omega_{1}$
$x^{2}$ $\frac{\rho}{3} + u_{0}^{2}$ $\omega_{0}$
$y^{2}$ $\frac{\rho}{3} + u_{1}^{2}$ $\omega_{0}$
$z^{2}$ $\frac{\rho}{3} + u_{2}^{2}$ $\omega_{0}$
$x y$ $u_{0} u_{1}$ $\omega_{0}$
$x z$ $u_{0} u_{2}$ $\omega_{0}$
$y z$ $u_{1} u_{2}$ $\omega_{0}$
$x^{2} y$ $\frac{u_{1}}{3}$ $\omega_{1}$
$x^{2} z$ $\frac{u_{2}}{3}$ $\omega_{1}$
$x y^{2}$ $\frac{u_{0}}{3}$ $\omega_{1}$
$x z^{2}$ $\frac{u_{0}}{3}$ $\omega_{1}$
$y^{2} z$ $\frac{u_{2}}{3}$ $\omega_{1}$
$y z^{2}$ $\frac{u_{1}}{3}$ $\omega_{1}$
$x^{2} y^{2}$ $\frac{\rho}{9} + \frac{u_{0}^{2}}{3} + \frac{u_{1}^{2}}{3}$ $\omega_{0}$
$x^{2} z^{2}$ $\frac{\rho}{9} + \frac{u_{0}^{2}}{3} + \frac{u_{2}^{2}}{3}$ $\omega_{0}$
$y^{2} z^{2}$ $\frac{\rho}{9} + \frac{u_{1}^{2}}{3} + \frac{u_{2}^{2}}{3}$ $\omega_{0}$

Next, the Chapman Enskog analysis object is created. This may take a while…

[3]:
analysis = ChapmanEnskogAnalysis(method)

This object now information about the method, e.g. the relation of relaxation rate to viscosities, if the method approximates the compressible or incompressible continuity equation …

[4]:
analysis.compressible
[4]:
False
[5]:
analysis.pressure_equation
[5]:
$$\left [ \frac{\rho}{3}\right ]$$
[6]:
analysis.get_kinematic_viscosity()
[6]:
$$- \frac{\omega_{0} - 2}{6 \omega_{0}}$$
[7]:
analysis.get_bulk_viscosity()
[7]:
$$- \frac{1}{9} + \frac{2}{9 \omega_{0}}$$

But also details of the analysis are available:

[8]:
sp.Matrix(analysis.get_macroscopic_equations())
[8]:
$$\left[\begin{matrix}{\partial_{t} \rho} + {\partial_{0} u_{0}} + {\partial_{1} u_{1}} + {\partial_{2} u_{2}}\\- \frac{\epsilon \omega_{0}}{2} {\partial_{2} {\Pi_{02}^{(1)}}} - \frac{\epsilon \omega_{0}}{2} {\partial_{1} {\Pi_{01}^{(1)}}} - \frac{\epsilon \omega_{0}}{2} {\partial_{0} {\Pi_{00}^{(1)}}} + \epsilon {\partial_{2} {\Pi_{02}^{(1)}}} + \epsilon {\partial_{1} {\Pi_{01}^{(1)}}} + \epsilon {\partial_{0} {\Pi_{00}^{(1)}}} + \frac{{\partial_{0} \rho}}{3} + {\partial_{t} u_{0}} + {\partial_{0} (u_{0}^{2}) } + {\partial_{1} (u_{0} u_{1}) } + {\partial_{2} (u_{0} u_{2}) }\\- \frac{\epsilon \omega_{0}}{2} {\partial_{2} {\Pi_{12}^{(1)}}} - \frac{\epsilon \omega_{0}}{2} {\partial_{1} {\Pi_{11}^{(1)}}} - \frac{\epsilon \omega_{0}}{2} {\partial_{0} {\Pi_{01}^{(1)}}} + \epsilon {\partial_{2} {\Pi_{12}^{(1)}}} + \epsilon {\partial_{1} {\Pi_{11}^{(1)}}} + \epsilon {\partial_{0} {\Pi_{01}^{(1)}}} + \frac{{\partial_{1} \rho}}{3} + {\partial_{t} u_{1}} + {\partial_{1} (u_{1}^{2}) } + {\partial_{0} (u_{0} u_{1}) } + {\partial_{2} (u_{1} u_{2}) }\\- \frac{\epsilon \omega_{0}}{2} {\partial_{2} {\Pi_{22}^{(1)}}} - \frac{\epsilon \omega_{0}}{2} {\partial_{1} {\Pi_{12}^{(1)}}} - \frac{\epsilon \omega_{0}}{2} {\partial_{0} {\Pi_{02}^{(1)}}} + \epsilon {\partial_{2} {\Pi_{22}^{(1)}}} + \epsilon {\partial_{1} {\Pi_{12}^{(1)}}} + \epsilon {\partial_{0} {\Pi_{02}^{(1)}}} + \frac{{\partial_{2} \rho}}{3} + {\partial_{t} u_{2}} + {\partial_{2} (u_{2}^{2}) } + {\partial_{0} (u_{0} u_{2}) } + {\partial_{1} (u_{1} u_{2}) }\end{matrix}\right]$$