LRR.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-2016 OpenFOAM Foundation
9 Copyright (C) 2019-2020 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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 "LRR.H"
30#include "fvOptions.H"
31#include "wallFvPatch.H"
32
33// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34
35namespace Foam
36{
37namespace RASModels
38{
39
40// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
41
42template<class BasicTurbulenceModel>
44{
45 this->nut_ = this->Cmu_*sqr(k_)/epsilon_;
46 this->nut_.correctBoundaryConditions();
47 fv::options::New(this->mesh_).correct(this->nut_);
48
49 BasicTurbulenceModel::correctNut();
50}
51
52
53// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
54
55template<class BasicTurbulenceModel>
57(
58 const alphaField& alpha,
59 const rhoField& rho,
60 const volVectorField& U,
61 const surfaceScalarField& alphaRhoPhi,
63 const transportModel& transport,
64 const word& propertiesName,
65 const word& type
66)
67:
68 ReynoldsStress<RASModel<BasicTurbulenceModel>>
69 (
70 type,
71 alpha,
72 rho,
73 U,
74 alphaRhoPhi,
75 phi,
76 transport,
77 propertiesName
78 ),
79
80 Cmu_
81 (
82 dimensioned<scalar>::getOrAddToDict
83 (
84 "Cmu",
85 this->coeffDict_,
86 0.09
87 )
88 ),
89 C1_
90 (
91 dimensioned<scalar>::getOrAddToDict
92 (
93 "C1",
94 this->coeffDict_,
95 1.8
96 )
97 ),
98 C2_
99 (
100 dimensioned<scalar>::getOrAddToDict
101 (
102 "C2",
103 this->coeffDict_,
104 0.6
105 )
106 ),
107 Ceps1_
108 (
109 dimensioned<scalar>::getOrAddToDict
110 (
111 "Ceps1",
112 this->coeffDict_,
113 1.44
114 )
115 ),
116 Ceps2_
117 (
118 dimensioned<scalar>::getOrAddToDict
119 (
120 "Ceps2",
121 this->coeffDict_,
122 1.92
123 )
124 ),
125 Cs_
126 (
127 dimensioned<scalar>::getOrAddToDict
128 (
129 "Cs",
130 this->coeffDict_,
131 0.25
132 )
133 ),
134 Ceps_
135 (
136 dimensioned<scalar>::getOrAddToDict
137 (
138 "Ceps",
139 this->coeffDict_,
140 0.15
141 )
142 ),
143
144 wallReflection_
145 (
146 Switch::getOrAddToDict
147 (
148 "wallReflection",
149 this->coeffDict_,
150 true
151 )
152 ),
153 kappa_
154 (
155 dimensioned<scalar>::getOrAddToDict
156 (
157 "kappa",
158 this->coeffDict_,
159 0.41
160 )
161 ),
162 Cref1_
163 (
164 dimensioned<scalar>::getOrAddToDict
165 (
166 "Cref1",
167 this->coeffDict_,
168 0.5
169 )
170 ),
171 Cref2_
172 (
173 dimensioned<scalar>::getOrAddToDict
174 (
175 "Cref2",
176 this->coeffDict_,
177 0.3
178 )
179 ),
180
181 k_
182 (
184 (
185 "k",
186 this->runTime_.timeName(),
187 this->mesh_,
188 IOobject::NO_READ,
189 IOobject::AUTO_WRITE
190 ),
191 0.5*tr(this->R_)
192 ),
193 epsilon_
194 (
196 (
197 "epsilon",
198 this->runTime_.timeName(),
199 this->mesh_,
200 IOobject::MUST_READ,
201 IOobject::AUTO_WRITE
202 ),
203 this->mesh_
204 )
205{
206 if (type == typeName)
207 {
208 this->printCoeffs(type);
209
210 this->boundNormalStress(this->R_);
211 bound(epsilon_, this->epsilonMin_);
212 k_ = 0.5*tr(this->R_);
213 }
214}
215
216
217// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
218
219template<class BasicTurbulenceModel>
221{
223 {
224 Cmu_.readIfPresent(this->coeffDict());
225 C1_.readIfPresent(this->coeffDict());
226 C2_.readIfPresent(this->coeffDict());
227 Ceps1_.readIfPresent(this->coeffDict());
228 Ceps2_.readIfPresent(this->coeffDict());
229 Cs_.readIfPresent(this->coeffDict());
230 Ceps_.readIfPresent(this->coeffDict());
231
232 wallReflection_.readIfPresent("wallReflection", this->coeffDict());
233 kappa_.readIfPresent(this->coeffDict());
234 Cref1_.readIfPresent(this->coeffDict());
235 Cref2_.readIfPresent(this->coeffDict());
236
237 return true;
238 }
239
240 return false;
241}
242
243
244template<class BasicTurbulenceModel>
246{
248 (
250 (
251 "DREff",
252 (Cs_*(this->k_/this->epsilon_))*this->R_ + I*this->nu()
253 )
254 );
255}
256
257
258template<class BasicTurbulenceModel>
260{
262 (
264 (
265 "DepsilonEff",
266 (Ceps_*(this->k_/this->epsilon_))*this->R_ + I*this->nu()
267 )
268 );
269}
270
271
272template<class BasicTurbulenceModel>
274{
275 if (!this->turbulence_)
276 {
277 return;
278 }
279
280 // Local references
281 const alphaField& alpha = this->alpha_;
282 const rhoField& rho = this->rho_;
283 const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
284 const volVectorField& U = this->U_;
285 volSymmTensorField& R = this->R_;
287
289
291 const volTensorField& gradU = tgradU();
292
293 volSymmTensorField P(-twoSymm(R & gradU));
294 volScalarField G(this->GName(), 0.5*mag(tr(P)));
295
296 // Update epsilon and G at the wall
297 epsilon_.boundaryFieldRef().updateCoeffs();
298
299 // Dissipation equation
301 (
302 fvm::ddt(alpha, rho, epsilon_)
303 + fvm::div(alphaRhoPhi, epsilon_)
304 - fvm::laplacian(alpha*rho*DepsilonEff(), epsilon_)
305 ==
306 Ceps1_*alpha*rho*G*epsilon_/k_
307 - fvm::Sp(Ceps2_*alpha*rho*epsilon_/k_, epsilon_)
308 + fvOptions(alpha, rho, epsilon_)
309 );
310
311 epsEqn.ref().relax();
312 fvOptions.constrain(epsEqn.ref());
313 epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
314 solve(epsEqn);
315 fvOptions.correct(epsilon_);
316 bound(epsilon_, this->epsilonMin_);
317
318
319 // Correct the trace of the tensorial production to be consistent
320 // with the near-wall generation from the wall-functions
321 const fvPatchList& patches = this->mesh_.boundary();
322
323 forAll(patches, patchi)
324 {
325 const fvPatch& curPatch = patches[patchi];
326
327 if (isA<wallFvPatch>(curPatch))
328 {
329 forAll(curPatch, facei)
330 {
331 label celli = curPatch.faceCells()[facei];
332 P[celli] *= min
333 (
334 G[celli]/(0.5*mag(tr(P[celli])) + SMALL),
335 1.0
336 );
337 }
338 }
339 }
340
341 // Reynolds stress equation
343 (
345 + fvm::div(alphaRhoPhi, R)
346 - fvm::laplacian(alpha*rho*DREff(), R)
347 + fvm::Sp(C1_*alpha*rho*epsilon_/k_, R)
348 ==
349 alpha*rho*P
350 - (2.0/3.0*(1 - C1_)*I)*alpha*rho*epsilon_
351 - C2_*alpha*rho*dev(P)
352 + fvOptions(alpha, rho, R)
353 );
354
355 // Optionally add wall-refection term
356 if (wallReflection_)
357 {
358 const volVectorField& n_(wallDist::New(this->mesh_).n());
359 const volScalarField& y_(wallDist::New(this->mesh_).y());
360
361 const volSymmTensorField reflect
362 (
363 Cref1_*R - ((Cref2_*C2_)*(k_/epsilon_))*dev(P)
364 );
365
366 REqn.ref() +=
367 ((3*pow(Cmu_, 0.75)/kappa_)*(alpha*rho*sqrt(k_)/y_))
368 *dev(symm((n_ & reflect)*n_));
369 }
370
371 REqn.ref().relax();
372 fvOptions.constrain(REqn.ref());
373 solve(REqn);
374 fvOptions.correct(R);
375
376 this->boundNormalStress(R);
377
378 k_ = 0.5*tr(R);
379
380 correctNut();
381
382 // Correct wall shear-stresses when applying wall-functions
383 this->correctWallShearStress(R);
384}
385
386
387// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
388
389} // End namespace RASModels
390} // End namespace Foam
391
392// ************************************************************************* //
scalar y
#define R(A, B, C, D, E, F, K, M)
label n
fv::options & fvOptions
surfaceScalarField & phi
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
Templated abstract base class for RAS turbulence models.
Definition: RASModel.H:55
virtual void printCoeffs(const word &type)
Print model coefficients.
Definition: RASModel.C:34
dimensionedScalar epsilonMin_
Lower limit of epsilon.
Definition: RASModel.H:77
Launder, Reece and Rodi Reynolds-stress turbulence model for incompressible and compressible flows.
Definition: LRR.H:106
volScalarField epsilon_
Definition: LRR.H:144
BasicTurbulenceModel::alphaField alphaField
Definition: LRR.H:155
tmp< volSymmTensorField > DepsilonEff() const
Return the effective diffusivity for epsilon.
Definition: LRR.C:259
BasicTurbulenceModel::rhoField rhoField
Definition: LRR.H:156
volScalarField k_
Definition: LRR.H:143
virtual void correct()
Solve the turbulence equations and correct eddy-Viscosity and.
Definition: LRR.C:273
tmp< volSymmTensorField > DREff() const
Return the effective diffusivity for R.
Definition: LRR.C:245
virtual void correctNut()
Update the eddy-viscosity.
Definition: LRR.C:43
BasicTurbulenceModel::transportModel transportModel
Definition: LRR.H:157
virtual bool read()
Read model coefficients if they have changed.
Definition: LRR.C:220
Reynolds-stress turbulence model base class.
void boundNormalStress(volSymmTensorField &R) const
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:78
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
Generic dimensioned Type class.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:71
virtual const labelUList & faceCells() const
Return faceCells.
Definition: fvPatch.C:113
Finite-volume options.
Definition: fvOptions.H:59
A class for managing temporary objects.
Definition: tmp.H:65
T & ref() const
Definition: tmpI.H:227
A class for handling words, derived from Foam::string.
Definition: word.H:68
U
Definition: pEqn.H:72
thermo correct()
const polyBoundaryMesh & patches
word timeName
Definition: getTimeIndex.H:3
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:54
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:48
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:48
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:48
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
Namespace for OpenFOAM.
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:108
dimensionedSymmTensor sqr(const dimensionedVector &dv)
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:35
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
static const Identity< scalar > I
Definition: Identity.H:94
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
volScalarField & nu
volScalarField & alpha
CEqn solve()
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333