ReynoldsStress.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) 2015-2017 OpenFOAM Foundation
9 Copyright (C) 2020-2022 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 "ReynoldsStress.H"
30#include "fvc.H"
31#include "fvm.H"
32#include "wallFvPatch.H"
33
34// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
35
36template<class BasicTurbulenceModel>
38(
40) const
41{
42 scalar kMin = this->kMin_.value();
43
44 R.max
45 (
47 (
48 "zero",
49 R.dimensions(),
51 (
52 kMin, -GREAT, -GREAT,
53 kMin, -GREAT,
54 kMin
55 )
56 )
57 );
58}
59
60
61template<class BasicTurbulenceModel>
63(
65) const
66{
67 const fvPatchList& patches = this->mesh_.boundary();
68
69 volSymmTensorField::Boundary& RBf = R.boundaryFieldRef();
70
71 forAll(patches, patchi)
72 {
73 const fvPatch& curPatch = patches[patchi];
74
75 if (isA<wallFvPatch>(curPatch))
76 {
77 symmTensorField& Rw = RBf[patchi];
78
79 const scalarField& nutw = this->nut_.boundaryField()[patchi];
80
81 const vectorField snGradU
82 (
83 this->U_.boundaryField()[patchi].snGrad()
84 );
86 const vectorField& faceAreas
87 = this->mesh_.Sf().boundaryField()[patchi];
88
89 const scalarField& magFaceAreas
90 = this->mesh_.magSf().boundaryField()[patchi];
91
92 forAll(curPatch, facei)
93 {
94 // Calculate near-wall velocity gradient
95 const tensor gradUw
96 = (faceAreas[facei]/magFaceAreas[facei])*snGradU[facei];
97
98 // Set the wall Reynolds-stress to the near-wall shear-stress
99 // Note: the spherical part of the normal stress is included in
100 // the pressure
101 Rw[facei] = -nutw[facei]*2*dev(symm(gradUw));
102 }
103 }
104 }
105}
106
107
108template<class BasicTurbulenceModel>
110(
111 const volSymmTensorField& R
112) const
113{
114 const label maxDiffs = 5;
115 label nDiffs = 0;
116
117 // (S:Eq. 4a-4c)
118 forAll(R, celli)
119 {
120 bool diff = false;
121
122 if (maxDiffs < nDiffs)
123 {
124 Info<< "More than " << maxDiffs << " times"
125 << " Reynolds-stress realizability checks failed."
126 << " Skipping further comparisons." << endl;
127 return;
128 }
130 const symmTensor& r = R[celli];
131
132 if (r.xx() < 0)
133 {
135 << "Reynolds stress " << r << " at cell " << celli
136 << " does not obey the constraint: Rxx >= 0"
137 << endl;
138 diff = true;
139 }
140
141 if (r.xx()*r.yy() - sqr(r.xy()) < 0)
142 {
144 << "Reynolds stress " << r << " at cell " << celli
145 << " does not obey the constraint: Rxx*Ryy - sqr(Rxy) >= 0"
146 << endl;
147 diff = true;
148 }
149
150 if (det(r) < 0)
151 {
153 << "Reynolds stress " << r << " at cell " << celli
154 << " does not obey the constraint: det(R) >= 0"
155 << endl;
156 diff = true;
157 }
158
159 if (diff)
160 {
161 ++nDiffs;
163 }
164}
165
166
167// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
168
169template<class BasicTurbulenceModel>
171(
172 const word& modelName,
174 const rhoField& rho,
175 const volVectorField& U,
176 const surfaceScalarField& alphaRhoPhi,
177 const surfaceScalarField& phi,
178 const transportModel& transport,
179 const word& propertiesName
180)
181:
182 BasicTurbulenceModel
183 (
184 modelName,
185 alpha,
186 rho,
187 U,
188 alphaRhoPhi,
189 phi,
190 transport,
191 propertiesName
192 ),
193
194 couplingFactor_
195 (
196 dimensioned<scalar>::getOrAddToDict
197 (
198 "couplingFactor",
199 this->coeffDict_,
200 0.0
201 )
202 ),
203
204 R_
205 (
207 (
208 IOobject::groupName("R", alphaRhoPhi.group()),
209 this->runTime_.timeName(),
210 this->mesh_,
211 IOobject::MUST_READ,
212 IOobject::AUTO_WRITE
213 ),
214 this->mesh_
215 ),
216
217 nut_
218 (
220 (
221 IOobject::groupName("nut", alphaRhoPhi.group()),
222 this->runTime_.timeName(),
223 this->mesh_,
224 IOobject::MUST_READ,
225 IOobject::AUTO_WRITE
226 ),
227 this->mesh_
228 )
229{
230 if (couplingFactor_.value() < 0.0 || couplingFactor_.value() > 1.0)
231 {
233 << "couplingFactor = " << couplingFactor_
234 << " is not in range 0 - 1" << nl
235 << exit(FatalError);
236 }
237}
238
239
240// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
241
242template<class BasicTurbulenceModel>
244{
246}
247
248
249template<class BasicTurbulenceModel>
252{
253 return R_;
254}
255
256
257template<class BasicTurbulenceModel>
260{
261 tmp<Foam::volScalarField> tk(0.5*tr(R_));
262 tk.ref().rename("k");
263 return tk;
264}
265
266
267template<class BasicTurbulenceModel>
270{
271 return devRhoReff(this->U_);
272}
273
274
275template<class BasicTurbulenceModel>
278(
279 const volVectorField& U
280) const
281{
283 (
285 (
287 (
288 IOobject::groupName("devRhoReff", this->alphaRhoPhi_.group()),
289 this->runTime_.timeName(),
290 this->mesh_,
291 IOobject::NO_READ,
292 IOobject::NO_WRITE
293 ),
294 this->alpha_*this->rho_*R_
295 - (this->alpha_*this->rho_*this->nu())
296 *dev(twoSymm(fvc::grad(U)))
297 )
298 );
299}
300
301
302template<class BasicTurbulenceModel>
303template<class RhoFieldType>
306(
307 const RhoFieldType& rho,
309) const
310{
311 if (couplingFactor_.value() > 0.0)
312 {
313 return
314 (
315 fvc::laplacian
316 (
317 (1.0 - couplingFactor_)*this->alpha_*rho*this->nut(),
318 U,
319 "laplacian(nuEff,U)"
320 )
321 + fvc::div
322 (
323 this->alpha_*rho*R_
324 + couplingFactor_
325 *this->alpha_*rho*this->nut()*fvc::grad(U),
326 "div(devRhoReff)"
327 )
328 - fvc::div(this->alpha_*rho*this->nu()*dev2(T(fvc::grad(U))))
329 - fvm::laplacian(this->alpha_*rho*this->nuEff(), U)
330 );
331 }
332 else
333 {
334 return
335 (
336 fvc::laplacian
337 (
338 this->alpha_*rho*this->nut(),
339 U,
340 "laplacian(nuEff,U)"
341 )
342 + fvc::div(this->alpha_*rho*R_)
343 - fvc::div(this->alpha_*rho*this->nu()*dev2(T(fvc::grad(U))))
344 - fvm::laplacian(this->alpha_*rho*this->nuEff(), U)
345 );
346 }
347}
348
349
350template<class BasicTurbulenceModel>
353(
355) const
356{
357 return DivDevRhoReff(this->rho_, U);
358}
359
360
361template<class BasicTurbulenceModel>
364(
365 const volScalarField& rho,
367) const
368{
369 return DivDevRhoReff(rho, U);
370}
371
372
373template<class BasicTurbulenceModel>
375{
376 correctNut();
377}
378
379
380template<class BasicTurbulenceModel>
382{
384}
385
386
387// ************************************************************************* //
#define R(A, B, C, D, E, F, K, M)
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
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
virtual bool read()
Re-read model coefficients if they have changed.
Reynolds-stress turbulence model base class.
virtual tmp< volSymmTensorField > R() const
Return the Reynolds stress tensor.
void boundNormalStress(volSymmTensorField &R) const
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
virtual void validate()
Validate the turbulence fields after construction.
virtual bool read()=0
Re-read model coefficients if they have changed.
tmp< fvVectorMatrix > DivDevRhoReff(const RhoFieldType &rho, volVectorField &U) const
Return the source term for the momentum equation.
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
BasicTurbulenceModel::transportModel transportModel
virtual tmp< volSymmTensorField > devRhoReff() const
Return the effective stress tensor.
void correctWallShearStress(volSymmTensorField &R) const
void checkRealizabilityConditions(const volSymmTensorField &R) const
const Cmpt & xx() const
Definition: SymmTensorI.H:97
const Cmpt & xy() const
Definition: SymmTensorI.H:103
const Cmpt & yy() const
Definition: SymmTensorI.H:121
Generic dimensioned Type class.
tmp< volVectorField > divDevRhoReff()
Return the effective viscous stress (laminar + turbulent).
Definition: momentumError.C:52
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:71
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
const volScalarField & T
const polyBoundaryMesh & patches
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
scalar nut
word timeName
Definition: getTimeIndex.H:3
#define WarningInFunction
Report a warning using Foam::Warning.
dimensionedSymmTensor dev2(const dimensionedSymmTensor &dt)
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
dimensionedScalar det(const dimensionedSphericalTensor &dt)
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
scalar diff(const triad &A, const triad &B)
Return a quantity of the difference between two triads.
Definition: triad.C:378
error FatalError
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
volScalarField & nu
volScalarField & alpha
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333