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-2020 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 "fvOptions.H"
28 #include "wallDist.H"
29 
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace RASModels
35 {
36 
37 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
38 
39 template<class BasicMomentumTransportModel>
41 {
42  this->nut_ = this->Cmu_*sqr(k_)/epsilon_;
43  this->nut_.correctBoundaryConditions();
44  fv::options::New(this->mesh_).correct(this->nut_);
45 }
46 
47 
48 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
49 
50 template<class BasicMomentumTransportModel>
52 (
53  const alphaField& alpha,
54  const rhoField& rho,
55  const volVectorField& U,
56  const surfaceScalarField& alphaRhoPhi,
57  const surfaceScalarField& phi,
58  const transportModel& transport,
59  const word& type
60 )
61 :
63  (
64  type,
65  alpha,
66  rho,
67  U,
68  alphaRhoPhi,
69  phi,
70  transport
71  ),
72 
73  Cmu_
74  (
76  (
77  "Cmu",
78  this->coeffDict_,
79  0.09
80  )
81  ),
82  C1_
83  (
85  (
86  "C1",
87  this->coeffDict_,
88  1.8
89  )
90  ),
91  C2_
92  (
94  (
95  "C2",
96  this->coeffDict_,
97  0.6
98  )
99  ),
100  Ceps1_
101  (
103  (
104  "Ceps1",
105  this->coeffDict_,
106  1.44
107  )
108  ),
109  Ceps2_
110  (
112  (
113  "Ceps2",
114  this->coeffDict_,
115  1.92
116  )
117  ),
118  Cs_
119  (
121  (
122  "Cs",
123  this->coeffDict_,
124  0.25
125  )
126  ),
127  Ceps_
128  (
130  (
131  "Ceps",
132  this->coeffDict_,
133  0.15
134  )
135  ),
136 
137  wallReflection_
138  (
140  (
141  "wallReflection",
142  this->coeffDict_,
143  true
144  )
145  ),
146  kappa_
147  (
149  (
150  "kappa",
151  this->coeffDict_,
152  0.41
153  )
154  ),
155  Cref1_
156  (
158  (
159  "Cref1",
160  this->coeffDict_,
161  0.5
162  )
163  ),
164  Cref2_
165  (
167  (
168  "Cref2",
169  this->coeffDict_,
170  0.3
171  )
172  ),
173 
174  k_
175  (
176  IOobject
177  (
178  "k",
179  this->runTime_.timeName(),
180  this->mesh_,
183  ),
184  0.5*tr(this->R_)
185  ),
186  epsilon_
187  (
188  IOobject
189  (
190  "epsilon",
191  this->runTime_.timeName(),
192  this->mesh_,
195  ),
196  this->mesh_
197  )
198 {
199  if (type == typeName)
200  {
201  this->printCoeffs(type);
202 
203  this->boundNormalStress(this->R_);
204  bound(epsilon_, this->epsilonMin_);
205  k_ = 0.5*tr(this->R_);
206  }
207 }
208 
209 
210 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
211 
212 template<class BasicMomentumTransportModel>
214 {
216  {
217  Cmu_.readIfPresent(this->coeffDict());
218  C1_.readIfPresent(this->coeffDict());
219  C2_.readIfPresent(this->coeffDict());
220  Ceps1_.readIfPresent(this->coeffDict());
221  Ceps2_.readIfPresent(this->coeffDict());
222  Cs_.readIfPresent(this->coeffDict());
223  Ceps_.readIfPresent(this->coeffDict());
224 
225  wallReflection_.readIfPresent("wallReflection", this->coeffDict());
226  kappa_.readIfPresent(this->coeffDict());
227  Cref1_.readIfPresent(this->coeffDict());
228  Cref2_.readIfPresent(this->coeffDict());
229 
230  return true;
231  }
232  else
233  {
234  return false;
235  }
236 }
237 
238 
239 template<class BasicMomentumTransportModel>
241 {
243  (
244  "DREff",
245  (Cs_*(this->k_/this->epsilon_))*this->R_ + I*this->nu()
246  );
247 }
248 
249 
250 template<class BasicMomentumTransportModel>
252 {
254  (
255  "DepsilonEff",
256  (Ceps_*(this->k_/this->epsilon_))*this->R_ + I*this->nu()
257  );
258 }
259 
260 
261 template<class BasicMomentumTransportModel>
263 {
264  if (!this->turbulence_)
265  {
266  return;
267  }
268 
269  // Local references
270  const alphaField& alpha = this->alpha_;
271  const rhoField& rho = this->rho_;
272  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
273  const volVectorField& U = this->U_;
274  volSymmTensorField& R = this->R_;
275  fv::options& fvOptions(fv::options::New(this->mesh_));
276 
278 
279  tmp<volTensorField> tgradU(fvc::grad(U));
280  const volTensorField& gradU = tgradU();
281 
282  volSymmTensorField P(-twoSymm(R & gradU));
283  volScalarField G(this->GName(), 0.5*mag(tr(P)));
284 
285  // Update epsilon and G at the wall
286  epsilon_.boundaryFieldRef().updateCoeffs();
287 
288  // Dissipation equation
289  tmp<fvScalarMatrix> epsEqn
290  (
291  fvm::ddt(alpha, rho, epsilon_)
292  + fvm::div(alphaRhoPhi, epsilon_)
293  - fvm::laplacian(alpha*rho*DepsilonEff(), epsilon_)
294  ==
295  Ceps1_*alpha*rho*G*epsilon_/k_
296  - fvm::Sp(Ceps2_*alpha*rho*epsilon_/k_, epsilon_)
297  + fvOptions(alpha, rho, epsilon_)
298  );
299 
300  epsEqn.ref().relax();
301  fvOptions.constrain(epsEqn.ref());
302  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
303  solve(epsEqn);
304  fvOptions.correct(epsilon_);
305  bound(epsilon_, this->epsilonMin_);
306 
307 
308  // Correct the trace of the tensorial production to be consistent
309  // with the near-wall generation from the wall-functions
310  const fvPatchList& patches = this->mesh_.boundary();
311 
312  forAll(patches, patchi)
313  {
314  const fvPatch& curPatch = patches[patchi];
315 
316  if (isA<wallFvPatch>(curPatch))
317  {
318  forAll(curPatch, facei)
319  {
320  label celli = curPatch.faceCells()[facei];
321  P[celli] *= min
322  (
323  G[celli]/(0.5*mag(tr(P[celli])) + small),
324  1.0
325  );
326  }
327  }
328  }
329 
330  // Reynolds stress equation
332  (
333  fvm::ddt(alpha, rho, R)
334  + fvm::div(alphaRhoPhi, R)
335  - fvm::laplacian(alpha*rho*DREff(), R)
336  + fvm::Sp(C1_*alpha*rho*epsilon_/k_, R)
337  ==
338  alpha*rho*P
339  - (2.0/3.0*(1 - C1_)*I)*alpha*rho*epsilon_
340  - C2_*alpha*rho*dev(P)
341  + fvOptions(alpha, rho, R)
342  );
343 
344  // Optionally add wall-refection term
345  if (wallReflection_)
346  {
347  const volVectorField& n_(wallDist::New(this->mesh_).n());
348  const volScalarField& y_(wallDist::New(this->mesh_).y());
349 
350  const volSymmTensorField reflect
351  (
352  Cref1_*R - ((Cref2_*C2_)*(k_/epsilon_))*dev(P)
353  );
354 
355  REqn.ref() +=
356  ((3*pow(Cmu_, 0.75)/kappa_)*(alpha*rho*sqrt(k_)/y_))
357  *dev(symm((n_ & reflect)*n_));
358  }
359 
360  REqn.ref().relax();
361  fvOptions.constrain(REqn.ref());
362  solve(REqn);
363  fvOptions.correct(R);
364 
365  this->boundNormalStress(R);
366 
367  k_ = 0.5*tr(R);
368 
369  correctNut();
370 
371  // Correct wall shear-stresses when applying wall-functions
372  this->correctWallShearStress(R);
373 }
374 
375 
376 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
377 
378 } // End namespace RASModels
379 } // End namespace Foam
380 
381 // ************************************************************************* //
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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
fv::options & fvOptions
tmp< volSymmTensorField > DREff() const
Return the effective diffusivity for R.
Definition: LRR.C:240
void correct(GeometricField< Type, fvPatchField, volMesh > &field)
Apply correction to field.
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
dimensionedSymmTensor sqr(const dimensionedVector &dv)
static tmp< GeometricField< symmTensor, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< symmTensor >> &)
Return a temporary field constructed from name,.
dimensionedScalar sqrt(const dimensionedScalar &ds)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
Generic dimensioned Type class.
patches[0]
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
virtual bool read()
Read model coefficients if they have changed.
Definition: LRR.C:213
Finite-volume options.
Definition: fvOptions.H:52
BasicMomentumTransportModel::alphaField alphaField
Definition: LRR.H:140
static Switch lookupOrAddToDict(const word &, dictionary &, const Switch &defaultValue=false)
Construct from dictionary, supplying default value so that if the.
Definition: Switch.C:111
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
tmp< volSymmTensorField > DepsilonEff() const
Return the effective diffusivity for epsilon.
Definition: LRR.C:251
phi
Definition: pEqn.H:104
static const wallDist & New(const fvMesh &mesh)
Definition: MeshObject.C:44
Templated abstract base class for RAS turbulence models.
Definition: RASModel.H:49
scalar y
bool read(const char *, int32_t &)
Definition: int32IO.C:85
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
static const Identity< scalar > I
Definition: Identity.H:93
BasicMomentumTransportModel::transportModel transportModel
Definition: LRR.H:142
virtual void correctNut()
Update the eddy-viscosity.
Definition: LRR.C:40
A class for handling words, derived from string.
Definition: word.H:59
virtual const labelUList & faceCells() const
Return faceCells.
Definition: fvPatch.C:93
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
void constrain(fvMatrix< Type > &eqn)
Apply constraints to equation.
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:46
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(rho0+psi *p, rhoMin);# 1 "/home/ubuntu/OpenFOAM-8/applications/solvers/multiphase/cavitatingFoam/alphavPsi.H" 1{ alphav=max(min((rho - rholSat)/(rhovSat - rholSat), scalar(1)), scalar(0));alphal=1.0 - alphav;Info<< "max-min alphav: "<< max(alphav).value()<< " "<< min(alphav).value()<< endl;psiModel-> correct()
Definition: pEqn.H:68
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
BasicMomentumTransportModel::rhoField rhoField
Definition: LRR.H:141
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
virtual void correct()
Solve the turbulence equations and correct eddy-Viscosity and.
Definition: LRR.C:262
label patchi
U
Definition: pEqn.H:72
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:33
#define R(A, B, C, D, E, F, K, M)
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:46
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:70
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.
dimensioned< scalar > mag(const dimensioned< Type > &)
label n
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
LRR(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &type=typeName)
Construct from components.
Definition: LRR.C:52
A class for managing temporary objects.
Definition: PtrList.H:53
static options & New(const fvMesh &mesh)
Construct fvOptions and register to datbase if not present.
Definition: fvOptions.C:101
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
const dimensionedScalar & G
Newtonian constant of gravitation.
Namespace for OpenFOAM.
Reynolds-stress turbulence model base class.