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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2019-2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "DeardorffDiffStress.H"
30 #include "fvOptions.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace LESModels
37 {
38 
39 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
40 
41 template<class BasicTurbulenceModel>
43 {
44  this->nut_ = Ck_*sqrt(this->k())*this->delta();
45  this->nut_.correctBoundaryConditions();
46  fv::options::New(this->mesh_).correct(this->nut_);
47 
48  BasicTurbulenceModel::correctNut();
49 }
50 
51 
52 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
53 
54 template<class BasicTurbulenceModel>
56 (
57  const alphaField& alpha,
58  const rhoField& rho,
59  const volVectorField& U,
60  const surfaceScalarField& alphaRhoPhi,
61  const surfaceScalarField& phi,
62  const transportModel& transport,
63  const word& propertiesName,
64  const word& type
65 )
66 :
68  (
69  type,
70  alpha,
71  rho,
72  U,
73  alphaRhoPhi,
74  phi,
75  transport,
76  propertiesName
77  ),
78 
79  Ck_
80  (
82  (
83  "Ck",
84  this->coeffDict_,
85  0.094
86  )
87  ),
88  Cm_
89  (
91  (
92  "Cm",
93  this->coeffDict_,
94  4.13
95  )
96  ),
97  Ce_
98  (
100  (
101  "Ce",
102  this->coeffDict_,
103  1.05
104  )
105  ),
106  Cs_
107  (
109  (
110  "Cs",
111  this->coeffDict_,
112  0.25
113  )
114  )
115 {
116  if (type == typeName)
117  {
118  this->printCoeffs(type);
119  this->boundNormalStress(this->R_);
120  }
121 }
122 
123 
124 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
125 
126 template<class BasicTurbulenceModel>
128 {
130  {
131  Ck_.readIfPresent(this->coeffDict());
132  Cm_.readIfPresent(this->coeffDict());
133  Ce_.readIfPresent(this->coeffDict());
134  Cs_.readIfPresent(this->coeffDict());
135 
136  return true;
137  }
138 
139  return false;
140 }
141 
142 
143 template<class BasicTurbulenceModel>
145 {
146  volScalarField k(this->k());
147 
148  return tmp<volScalarField>
149  (
150  new volScalarField
151  (
152  IOobject
153  (
154  IOobject::groupName("epsilon", this->alphaRhoPhi_.group()),
155  this->runTime_.timeName(),
156  this->mesh_,
159  ),
160  this->Ce_*k*sqrt(k)/this->delta()
161  )
162  );
163 }
164 
165 
166 template<class BasicTurbulenceModel>
168 {
169  volScalarField k(this->k());
170  volScalarField epsilon(this->Ce_*k*sqrt(k)/this->delta());
171 
173  (
174  IOobject
175  (
176  IOobject::groupName("omega", this->alphaRhoPhi_.group()),
177  this->runTime_.timeName(),
178  this->mesh_
179  ),
180  epsilon/(0.09*k)
181  );
182 }
183 
184 
185 template<class BasicTurbulenceModel>
187 {
188  if (!this->turbulence_)
189  {
190  return;
191  }
192 
193  // Local references
194  const alphaField& alpha = this->alpha_;
195  const rhoField& rho = this->rho_;
196  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
197  const volVectorField& U = this->U_;
198  volSymmTensorField& R = this->R_;
199  fv::options& fvOptions(fv::options::New(this->mesh_));
200 
202 
204  const volTensorField& gradU = tgradU();
205 
206  volSymmTensorField D(symm(gradU));
207 
208  volSymmTensorField P(-twoSymm(R & gradU));
209 
210  volScalarField k(this->k());
211 
213  (
214  fvm::ddt(alpha, rho, R)
215  + fvm::div(alphaRhoPhi, R)
216  - fvm::laplacian(I*this->nu() + Cs_*(k/this->epsilon())*R, R)
217  + fvm::Sp(Cm_*alpha*rho*sqrt(k)/this->delta(), R)
218  ==
219  alpha*rho*P
220  + (4.0/5.0)*alpha*rho*k*D
221  - ((2.0/3.0)*(1.0 - Cm_/this->Ce_)*I)*(alpha*rho*this->epsilon())
222  + fvOptions(alpha, rho, R)
223  );
224 
225  REqn.ref().relax();
226  fvOptions.constrain(REqn.ref());
227  REqn.ref().solve();
228  fvOptions.correct(R);
229  this->boundNormalStress(R);
230 
231  correctNut();
232 }
233 
234 
235 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
236 
237 } // End namespace LESModels
238 } // End namespace Foam
239 
240 // ************************************************************************* //
Foam::IOobject::NO_WRITE
Definition: IOobject.H:130
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::LESModels::DeardorffDiffStress
Differential SGS Stress Equation Model for incompressible and compressible flows.
Definition: DeardorffDiffStress.H:85
Foam::symm
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:84
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:54
Foam::LESModels::DeardorffDiffStress::omega
virtual tmp< volScalarField > omega() const
Return the specific dissipation rate.
Definition: DeardorffDiffStress.C:167
fvOptions.H
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::fv::optionList::correct
void correct(GeometricField< Type, fvPatchField, volMesh > &field)
Apply correction to field.
Definition: fvOptionListTemplates.C:323
Foam::LESModels::DeardorffDiffStress::correctNut
virtual void correctNut()
Update the eddy-viscosity.
Definition: DeardorffDiffStress.C:42
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::fv::options::New
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present.
Definition: fvOptions.C:104
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
Foam::read
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:108
Foam::LESModels::DeardorffDiffStress::read
virtual bool read()
Read model coefficients if they have changed.
Definition: DeardorffDiffStress.C:127
rho
rho
Definition: readInitialConditions.H:88
nu
volScalarField & nu
Definition: readMechanicalProperties.H:176
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:144
Foam::ReynoldsStress
Reynolds-stress turbulence model base class.
Definition: ReynoldsStress.H:53
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:228
Foam::fvm::Sp
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
fvOptions
fv::options & fvOptions
Definition: setRegionFluidFields.H:23
Foam::fvm::laplacian
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:48
correct
fvOptions correct(rho)
delta
scalar delta
Definition: LISASMDCalcMethod2.H:8
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::fv::options
Finite-volume options.
Definition: fvOptions.H:55
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:43
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::LESModels::DeardorffDiffStress::alphaField
BasicTurbulenceModel::alphaField alphaField
Definition: DeardorffDiffStress.H:118
U
U
Definition: pEqn.H:72
Foam::LESModels::DeardorffDiffStress::transportModel
BasicTurbulenceModel::transportModel transportModel
Definition: DeardorffDiffStress.H:120
Foam::fvm::ddt
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:48
Foam::LESModels::DeardorffDiffStress::rhoField
BasicTurbulenceModel::rhoField rhoField
Definition: DeardorffDiffStress.H:119
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
k
label k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
Foam::tmp::New
static tmp< T > New(Args &&... args)
Construct tmp of T with forwarding arguments.
D
const dimensionedScalar & D
Definition: solveBulkSurfactant.H:4
Foam::LESModel
Templated abstract base class for LES SGS models.
Definition: LESModel.H:58
epsilon
scalar epsilon
Definition: evaluateNearWall.H:7
Foam::IOobject::groupName
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
Foam::GeometricField< vector, fvPatchField, volMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:123
Foam::fvm::div
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:48
Foam::twoSymm
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:95
Foam::LESModels::DeardorffDiffStress::correct
virtual void correct()
Correct sub-grid stress, eddy-Viscosity and related properties.
Definition: DeardorffDiffStress.C:186
Foam::I
static const Identity< scalar > I
Definition: Identity.H:95