All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
LRR.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2023 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 "LRR.H"
27 #include "fvModels.H"
28 #include "fvConstraints.H"
29 #include "wallDist.H"
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace RASModels
36 {
37 
38 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
39 
40 template<class BasicMomentumTransportModel>
42 {
43  tmp<volScalarField> tCmuk2(Cmu_*sqr(k_));
44  epsilon_ = max(epsilon_, tCmuk2()/(this->nutMaxCoeff_*this->nu()));
45  return tCmuk2;
46 }
47 
48 
49 template<class BasicMomentumTransportModel>
51 {
52  this->nut_ = boundEpsilon()/epsilon_;
53  this->nut_.correctBoundaryConditions();
54  fvConstraints::New(this->mesh_).constrain(this->nut_);
55 }
56 
57 
58 template<class BasicMomentumTransportModel>
60 {
61  return tmp<fvScalarMatrix>
62  (
63  new fvScalarMatrix
64  (
65  epsilon_,
66  dimVolume*this->rho_.dimensions()*epsilon_.dimensions()/dimTime
67  )
68  );
69 }
70 
71 
72 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
73 
74 template<class BasicMomentumTransportModel>
76 (
77  const alphaField& alpha,
78  const rhoField& rho,
79  const volVectorField& U,
80  const surfaceScalarField& alphaRhoPhi,
81  const surfaceScalarField& phi,
82  const viscosity& viscosity,
83  const word& type
84 )
85 :
86  ReynoldsStress<RASModel<BasicMomentumTransportModel>>
87  (
88  type,
89  alpha,
90  rho,
91  U,
92  alphaRhoPhi,
93  phi,
94  viscosity
95  ),
96 
97  Cmu_
98  (
99  dimensioned<scalar>::lookupOrAddToDict
100  (
101  "Cmu",
102  this->coeffDict_,
103  0.09
104  )
105  ),
106  C1_
107  (
108  dimensioned<scalar>::lookupOrAddToDict
109  (
110  "C1",
111  this->coeffDict_,
112  1.8
113  )
114  ),
115  C2_
116  (
117  dimensioned<scalar>::lookupOrAddToDict
118  (
119  "C2",
120  this->coeffDict_,
121  0.6
122  )
123  ),
124  Ceps1_
125  (
126  dimensioned<scalar>::lookupOrAddToDict
127  (
128  "Ceps1",
129  this->coeffDict_,
130  1.44
131  )
132  ),
133  Ceps2_
134  (
135  dimensioned<scalar>::lookupOrAddToDict
136  (
137  "Ceps2",
138  this->coeffDict_,
139  1.92
140  )
141  ),
142  Cs_
143  (
144  dimensioned<scalar>::lookupOrAddToDict
145  (
146  "Cs",
147  this->coeffDict_,
148  0.25
149  )
150  ),
151  Ceps_
152  (
153  dimensioned<scalar>::lookupOrAddToDict
154  (
155  "Ceps",
156  this->coeffDict_,
157  0.15
158  )
159  ),
160 
161  wallReflection_
162  (
163  Switch::lookupOrAddToDict
164  (
165  "wallReflection",
166  this->coeffDict_,
167  true
168  )
169  ),
170  kappa_
171  (
172  dimensioned<scalar>::lookupOrAddToDict
173  (
174  "kappa",
175  this->coeffDict_,
176  0.41
177  )
178  ),
179  Cref1_
180  (
181  dimensioned<scalar>::lookupOrAddToDict
182  (
183  "Cref1",
184  this->coeffDict_,
185  0.5
186  )
187  ),
188  Cref2_
189  (
190  dimensioned<scalar>::lookupOrAddToDict
191  (
192  "Cref2",
193  this->coeffDict_,
194  0.3
195  )
196  ),
197 
198  k_
199  (
200  IOobject
201  (
202  this->groupName("k"),
203  this->runTime_.name(),
204  this->mesh_,
205  IOobject::NO_READ,
206  IOobject::AUTO_WRITE
207  ),
208  0.5*tr(this->R_)
209  ),
210  epsilon_
211  (
212  IOobject
213  (
214  this->groupName("epsilon"),
215  this->runTime_.name(),
216  this->mesh_,
217  IOobject::MUST_READ,
218  IOobject::AUTO_WRITE
219  ),
220  this->mesh_
221  )
222 {
223  if (type == typeName)
224  {
225  this->printCoeffs(type);
226 
227  this->boundNormalStress(this->R_);
228  boundEpsilon();
229  k_ = 0.5*tr(this->R_);
230  }
231 }
232 
233 
234 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
235 
236 template<class BasicMomentumTransportModel>
238 {
240  {
241  Cmu_.readIfPresent(this->coeffDict());
242  C1_.readIfPresent(this->coeffDict());
243  C2_.readIfPresent(this->coeffDict());
244  Ceps1_.readIfPresent(this->coeffDict());
245  Ceps2_.readIfPresent(this->coeffDict());
246  Cs_.readIfPresent(this->coeffDict());
247  Ceps_.readIfPresent(this->coeffDict());
248 
249  wallReflection_.readIfPresent("wallReflection", this->coeffDict());
250  kappa_.readIfPresent(this->coeffDict());
251  Cref1_.readIfPresent(this->coeffDict());
252  Cref2_.readIfPresent(this->coeffDict());
253 
254  return true;
255  }
256  else
257  {
258  return false;
259  }
260 }
261 
262 
263 template<class BasicMomentumTransportModel>
265 {
267  (
268  "DREff",
269  (Cs_*(this->k_/this->epsilon_))*this->R_ + I*this->nu()
270  );
271 }
272 
273 
274 template<class BasicMomentumTransportModel>
276 {
278  (
279  "DepsilonEff",
280  (Ceps_*(this->k_/this->epsilon_))*this->R_ + I*this->nu()
281  );
282 }
283 
284 
285 template<class BasicMomentumTransportModel>
287 {
288  if (!this->turbulence_)
289  {
290  return;
291  }
292 
293  // Local references
294  const alphaField& alpha = this->alpha_;
295  const rhoField& rho = this->rho_;
296  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
297  const volVectorField& U = this->U_;
298  volSymmTensorField& R = this->R_;
299  const Foam::fvModels& fvModels(Foam::fvModels::New(this->mesh_));
301  (
302  Foam::fvConstraints::New(this->mesh_)
303  );
304 
306 
308  const volTensorField& gradU = tgradU();
309 
310  volSymmTensorField P(-twoSymm(R & gradU));
311  volScalarField G(this->GName(), 0.5*mag(tr(P)));
312 
313  // Update epsilon and G at the wall
314  epsilon_.boundaryFieldRef().updateCoeffs();
315 
316  // Dissipation equation
317  tmp<fvScalarMatrix> epsEqn
318  (
319  fvm::ddt(alpha, rho, epsilon_)
320  + fvm::div(alphaRhoPhi, epsilon_)
321  - fvm::laplacian(alpha*rho*DepsilonEff(), epsilon_)
322  ==
323  Ceps1_*alpha*rho*G*epsilon_/k_
324  - fvm::Sp(Ceps2_*alpha*rho*epsilon_/k_, epsilon_)
325  + epsilonSource()
326  + fvModels.source(alpha, rho, epsilon_)
327  );
328 
329  epsEqn.ref().relax();
330  fvConstraints.constrain(epsEqn.ref());
331  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
332  solve(epsEqn);
333  fvConstraints.constrain(epsilon_);
334  boundEpsilon();
335 
336 
337  // Correct the trace of the tensorial production to be consistent
338  // with the near-wall generation from the wall-functions
339  const fvPatchList& patches = this->mesh_.boundary();
340 
342  {
343  const fvPatch& curPatch = patches[patchi];
344 
345  if (isA<wallFvPatch>(curPatch))
346  {
347  forAll(curPatch, facei)
348  {
349  label celli = curPatch.faceCells()[facei];
350  P[celli] *= min
351  (
352  G[celli]/(0.5*mag(tr(P[celli])) + small),
353  1.0
354  );
355  }
356  }
357  }
358 
359  // Reynolds stress equation
361  (
362  fvm::ddt(alpha, rho, R)
363  + fvm::div(alphaRhoPhi, R)
364  - fvm::laplacian(alpha*rho*DREff(), R)
365  + fvm::Sp(C1_*alpha*rho*epsilon_/k_, R)
366  ==
367  alpha*rho*P
368  - (2.0/3.0*(1 - C1_)*I)*alpha*rho*epsilon_
369  - C2_*alpha*rho*dev(P)
370  + this->RSource()
371  + fvModels.source(alpha, rho, R)
372  );
373 
374  // Optionally add wall-refection term
375  if (wallReflection_)
376  {
377  const volVectorField& n(wallDist::New(this->mesh_).n());
378  const volScalarField& y(wallDist::New(this->mesh_).y());
379 
380  const volSymmTensorField reflect
381  (
382  Cref1_*R - ((Cref2_*C2_)*(k_/epsilon_))*dev(P)
383  );
384 
385  REqn.ref() +=
386  ((3*pow(Cmu_, 0.75)/kappa_)*(alpha*rho*sqrt(k_)/y))
387  *dev(symm((n & reflect)*n));
388  }
389 
390  REqn.ref().relax();
391  fvConstraints.constrain(REqn.ref());
392  solve(REqn);
394 
395  this->boundNormalStress(R);
396 
397  k_ = 0.5*tr(R);
398 
399  correctNut();
400 
401  // Correct wall shear-stresses when applying wall-functions
402  this->correctWallShearStress(R);
403 }
404 
405 
406 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
407 
408 } // End namespace RASModels
409 } // End namespace Foam
410 
411 // ************************************************************************* //
scalar y
label n
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
static wallDist & New(const word &name, const fvMesh &mesh)
Construct and return the named DemandDrivenMeshObject.
Generic GeometricField class.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
Templated abstract base class for RAS turbulence models.
Definition: RASModel.H:56
virtual void printCoeffs(const word &type)
Print model coefficients.
Definition: RASModel.C:32
BasicMomentumTransportModel::alphaField alphaField
Definition: LRR.H:146
tmp< volSymmTensorField > DepsilonEff() const
Return the effective diffusivity for epsilon.
Definition: LRR.C:275
virtual tmp< fvScalarMatrix > epsilonSource() const
Source term for the epsilon equation.
Definition: LRR.C:59
volScalarField k_
Definition: LRR.H:128
LRR(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity, const word &type=typeName)
Construct from components.
Definition: LRR.C:76
virtual void correct()
Solve the turbulence equations and correct eddy-Viscosity and.
Definition: LRR.C:286
tmp< volScalarField > boundEpsilon()
Bound epsilon and return Cmu*sqr(k) for nut.
Definition: LRR.C:41
tmp< volSymmTensorField > DREff() const
Return the effective diffusivity for R.
Definition: LRR.C:264
virtual void correctNut()
Correct the eddy-viscosity nut.
Definition: LRR.C:50
virtual bool read()
Read model coefficients if they have changed.
Definition: LRR.C:237
BasicMomentumTransportModel::rhoField rhoField
Definition: LRR.H:147
Reynolds-stress turbulence model base class.
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:61
static autoPtr< dictionary > New(Istream &)
Construct top-level dictionary on freestore from Istream.
Definition: dictionaryIO.C:100
Generic dimensioned Type class.
Finite volume constraints.
Definition: fvConstraints.H:67
bool constrain(fvMatrix< Type > &eqn) const
Apply constraints to an equation.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
Finite volume models.
Definition: fvModels.H:65
tmp< fvMatrix< Type > > source(const VolField< Type > &field) const
Return source for an equation.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
virtual const labelUList & faceCells() const
Return faceCells.
Definition: fvPatch.C:99
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
Abstract base class for all fluid physical properties.
Definition: viscosity.H:50
A class for handling words, derived from string.
Definition: word.H:62
Foam::fvConstraints & fvConstraints(Foam::fvConstraints::New(mesh))
Foam::fvModels & fvModels(Foam::fvModels::New(mesh))
label patchi
const fvPatchList & patches
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
void correct(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiCorr, const SpType &Sp, const SuType &Su)
const dimensionedScalar G
Newtonian constant of gravitation.
tmp< VolField< typename outerProduct< vector, Type >::type > > grad(const SurfaceField< Type > &ssf)
Definition: fvcGrad.C:46
tmp< fvMatrix< Type > > laplacian(const VolField< Type > &vf, const word &name)
Definition: fvmLaplacian.C:47
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const VolField< Type > &vf, const word &name)
Definition: fvmDiv.C:48
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const VolField< Type > &)
tmp< fvMatrix< Type > > ddt(const VolField< Type > &vf)
Definition: fvmDdt.C:46
Namespace for OpenFOAM.
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
bool read(const char *, int32_t &)
Definition: int32IO.C:85
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
dimensionedSymmTensor sqr(const dimensionedVector &dv)
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
static const Identity< scalar > I
Definition: Identity.H:93
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
const dimensionSet dimTime
dimensionedScalar sqrt(const dimensionedScalar &ds)
const dimensionSet dimVolume
dimensioned< scalar > mag(const dimensioned< Type > &)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
static scalar R(const scalar a, const scalar x)
Definition: invIncGamma.C:102
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
SolverPerformance< Type > solve(fvMatrix< Type > &, const word &)
Solve returning the solution statistics given convergence tolerance.