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-2023 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 "fvModels.H"
28 #include "fvConstraints.H"
29 #include "porosityForce.H"
30 #include "bound.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace RASModels
37 {
38 
39 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
40 
41 template<class BasicMomentumTransportModel>
43 (
46 )
47 {
48  if (pm.dict().found(C.name()))
49  {
50  const scalar Cpm = pm.dict().lookup<scalar>(C.name());
51 
52  const labelList& cells = this->mesh_.cellZones()[pm.zoneName()];
53 
54  forAll(cells, i)
55  {
56  const label celli = cells[i];
57  C[celli] = Cpm;
58  }
59  }
60 }
61 
62 
63 template<class BasicMomentumTransportModel>
65 (
68 )
69 {
70  if (pm.dict().found(C.name()))
71  {
72  const scalarField& Av = pm.Av();
73 
74  const scalar Cpm = pm.dict().lookup<scalar>(C.name());
75 
76  const labelList& cells = this->mesh_.cellZones()[pm.zoneName()];
77 
78  forAll(cells, i)
79  {
80  const label celli = cells[i];
81  C[celli] = Cpm*Av[celli];
82  }
83  }
84 }
85 
86 
87 template<class BasicMomentumTransportModel>
90 {
92  (
93  Foam::fvModels::New(this->mesh_)
94  );
95 
96  forAll(fvModels, i)
97  {
98  if (isA<fv::porosityForce>(fvModels[i]))
99  {
100  const fv::porosityForce& eps =
101  refCast<const fv::porosityForce>(fvModels[i]);
102 
103  if (isA<porosityModels::powerLawLopesdaCosta>(eps.model()))
104  {
106  refCast<const porosityModels::powerLawLopesdaCosta>
107  (
108  eps.model()
109  );
110 
111  setPorosityCoefficient(Cmu_, pm);
112  setPorosityCoefficient(C1_, pm);
113  setPorosityCoefficient(C2_, pm);
114  setPorosityCoefficient(sigmak_, pm);
115  setPorosityCoefficient(sigmaEps_, pm);
116 
117  setCdAv(CdAv_, pm);
118  setPorosityCoefficient(betap_, pm);
119  setPorosityCoefficient(betad_, pm);
120  setPorosityCoefficient(C4_, pm);
121  setPorosityCoefficient(C5_, pm);
122  }
123  }
124  }
125 }
126 
127 
128 template<class BasicMomentumTransportModel>
131 {
132  tmp<volScalarField> tCmuk2(Cmu_*sqr(k_));
133  epsilon_ = max(epsilon_, tCmuk2()/(this->nutMaxCoeff_*this->nu()));
134  return tCmuk2;
135 }
136 
137 
138 template<class BasicMomentumTransportModel>
140 {
141  this->nut_ = boundEpsilon()/epsilon_;
142  this->nut_.correctBoundaryConditions();
143  fvConstraints::New(this->mesh_).constrain(this->nut_);
144 }
145 
146 
147 template<class BasicMomentumTransportModel>
149 (
150  const volScalarField::Internal& magU,
151  const volScalarField::Internal& magU3
152 ) const
153 {
154  return fvm::Su(CdAv_*(betap_*magU3 - betad_*magU*k_()), k_);
155 }
156 
157 
158 template<class BasicMomentumTransportModel>
161 (
162  const volScalarField::Internal& magU,
163  const volScalarField::Internal& magU3
164 ) const
165 {
166  return fvm::Su
167  (
168  CdAv_
169  *(C4_*betap_*epsilon_()/k_()*magU3 - C5_*betad_*magU*epsilon_()),
170  epsilon_
171  );
172 }
173 
174 
175 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
176 
177 template<class BasicMomentumTransportModel>
179 (
180  const alphaField& alpha,
181  const rhoField& rho,
182  const volVectorField& U,
183  const surfaceScalarField& alphaRhoPhi,
184  const surfaceScalarField& phi,
185  const viscosity& viscosity,
186  const word& type
187 )
188 :
189  eddyViscosity<RASModel<BasicMomentumTransportModel>>
190  (
191  type,
192  alpha,
193  rho,
194  U,
195  alphaRhoPhi,
196  phi,
197  viscosity
198  ),
199 
200  Cmu_
201  (
202  IOobject
203  (
204  "Cmu",
205  this->runTime_.name(),
206  this->mesh_
207  ),
208  this->mesh_,
209  dimensioned<scalar>::lookupOrAddToDict
210  (
211  "Cmu",
212  this->coeffDict_,
213  0.09
214  )
215  ),
216  C1_
217  (
218  IOobject
219  (
220  "C1",
221  this->runTime_.name(),
222  this->mesh_
223  ),
224  this->mesh_,
225  dimensioned<scalar>::lookupOrAddToDict
226  (
227  "C1",
228  this->coeffDict_,
229  1.44
230  )
231  ),
232  C2_
233  (
234  IOobject
235  (
236  "C2",
237  this->runTime_.name(),
238  this->mesh_
239  ),
240  this->mesh_,
241  dimensioned<scalar>::lookupOrAddToDict
242  (
243  "C2",
244  this->coeffDict_,
245  1.92
246  )
247  ),
248  sigmak_
249  (
250  IOobject
251  (
252  "sigmak",
253  this->runTime_.name(),
254  this->mesh_
255  ),
256  this->mesh_,
257  dimensioned<scalar>::lookupOrAddToDict
258  (
259  "sigmak",
260  this->coeffDict_,
261  1.0
262  )
263  ),
264  sigmaEps_
265  (
266  IOobject
267  (
268  "sigmaEps",
269  this->runTime_.name(),
270  this->mesh_
271  ),
272  this->mesh_,
273  dimensioned<scalar>::lookupOrAddToDict
274  (
275  "sigmaEps",
276  this->coeffDict_,
277  1.3
278  )
279  ),
280 
281  CdAv_
282  (
283  IOobject
284  (
285  "CdAv",
286  this->runTime_.name(),
287  this->mesh_
288  ),
289  this->mesh_,
291  ),
292  betap_
293  (
294  IOobject
295  (
296  "betap",
297  this->runTime_.name(),
298  this->mesh_
299  ),
300  this->mesh_,
302  ),
303  betad_
304  (
305  IOobject
306  (
307  "betad",
308  this->runTime_.name(),
309  this->mesh_
310  ),
311  this->mesh_,
313  ),
314  C4_
315  (
316  IOobject
317  (
318  "C4",
319  this->runTime_.name(),
320  this->mesh_
321  ),
322  this->mesh_,
324  ),
325  C5_
326  (
327  IOobject
328  (
329  "C5",
330  this->runTime_.name(),
331  this->mesh_
332  ),
333  this->mesh_,
335  ),
336 
337  k_
338  (
339  IOobject
340  (
341  IOobject::groupName("k", alphaRhoPhi.group()),
342  this->runTime_.name(),
343  this->mesh_,
344  IOobject::MUST_READ,
345  IOobject::AUTO_WRITE
346  ),
347  this->mesh_
348  ),
349  epsilon_
350  (
351  IOobject
352  (
353  IOobject::groupName("epsilon", alphaRhoPhi.group()),
354  this->runTime_.name(),
355  this->mesh_,
356  IOobject::MUST_READ,
357  IOobject::AUTO_WRITE
358  ),
359  this->mesh_
360  )
361 {
362  bound(k_, this->kMin_);
363  boundEpsilon();
364 
365  if (type == typeName)
366  {
367  this->printCoeffs(type);
368  }
369 
371 }
372 
373 
374 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
375 
376 template<class BasicMomentumTransportModel>
378 {
380  {
381  return true;
382  }
383  else
384  {
385  return false;
386  }
387 }
388 
389 
390 template<class BasicMomentumTransportModel>
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  const Foam::fvModels& fvModels(Foam::fvModels::New(this->mesh_));
406  (
407  Foam::fvConstraints::New(this->mesh_)
408  );
409 
411 
413  (
414  fvc::div(fvc::absolute(this->phi(), U))().v()
415  );
416 
417  tmp<volTensorField> tgradU = fvc::grad(U);
419  (
420  this->GName(),
421  nut.v()*(dev(twoSymm(tgradU().v())) && tgradU().v())
422  );
423  tgradU.clear();
424 
425  // Update epsilon and G at the wall
426  epsilon_.boundaryFieldRef().updateCoeffs();
427 
429  volScalarField::Internal magU3(pow3(magU));
430 
431  // Dissipation equation
432  tmp<fvScalarMatrix> epsEqn
433  (
434  fvm::ddt(alpha, rho, epsilon_)
435  + fvm::div(alphaRhoPhi, epsilon_)
436  - fvm::laplacian(alpha*rho*DepsilonEff(), epsilon_)
437  ==
438  C1_*alpha()*rho()*G*epsilon_()/k_()
439  - fvm::SuSp(((2.0/3.0)*C1_)*alpha()*rho()*divU, epsilon_)
440  - fvm::Sp(C2_*alpha()*rho()*epsilon_()/k_(), epsilon_)
441  + epsilonSource(magU, magU3)
442  + fvModels.source(alpha, rho, epsilon_)
443  );
444 
445  epsEqn.ref().relax();
446  fvConstraints.constrain(epsEqn.ref());
447  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
448  solve(epsEqn);
449  fvConstraints.constrain(epsilon_);
450  boundEpsilon();
451 
452  // Turbulent kinetic energy equation
454  (
455  fvm::ddt(alpha, rho, k_)
456  + fvm::div(alphaRhoPhi, k_)
457  - fvm::laplacian(alpha*rho*DkEff(), k_)
458  ==
459  alpha()*rho()*G
460  - fvm::SuSp((2.0/3.0)*alpha()*rho()*divU, k_)
461  - fvm::Sp(alpha()*rho()*epsilon_()/k_(), k_)
462  + kSource(magU, magU3)
463  + fvModels.source(alpha, rho, k_)
464  );
465 
466  kEqn.ref().relax();
467  fvConstraints.constrain(kEqn.ref());
468  solve(kEqn);
470  bound(k_, this->kMin_);
471 
472  correctNut();
473 }
474 
475 
476 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
477 
478 } // End namespace RASModels
479 } // End namespace Foam
480 
481 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Bound the given scalar field where it is below the specified minimum.
Graphite solid properties.
Definition: C.H:51
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic GeometricField class.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
Templated abstract base class for RAS turbulence models.
Definition: RASModel.H:56
BasicMomentumTransportModel::alphaField alphaField
virtual tmp< fvScalarMatrix > kSource(const volScalarField::Internal &magU, const volScalarField::Internal &magU3) const
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
tmp< volScalarField > boundEpsilon()
Bound epsilon and return Cmu*sqr(k) for nut.
kEpsilonLopesdaCosta(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity, const word &type=typeName)
Construct from components.
void setCdAv(volScalarField::Internal &C, const porosityModels::powerLawLopesdaCosta &pm)
virtual void correctNut()
Correct the eddy-viscosity nut.
virtual tmp< fvScalarMatrix > epsilonSource(const volScalarField::Internal &magU, const volScalarField::Internal &magU3) const
void setPorosityCoefficient(volScalarField::Internal &C, const porosityModels::powerLawLopesdaCosta &pm)
virtual bool read()
Re-read model coefficients if they have changed.
BasicMomentumTransportModel::rhoField rhoField
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:710
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:509
static autoPtr< dictionary > New(Istream &)
Construct top-level dictionary on freestore from Istream.
Definition: dictionaryIO.C:100
Generic dimensioned Type class.
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:52
Finite volume constraints.
Definition: fvConstraints.H:67
bool constrain(fvMatrix< Type > &eqn) const
Apply constraints to an equation.
Finite volume models.
Definition: fvModels.H:65
tmp< fvMatrix< Type > > source(const VolField< Type > &field) const
Return source for an equation.
This model applies the force exerted on the fluid by a porous media.
Definition: porosityForce.H:86
const porosityModel & model() const
Return the porosity model.
const word & zoneName() const
Return const access to the cell zone name.
const dictionary & dict() const
Return dictionary used for model construction.
const scalarField & Av() const
Return the porosity surface area per unit volume zone field.
Variant of the power law porosity model with spatially varying drag coefficient.
A class for managing temporary objects.
Definition: tmp.H:55
void clear() const
If object pointer points to valid object:
Definition: tmpI.H:237
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
Abstract base class for all fluid physical properties.
Definition: viscosity.H:50
A class for handling words, derived from string.
Definition: word.H:62
Foam::fvConstraints & fvConstraints(Foam::fvConstraints::New(mesh))
Foam::fvModels & fvModels(Foam::fvModels::New(mesh))
const scalar nut
const cellShapeList & cells
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
void correct(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiCorr, const SpType &Sp, const SuType &Su)
const char *const group
Group name for atomic constants.
const dimensionedScalar G
Newtonian constant of gravitation.
tmp< VolField< typename outerProduct< vector, Type >::type > > grad(const SurfaceField< Type > &ssf)
Definition: fvcGrad.C:46
tmp< VolField< Type > > div(const SurfaceField< Type > &ssf)
Definition: fvcDiv.C:47
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:202
tmp< fvMatrix< Type > > Su(const DimensionedField< Type, volMesh > &, const VolField< Type > &)
tmp< fvMatrix< Type > > laplacian(const VolField< Type > &vf, const word &name)
Definition: fvmLaplacian.C:47
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const VolField< Type > &vf, const word &name)
Definition: fvmDiv.C:48
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const VolField< Type > &)
tmp< fvMatrix< Type > > SuSp(const volScalarField::Internal &, const VolField< Type > &)
tmp< fvMatrix< Type > > ddt(const VolField< Type > &vf)
Definition: fvmDdt.C:46
Namespace for OpenFOAM.
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
bool read(const char *, int32_t &)
Definition: int32IO.C:85
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
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar pow3(const dimensionedScalar &ds)
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
const dimensionSet dimless
const dimensionSet dimLength
dimensioned< scalar > mag(const dimensioned< Type > &)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
bool bound(volScalarField &, const dimensionedScalar &min)
Bound the given scalar field where it is below the specified min value.
Definition: bound.C:31
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.