DeardorffDiffStress.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) 2011-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 "DeardorffDiffStress.H"
27 
28 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
29 
30 namespace Foam
31 {
32 namespace LESModels
33 {
34 
35 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
36 
37 template<class BasicTurbulenceModel>
39 {
40  this->nut_ = Ck_*sqrt(this->k())*this->delta();
41  this->nut_.correctBoundaryConditions();
42 
43  BasicTurbulenceModel::correctNut();
44 }
45 
46 
47 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48 
49 template<class BasicTurbulenceModel>
51 (
52  const alphaField& alpha,
53  const rhoField& rho,
54  const volVectorField& U,
55  const surfaceScalarField& alphaRhoPhi,
56  const surfaceScalarField& phi,
57  const transportModel& transport,
58  const word& propertiesName,
59  const word& type
60 )
61 :
63  (
64  type,
65  alpha,
66  rho,
67  U,
68  alphaRhoPhi,
69  phi,
70  transport,
71  propertiesName
72  ),
73 
74  Ck_
75  (
77  (
78  "Ck",
79  this->coeffDict_,
80  0.094
81  )
82  ),
83  Cm_
84  (
86  (
87  "Cm",
88  this->coeffDict_,
89  4.13
90  )
91  ),
92  Ce_
93  (
95  (
96  "Ce",
97  this->coeffDict_,
98  1.05
99  )
100  ),
101  Cs_
102  (
104  (
105  "Cs",
106  this->coeffDict_,
107  0.25
108  )
109  )
110 {
111  if (type == typeName)
112  {
113  this->printCoeffs(type);
114  this->boundNormalStress(this->R_);
115  }
116 }
117 
118 
119 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
120 
121 template<class BasicTurbulenceModel>
123 {
125  {
126  Ck_.readIfPresent(this->coeffDict());
127  Cm_.readIfPresent(this->coeffDict());
128  Ce_.readIfPresent(this->coeffDict());
129  Cs_.readIfPresent(this->coeffDict());
130 
131  return true;
132  }
133  else
134  {
135  return false;
136  }
137 }
138 
139 
140 template<class BasicTurbulenceModel>
142 {
143  volScalarField k(this->k());
144 
145  return tmp<volScalarField>
146  (
147  new volScalarField
148  (
149  IOobject
150  (
151  IOobject::groupName("epsilon", this->U_.group()),
152  this->runTime_.timeName(),
153  this->mesh_,
156  ),
157  this->Ce_*k*sqrt(k)/this->delta()
158  )
159  );
160 }
161 
162 
163 template<class BasicTurbulenceModel>
165 {
166  if (!this->turbulence_)
167  {
168  return;
169  }
170 
171  // Local references
172  const alphaField& alpha = this->alpha_;
173  const rhoField& rho = this->rho_;
174  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
175  const volVectorField& U = this->U_;
176  volSymmTensorField& R = this->R_;
177 
179 
181  const volTensorField& gradU = tgradU();
182 
183  volSymmTensorField D(symm(gradU));
184 
185  volSymmTensorField P(-twoSymm(R & gradU));
186 
187  volScalarField k(this->k());
188 
190  (
191  fvm::ddt(alpha, rho, R)
192  + fvm::div(alphaRhoPhi, R)
193  - fvm::laplacian(I*this->nu() + Cs_*(k/this->epsilon())*R, R)
194  + fvm::Sp(Cm_*alpha*rho*sqrt(k)/this->delta(), R)
195  ==
196  alpha*rho*P
197  + (4.0/5.0)*alpha*rho*k*D
198  - ((2.0/3.0)*(1.0 - Cm_/this->Ce_)*I)*(alpha*rho*this->epsilon())
199  );
200 
201  REqn().relax();
202  REqn().solve();
203 
204  this->boundNormalStress(R);
205  correctNut();
206 }
207 
208 
209 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
210 
211 } // End namespace LESModels
212 } // End namespace Foam
213 
214 // ************************************************************************* //
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::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::LESModels::DeardorffDiffStress::correctNut
virtual void correctNut()
Update the eddy-viscosity.
Definition: DeardorffDiffStress.C:38
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
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::read
bool read(const char *, int32_t &)
Definition: int32IO.C:87
Foam::LESModels::DeardorffDiffStress::DeardorffDiffStress
DeardorffDiffStress(const DeardorffDiffStress &)
Foam::LESModels::DeardorffDiffStress::read
virtual bool read()
Read model coefficients if they have changed.
Definition: DeardorffDiffStress.C:122
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
U
U
Definition: pEqn.H:46
Foam::fvm::Sp
tmp< fvMatrix< Type > > Sp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
nu
volScalarField & nu
Definition: readMechanicalProperties.H:179
DeardorffDiffStress.H
R
#define R(A, B, C, D, E, F, K, M)
Foam::LESModels::DeardorffDiffStress::epsilon
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
Definition: DeardorffDiffStress.C:141
Foam::ReynoldsStress
Reynolds-stress turbulence model base class.
Definition: ReynoldsStress.H:50
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
Foam::fvm::laplacian
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:46
correct
fvOptions correct(rho)
delta
scalar delta
Definition: LISASMDCalcMethod2.H:8
Foam::I
static const sphericalTensor I(1)
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
rho
rho
Definition: pEqn.H:3
Foam::LESModels::DeardorffDiffStress::alphaField
BasicTurbulenceModel::alphaField alphaField
Definition: DeardorffDiffStress.H:113
Foam::LESModels::DeardorffDiffStress::transportModel
BasicTurbulenceModel::transportModel transportModel
Definition: DeardorffDiffStress.H:115
epsilon
epsilon
Definition: createTDFields.H:56
Foam::fvm::ddt
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
Foam::LESModels::DeardorffDiffStress::rhoField
BasicTurbulenceModel::rhoField rhoField
Definition: DeardorffDiffStress.H:114
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:142
k
label k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
Foam::LESModel
Templated abstract base class for LES SGS models.
Definition: LESModel.H:56
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::fvm::div
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:46
Foam::twoSymm
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:93
Foam::LESModels::DeardorffDiffStress::correct
virtual void correct()
Correct sub-grid stress, eddy-Viscosity and related properties.
Definition: DeardorffDiffStress.C:164