WALE.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2015 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 
28 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
29 
30 namespace Foam
31 {
32 namespace LESModels
33 {
34 
35 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
36 
37 template<class BasicTurbulenceModel>
38 tmp<volSymmTensorField> WALE<BasicTurbulenceModel>::Sd
39 (
40  const volTensorField& gradU
41 ) const
42 {
43  return dev(symm(gradU & gradU));
44 }
45 
46 
47 template<class BasicTurbulenceModel>
49 (
50  const volTensorField& gradU
51 ) const
52 {
53  volScalarField magSqrSd(magSqr(Sd(gradU)));
54 
55  return tmp<volScalarField>
56  (
57  new volScalarField
58  (
59  IOobject
60  (
61  IOobject::groupName("k", this->U_.group()),
62  this->runTime_.timeName(),
63  this->mesh_
64  ),
65  sqr(sqr(Cw_)*this->delta()/Ck_)*
66  (
67  pow3(magSqrSd)
68  /(
69  sqr
70  (
71  pow(magSqr(symm(gradU)), 5.0/2.0)
72  + pow(magSqrSd, 5.0/4.0)
73  )
75  (
76  "SMALL",
77  dimensionSet(0, 0, -10, 0, 0),
78  SMALL
79  )
80  )
81  )
82  )
83  );
84 }
85 
86 
87 template<class BasicTurbulenceModel>
89 {
90  this->nut_ = Ck_*this->delta()*sqrt(this->k(fvc::grad(this->U_)));
91  this->nut_.correctBoundaryConditions();
92 
93  BasicTurbulenceModel::correctNut();
94 }
95 
96 
97 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
98 
99 template<class BasicTurbulenceModel>
101 (
102  const alphaField& alpha,
103  const rhoField& rho,
104  const volVectorField& U,
105  const surfaceScalarField& alphaRhoPhi,
106  const surfaceScalarField& phi,
107  const transportModel& transport,
108  const word& propertiesName,
109  const word& type
110 )
111 :
113  (
114  type,
115  alpha,
116  rho,
117  U,
118  alphaRhoPhi,
119  phi,
120  transport,
121  propertiesName
122  ),
123 
124  Ck_
125  (
127  (
128  "Ck",
129  this->coeffDict_,
130  0.094
131  )
132  ),
133 
134  Cw_
135  (
137  (
138  "Cw",
139  this->coeffDict_,
140  0.325
141  )
142  )
143 {
144  if (type == typeName)
145  {
146  this->printCoeffs(type);
147  }
148 }
149 
150 
151 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
152 
153 template<class BasicTurbulenceModel>
155 {
157  {
158  Ck_.readIfPresent(this->coeffDict());
159  Cw_.readIfPresent(this->coeffDict());
160 
161  return true;
162  }
163  else
164  {
165  return false;
166  }
167 }
168 
169 
170 template<class BasicTurbulenceModel>
172 {
173  volScalarField k(this->k(fvc::grad(this->U_)));
174 
175  return tmp<volScalarField>
176  (
177  new volScalarField
178  (
179  IOobject
180  (
181  IOobject::groupName("epsilon", this->U_.group()),
182  this->runTime_.timeName(),
183  this->mesh_,
186  ),
187  this->Ce_*k*sqrt(k)/this->delta()
188  )
189  );
190 }
191 
192 
193 template<class BasicTurbulenceModel>
195 {
197  correctNut();
198 }
199 
200 
201 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
202 
203 } // End namespace LESModels
204 } // End namespace Foam
205 
206 // ************************************************************************* //
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::IOobject::groupName
static word groupName(Name name, const word &group)
Foam::symm
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:82
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
Foam::LESModels::WALE::WALE
WALE(const WALE &)
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::LESModels::WALE::epsilon
virtual tmp< volScalarField > epsilon() const
Return sub-grid disipation rate.
Definition: WALE.C:171
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::LESModels::WALE::correctNut
virtual void correctNut()
Update the SGS eddy-viscosity.
Definition: WALE.C:88
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:216
Foam::dimensionSet
Dimension set for the base types.
Definition: dimensionSet.H:116
Foam::LESModels::LESeddyViscosity::transportModel
BasicTurbulenceModel::transportModel transportModel
Definition: LESeddyViscosity.H:77
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
U
U
Definition: pEqn.H:46
Foam::LESModels::WALE::Sd
tmp< volSymmTensorField > Sd(const volTensorField &gradU) const
Return the deviatoric symmetric part of the square of the given.
Definition: WALE.C:39
Foam::pow3
dimensionedScalar pow3(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:87
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
delta
scalar delta
Definition: LISASMDCalcMethod2.H:8
Foam::LESModels::WALE::k
virtual tmp< volScalarField > k() const
Return SGS kinetic energy.
Definition: WALE.H:147
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:73
Foam::LESModels::LESeddyViscosity
Eddy viscosity LES SGS model base class.
Definition: LESeddyViscosity.H:55
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::LESModels::LESeddyViscosity::rhoField
BasicTurbulenceModel::rhoField rhoField
Definition: LESeddyViscosity.H:76
rho
rho
Definition: pEqn.H:3
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:49
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:142
Foam::LESModels::WALE::read
virtual bool read()
Read model coefficients if they have changed.
Definition: WALE.C:154
k
label k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
Foam::LESModels::WALE::correct
virtual void correct()
Correct Eddy-Viscosity and related properties.
Definition: WALE.C:194
WALE.H
Foam::eddyViscosity< LESModel< BasicTurbulenceModel > >::correct
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
Definition: eddyViscosity.C:132
Foam::type
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:588
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::magSqr
dimensioned< scalar > magSqr(const dimensioned< Type > &)
Foam::LESModels::LESeddyViscosity::alphaField
BasicTurbulenceModel::alphaField alphaField
Definition: LESeddyViscosity.H:75
Foam::dev
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:104