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-2024 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_("Cmu", this->coeffDict(), 0.09),
98  C1_("C1", this->coeffDict(), 1.8),
99  C2_("C2", this->coeffDict(), 0.6),
100  Ceps1_("Ceps1", this->coeffDict(), 1.44),
101  Ceps2_("Ceps2", this->coeffDict(), 1.92),
102  Cs_("Cs", this->coeffDict(), 0.25),
103  Ceps_("Ceps", this->coeffDict(), 0.15),
104 
105  wallReflection_
106  (
107  this->coeffDict().template lookupOrDefault<Switch>
108  (
109  "wallReflection",
110  true
111  )
112  ),
113  kappa_("kappa", this->coeffDict(), 0.41),
114  Cref1_("Cref1", this->coeffDict(), 0.5),
115  Cref2_("Cref2", this->coeffDict(), 0.3),
116 
117  k_
118  (
119  IOobject
120  (
121  this->groupName("k"),
122  this->runTime_.name(),
123  this->mesh_,
124  IOobject::NO_READ,
125  IOobject::AUTO_WRITE
126  ),
127  0.5*tr(this->R_)
128  ),
129  epsilon_
130  (
131  IOobject
132  (
133  this->groupName("epsilon"),
134  this->runTime_.name(),
135  this->mesh_,
136  IOobject::MUST_READ,
137  IOobject::AUTO_WRITE
138  ),
139  this->mesh_
140  )
141 {
142  if (type == typeName)
143  {
144  this->boundNormalStress(this->R_);
145  boundEpsilon();
146  k_ = 0.5*tr(this->R_);
147  }
148 }
149 
150 
151 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
152 
153 template<class BasicMomentumTransportModel>
155 {
157  {
158  Cmu_.readIfPresent(this->coeffDict());
159  C1_.readIfPresent(this->coeffDict());
160  C2_.readIfPresent(this->coeffDict());
161  Ceps1_.readIfPresent(this->coeffDict());
162  Ceps2_.readIfPresent(this->coeffDict());
163  Cs_.readIfPresent(this->coeffDict());
164  Ceps_.readIfPresent(this->coeffDict());
165 
166  wallReflection_.readIfPresent("wallReflection", this->coeffDict());
167  kappa_.readIfPresent(this->coeffDict());
168  Cref1_.readIfPresent(this->coeffDict());
169  Cref2_.readIfPresent(this->coeffDict());
170 
171  return true;
172  }
173  else
174  {
175  return false;
176  }
177 }
178 
179 
180 template<class BasicMomentumTransportModel>
182 {
184  (
185  "DREff",
186  (Cs_*(this->k_/this->epsilon_))*this->R_ + I*this->nu()
187  );
188 }
189 
190 
191 template<class BasicMomentumTransportModel>
193 {
195  (
196  "DepsilonEff",
197  (Ceps_*(this->k_/this->epsilon_))*this->R_ + I*this->nu()
198  );
199 }
200 
201 
202 template<class BasicMomentumTransportModel>
204 {
205  if (!this->turbulence_)
206  {
207  return;
208  }
209 
210  // Local references
211  const alphaField& alpha = this->alpha_;
212  const rhoField& rho = this->rho_;
213  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
214  const volVectorField& U = this->U_;
215  volSymmTensorField& R = this->R_;
216  const Foam::fvModels& fvModels(Foam::fvModels::New(this->mesh_));
218  (
219  Foam::fvConstraints::New(this->mesh_)
220  );
221 
223 
225  const volTensorField& gradU = tgradU();
226 
227  volSymmTensorField P(-twoSymm(R & gradU));
228  volScalarField G(this->GName(), 0.5*mag(tr(P)));
229 
230  // Update epsilon and G at the wall
231  epsilon_.boundaryFieldRef().updateCoeffs();
232 
233  // Dissipation equation
234  tmp<fvScalarMatrix> epsEqn
235  (
236  fvm::ddt(alpha, rho, epsilon_)
237  + fvm::div(alphaRhoPhi, epsilon_)
238  - fvm::laplacian(alpha*rho*DepsilonEff(), epsilon_)
239  ==
240  Ceps1_*alpha*rho*G*epsilon_/k_
241  - fvm::Sp(Ceps2_*alpha*rho*epsilon_/k_, epsilon_)
242  + epsilonSource()
243  + fvModels.source(alpha, rho, epsilon_)
244  );
245 
246  epsEqn.ref().relax();
247  fvConstraints.constrain(epsEqn.ref());
248  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
249  solve(epsEqn);
250  fvConstraints.constrain(epsilon_);
251  boundEpsilon();
252 
253 
254  // Correct the trace of the tensorial production to be consistent
255  // with the near-wall generation from the wall-functions
256  const fvPatchList& patches = this->mesh_.boundary();
257 
259  {
260  const fvPatch& curPatch = patches[patchi];
261 
262  if (isA<wallFvPatch>(curPatch))
263  {
264  forAll(curPatch, facei)
265  {
266  label celli = curPatch.faceCells()[facei];
267  P[celli] *= min
268  (
269  G[celli]/(0.5*mag(tr(P[celli])) + small),
270  1.0
271  );
272  }
273  }
274  }
275 
276  // Reynolds stress equation
278  (
279  fvm::ddt(alpha, rho, R)
280  + fvm::div(alphaRhoPhi, R)
281  - fvm::laplacian(alpha*rho*DREff(), R)
282  + fvm::Sp(C1_*alpha*rho*epsilon_/k_, R)
283  ==
284  alpha*rho*P
285  - (2.0/3.0*(1 - C1_)*I)*alpha*rho*epsilon_
286  - C2_*alpha*rho*dev(P)
287  + this->RSource()
288  + fvModels.source(alpha, rho, R)
289  );
290 
291  // Optionally add wall-refection term
292  if (wallReflection_)
293  {
294  const volVectorField& n(wallDist::New(this->mesh_).n());
295  const volScalarField& y(wallDist::New(this->mesh_).y());
296 
298  (
299  Cref1_*R - ((Cref2_*C2_)*(k_/epsilon_))*dev(P)
300  );
301 
302  REqn.ref() +=
303  ((3*pow(Cmu_, 0.75)/kappa_)*(alpha*rho*sqrt(k_)/y))
304  *dev(symm((n & reflect)*n));
305  }
306 
307  REqn.ref().relax();
308  fvConstraints.constrain(REqn.ref());
309  solve(REqn);
311 
312  this->boundNormalStress(R);
313 
314  k_ = 0.5*tr(R);
315 
316  correctNut();
317 
318  // Correct wall shear-stresses when applying wall-functions
319  this->correctWallShearStress(R);
320 }
321 
322 
323 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
324 
325 } // End namespace RASModels
326 } // End namespace Foam
327 
328 // ************************************************************************* //
scalar y
label n
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
static wallDist & New(const word &name, const fvMesh &mesh)
Construct and return the named DemandDrivenMeshObject.
Generic GeometricField class.
static tmp< GeometricField< Type, GeoMesh, PrimitiveField > > New(const word &name, const Internal &, const PtrList< Patch > &, 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
BasicMomentumTransportModel::alphaField alphaField
Definition: LRR.H:146
tmp< volSymmTensorField > DepsilonEff() const
Return the effective diffusivity for epsilon.
Definition: LRR.C:192
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:203
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:181
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:154
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:103
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:197
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 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
void reflect(barycentric &coordinates)
Reflection transform. Corrects the coordinates when the track moves.
Definition: tracking.C:1028
Namespace for OpenFOAM.
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
void dev(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
void twoSymm(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
static const Identity< scalar > I
Definition: Identity.H:93
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
const dimensionSet dimTime
void symm(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
void pow(LagrangianPatchField< typename powProduct< Type, r >::type > &f, const LagrangianPatchField< Type > &f1)
const dimensionSet dimVolume
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
static scalar R(const scalar a, const scalar x)
Definition: invIncGamma.C:102
void sqr(LagrangianPatchField< typename outerProduct< Type, Type >::type > &f, const LagrangianPatchField< Type > &f1)
void tr(LagrangianPatchField< scalar > &f, const LagrangianPatchField< tensor > &f1)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
void sqrt(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
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.