kEpsilonLopesdaCosta.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) 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 "kEpsilonLopesdaCosta.H"
27 #include "fvOptions.H"
28 #include "explicitPorositySource.H"
29 #include "bound.H"
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace RASModels
36 {
37 
38 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
39 
40 template<class BasicTurbulenceModel>
42 (
45 )
46 {
47  if (pm.dict().found(C.name()))
48  {
49  const labelList& cellZoneIDs = pm.cellZoneIDs();
50 
51  const scalar Cpm = readScalar(pm.dict().lookup(C.name()));
52 
53  forAll(cellZoneIDs, zonei)
54  {
55  const labelList& cells =
56  this->mesh_.cellZones()[cellZoneIDs[zonei]];
57 
58  forAll(cells, i)
59  {
60  const label celli = cells[i];
61  C[celli] = Cpm;
62  }
63  }
64  }
65 }
66 
67 
68 template<class BasicTurbulenceModel>
70 (
73 )
74 {
75  if (pm.dict().found(C.name()))
76  {
77  const labelList& cellZoneIDs = pm.cellZoneIDs();
78  const scalarField& Sigma = pm.Sigma();
79 
80  const scalar Cpm = readScalar(pm.dict().lookup(C.name()));
81 
82  forAll(cellZoneIDs, zonei)
83  {
84  const labelList& cells =
85  this->mesh_.cellZones()[cellZoneIDs[zonei]];
86 
87  forAll(cells, i)
88  {
89  const label celli = cells[i];
90  C[celli] = Cpm*Sigma[celli];
91  }
92  }
93  }
94 }
95 
96 
97 template<class BasicTurbulenceModel>
99 {
101 
102  forAll(fvOptions, i)
103  {
104  if (isA<fv::explicitPorositySource>(fvOptions[i]))
105  {
106  const fv::explicitPorositySource& eps =
107  refCast<const fv::explicitPorositySource>(fvOptions[i]);
108 
109  if (isA<porosityModels::powerLawLopesdaCosta>(eps.model()))
110  {
112  refCast<const porosityModels::powerLawLopesdaCosta>
113  (
114  eps.model()
115  );
116 
117  setPorosityCoefficient(Cmu_, pm);
118  setPorosityCoefficient(C1_, pm);
119  setPorosityCoefficient(C2_, pm);
120  setPorosityCoefficient(sigmak_, pm);
121  setPorosityCoefficient(sigmaEps_, pm);
122 
123  setCdSigma(CdSigma_, pm);
124  setPorosityCoefficient(betap_, pm);
125  setPorosityCoefficient(betad_, pm);
126  setPorosityCoefficient(C4_, pm);
127  setPorosityCoefficient(C5_, pm);
128  }
129  }
130  }
131 }
132 
133 
134 template<class BasicTurbulenceModel>
136 {
137  this->nut_ = Cmu_*sqr(k_)/epsilon_;
138  this->nut_.correctBoundaryConditions();
139  fv::options::New(this->mesh_).correct(this->nut_);
140 
141  BasicTurbulenceModel::correctNut();
142 }
143 
144 
145 template<class BasicTurbulenceModel>
147 (
148  const volScalarField::Internal& magU,
149  const volScalarField::Internal& magU3
150 ) const
151 {
152  return fvm::Su(CdSigma_*(betap_*magU3 - betad_*magU*k_()), k_);
153 }
154 
155 
156 template<class BasicTurbulenceModel>
159 (
160  const volScalarField::Internal& magU,
161  const volScalarField::Internal& magU3
162 ) const
163 {
164  return fvm::Su
165  (
166  CdSigma_
167  *(C4_*betap_*epsilon_()/k_()*magU3 - C5_*betad_*magU*epsilon_()),
168  epsilon_
169  );
170 }
171 
172 
173 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
174 
175 template<class BasicTurbulenceModel>
177 (
178  const alphaField& alpha,
179  const rhoField& rho,
180  const volVectorField& U,
181  const surfaceScalarField& alphaRhoPhi,
182  const surfaceScalarField& phi,
183  const transportModel& transport,
184  const word& propertiesName,
185  const word& type
186 )
187 :
189  (
190  type,
191  alpha,
192  rho,
193  U,
194  alphaRhoPhi,
195  phi,
196  transport,
197  propertiesName
198  ),
199 
200  Cmu_
201  (
202  IOobject
203  (
204  "Cmu",
205  this->runTime_.timeName(),
206  this->mesh_
207  ),
208  this->mesh_,
210  (
211  "Cmu",
212  this->coeffDict_,
213  0.09
214  )
215  ),
216  C1_
217  (
218  IOobject
219  (
220  "C1",
221  this->runTime_.timeName(),
222  this->mesh_
223  ),
224  this->mesh_,
226  (
227  "C1",
228  this->coeffDict_,
229  1.44
230  )
231  ),
232  C2_
233  (
234  IOobject
235  (
236  "C2",
237  this->runTime_.timeName(),
238  this->mesh_
239  ),
240  this->mesh_,
242  (
243  "C2",
244  this->coeffDict_,
245  1.92
246  )
247  ),
248  sigmak_
249  (
250  IOobject
251  (
252  "sigmak",
253  this->runTime_.timeName(),
254  this->mesh_
255  ),
256  this->mesh_,
258  (
259  "sigmak",
260  this->coeffDict_,
261  1.0
262  )
263  ),
264  sigmaEps_
265  (
266  IOobject
267  (
268  "sigmaEps",
269  this->runTime_.timeName(),
270  this->mesh_
271  ),
272  this->mesh_,
274  (
275  "sigmaEps",
276  this->coeffDict_,
277  1.3
278  )
279  ),
280 
281  CdSigma_
282  (
283  IOobject
284  (
285  "CdSigma",
286  this->runTime_.timeName(),
287  this->mesh_
288  ),
289  this->mesh_,
291  ),
292  betap_
293  (
294  IOobject
295  (
296  "betap",
297  this->runTime_.timeName(),
298  this->mesh_
299  ),
300  this->mesh_,
302  ),
303  betad_
304  (
305  IOobject
306  (
307  "betad",
308  this->runTime_.timeName(),
309  this->mesh_
310  ),
311  this->mesh_,
313  ),
314  C4_
315  (
316  IOobject
317  (
318  "C4",
319  this->runTime_.timeName(),
320  this->mesh_
321  ),
322  this->mesh_,
324  ),
325  C5_
326  (
327  IOobject
328  (
329  "C5",
330  this->runTime_.timeName(),
331  this->mesh_
332  ),
333  this->mesh_,
335  ),
336 
337  k_
338  (
339  IOobject
340  (
341  IOobject::groupName("k", alphaRhoPhi.group()),
342  this->runTime_.timeName(),
343  this->mesh_,
346  ),
347  this->mesh_
348  ),
349  epsilon_
350  (
351  IOobject
352  (
353  IOobject::groupName("epsilon", alphaRhoPhi.group()),
354  this->runTime_.timeName(),
355  this->mesh_,
358  ),
359  this->mesh_
360  )
361 {
362  bound(k_, this->kMin_);
363  bound(epsilon_, this->epsilonMin_);
364 
365  if (type == typeName)
366  {
367  this->printCoeffs(type);
368  }
369 
370  setPorosityCoefficients();
371 }
372 
373 
374 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
375 
376 template<class BasicTurbulenceModel>
378 {
380  {
381  return true;
382  }
383  else
384  {
385  return false;
386  }
387 }
388 
389 
390 template<class BasicTurbulenceModel>
392 {
393  if (!this->turbulence_)
394  {
395  return;
396  }
397 
398  // Local references
399  const alphaField& alpha = this->alpha_;
400  const rhoField& rho = this->rho_;
401  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
402  const volVectorField& U = this->U_;
403  volScalarField& nut = this->nut_;
404  fv::options& fvOptions(fv::options::New(this->mesh_));
405 
407 
409  (
410  fvc::div(fvc::absolute(this->phi(), U))().v()
411  );
412 
413  tmp<volTensorField> tgradU = fvc::grad(U);
415  (
416  this->GName(),
417  nut.v()*(dev(twoSymm(tgradU().v())) && tgradU().v())
418  );
419  tgradU.clear();
420 
421  // Update epsilon and G at the wall
422  epsilon_.boundaryFieldRef().updateCoeffs();
423 
424  volScalarField::Internal magU(mag(U));
425  volScalarField::Internal magU3(pow3(magU));
426 
427  // Dissipation equation
428  tmp<fvScalarMatrix> epsEqn
429  (
430  fvm::ddt(alpha, rho, epsilon_)
431  + fvm::div(alphaRhoPhi, epsilon_)
432  - fvm::laplacian(alpha*rho*DepsilonEff(), epsilon_)
433  ==
434  C1_*alpha()*rho()*G*epsilon_()/k_()
435  - fvm::SuSp(((2.0/3.0)*C1_)*alpha()*rho()*divU, epsilon_)
436  - fvm::Sp(C2_*alpha()*rho()*epsilon_()/k_(), epsilon_)
437  + epsilonSource(magU, magU3)
438  + fvOptions(alpha, rho, epsilon_)
439  );
440 
441  epsEqn.ref().relax();
442  fvOptions.constrain(epsEqn.ref());
443  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
444  solve(epsEqn);
445  fvOptions.correct(epsilon_);
446  bound(epsilon_, this->epsilonMin_);
447 
448  // Turbulent kinetic energy equation
450  (
451  fvm::ddt(alpha, rho, k_)
452  + fvm::div(alphaRhoPhi, k_)
453  - fvm::laplacian(alpha*rho*DkEff(), k_)
454  ==
455  alpha()*rho()*G
456  - fvm::SuSp((2.0/3.0)*alpha()*rho()*divU, k_)
457  - fvm::Sp(alpha()*rho()*epsilon_()/k_(), k_)
458  + kSource(magU, magU3)
459  + fvOptions(alpha, rho, k_)
460  );
461 
462  kEqn.ref().relax();
463  fvOptions.constrain(kEqn.ref());
464  solve(kEqn);
465  fvOptions.correct(k_);
466  bound(k_, this->kMin_);
467 
468  correctNut();
469 }
470 
471 
472 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
473 
474 } // End namespace RASModels
475 } // End namespace Foam
476 
477 // ************************************************************************* //
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
void setPorosityCoefficient(volScalarField::Internal &C, const porosityModels::powerLawLopesdaCosta &pm)
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:438
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
tmp< fvMatrix< Type > > SuSp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
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 clear() const
If object pointer points to valid object:
Definition: tmpI.H:230
fv::options & fvOptions
surfaceScalarField & phi
Variant of the power law porosity model with spatially varying drag coefficient.
const dictionary & dict() const
Return dictionary used for model construction.
const scalarField & Sigma() const
Return the porosity surface area per unit volume zone field.
void correct(GeometricField< Type, fvPatchField, volMesh > &field)
Apply correction to field.
const dimensionedScalar G
Newtonian constant of gravitation.
tmp< fvMatrix< Type > > Su(const DimensionedField< Type, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
dimensionedSymmTensor sqr(const dimensionedVector &dv)
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
BasicTurbulenceModel::rhoField rhoField
virtual tmp< fvScalarMatrix > kSource(const volScalarField::Internal &magU, const volScalarField::Internal &magU3) const
Generic dimensioned Type class.
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
virtual bool read()
Re-read model coefficients if they have changed.
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:49
Finite-volume options.
Definition: fvOptions.H:52
DimensionedField< scalar, volMesh > Internal
Type of the internal field from which this GeometricField is derived.
const labelList & cellZoneIDs() const
Return const access to the cell zone IDs.
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
Templated abstract base class for RAS turbulence models.
Definition: RASModel.H:49
const porosityModel & model() const
virtual tmp< fvScalarMatrix > epsilonSource(const volScalarField::Internal &magU, const volScalarField::Internal &magU3) const
bool read(const char *, int32_t &)
Definition: int32IO.C:85
const Internal & v() const
Return a const-reference to the dimensioned internal field.
Definition: volFieldsI.H:31
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
const cellShapeList & cells
kEpsilonLopesdaCosta(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.
A class for handling words, derived from string.
Definition: word.H:59
static word groupName(Name name, const word &group)
optionList(const fvMesh &mesh)
Construct null.
Definition: fvOptionList.C:98
void setCdSigma(volScalarField::Internal &C, const porosityModels::powerLawLopesdaCosta &pm)
BasicTurbulenceModel::transportModel transportModel
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if successful.
Definition: doubleScalar.H:68
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
Bound the given scalar field if it has gone unbounded.
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
dimensionedScalar pow3(const dimensionedScalar &ds)
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:188
U
Definition: pEqn.H:72
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:33
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:46
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 > &)
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
zeroField divU
Definition: alphaSuSp.H:3
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
BasicTurbulenceModel::alphaField alphaField
scalar nut
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:583