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