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-2018 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 BasicTurbulenceModel>
39 tmp<volSymmTensorField> WALE<BasicTurbulenceModel>::Sd
40 (
41  const volTensorField& gradU
42 ) const
43 {
44  return dev(symm(gradU & gradU));
45 }
46 
47 
48 template<class BasicTurbulenceModel>
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 BasicTurbulenceModel>
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  BasicTurbulenceModel::correctNut();
88 }
89 
90 
91 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
92 
93 template<class BasicTurbulenceModel>
95 (
96  const alphaField& alpha,
97  const rhoField& rho,
98  const volVectorField& U,
99  const surfaceScalarField& alphaRhoPhi,
100  const surfaceScalarField& phi,
101  const transportModel& transport,
102  const word& propertiesName,
103  const word& type
104 )
105 :
107  (
108  type,
109  alpha,
110  rho,
111  U,
112  alphaRhoPhi,
113  phi,
114  transport,
115  propertiesName
116  ),
117 
118  Ck_
119  (
121  (
122  "Ck",
123  this->coeffDict_,
124  0.094
125  )
126  ),
127 
128  Cw_
129  (
131  (
132  "Cw",
133  this->coeffDict_,
134  0.325
135  )
136  )
137 {
138  if (type == typeName)
139  {
140  this->printCoeffs(type);
141  }
142 }
143 
144 
145 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
146 
147 template<class BasicTurbulenceModel>
149 {
151  {
152  Ck_.readIfPresent(this->coeffDict());
153  Cw_.readIfPresent(this->coeffDict());
154 
155  return true;
156  }
157  else
158  {
159  return false;
160  }
161 }
162 
163 
164 template<class BasicTurbulenceModel>
166 {
167  volScalarField k(this->k(fvc::grad(this->U_)));
168 
169  return volScalarField::New
170  (
171  IOobject::groupName("epsilon", this->alphaRhoPhi_.group()),
172  this->Ce_*k*sqrt(k)/this->delta()
173  );
174 }
175 
176 
177 template<class BasicTurbulenceModel>
179 {
181  correctNut();
182 }
183 
184 
185 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
186 
187 } // End namespace LESModels
188 } // End namespace Foam
189 
190 // ************************************************************************* //
virtual bool read()
Read model coefficients if they have changed.
Definition: WALE.C:148
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
surfaceScalarField & phi
BasicTurbulenceModel::transportModel transportModel
void correct(GeometricField< Type, fvPatchField, volMesh > &field)
Apply correction to field.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
virtual tmp< volScalarField > epsilon() const
Return sub-grid disipation rate.
Definition: WALE.C:165
dimensionedScalar sqrt(const dimensionedScalar &ds)
label k
Boltzmann constant.
Generic dimensioned Type class.
BasicTurbulenceModel::alphaField alphaField
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Mesh &, const dimensionSet &, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
Return a temporary field constructed from name, mesh, dimensionSet.
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 > &)
WALE(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName=turbulenceModel::propertiesName, const word &type=typeName)
Construct from components.
Definition: WALE.C:95
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar pow3(const dimensionedScalar &ds)
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
BasicTurbulenceModel::rhoField rhoField
virtual void correct()
Correct Eddy-Viscosity and related properties.
Definition: WALE.C:178
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:140
Namespace for OpenFOAM.