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-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 "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 BasicMomentumTransportModel>
42 (
45 )
46 {
47  if (pm.dict().found(C.name()))
48  {
49  const labelList& cellZoneIDs = pm.cellZoneIDs();
50 
51  const scalar Cpm = pm.dict().lookup<scalar>(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 BasicMomentumTransportModel>
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 = pm.dict().lookup<scalar>(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 BasicMomentumTransportModel>
100 {
102 
103  forAll(fvOptions, i)
104  {
105  if (isA<fv::explicitPorositySource>(fvOptions[i]))
106  {
107  const fv::explicitPorositySource& eps =
108  refCast<const fv::explicitPorositySource>(fvOptions[i]);
109 
110  if (isA<porosityModels::powerLawLopesdaCosta>(eps.model()))
111  {
113  refCast<const porosityModels::powerLawLopesdaCosta>
114  (
115  eps.model()
116  );
117 
118  setPorosityCoefficient(Cmu_, pm);
119  setPorosityCoefficient(C1_, pm);
120  setPorosityCoefficient(C2_, pm);
121  setPorosityCoefficient(sigmak_, pm);
122  setPorosityCoefficient(sigmaEps_, pm);
123 
124  setCdSigma(CdSigma_, pm);
125  setPorosityCoefficient(betap_, pm);
126  setPorosityCoefficient(betad_, pm);
127  setPorosityCoefficient(C4_, pm);
128  setPorosityCoefficient(C5_, pm);
129  }
130  }
131  }
132 }
133 
134 
135 template<class BasicMomentumTransportModel>
137 {
138  this->nut_ = Cmu_*sqr(k_)/epsilon_;
139  this->nut_.correctBoundaryConditions();
140  fv::options::New(this->mesh_).correct(this->nut_);
141 }
142 
143 
144 template<class BasicMomentumTransportModel>
146 (
147  const volScalarField::Internal& magU,
148  const volScalarField::Internal& magU3
149 ) const
150 {
151  return fvm::Su(CdSigma_*(betap_*magU3 - betad_*magU*k_()), k_);
152 }
153 
154 
155 template<class BasicMomentumTransportModel>
158 (
159  const volScalarField::Internal& magU,
160  const volScalarField::Internal& magU3
161 ) const
162 {
163  return fvm::Su
164  (
165  CdSigma_
166  *(C4_*betap_*epsilon_()/k_()*magU3 - C5_*betad_*magU*epsilon_()),
167  epsilon_
168  );
169 }
170 
171 
172 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
173 
174 template<class BasicMomentumTransportModel>
176 (
177  const alphaField& alpha,
178  const rhoField& rho,
179  const volVectorField& U,
180  const surfaceScalarField& alphaRhoPhi,
181  const surfaceScalarField& phi,
182  const transportModel& transport,
183  const word& type
184 )
185 :
187  (
188  type,
189  alpha,
190  rho,
191  U,
192  alphaRhoPhi,
193  phi,
194  transport
195  ),
196 
197  Cmu_
198  (
199  IOobject
200  (
201  "Cmu",
202  this->runTime_.timeName(),
203  this->mesh_
204  ),
205  this->mesh_,
207  (
208  "Cmu",
209  this->coeffDict_,
210  0.09
211  )
212  ),
213  C1_
214  (
215  IOobject
216  (
217  "C1",
218  this->runTime_.timeName(),
219  this->mesh_
220  ),
221  this->mesh_,
223  (
224  "C1",
225  this->coeffDict_,
226  1.44
227  )
228  ),
229  C2_
230  (
231  IOobject
232  (
233  "C2",
234  this->runTime_.timeName(),
235  this->mesh_
236  ),
237  this->mesh_,
239  (
240  "C2",
241  this->coeffDict_,
242  1.92
243  )
244  ),
245  sigmak_
246  (
247  IOobject
248  (
249  "sigmak",
250  this->runTime_.timeName(),
251  this->mesh_
252  ),
253  this->mesh_,
255  (
256  "sigmak",
257  this->coeffDict_,
258  1.0
259  )
260  ),
261  sigmaEps_
262  (
263  IOobject
264  (
265  "sigmaEps",
266  this->runTime_.timeName(),
267  this->mesh_
268  ),
269  this->mesh_,
271  (
272  "sigmaEps",
273  this->coeffDict_,
274  1.3
275  )
276  ),
277 
278  CdSigma_
279  (
280  IOobject
281  (
282  "CdSigma",
283  this->runTime_.timeName(),
284  this->mesh_
285  ),
286  this->mesh_,
288  ),
289  betap_
290  (
291  IOobject
292  (
293  "betap",
294  this->runTime_.timeName(),
295  this->mesh_
296  ),
297  this->mesh_,
299  ),
300  betad_
301  (
302  IOobject
303  (
304  "betad",
305  this->runTime_.timeName(),
306  this->mesh_
307  ),
308  this->mesh_,
310  ),
311  C4_
312  (
313  IOobject
314  (
315  "C4",
316  this->runTime_.timeName(),
317  this->mesh_
318  ),
319  this->mesh_,
321  ),
322  C5_
323  (
324  IOobject
325  (
326  "C5",
327  this->runTime_.timeName(),
328  this->mesh_
329  ),
330  this->mesh_,
332  ),
333 
334  k_
335  (
336  IOobject
337  (
338  IOobject::groupName("k", alphaRhoPhi.group()),
339  this->runTime_.timeName(),
340  this->mesh_,
343  ),
344  this->mesh_
345  ),
346  epsilon_
347  (
348  IOobject
349  (
350  IOobject::groupName("epsilon", alphaRhoPhi.group()),
351  this->runTime_.timeName(),
352  this->mesh_,
355  ),
356  this->mesh_
357  )
358 {
359  bound(k_, this->kMin_);
360  bound(epsilon_, this->epsilonMin_);
361 
362  if (type == typeName)
363  {
364  this->printCoeffs(type);
365  }
366 
367  setPorosityCoefficients();
368 }
369 
370 
371 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
372 
373 template<class BasicMomentumTransportModel>
375 {
377  {
378  return true;
379  }
380  else
381  {
382  return false;
383  }
384 }
385 
386 
387 template<class BasicMomentumTransportModel>
389 {
390  if (!this->turbulence_)
391  {
392  return;
393  }
394 
395  // Local references
396  const alphaField& alpha = this->alpha_;
397  const rhoField& rho = this->rho_;
398  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
399  const volVectorField& U = this->U_;
400  volScalarField& nut = this->nut_;
401  fv::options& fvOptions(fv::options::New(this->mesh_));
402 
404 
406  (
407  fvc::div(fvc::absolute(this->phi(), U))().v()
408  );
409 
410  tmp<volTensorField> tgradU = fvc::grad(U);
412  (
413  this->GName(),
414  nut.v()*(dev(twoSymm(tgradU().v())) && tgradU().v())
415  );
416  tgradU.clear();
417 
418  // Update epsilon and G at the wall
419  epsilon_.boundaryFieldRef().updateCoeffs();
420 
421  volScalarField::Internal magU(mag(U));
422  volScalarField::Internal magU3(pow3(magU));
423 
424  // Dissipation equation
425  tmp<fvScalarMatrix> epsEqn
426  (
427  fvm::ddt(alpha, rho, epsilon_)
428  + fvm::div(alphaRhoPhi, epsilon_)
429  - fvm::laplacian(alpha*rho*DepsilonEff(), epsilon_)
430  ==
431  C1_*alpha()*rho()*G*epsilon_()/k_()
432  - fvm::SuSp(((2.0/3.0)*C1_)*alpha()*rho()*divU, epsilon_)
433  - fvm::Sp(C2_*alpha()*rho()*epsilon_()/k_(), epsilon_)
434  + epsilonSource(magU, magU3)
435  + fvOptions(alpha, rho, epsilon_)
436  );
437 
438  epsEqn.ref().relax();
439  fvOptions.constrain(epsEqn.ref());
440  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
441  solve(epsEqn);
442  fvOptions.correct(epsilon_);
443  bound(epsilon_, this->epsilonMin_);
444 
445  // Turbulent kinetic energy equation
447  (
448  fvm::ddt(alpha, rho, k_)
449  + fvm::div(alphaRhoPhi, k_)
450  - fvm::laplacian(alpha*rho*DkEff(), k_)
451  ==
452  alpha()*rho()*G
453  - fvm::SuSp((2.0/3.0)*alpha()*rho()*divU, k_)
454  - fvm::Sp(alpha()*rho()*epsilon_()/k_(), k_)
455  + kSource(magU, magU3)
456  + fvOptions(alpha, rho, k_)
457  );
458 
459  kEqn.ref().relax();
460  fvOptions.constrain(kEqn.ref());
461  solve(kEqn);
462  fvOptions.correct(k_);
463  bound(k_, this->kMin_);
464 
465  correctNut();
466 }
467 
468 
469 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
470 
471 } // End namespace RASModels
472 } // End namespace Foam
473 
474 // ************************************************************************* //
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:667
#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
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.
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
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)
phi
Definition: pEqn.H:104
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
BasicMomentumTransportModel::transportModel transportModel
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)
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
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-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
BasicMomentumTransportModel::rhoField rhoField
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
kEpsilonLopesdaCosta(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.
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
const dimensionedScalar & G
Newtonian constant of gravitation.
BasicMomentumTransportModel::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:812