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-2023 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 "ReynoldsStress.H"
30 #include "fvc.H"
31 #include "fvm.H"
32 #include "wallFvPatch.H"
33 
34 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
35 
36 template<class BasicTurbulenceModel>
38 (
40 ) const
41 {
42  const scalar kMin = this->kMin_.value();
43 
44  R.clamp_min
45  (
47  (
48  kMin, -GREAT, -GREAT,
49  kMin, -GREAT,
50  kMin
51  )
52  );
53 }
54 
55 
56 template<class BasicTurbulenceModel>
58 (
60 ) const
61 {
62  const fvPatchList& patches = this->mesh_.boundary();
63 
64  volSymmTensorField::Boundary& RBf = R.boundaryFieldRef();
65 
66  forAll(patches, patchi)
67  {
68  const fvPatch& curPatch = patches[patchi];
69 
70  if (isA<wallFvPatch>(curPatch))
71  {
72  symmTensorField& Rw = RBf[patchi];
73 
74  const scalarField& nutw = this->nut_.boundaryField()[patchi];
75 
76  const vectorField snGradU
77  (
78  this->U_.boundaryField()[patchi].snGrad()
79  );
80 
81  const vectorField& faceAreas
82  = this->mesh_.Sf().boundaryField()[patchi];
83 
84  const scalarField& magFaceAreas
85  = this->mesh_.magSf().boundaryField()[patchi];
86 
87  forAll(curPatch, facei)
88  {
89  // Calculate near-wall velocity gradient
90  const tensor gradUw
91  = (faceAreas[facei]/magFaceAreas[facei])*snGradU[facei];
92 
93  // Set the wall Reynolds-stress to the near-wall shear-stress
94  // Note: the spherical part of the normal stress is included in
95  // the pressure
96  Rw[facei] = -nutw[facei]*2*devSymm(gradUw);
97  }
98  }
99  }
100 }
101 
102 
103 template<class BasicTurbulenceModel>
105 (
106  const volSymmTensorField& R
107 ) const
108 {
109  const label maxDiffs = 5;
110  label nDiffs = 0;
111 
112  // (S:Eq. 4a-4c)
113  forAll(R, celli)
114  {
115  bool diff = false;
116 
117  if (maxDiffs < nDiffs)
118  {
119  Info<< "More than " << maxDiffs << " times"
120  << " Reynolds-stress realizability checks failed."
121  << " Skipping further comparisons." << endl;
122  return;
123  }
124 
125  const symmTensor& r = R[celli];
126 
127  if (r.xx() < 0)
128  {
130  << "Reynolds stress " << r << " at cell " << celli
131  << " does not obey the constraint: Rxx >= 0"
132  << endl;
133  diff = true;
134  }
135 
136  if (r.xx()*r.yy() - sqr(r.xy()) < 0)
137  {
139  << "Reynolds stress " << r << " at cell " << celli
140  << " does not obey the constraint: Rxx*Ryy - sqr(Rxy) >= 0"
141  << endl;
142  diff = true;
143  }
144 
145  if (det(r) < 0)
146  {
148  << "Reynolds stress " << r << " at cell " << celli
149  << " does not obey the constraint: det(R) >= 0"
150  << endl;
151  diff = true;
152  }
153 
154  if (diff)
155  {
156  ++nDiffs;
157  }
158  }
159 }
160 
161 
162 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
163 
164 template<class BasicTurbulenceModel>
166 (
167  const word& modelName,
168  const alphaField& alpha,
169  const rhoField& rho,
170  const volVectorField& U,
171  const surfaceScalarField& alphaRhoPhi,
172  const surfaceScalarField& phi,
173  const transportModel& transport,
174  const word& propertiesName
175 )
176 :
177  BasicTurbulenceModel
178  (
179  modelName,
180  alpha,
181  rho,
182  U,
183  alphaRhoPhi,
184  phi,
185  transport,
186  propertiesName
187  ),
188 
189  couplingFactor_
190  (
191  dimensioned<scalar>::getOrAddToDict
192  (
193  "couplingFactor",
194  this->coeffDict_,
195  0.0
196  )
197  ),
198 
199  R_
200  (
201  IOobject
202  (
203  IOobject::groupName("R", alphaRhoPhi.group()),
204  this->runTime_.timeName(),
205  this->mesh_,
206  IOobject::MUST_READ,
207  IOobject::AUTO_WRITE
208  ),
209  this->mesh_
210  ),
211 
212  nut_
213  (
214  IOobject
215  (
216  IOobject::groupName("nut", alphaRhoPhi.group()),
217  this->runTime_.timeName(),
218  this->mesh_,
219  IOobject::MUST_READ,
220  IOobject::AUTO_WRITE
221  ),
222  this->mesh_
223  )
224 {
225  if (couplingFactor_.value() < 0.0 || couplingFactor_.value() > 1.0)
226  {
228  << "couplingFactor = " << couplingFactor_
229  << " is not in range 0 - 1" << nl
230  << exit(FatalError);
231  }
232 }
233 
234 
235 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
236 
237 template<class BasicTurbulenceModel>
239 {
241 }
242 
243 
244 template<class BasicTurbulenceModel>
247 {
248  return R_;
249 }
250 
251 
252 template<class BasicTurbulenceModel>
255 {
256  tmp<Foam::volScalarField> tk(0.5*tr(R_));
257  tk.ref().rename("k");
258  return tk;
259 }
260 
261 
262 template<class BasicTurbulenceModel>
265 {
266  return devRhoReff(this->U_);
267 }
268 
269 
270 template<class BasicTurbulenceModel>
273 (
274  const volVectorField& U
275 ) const
276 {
278  (
280  (
281  IOobject
282  (
283  IOobject::groupName("devRhoReff", this->alphaRhoPhi_.group()),
284  this->runTime_.timeName(),
285  this->mesh_,
286  IOobject::NO_READ,
287  IOobject::NO_WRITE
288  ),
289  this->alpha_*this->rho_*R_
290  - (this->alpha_*this->rho_*this->nu())
291  *devTwoSymm(fvc::grad(U))
292  )
293  );
294 }
295 
296 
297 template<class BasicTurbulenceModel>
298 template<class RhoFieldType>
301 (
302  const RhoFieldType& rho,
304 ) const
305 {
306  if (couplingFactor_.value() > 0.0)
307  {
308  return
309  (
311  (
312  (1.0 - couplingFactor_)*this->alpha_*rho*this->nut(),
313  U,
314  "laplacian(nuEff,U)"
315  )
316  + fvc::div
317  (
318  this->alpha_*rho*R_
319  + couplingFactor_
320  *this->alpha_*rho*this->nut()*fvc::grad(U),
321  "div(devRhoReff)"
322  )
323  - fvc::div(this->alpha_*rho*this->nu()*dev2(T(fvc::grad(U))))
324  - fvm::laplacian(this->alpha_*rho*this->nuEff(), U)
325  );
326  }
327  else
328  {
329  return
330  (
332  (
333  this->alpha_*rho*this->nut(),
334  U,
335  "laplacian(nuEff,U)"
336  )
337  + fvc::div(this->alpha_*rho*R_)
338  - fvc::div(this->alpha_*rho*this->nu()*dev2(T(fvc::grad(U))))
339  - fvm::laplacian(this->alpha_*rho*this->nuEff(), U)
340  );
341  }
342 }
343 
344 
345 template<class BasicTurbulenceModel>
348 (
350 ) const
351 {
352  return DivDevRhoReff(this->rho_, U);
353 }
354 
355 
356 template<class BasicTurbulenceModel>
359 (
360  const volScalarField& rho,
362 ) const
363 {
364  return DivDevRhoReff(rho, U);
365 }
366 
367 
368 template<class BasicTurbulenceModel>
370 {
371  correctNut();
372 }
373 
374 
375 template<class BasicTurbulenceModel>
377 {
379 }
380 
381 
382 // ************************************************************************* //
void correctWallShearStress(volSymmTensorField &R) const
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
Definition: volFieldsFwd.H:85
scalar diff(const triad &A, const triad &B)
Return a quantity of the difference between two triads.
Definition: triad.C:373
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
tmp< fvVectorMatrix > DivDevRhoReff(const RhoFieldType &rho, volVectorField &U) const
Return the source term for the momentum equation.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
tmp< GeometricField< Type, faPatchField, areaMesh > > div(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Definition: facDiv.C:43
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
virtual tmp< volSymmTensorField > R() const
Return the Reynolds stress tensor.
SymmTensor< Cmpt > devSymm(const SymmTensor< Cmpt > &st)
Return the deviatoric part of the symmetric part of a SymmTensor.
Definition: SymmTensorI.H:481
virtual tmp< fvVectorMatrix > divDevRhoReff(volVectorField &U) const
Return the source term for the momentum equation.
Generic dimensioned Type class.
dimensionedScalar det(const dimensionedSphericalTensor &dt)
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:82
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
constexpr const char *const group
Group name for atomic constants.
void checkRealizabilityConditions(const volSymmTensorField &R) const
word timeName
Definition: getTimeIndex.H:3
virtual bool read()=0
Re-read model coefficients if they have changed.
SymmTensor< scalar > symmTensor
SymmTensor of scalars, i.e. SymmTensor<scalar>.
Definition: symmTensor.H:55
A class for handling words, derived from Foam::string.
Definition: word.H:63
tmp< GeometricField< Type, faPatchField, areaMesh > > laplacian(const GeometricField< Type, faPatchField, areaMesh > &vf, const word &name)
Definition: facLaplacian.C:40
virtual tmp< volSymmTensorField > devRhoReff() const
Return the effective stress tensor.
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
virtual void validate()
Validate the turbulence fields after construction.
void boundNormalStress(volSymmTensorField &R) const
void read(Istream &, label &val, const dictionary &)
In-place read with dictionary lookup.
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(psi *p+alphal *rhol0+((alphav *psiv+alphal *psil) - psi) *pSat, rhoMin);# 1 "/home/chef2/andy/OpenFOAM/release/v2312/OpenFOAM-v2312/applications/solvers/multiphase/cavitatingFoam/alphavPsi.H" 1{ alphav=clamp((rho - rholSat)/(rhovSat - rholSat), zero_one{});alphal=1.0 - alphav;Info<< "max-min alphav: "<< max(alphav).value()<< " "<< min(alphav).value()<< endl;psiModel-> correct()
Definition: pEqn.H:63
const volScalarField & T
ReynoldsStress(const word &modelName, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName)
Construct from components.
dimensionedSymmTensor dev2(const dimensionedSymmTensor &dt)
U
Definition: pEqn.H:72
#define R(A, B, C, D, E, F, K, M)
#define WarningInFunction
Report a warning using Foam::Warning.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
Base-class for all transport models used by the incompressible turbulence models. ...
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
const polyBoundaryMesh & patches
messageStream Info
Information stream (stdout output on master, null elsewhere)
A class for managing temporary objects.
Definition: HashPtrTable.H:50
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Tensor of scalars, i.e. Tensor<scalar>.
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
volScalarField & nu
scalar nut
tmp< GeometricField< typename outerProduct< vector, Type >::type, faPatchField, areaMesh >> grad(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Definition: facGrad.C:51
SymmTensor< Cmpt > devTwoSymm(const SymmTensor< Cmpt > &st)
Return the deviatoric part of twice the symmetric part of a SymmTensor.
Definition: SymmTensorI.H:491