4th Symposium on Fluid-Structure-Sound Interactions and Control
Month: December 2016
32th TSFD Symposium
第32回生研TSFD (Turbulence Simulation and Flow Design) シンポジウム 「乱流シミュレーションと流れの設計 -数値計算による熱流体現象の理解、予測、制御」 が,2017年3月22日(水) に東京大学生産技術研究所 (駒場IIリサーチキャンパス) において開催されるそうです.プログラム等はまだ未確定ですが,ぜひ勉強してきたいと思います.
WALE SGS model in OpenFOAM
Last Updated: May 6, 2019
The WALE (Wall-Adapting Local Eddy-viscosity) SGS model is one of the major SGS models. It is an algebraic eddy viscosity model (0-equation model) as with the Smagorinsky SGS model, but it has some excellent features that the Smagorinsky model does not have.
Implementation in OpenFOAM
In the WALE model, the subgrid scale eddy viscosity \(\nu_{sgs}\) is evaluated as
\begin{equation}
\nu_{sgs} = C_{k} \Delta \sqrt{k_{sgs}}, \tag{1} \label{eq:nusgs}
\end{equation}
where \(C_{k}\) is a model constant and \(k_{sgs}\) is the subgrid scale kinetic energy. You can find its definition in this post.
1 2 3 4 5 6 7 8 9 |
template<class BasicTurbulenceModel> void WALE<BasicTurbulenceModel>::correctNut() { this->nut_ = Ck_*this->delta()*sqrt(this->k(fvc::grad(this->U_))); this->nut_.correctBoundaryConditions(); fv::options::New(this->mesh_).correct(this->nut_); BasicTurbulenceModel::correctNut(); } |
The traceless symmetric part of the square of the velocity gradient tensor \(S^d\) is calculated in the following function.
\begin{eqnarray}
S_{ij}^d = \frac{1}{2} \left( \frac{\partial \overline{u}_k}{\partial x_i}\frac{\partial \overline{u}_j}{\partial x_k} + \frac{\partial \overline{u}_k}{\partial x_j}\frac{\partial \overline{u}_i}{\partial x_k} \right) – \frac{1}{3} \delta_{ij} \frac{\partial \overline{u}_k}{\partial x_l}\frac{\partial \overline{u}_l}{\partial x_k}, \tag{2} \label{eq:sd}
\end{eqnarray}
where \(\delta_{ij}\) is the Kronecker delta.
1 2 3 4 5 6 7 8 |
template<class BasicTurbulenceModel> tmp<volSymmTensorField> WALE<BasicTurbulenceModel>::Sd ( const volTensorField& gradU ) const { return dev(symm(gradU & gradU)); } |
These tensor operations used above are summarized in the following post:
The subgrid scale kinetic energy \(k_{sgs}\) is
\begin{equation}
k_{sgs} = \left( \frac{C_w^2 \Delta}{C_k} \right)^2 \frac{\left( S_{ij}^d S_{ij}^d \right)^3}{\left( \left( \overline{S}_{ij} \overline{S}_{ij} \right)^{5/2} + \left( S_{ij}^d S_{ij}^d \right)^{5/4} \right)^2}, \tag{3} \label{eq:ksgs}
\end{equation}
where
\begin{equation}
\overline{S}_{ij} = \frac{1}{2} \left( \frac{\partial \overline{u}_{i}}{\partial x_{j}} + \frac{\partial \overline{u}_{j}}{\partial x_{i}}\right), \tag{4} \label{eq:s}
\end{equation}
is the resolved-scale strain rate tensor.
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 |
template<class BasicTurbulenceModel> tmp<volScalarField> WALE<BasicTurbulenceModel>::k ( const volTensorField& gradU ) const { volScalarField magSqrSd(magSqr(Sd(gradU))); return volScalarField::New ( IOobject::groupName("k", this->alphaRhoPhi_.group()), sqr(sqr(Cw_)*this->delta()/Ck_)* ( pow3(magSqrSd) /( sqr ( pow(magSqr(symm(gradU)), 5.0/2.0) + pow(magSqrSd, 5.0/4.0) ) + dimensionedScalar ( "small", dimensionSet(0, 0, -10, 0, 0), small ) ) ) ); } |
Finally, we can get the following expression by substituting Eq. \eqref{eq:ksgs} into Eq. \eqref{eq:nusgs}:
\begin{eqnarray}
\nu_{sgs} = \left( C_w \Delta \right)^2 \frac{\left( S_{ij}^d S_{ij}^d \right)^{3/2}}{\left( \overline{S}_{ij} \overline{S}_{ij} \right)^{5/2} + \left( S_{ij}^d S_{ij}^d \right)^{5/4}}, \tag{5} \label{eq:nusgs2}
\end{eqnarray}
It is the same as Eq. (13) in [1].
Features
- Algebraic eddy viscosity model (0-equation model)
- The rotation rate is taken into account in the calculation of \(\nu_{sgs}\)
- To be able to handle transition
- Damping is Not necessary for \(\nu_{sgs}\) in the near-wall region
References
Final expression of the SGS eddy viscosity is relatively simple and the implementation into OpenFOAM is not so complicated but the process of deriving the expression described in [1] is not easy to understand. I want to develop my SGS model someday 🙂
Thanks 2016!
Thank you very much for visiting my page and sending encouragement messages.
I launched this blog this year and I’m glad that many visitors all over the world have read my posts.
I’m going to continue this blog next year. Additionally, I want to collaborate with many contributors all over the world to create a synergy effect and provide more useful, reliable and accurate information.
Have a happy New Year!
Fumiya
Some Topics of My Interest on (Computational) Fluid Dynamics
I will accumulate more and more information of some topics of interest in the field of fluid dynamics. When I get a better understanding about the following topics, I will prepare dedicated pages.
Physical mechanism of energy cascade |
I want to understand the physical mechanism of the energy cascade in turbulent flows and have found some useful references:
- S. Goto and G. Kawahara, GENERATION MECHANISM OF HIERARCHY OF COHERENT VORTICES IN TURBULENCE SUSTAINED BY STEADY FORCE
- H. Tennekes and J. L. Lumley, A FIRST COURSE IN TURBULENCE
If you have any other information on this subject, could you inform me of the current status of research?
Boundary layer flow |
I found a very informative and beautiful video of a spatially developing turbulent boundary layer on a flat plate.
A video is worth ten thousand words!
Video credit: J. H. Lee, Y. S. Kwon, N. Hutchins and J. P. Monty
As we can see, the turbulent boundary layer becomes thinner and the eddies are smaller in the higher Reynolds number \(Re_{\tau}\) case. The Reynolds number \(Re_{\tau}\) based on friction velocity \(u_{\tau}\) is often used in studies of wall-bounded turbulence, e.g., the turbulent channel flow.
Wall-Modeled Large Eddy Simulation |
TTT Technical Highlight: Wall-Modeled Large Eddy Simulation of High Reynolds Number Flows
I will add topics and update information regularly!
Tensor Operations in OpenFOAM
Last Updated: May 5, 2019
The gradient of the velocity field is a strain-rate tensor field, that is, a second rank tensor field. It appears in the diffusion term of the Navier-Stokes equation.
A second rank tensor has nine components and can be expressed as a 3×3 matrix as shown in the above image.
In this blog post, I will pick out some typical tensor operations and give brief explanations of them with some usage examples in OpenFOAM.
Keywords
strain rate tensor, vorticity tensor, Q-criterion, Hodge dual
Gradient of a Vector Field | fvc::grad(u)
The gradient of a velocity vector \(\boldsymbol{u}\) returns a velocity gradient tensor (second rank tensor).
\begin{eqnarray}
\nabla \boldsymbol{u} &\equiv& \partial_{i} u_j \\
&=& \left(
\begin{matrix}
\partial u_1/\partial x_1 & \partial u_2/\partial x_1 & \partial u_3/\partial x_1 \\
\partial u_1/\partial x_2 & \partial u_2/\partial x_2 & \partial u_3/\partial x_2 \\
\partial u_1/\partial x_3 & \partial u_2/\partial x_3 & \partial u_3/\partial x_3
\end{matrix} \right)
\end{eqnarray}
Symmetric Part of a Second Rank Tensor | symm(T)
\begin{eqnarray}
{\rm symm}(\boldsymbol{T}) &\equiv& \frac{1}{2} (\boldsymbol{T} + \boldsymbol{T}^T) \\
&=& \frac{1}{2} \left(
\begin{matrix}
2T_{11} & T_{12} + T_{21} & T_{13} + T_{31} \\
T_{21} + T_{12} & 2T_{22} & T_{23} + T_{32} \\
T_{31} + T_{13} & T_{32} + T_{23} & 2T_{33}
\end{matrix} \right)
\end{eqnarray}
1 2 3 4 5 6 7 8 9 10 11 |
//- Return the symmetric part of a tensor template<class Cmpt> inline SymmTensor<Cmpt> symm(const Tensor<Cmpt>& t) { return SymmTensor<Cmpt> ( t.xx(), 0.5*(t.xy() + t.yx()), 0.5*(t.xz() + t.zx()), t.yy(), 0.5*(t.yz() + t.zy()), t.zz() ); } |
Twice the Symmetric Part of a Second Rank Tensor | twoSymm(T)
\begin{eqnarray}
{\rm twoSymm}(\boldsymbol{T}) &\equiv& \boldsymbol{T} + \boldsymbol{T}^T \\
&=& \left(
\begin{matrix}
2T_{11} & T_{12} + T_{21} & T_{13} + T_{31} \\
T_{21} + T_{12} & 2T_{22} & T_{23} + T_{32} \\
T_{31} + T_{13} & T_{32} + T_{23} & 2T_{33}
\end{matrix} \right)
\end{eqnarray}
1 2 3 4 5 6 7 8 9 10 11 |
//- Return twice the symmetric part of a tensor template<class Cmpt> inline SymmTensor<Cmpt> twoSymm(const Tensor<Cmpt>& t) { return SymmTensor<Cmpt> ( 2*t.xx(), (t.xy() + t.yx()), (t.xz() + t.zx()), 2*t.yy(), (t.yz() + t.zy()), 2*t.zz() ); } |
Skew-symmetric Part of a Second Rank Tensor | skew(T)
\begin{eqnarray}
{\rm skew}(\boldsymbol{T}) &\equiv& \frac{1}{2} (\boldsymbol{T} – \boldsymbol{T}^T) \\
&=& \frac{1}{2} \left(
\begin{matrix}
0 & T_{12} – T_{21} & T_{13} – T_{31} \\
T_{21} – T_{12} & 0 & T_{23} – T_{32} \\
T_{31} – T_{13} & T_{32} – T_{23} & 0
\end{matrix} \right)
\end{eqnarray}
1 2 3 4 5 6 7 8 9 10 11 |
//- Return the skew-symmetric part of a tensor template<class Cmpt> inline Tensor<Cmpt> skew(const Tensor<Cmpt>& t) { return Tensor<Cmpt> ( 0.0, 0.5*(t.xy() - t.yx()), 0.5*(t.xz() - t.zx()), 0.5*(t.yx() - t.xy()), 0.0, 0.5*(t.yz() - t.zy()), 0.5*(t.zx() - t.xz()), 0.5*(t.zy() - t.yz()), 0.0 ); } |
The symm and skew operations of the velocity gradient tensor field frequently appear in the source code. The following is a typical example.
1 2 3 4 5 6 |
tmp<volTensorField> tgradU(fvc::grad(U)); const volTensorField& gradU = tgradU(); volSymmTensorField b(dev(R)/(2*k_)); volSymmTensorField S(symm(gradU)); volTensorField Omega(skew(gradU)); |
In the above code, the symmetric and antisymmetric parts of the velocity gradient tensor \(\partial u_j/\partial x_i\) are defined as follows:
\begin{eqnarray}
S_{ij} &=& \frac{1}{2} \left( \frac{\partial u_j}{\partial x_i} + \frac{\partial u_i}{\partial x_j} \right), \\
\Omega_{ij} &=& \frac{1}{2} \left( \frac{\partial u_j}{\partial x_i} – \frac{\partial u_i}{\partial x_j} \right),
\end{eqnarray}
where \(S_{ij}\) is the strain rate tensor and \(-\Omega_{ij}\) is the vorticity (spin) tensor.
Hodge Dual | *T
\begin{equation}
*T = \left( T_{23},\;-T_{13},\;T_{12} \right)
\end{equation}
1 2 3 4 5 |
template<class Cmpt> inline Vector<Cmpt> operator*(const Tensor<Cmpt>& t) { return Vector<Cmpt>(t.yz(), -t.xz(), t.xy()); } |
The vorticity vector \(\boldsymbol{\omega}\) is calculated as the Hodge dual of the skew-symmetric part of the velocity gradient tensor.
\begin{eqnarray}
\boldsymbol{\omega} &=& 2 \times \left( * \Omega_{ij} \right) \\
&=& 2 \times \left( \Omega_{23},\;-\Omega_{13},\;\Omega_{12} \right) \\
&=& \left( \frac{\partial u_3}{\partial x_2}-\frac{\partial u_2}{\partial x_3},\;\frac{\partial u_1}{\partial x_3}-\frac{\partial u_3}{\partial x_1},\;\frac{\partial u_2}{\partial x_1}-\frac{\partial u_1}{\partial x_2} \right)
\end{eqnarray}
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 |
template<class Type> tmp<GeometricField<Type, fvPatchField, volMesh>> curl ( const GeometricField<Type, fvPatchField, volMesh>& vf ) { word nameCurlVf = "curl(" + vf.name() + ')'; // Gausses theorem curl // tmp<GeometricField<Type, fvPatchField, volMesh>> tcurlVf = // fvc::surfaceIntegrate(vf.mesh().Sf() ^ fvc::interpolate(vf)); // Calculate curl as the Hodge dual of the skew-symmetric part of grad tmp<GeometricField<Type, fvPatchField, volMesh>> tcurlVf = 2.0*(*skew(fvc::grad(vf, nameCurlVf))); tcurlVf.ref().rename(nameCurlVf); return tcurlVf; } |
The operator * in front of the skew is the Hodge dual operator.
Inner Product of Two Second Rank Tensors | T & S
\begin{equation}
P_{ij} = T_{ik} S_{kj}
\end{equation}
1 2 3 4 5 6 7 8 9 10 11 |
tmp<volSymmTensorField> tS = dev(symm(gradU)); const volSymmTensorField& S = tS(); volScalarField W ( (2*sqrt(2.0))*((S&S)&&S) /( magS*S2 + dimensionedScalar("small", dimensionSet(0, 0, -3, 0, 0), SMALL) ) ); |
Double Inner Product of Two Second Rank Tensors | T && S
\begin{align}
s = T_{ij}S_{ij} &= T_{11}S_{11} + T_{12}S_{12} + T_{13}S_{13} \\
&+ T_{21}S_{21} + T_{22}S_{22} + T_{23}S_{23} \\
&+ T_{31}S_{31} + T_{32}S_{32} + T_{33}S_{33}
\end{align}
1 2 3 4 5 6 |
tmp<volTensorField> tgradU = fvc::grad(U); volScalarField::Internal G ( this->GName(), nut.v()*(dev(twoSymm(tgradU().v())) && tgradU().v()) ); |
Trace of a Second Rank Tensor | tr(T)
\begin{equation}
{\rm tr}(\boldsymbol{T}) \equiv T_{11} + T_{22} + T_{33} = \boldsymbol{T} {\rm \&}{\rm \&} \boldsymbol{I}
\end{equation}
1 2 3 4 5 6 |
//- Return the trace of a tensor template<class Cmpt> inline Cmpt tr(const Tensor<Cmpt>& t) { return t.xx() + t.yy() + t.zz(); } |
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 |
bool Foam::functionObjects::Q::calc() { if (foundObject<volVectorField>(fieldName_)) { const volVectorField& U = lookupObject<volVectorField>(fieldName_); const tmp<volTensorField> tgradU(fvc::grad(U)); const volTensorField& gradU = tgradU(); return store ( resultName_, 0.5*(sqr(tr(gradU)) - tr(((gradU) & (gradU)))) ); } else { return false; } } |
The Q-criterion can be used to identify vortex cores:
\begin{align}
Q &= \frac{1}{2} \left[ (\nabla \cdot \boldsymbol{u})^2 – {\rm tr}(\nabla \boldsymbol{u}^2)\right] \\
&= \frac{1}{2} \left[ (\nabla \cdot \boldsymbol{u})^2 + \Omega_{ij}\Omega_{ij} – S_{ij}S_{ij} \right].
\end{align}
For incompressible flows, it can be simplified as follows:
\begin{equation}
Q = \frac{1}{2} \left[ \Omega_{ij}\Omega_{ij} – S_{ij}S_{ij} \right].
\end{equation}
References
[1] OpenFOAM Programmer’s Guide
[2] CFD Direct | Tensor Mathematics
National Committee for Fluid Mechanics Films
We can find many informative films on fluid mechanics in the National Committee for Fluid Mechanics Films (NCFMF) site.
Notes are also available for most of the movies.
Turbulence: An Introduction for Scientists and Engineers by P. A. Davidson
“Turbulence: An Introduction for Scientists and Engineers” by P. A. Davidson has been translated into Japanese and is going to be published soon. I’m going to read it during year-end through New Year holidays.
One equation eddy-viscosity SGS model in OpenFOAM
OpenFOAM Version: OpenFOAM-dev, OpenFOAM-1606+
Implementation in OpenFOAM |
Assumption 1: Eddy viscosity approximation
As in the case of the Smagorinsky SGS model, the one equation eddy viscosity SGS model uses the eddy viscosity approximation as the name suggests, so the subgrid scale stress tensor \(\tau_{ij}\) is approximated as follows:
\begin{equation}
\tau_{ij} \approx \frac{2}{3} k_{sgs} \delta_{ij} – 2 \nu_{sgs} dev(\overline{D})_{ij}, \tag{1} \label{eq:tauij}
\end{equation}
where \(\nu_{sgs}\) is the subgrid scale eddy viscosity, \(\overline{D}\) is the resolved-scale strain rate tensor defined as
\begin{equation}
\overline{D}_{ij} = \frac{1}{2} \left( \frac{\partial \overline{u}_{i}}{\partial x_{j}} + \frac{\partial \overline{u}_{j}}{\partial x_{i}}\right), \tag{2} \label{eq:Dij}
\end{equation}
and the subgrid scale kinetic energy \(k_{sgs}\) is
\begin{equation}
k_{sgs} = \frac{1}{2} \tau_{kk} = \frac{1}{2} \left( \overline{u_{k}u_{k}} – \overline{u}_{k}\overline{u}_{k} \right). \tag{3} \label{eq:ksgs}
\end{equation}
The subgrid scale eddy viscosity \(\nu_{sgs}\) is computed using \(k_{sgs}\)
\begin{equation}
\nu_{sgs} = C_{k} \sqrt{k_{sgs}} \Delta, \tag{4} \label{eq:nusgs}
\end{equation}
where \(C_{k}\) is a model constant whose default value is \(0.094\).
1 2 3 4 5 6 7 8 9 |
template<class BasicTurbulenceModel> void kEqn<BasicTurbulenceModel>::correctNut() { this->nut_ = Ck_*sqrt(k_)*this->delta(); this->nut_.correctBoundaryConditions(); fv::options::New(this->mesh_).correct(this->nut_); BasicTurbulenceModel::correctNut(); } |
The procedure so far is the same as the Smagorinsky SGS model but there is a difference in what follows. These models are different in terms of how to compute \(k_{sgs}\). The Smagorinsky model assumes the local equilibrium to compute \(k_{sgs}\) but the one equation eddy viscosity model solves a transport equation of \(k_{sgs}\) [1].
The second category of SGS model is one-equation eddy viscosity models. The main reason to develop the one-equation SGS models is to overcome the deficiency of local balance assumption between the SGS energy production and dissipation adopted in algebraic eddy viscosity models. Such a phenomenon may occur in high Reynolds number flows and/or in the cases of coarse grid resolution. The first one-equation eddy viscosity SGS model was theoretically derived by Yoshizawa and Horiuti [2], in which the SGS kinetic energy is defined as \(k_{sgs} = \frac{1}{2}(\overline{u_k^2} – \overline{u}_k^2)\), and the SGS eddy viscosity, \(\nu_{sgs}\), is computed using \(k_{sgs}\) as \(\nu_{sgs} = C_{k} k_{sgs}^{1/2} \Delta\). A transportation equation is derived to account for the historic effect of \(k_{sgs}\) due to production, dissipation and diffusion:
Transport Equation of \(k_{sgs}\)
\begin{align}
&\frac{\partial (\rho k_{sgs})}{\partial t} + \frac{\partial (\rho \overline{u}_j k_{sgs})}{\partial x_{j}} – \frac{\partial}{\partial x_{j}} \left[ \rho \left(\nu + \nu_{sgs}\right) \frac{\partial k_{sgs}}{\partial x_{j}} \right] \\
&= \; – \rho \tau_{ij} : \overline{D}_{ij} \; – \; C_{\epsilon} \frac{\rho k_{sgs}^{3/2}}{\Delta}, \tag{5} \label{eq:ksgsEqn}
\end{align}
where the operator : is a double inner product, \(\rho\) is the density, \(\nu\) is the kinetic viscosity and \(C_{\epsilon}\) is another model constant. The terms in Eq. \eqref{eq:ksgsEqn} are, starting from the left, the time derivative term, convective term, diffusion term, production term and dissipation term. In the case of the Smagorinsky SGS model, only the production and dissipation terms are taken into account with the assumption of the local equilibrium.
The production term in Eq. \eqref{eq:ksgsEqn} can be rearranged to yield the expression used in the source code as follows:
\begin{align}
– \rho \tau_{ij} : \overline{D}_{ij} &= \left( – \frac{2}{3} \rho k_{sgs} \delta_{ij} + 2 \rho \nu_{sgs} dev(\overline{D})_{ij} \right) : \overline{D}_{ij} \\
&=\;- \frac{2}{3} \rho k_{sgs} \frac{\partial \overline{u}_{k}}{\partial x_{k}} \\
&\;\;\;+ \rho \nu_{sgs} \frac{\partial \overline{u}_i}{\partial x_{j}} \left( 2\overline{D}_{ij} – \frac{1}{3} tr(2\overline{D}) \delta_{ij} \right). \tag{6} \label{eq:production}
\end{align}
Please see the appendix below for more detailed process of this deformation.
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 |
template<class BasicTurbulenceModel> void kEqn<BasicTurbulenceModel>::correct() { if (!this->turbulence_) { return; } // Local references const alphaField& alpha = this->alpha_; const rhoField& rho = this->rho_; const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_; const volVectorField& U = this->U_; volScalarField& nut = this->nut_; fv::options& fvOptions(fv::options::New(this->mesh_)); LESeddyViscosity<BasicTurbulenceModel>::correct(); volScalarField divU(fvc::div(fvc::absolute(this->phi(), U))); tmp<volTensorField> tgradU(fvc::grad(U)); volScalarField G(this->GName(), nut*(tgradU() && dev(twoSymm(tgradU())))); tgradU.clear(); tmp<fvScalarMatrix> kEqn ( fvm::ddt(alpha, rho, k_) + fvm::div(alphaRhoPhi, k_) - fvm::laplacian(alpha*rho*DkEff(), k_) == alpha*rho*G - fvm::SuSp((2.0/3.0)*alpha*rho*divU, k_) - fvm::Sp(this->Ce_*alpha*rho*sqrt(k_)/this->delta(), k_) + kSource() + fvOptions(alpha, rho, k_) ); kEqn.ref().relax(); fvOptions.constrain(kEqn.ref()); solve(kEqn); fvOptions.correct(k_); bound(k_, this->kMin_); correctNut(); } |
1 2 3 4 5 6 7 8 |
//- Return the effective diffusivity for k tmp<volScalarField> DkEff() const { return tmp<volScalarField> ( new volScalarField("DkEff", this->nut_ + this->nu()) ); } |
References |
[1] S. Huang and Q. S. Li, A new dynamic one-equation subgrid-scale model for large eddy simulations. INTERNATIONAL JOURNAL FOR NUMERICAL METHODS IN ENGINEERING, 81, 835-865, 2010.
[2] A. Yoshizawa and K. Horiuti, A Statistically-Derived Subgrid-Scale Kinetic Energy Model for the Large-Eddy Simulation of Turbulent Flows. Journal of the Physical Society of Japan, 54(8), 2834-2839, 1985.
Appendix: Deformation of the Production Term |
\begin{align}
&\; 2 \rho \nu_{sgs} dev(\overline{D})_{ij} : \overline{D}_{ij} \\
=&\; \rho \nu_{sgs} \frac{1}{2} \left( \frac{\partial \overline{u}_i}{\partial x_j} + \frac{\partial \overline{u}_j}{\partial x_i} – \frac{2}{3} tr(\overline{D}) \delta_{ij} \right) \left( \frac{\partial \overline{u}_i}{\partial x_j} + \frac{\partial \overline{u}_j}{\partial x_i} \right) \\
=&\; \rho \nu_{sgs} \frac{1}{2} \left( 2 \frac{\partial \overline{u}_i}{\partial x_j}\frac{\partial \overline{u}_i}{\partial x_j} + 2 \frac{\partial \overline{u}_i}{\partial x_j}\frac{\partial \overline{u}_j}{\partial x_i} – \frac{4}{3} tr(\overline{D}) \delta_{ij} \frac{\partial \overline{u}_i}{\partial x_j} \right) \\
=&\; \rho \nu_{sgs} \frac{\partial \overline{u}_i}{\partial x_j} \left( \frac{\partial \overline{u}_i}{\partial x_j} + \frac{\partial \overline{u}_j}{\partial x_i} – \frac{2}{3} tr(\overline{D}) \delta_{ij} \right) \\
=&\; \rho \nu_{sgs} \frac{\partial \overline{u}_i}{\partial x_j} \left( 2 \overline{D}_{ij} \; – \frac{1}{3} tr(\overline{2D}) \delta_{ij} \right) \\
=&\; \rho G \tag{7} \label{eq:productionTerm}
\end{align}
Smagorinsky SGS model in OpenFOAM
Last Updated: May 4, 2019
The Smagorinsky subgrid scale (SGS) model was developed by Joseph Smagorinsky in the meteorological community in the 1960s. It is based on the eddy viscosity assumption, which postulates a linear relationship between the SGS shear stress and the resolved rate of strain tensor. This model serves as a base for other SGS models.
OpenFOAM Version: OpenFOAM-dev, OpenFOAM-1606+
Implementation in OpenFOAM
The Smagorinsky SGS model is the oldest and best known subgrid scale model. The subgrid scale stress tensor \(\tau_{ij}\) is
\begin{align}
\tau_{ij} &= \overline{u_{i}u_{j}} – \overline{u}_{i}\overline{u}_{j} \tag{1a}\\
&= \frac{1}{3} \tau_{kk} \delta_{ij} + \left( \tau_{ij} \; – \frac{1}{3} \tau_{kk} \delta_{ij} \right) \tag{1b} \label{eq:smg2}\\
&\approx \frac{1}{3} \tau_{kk} \delta_{ij} – 2 \nu_{sgs} dev(\overline{D})_{ij} \tag{1c} \label{eq:smg3}\\
&= \frac{2}{3} k_{sgs} \delta_{ij} – 2 \nu_{sgs} dev(\overline{D})_{ij}, \tag{1d} \label{eq:smg4}
\end{align}
where \(\nu_{sgs}\) is the subgrid scale eddy viscosity and the resolved-scale strain rate tensor \(\overline{D}_{ij}\) is defined as
\begin{equation}
\overline{D}_{ij} = \frac{1}{2} \left( \frac{\partial \overline{u}_{i}}{\partial x_{j}} + \frac{\partial \overline{u}_{j}}{\partial x_{i}}\right), \tag{2} \label{eq:D}
\end{equation}
and the subgrid scale kinetic energy \(k_{sgs}\) is
\begin{equation}
k_{sgs} = \frac{1}{2} \tau_{kk} = \frac{1}{2} \left( \overline{u_{k}u_{k}} – \overline{u}_{k}\overline{u}_{k} \right). \tag{3} \label{eq:ksgs}
\end{equation}
The subgrid scale stress tensor \(\tau_{ij}\) is split into an isotropic part \(\frac{1}{3} \tau_{kk} \delta_{ij}\) and an anisotropic part \(\tau_{ij} \; – \frac{1}{3} \tau_{kk} \delta_{ij}\) in Eq. \eqref{eq:smg2}. I think \(dev(\overline{D})\) is used instead of \(\overline{D}\) because the anisotropic part is a traceless tensor in Eq. \eqref{eq:smg3}.
Assumption 1: Eddy viscosity approximation
In analogy with the molecular viscous stress in laminar flows, the anisotropic part is approximated by relating it to the resolved rate of strain tensor \(\overline{D}_{ij}\) as in Eq. \eqref{eq:smg3}
\begin{equation}
\tau_{ij} \; – \frac{1}{3} \tau_{kk} \delta_{ij} \approx \; – 2 \nu_{sgs} dev(\overline{D})_{ij}. \tag{1c}
\end{equation}
In OpenFOAM, the subgrid scale viscosity is computed as
\begin{equation}
\nu_{sgs} = C_{k} \Delta \sqrt{k_{sgs}}, \tag{4} \label{eq:nusgs}
\end{equation}
where \(C_{k}\) is a model constant whose default value is \(0.094\) and \(\Delta\) is the grid size that defines the subgrid length scale. The method for calculating \(\Delta\) is specified in the turbulenceProperties file and the available options are as follows:
- cubeRootVol
- maxDeltaxyz
- maxDeltaxyzCubeRoot
- smooth
- vanDriest
- Prandtl
- IDDESDelta
For the Smagorinsky SGS model, the vanDriest option is the first choice.
1 2 3 4 5 6 |
template<class BasicTurbulenceModel> void Smagorinsky<BasicTurbulenceModel>::correct() { LESeddyViscosity<BasicTurbulenceModel>::correct(); correctNut(); } |
1 2 3 4 5 6 7 8 9 10 11 |
template<class BasicTurbulenceModel> void Smagorinsky<BasicTurbulenceModel>::correctNut() { volScalarField k(this->k(fvc::grad(this->U_))); this->nut_ = Ck_*this->delta()*sqrt(k); this->nut_.correctBoundaryConditions(); fv::options::New(this->mesh_).correct(this->nut_); BasicTurbulenceModel::correctNut(); } |
Our remaining task is to evaluate the distribution of subgrid scale kinetic energy \(k_{sgs}\).
Assumption 2: Local equilibrium
The SGS kinetic energy \(k_{sgs}\) is computed with the assumption of the balance between the subgrid scale energy production and dissipation (local equilibrium)
\begin{equation}
\overline{D} : \tau_{ij} + C_{\epsilon} \frac{k_{sgs}^{1.5}}{\Delta} = 0, \tag{5} \label{eq:equilibrium}
\end{equation}
where the operator : is a double inner product of two second-rank tensors that can be evaluated as the sum of the 9 products of the tensor components. We can compute \(k_{sgs}\) by solving Eq. \eqref{eq:equilibrium} as shown below:
\begin{align}
&\overline{D} : \left( \frac{2}{3} k_{sgs} I -2 \nu_{sgs} dev(\overline{D}) \right) + C_{\epsilon} \frac{k_{sgs}^{1.5}}{\Delta} = 0 \\
\Leftrightarrow &\quad \overline{D} : \left( \frac{2}{3} k_{sgs} I -2 C_{k} \Delta \sqrt{k_{sgs}} dev(\overline{D}) \right) + C_{\epsilon} \frac{k_{sgs}^{1.5}}{\Delta} = 0 \\
\Leftrightarrow &\quad \sqrt{k_{sgs}} \left( \frac{C_{\epsilon}}{\Delta} k_{sgs} + \frac{2}{3} tr(\overline{D}) \sqrt{k_{sgs}} -2 C_{k} \Delta \left( dev(\overline{D}) : \overline{D} \right) \right) = 0 \\
\Leftrightarrow &\quad a k_{sgs} + b \sqrt{k_{sgs}} – c = 0 \\
\Rightarrow &\quad k_{sgs} = \left( \frac{-b + \sqrt{b^2 + 4ac}}{2a} \right)^2, \tag{6} \label{eq:ksgs2}
\end{align}
where
\begin{eqnarray}
\left\{
\begin{array}{l}
a = \frac{C_{\epsilon}}{\Delta} \\
b = \frac{2}{3} tr(\overline{D}) \tag{7} \label{eq:coeffs}\\
c = 2 C_{k} \Delta \left( dev(\overline{D}) : \overline{D} \right).
\end{array}
\right.
\end{eqnarray}
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 |
template<class BasicTurbulenceModel> tmp<volScalarField> Smagorinsky<BasicTurbulenceModel>::k ( const tmp<volTensorField>& gradU ) const { volSymmTensorField D(symm(gradU)); volScalarField a(this->Ce_/this->delta()); volScalarField b((2.0/3.0)*tr(D)); volScalarField c(2*Ck_*this->delta()*(dev(D) && D)); return tmp<volScalarField> ( new volScalarField ( IOobject ( IOobject::groupName("k", this->U_.group()), this->runTime_.timeName(), this->mesh_ ), sqr((-b + sqrt(sqr(b) + 4*a*c))/(2*a)) ) ); } |
In the case of incompressible flows, Eq. \eqref{eq:coeffs} reduces to
\begin{eqnarray}
\left\{
\begin{array}{l}
b = \frac{2}{3} tr(\overline{D}) = 0 \\
c = 2 C_{k} \Delta \left( dev(\overline{D}) : \overline{D} \right) = C_{k} \Delta \vert \overline{D} \vert^2, \tag{8} \label{eq:coeffsInc}
\end{array}
\right.
\end{eqnarray}
where
\begin{equation}
\vert \overline{D} \vert = \sqrt{2 \overline{D} : \overline{D}}. \tag{9} \label{eq:abs}
\end{equation}
By substituting these relations into Eq. \eqref{eq:ksgs2}, we can get
\begin{equation}
k_{sgs} = \frac{c}{a} = \frac{C_{k} \Delta^2 \vert \overline{D} \vert^2}{C_{\epsilon}}. \tag{10} \label{eq:ksgsInc}
\end{equation}
We can get the following expression for the SGS eddy viscosity in the case of incompressible flows by substituting the above equation into Eq. \eqref{eq:nusgs}
\begin{equation}
\nu_{sgs} = C_{k} \sqrt{\frac{C_{k}}{C_{\epsilon}}} \Delta^2 \vert \overline{D} \vert. \tag{11} \label{eq:nusgsInc}
\end{equation}
By comparing it with the following common expression in the literature
\begin{equation}
\nu_{sgs} = \left( C_{s} \Delta \right)^2 \vert \overline{D} \vert, \tag{12} \label{eq:nusgstext}
\end{equation}
we can get the following relation for the Smagorinsky constant \(C_{s}\):
\begin{equation}
C_{s}^2 = C_{k} \sqrt{\frac{C_{k}}{C_{\epsilon}}}. \tag{13} \label{eq:smgConstInc}
\end{equation}
Features
- Algebraic eddy viscosity model (0-equation model)
- The rotation rate is not taken into account in the calculation of \(\nu_{sgs}\)
- Based on the local equilibrium assumption
- Appropriate value of the Smagorinsky constant is problem dependent
- Not to be able to handle transition
- Not to be able to deal with energy backscatter
- Van-Driest damping is necessary for the SGS eddy viscosity in the near-wall region
The last three items arise from the fact that \(\nu_{sgs} \ge 0\) in Eq. \eqref{eq:nusgs}.
References
[1] (Free Access) J. Smagorinsky, General circulation experiments with the primitive equations: I. The basic experiment*. Monthly weather review, 91(3), 99-164, 1963.