WALE.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2015-2020 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "WALE.H"
27 #include "fvOptions.H"
28 
29 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace LESModels
34 {
35 
36 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
37 
38 template<class BasicMomentumTransportModel>
39 tmp<volSymmTensorField> WALE<BasicMomentumTransportModel>::Sd
40 (
41  const volTensorField& gradU
42 ) const
43 {
44  return dev(symm(gradU & gradU));
45 }
46 
47 
48 template<class BasicMomentumTransportModel>
50 (
51  const volTensorField& gradU
52 ) const
53 {
54  volScalarField magSqrSd(magSqr(Sd(gradU)));
55 
56  return volScalarField::New
57  (
58  IOobject::groupName("k", this->alphaRhoPhi_.group()),
59  sqr(sqr(Cw_)*this->delta()/Ck_)*
60  (
61  pow3(magSqrSd)
62  /(
63  sqr
64  (
65  pow(magSqr(symm(gradU)), 5.0/2.0)
66  + pow(magSqrSd, 5.0/4.0)
67  )
69  (
70  "small",
71  dimensionSet(0, 0, -10, 0, 0),
72  small
73  )
74  )
75  )
76  );
77 }
78 
79 
80 template<class BasicMomentumTransportModel>
82 {
83  this->nut_ = Ck_*this->delta()*sqrt(this->k(fvc::grad(this->U_)));
84  this->nut_.correctBoundaryConditions();
85  fv::options::New(this->mesh_).correct(this->nut_);
86 }
87 
88 
89 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
90 
91 template<class BasicMomentumTransportModel>
93 (
94  const alphaField& alpha,
95  const rhoField& rho,
96  const volVectorField& U,
97  const surfaceScalarField& alphaRhoPhi,
98  const surfaceScalarField& phi,
99  const transportModel& transport,
100  const word& type
101 )
102 :
104  (
105  type,
106  alpha,
107  rho,
108  U,
109  alphaRhoPhi,
110  phi,
111  transport
112  ),
113 
114  Ck_
115  (
117  (
118  "Ck",
119  this->coeffDict_,
120  0.094
121  )
122  ),
123 
124  Cw_
125  (
127  (
128  "Cw",
129  this->coeffDict_,
130  0.325
131  )
132  )
133 {
134  if (type == typeName)
135  {
136  this->printCoeffs(type);
137  }
138 }
139 
140 
141 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
142 
143 template<class BasicMomentumTransportModel>
145 {
147  {
148  Ck_.readIfPresent(this->coeffDict());
149  Cw_.readIfPresent(this->coeffDict());
150 
151  return true;
152  }
153  else
154  {
155  return false;
156  }
157 }
158 
159 
160 template<class BasicMomentumTransportModel>
162 {
163  volScalarField k(this->k(fvc::grad(this->U_)));
164 
165  return volScalarField::New
166  (
167  IOobject::groupName("epsilon", this->alphaRhoPhi_.group()),
168  this->Ce_*k*sqrt(k)/this->delta()
169  );
170 }
171 
172 
173 template<class BasicMomentumTransportModel>
175 {
177  correctNut();
178 }
179 
180 
181 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
182 
183 } // End namespace LESModels
184 } // End namespace Foam
185 
186 // ************************************************************************* //
virtual bool read()
Read model coefficients if they have changed.
Definition: WALE.C:144
Eddy viscosity LES SGS model base class.
scalar delta
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
BasicMomentumTransportModel::rhoField rhoField
void correct(GeometricField< Type, fvPatchField, volMesh > &field)
Apply correction to field.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
virtual tmp< volScalarField > epsilon() const
Return sub-grid disipation rate.
Definition: WALE.C:161
dimensionedScalar sqrt(const dimensionedScalar &ds)
label k
Boltzmann constant.
Generic dimensioned Type class.
WALE(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &type=typeName)
Construct from components.
Definition: WALE.C:93
phi
Definition: pEqn.H:104
virtual void correctNut()
Update the SGS eddy-viscosity.
Definition: WALE.C:81
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
Dimension set for the base types.
Definition: dimensionSet.H:120
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
tmp< volSymmTensorField > Sd(const volTensorField &gradU) const
Return the deviatoric symmetric part of the square of the given.
Definition: WALE.C:40
A class for handling words, derived from string.
Definition: word.H:59
static word groupName(Name name, const word &group)
dimensioned< scalar > magSqr(const dimensioned< Type > &)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar pow3(const dimensionedScalar &ds)
BasicMomentumTransportModel::alphaField alphaField
U
Definition: pEqn.H:72
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
BasicMomentumTransportModel::transportModel transportModel
virtual void correct()
Correct Eddy-Viscosity and related properties.
Definition: WALE.C:174
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
A class for managing temporary objects.
Definition: PtrList.H:53
static options & New(const fvMesh &mesh)
Construct fvOptions and register to datbase if not present.
Definition: fvOptions.C:101
virtual tmp< volScalarField > k() const
Return SGS kinetic energy.
Definition: WALE.H:139
Namespace for OpenFOAM.