[1]:
from lbmpy.session import *
from lbmpy.chapman_enskog import ChapmanEnskogAnalysis
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]:
lb_config = LBMConfig(method=Method.TRT, stencil=Stencil.D3Q19, compressible=False, zero_centered=False)
method = create_lb_method(lbm_config=lb_config)
method
[2]:
Moment-Based Method | Stencil: D3Q19 | Zero-Centered Storage: ✗ | Force Model: None |
---|
Continuous Hydrodynamic Maxwellian Equilibrium | $f (\rho, \left( u_{0}, \ u_{1}, \ u_{2}\right), \left( v_{0}, \ v_{1}, \ v_{2}\right)) = \frac{3 \sqrt{6} \delta_{\rho} e^{- \frac{3 v_{0}^{2}}{2} - \frac{3 v_{1}^{2}}{2} - \frac{3 v_{2}^{2}}{2}}}{4 \pi^{\frac{3}{2}}} + \frac{3 \sqrt{6} e^{- \frac{3 \left(- u_{0} + v_{0}\right)^{2}}{2} - \frac{3 \left(- u_{1} + v_{1}\right)^{2}}{2} - \frac{3 \left(- u_{2} + v_{2}\right)^{2}}{2}}}{4 \pi^{\frac{3}{2}}}$ | ||
---|---|---|---|
Compressible: ✗ | Deviation Only: ✗ | Order: 2 |
Relaxation Info | ||
---|---|---|
Moment | Eq. Value | Relaxation Rate |
$1$ | $\rho$ | $\omega$ |
$x$ | $u_{0}$ | $\omega$ |
$y$ | $u_{1}$ | $\omega$ |
$z$ | $u_{2}$ | $\omega$ |
$x^{2}$ | $\frac{\rho}{3} + u_{0}^{2}$ | $\omega$ |
$y^{2}$ | $\frac{\rho}{3} + u_{1}^{2}$ | $\omega$ |
$z^{2}$ | $\frac{\rho}{3} + u_{2}^{2}$ | $\omega$ |
$x y$ | $u_{0} u_{1}$ | $\omega$ |
$x z$ | $u_{0} u_{2}$ | $\omega$ |
$y z$ | $u_{1} u_{2}$ | $\omega$ |
$x^{2} y$ | $\frac{u_{1}}{3}$ | $\omega$ |
$x^{2} z$ | $\frac{u_{2}}{3}$ | $\omega$ |
$x y^{2}$ | $\frac{u_{0}}{3}$ | $\omega$ |
$x z^{2}$ | $\frac{u_{0}}{3}$ | $\omega$ |
$y^{2} z$ | $\frac{u_{2}}{3}$ | $\omega$ |
$y z^{2}$ | $\frac{u_{1}}{3}$ | $\omega$ |
$x^{2} y^{2}$ | $\frac{\rho}{9} + \frac{u_{0}^{2}}{3} + \frac{u_{1}^{2}}{3}$ | $\omega$ |
$x^{2} z^{2}$ | $\frac{\rho}{9} + \frac{u_{0}^{2}}{3} + \frac{u_{2}^{2}}{3}$ | $\omega$ |
$y^{2} z^{2}$ | $\frac{\rho}{9} + \frac{u_{1}^{2}}{3} + \frac{u_{2}^{2}}{3}$ | $\omega$ |
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]:
$\displaystyle \left[ \frac{\rho}{3}\right]$
[6]:
analysis.get_kinematic_viscosity()
[6]:
$\displaystyle - \frac{\omega - 2}{6 \omega}$
[7]:
analysis.get_bulk_viscosity()
[7]:
$\displaystyle - \frac{1}{9} + \frac{2}{9 \omega}$
But also details of the analysis are available:
[8]:
sp.Matrix(analysis.get_macroscopic_equations())
[8]:
$\displaystyle \left[\begin{matrix}{\partial_{t} \rho} + {\partial_{0} u_{0}} + {\partial_{1} u_{1}} + {\partial_{2} u_{2}}\\- \frac{\epsilon \omega {\partial_{0} {\Pi_{00}^{(1)}}}}{2} - \frac{\epsilon \omega {\partial_{1} {\Pi_{01}^{(1)}}}}{2} - \frac{\epsilon \omega {\partial_{2} {\Pi_{02}^{(1)}}}}{2} + \epsilon {\partial_{0} {\Pi_{00}^{(1)}}} + \epsilon {\partial_{1} {\Pi_{01}^{(1)}}} + \epsilon {\partial_{2} {\Pi_{02}^{(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 {\partial_{0} {\Pi_{01}^{(1)}}}}{2} - \frac{\epsilon \omega {\partial_{1} {\Pi_{11}^{(1)}}}}{2} - \frac{\epsilon \omega {\partial_{2} {\Pi_{12}^{(1)}}}}{2} + \epsilon {\partial_{0} {\Pi_{01}^{(1)}}} + \epsilon {\partial_{1} {\Pi_{11}^{(1)}}} + \epsilon {\partial_{2} {\Pi_{12}^{(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 {\partial_{0} {\Pi_{02}^{(1)}}}}{2} - \frac{\epsilon \omega {\partial_{1} {\Pi_{12}^{(1)}}}}{2} - \frac{\epsilon \omega {\partial_{2} {\Pi_{22}^{(1)}}}}{2} + \epsilon {\partial_{0} {\Pi_{02}^{(1)}}} + \epsilon {\partial_{1} {\Pi_{12}^{(1)}}} + \epsilon {\partial_{2} {\Pi_{22}^{(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]$