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-2018 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 BasicTurbulenceModel>
41 {
42  this->nut_ = this->Cmu_*sqr(k_)/epsilon_;
43  this->nut_.correctBoundaryConditions();
44  fv::options::New(this->mesh_).correct(this->nut_);
45 
46  BasicTurbulenceModel::correctNut();
47 }
48 
49 
50 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
51 
52 template<class BasicTurbulenceModel>
54 (
55  const alphaField& alpha,
56  const rhoField& rho,
57  const volVectorField& U,
58  const surfaceScalarField& alphaRhoPhi,
59  const surfaceScalarField& phi,
60  const transportModel& transport,
61  const word& propertiesName,
62  const word& type
63 )
64 :
66  (
67  type,
68  alpha,
69  rho,
70  U,
71  alphaRhoPhi,
72  phi,
73  transport,
74  propertiesName
75  ),
76 
77  Cmu_
78  (
80  (
81  "Cmu",
82  this->coeffDict_,
83  0.09
84  )
85  ),
86  C1_
87  (
89  (
90  "C1",
91  this->coeffDict_,
92  1.8
93  )
94  ),
95  C2_
96  (
98  (
99  "C2",
100  this->coeffDict_,
101  0.6
102  )
103  ),
104  Ceps1_
105  (
107  (
108  "Ceps1",
109  this->coeffDict_,
110  1.44
111  )
112  ),
113  Ceps2_
114  (
116  (
117  "Ceps2",
118  this->coeffDict_,
119  1.92
120  )
121  ),
122  Cs_
123  (
125  (
126  "Cs",
127  this->coeffDict_,
128  0.25
129  )
130  ),
131  Ceps_
132  (
134  (
135  "Ceps",
136  this->coeffDict_,
137  0.15
138  )
139  ),
140 
141  wallReflection_
142  (
144  (
145  "wallReflection",
146  this->coeffDict_,
147  true
148  )
149  ),
150  kappa_
151  (
153  (
154  "kappa",
155  this->coeffDict_,
156  0.41
157  )
158  ),
159  Cref1_
160  (
162  (
163  "Cref1",
164  this->coeffDict_,
165  0.5
166  )
167  ),
168  Cref2_
169  (
171  (
172  "Cref2",
173  this->coeffDict_,
174  0.3
175  )
176  ),
177 
178  k_
179  (
180  IOobject
181  (
182  "k",
183  this->runTime_.timeName(),
184  this->mesh_,
187  ),
188  0.5*tr(this->R_)
189  ),
190  epsilon_
191  (
192  IOobject
193  (
194  "epsilon",
195  this->runTime_.timeName(),
196  this->mesh_,
199  ),
200  this->mesh_
201  )
202 {
203  if (type == typeName)
204  {
205  this->printCoeffs(type);
206 
207  this->boundNormalStress(this->R_);
208  bound(epsilon_, this->epsilonMin_);
209  k_ = 0.5*tr(this->R_);
210  }
211 }
212 
213 
214 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
215 
216 template<class BasicTurbulenceModel>
218 {
220  {
221  Cmu_.readIfPresent(this->coeffDict());
222  C1_.readIfPresent(this->coeffDict());
223  C2_.readIfPresent(this->coeffDict());
224  Ceps1_.readIfPresent(this->coeffDict());
225  Ceps2_.readIfPresent(this->coeffDict());
226  Cs_.readIfPresent(this->coeffDict());
227  Ceps_.readIfPresent(this->coeffDict());
228 
229  wallReflection_.readIfPresent("wallReflection", this->coeffDict());
230  kappa_.readIfPresent(this->coeffDict());
231  Cref1_.readIfPresent(this->coeffDict());
232  Cref2_.readIfPresent(this->coeffDict());
233 
234  return true;
235  }
236  else
237  {
238  return false;
239  }
240 }
241 
242 
243 template<class BasicTurbulenceModel>
245 {
247  (
248  "DREff",
249  (Cs_*(this->k_/this->epsilon_))*this->R_ + I*this->nu()
250  );
251 }
252 
253 
254 template<class BasicTurbulenceModel>
256 {
258  (
259  "DepsilonEff",
260  (Ceps_*(this->k_/this->epsilon_))*this->R_ + I*this->nu()
261  );
262 }
263 
264 
265 template<class BasicTurbulenceModel>
267 {
268  if (!this->turbulence_)
269  {
270  return;
271  }
272 
273  // Local references
274  const alphaField& alpha = this->alpha_;
275  const rhoField& rho = this->rho_;
276  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
277  const volVectorField& U = this->U_;
278  volSymmTensorField& R = this->R_;
279  fv::options& fvOptions(fv::options::New(this->mesh_));
280 
282 
283  tmp<volTensorField> tgradU(fvc::grad(U));
284  const volTensorField& gradU = tgradU();
285 
286  volSymmTensorField P(-twoSymm(R & gradU));
287  volScalarField G(this->GName(), 0.5*mag(tr(P)));
288 
289  // Update epsilon and G at the wall
290  epsilon_.boundaryFieldRef().updateCoeffs();
291 
292  // Dissipation equation
293  tmp<fvScalarMatrix> epsEqn
294  (
295  fvm::ddt(alpha, rho, epsilon_)
296  + fvm::div(alphaRhoPhi, epsilon_)
297  - fvm::laplacian(alpha*rho*DepsilonEff(), epsilon_)
298  ==
299  Ceps1_*alpha*rho*G*epsilon_/k_
300  - fvm::Sp(Ceps2_*alpha*rho*epsilon_/k_, epsilon_)
301  + fvOptions(alpha, rho, epsilon_)
302  );
303 
304  epsEqn.ref().relax();
305  fvOptions.constrain(epsEqn.ref());
306  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
307  solve(epsEqn);
308  fvOptions.correct(epsilon_);
309  bound(epsilon_, this->epsilonMin_);
310 
311 
312  // Correct the trace of the tensorial production to be consistent
313  // with the near-wall generation from the wall-functions
314  const fvPatchList& patches = this->mesh_.boundary();
315 
316  forAll(patches, patchi)
317  {
318  const fvPatch& curPatch = patches[patchi];
319 
320  if (isA<wallFvPatch>(curPatch))
321  {
322  forAll(curPatch, facei)
323  {
324  label celli = curPatch.faceCells()[facei];
325  P[celli] *= min
326  (
327  G[celli]/(0.5*mag(tr(P[celli])) + small),
328  1.0
329  );
330  }
331  }
332  }
333 
334  // Reynolds stress equation
336  (
337  fvm::ddt(alpha, rho, R)
338  + fvm::div(alphaRhoPhi, R)
339  - fvm::laplacian(alpha*rho*DREff(), R)
340  + fvm::Sp(C1_*alpha*rho*epsilon_/k_, R)
341  ==
342  alpha*rho*P
343  - (2.0/3.0*(1 - C1_)*I)*alpha*rho*epsilon_
344  - C2_*alpha*rho*dev(P)
345  + fvOptions(alpha, rho, R)
346  );
347 
348  // Optionally add wall-refection term
349  if (wallReflection_)
350  {
351  const volVectorField& n_(wallDist::New(this->mesh_).n());
352  const volScalarField& y_(wallDist::New(this->mesh_).y());
353 
354  const volSymmTensorField reflect
355  (
356  Cref1_*R - ((Cref2_*C2_)*(k_/epsilon_))*dev(P)
357  );
358 
359  REqn.ref() +=
360  ((3*pow(Cmu_, 0.75)/kappa_)*(alpha*rho*sqrt(k_)/y_))
361  *dev(symm((n_ & reflect)*n_));
362  }
363 
364  REqn.ref().relax();
365  fvOptions.constrain(REqn.ref());
366  solve(REqn);
367  fvOptions.correct(R);
368 
369  this->boundNormalStress(R);
370 
371  k_ = 0.5*tr(R);
372 
373  correctNut();
374 
375  // Correct wall shear-stresses when applying wall-functions
376  this->correctWallShearStress(R);
377 }
378 
379 
380 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
381 
382 } // End namespace RASModels
383 } // End namespace Foam
384 
385 // ************************************************************************* //
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
surfaceScalarField & phi
tmp< volSymmTensorField > DREff() const
Return the effective diffusivity for R.
Definition: LRR.C:244
void correct(GeometricField< Type, fvPatchField, volMesh > &field)
Apply correction to field.
const dimensionedScalar G
Newtonian constant of gravitation.
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
dimensionedSymmTensor sqr(const dimensionedVector &dv)
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:217
Finite-volume options.
Definition: fvOptions.H:52
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:255
static const wallDist & New(const fvMesh &mesh)
Definition: MeshObject.C:44
static tmp< GeometricField< symmTensor, fvPatchField, volMesh > > New(const word &name, const Mesh &, const dimensionSet &, const word &patchFieldType=fvPatchField< symmTensor >::calculatedType())
Return a temporary field constructed from name, mesh, dimensionSet.
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
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
BasicTurbulenceModel::transportModel transportModel
Definition: LRR.H:142
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
LRR(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName=turbulenceModel::propertiesName, const word &type=typeName)
Construct from components.
Definition: LRR.C:54
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-7/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 > &)
BasicTurbulenceModel::alphaField alphaField
Definition: LRR.H:140
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
virtual void correct()
Solve the turbulence equations and correct eddy-Viscosity and.
Definition: LRR.C:266
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
BasicTurbulenceModel::rhoField rhoField
Definition: LRR.H:141
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
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
volScalarField & nu
Namespace for OpenFOAM.
Reynolds-stress turbulence model base class.