.. DO NOT EDIT. .. THIS FILE WAS AUTOMATICALLY GENERATED BY SPHINX-GALLERY. .. TO MAKE CHANGES, EDIT THE SOURCE PYTHON FILE: .. "auto_how_to/example_turbulence.py" .. LINE NUMBERS ARE GIVEN BELOW. .. only:: html .. note:: :class: sphx-glr-download-link-note :ref:`Go to the end ` to download the full example code. .. rst-class:: sphx-glr-example-title .. _sphx_glr_auto_how_to_example_turbulence.py: Build an incompressible turbulence model ======================================== The :mod:`pybFoam.turbulence` submodule wraps OpenFOAM's incompressible RAS / LES factory. This how-to reads the shipped ``pitzDaily`` case (backward-facing step, k-epsilon RAS), constructs the turbulence model from Python, and inspects its bound state — ``nut``, ``k``, ``epsilon``, ``correct()``, ``divDevReff()``. Prerequisite: finished :doc:`/auto_tutorials/example_01_dictionaries`. .. GENERATED FROM PYTHON SOURCE LINES 16-18 Clone the case -------------- .. GENERATED FROM PYTHON SOURCE LINES 18-23 .. code-block:: Python from pybFoam import clone_example case = clone_example("pitzDaily") .. GENERATED FROM PYTHON SOURCE LINES 24-26 Build mesh and read U --------------------- .. GENERATED FROM PYTHON SOURCE LINES 26-44 .. code-block:: Python from pybFoam import ( Time, argList, createPhi, dictionary, fvMesh, volVectorField, ) from pybFoam.meshing import generate_blockmesh time = Time(argList([str(case), "-case", str(case)])) generate_blockmesh(time, dictionary.read(str(case / "system" / "blockMeshDict"))) mesh = fvMesh(time) U = volVectorField.read_field(mesh, "U") phi = createPhi(U) .. GENERATED FROM PYTHON SOURCE LINES 45-51 Construct the transport and turbulence models --------------------------------------------- ``singlePhaseTransportModel(U, phi)`` reads the Newtonian ``nu`` from ``constant/transportProperties``. ``incompressibleTurbulenceModel.New`` dispatches on ``simulationType`` and ``RASModel`` from ``turbulenceProperties`` to instantiate kEpsilon, kOmegaSST, etc. .. GENERATED FROM PYTHON SOURCE LINES 51-60 .. code-block:: Python from pybFoam.turbulence import ( incompressibleTurbulenceModel, singlePhaseTransportModel, ) transport = singlePhaseTransportModel(U, phi) model = incompressibleTurbulenceModel.New(U, phi, transport) .. GENERATED FROM PYTHON SOURCE LINES 61-65 Inspect the model ----------------- Each accessor returns a ``tmp``; call ``()`` to get the underlying field, then ``["internalField"]`` for the NumPy view. .. GENERATED FROM PYTHON SOURCE LINES 65-76 .. code-block:: Python import numpy as np nut = np.asarray(model.nut()()["internalField"]) k = np.asarray(model.k()()["internalField"]) eps = np.asarray(model.epsilon()()["internalField"]) print(f"nut : shape={nut.shape} min={nut.min():.3e} max={nut.max():.3e}") print(f"k : shape={k.shape} min={k.min():.3e} max={k.max():.3e}") print(f"epsilon : shape={eps.shape} min={eps.min():.3e} max={eps.max():.3e}") .. rst-class:: sphx-glr-script-out .. code-block:: none nut : shape=(12225,) min=0.000e+00 max=0.000e+00 k : shape=(12225,) min=3.750e-01 max=3.750e-01 epsilon : shape=(12225,) min=1.486e+01 max=1.486e+01 .. GENERATED FROM PYTHON SOURCE LINES 77-82 Update the eddy viscosity ------------------------- ``correct()`` is what a PISO/PIMPLE solver calls at the end of every outer iteration. For k-epsilon it solves the k and epsilon transport equations and recomputes ``nut = Cmu * k^2 / epsilon``. .. GENERATED FROM PYTHON SOURCE LINES 82-88 .. code-block:: Python model.correct() nut_after = np.asarray(model.nut()()["internalField"]) print(f"nut after correct() : min={nut_after.min():.3e} max={nut_after.max():.3e}") .. rst-class:: sphx-glr-script-out .. code-block:: none nut after correct() : min=1.938e-07 max=2.148e-03 .. GENERATED FROM PYTHON SOURCE LINES 89-97 divDevReff in a momentum equation --------------------------------- In a turbulent solver the laminar diffusion ``laplacian(nu, U)`` is replaced by the model's stress divergence: .. code-block:: python UEqn = fvVectorMatrix(fvm.ddt(U) + fvm.div(phi, U) + model.divDevReff(U)) .. GENERATED FROM PYTHON SOURCE LINES 97-101 .. code-block:: Python UEqn_term = model.divDevReff(U) print(f"divDevReff term : {type(UEqn_term).__name__}") .. rst-class:: sphx-glr-script-out .. code-block:: none divDevReff term : tmp_fvVectorMatrix .. GENERATED FROM PYTHON SOURCE LINES 102-109 See also -------- - :doc:`example_modify_boundary_conditions` — change a patch field (e.g. inlet velocity) before constructing the turbulence model. - :doc:`/auto_tutorials/example_04_run_and_sample` — drive a solver loop; add ``model.correct()`` inside it for a turbulent run. .. rst-class:: sphx-glr-timing **Total running time of the script:** (0 minutes 0.222 seconds) .. _sphx_glr_download_auto_how_to_example_turbulence.py: .. only:: html .. container:: sphx-glr-footer sphx-glr-footer-example .. container:: sphx-glr-download sphx-glr-download-jupyter :download:`Download Jupyter notebook: example_turbulence.ipynb ` .. container:: sphx-glr-download sphx-glr-download-python :download:`Download Python source code: example_turbulence.py ` .. container:: sphx-glr-download sphx-glr-download-zip :download:`Download zipped: example_turbulence.zip ` .. only:: html .. rst-class:: sphx-glr-signature `Gallery generated by Sphinx-Gallery `_