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-2024 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>("Cmu", this->coeffDict(), 0.09)
210  ),
211  C1_
212  (
213  IOobject
214  (
215  "C1",
216  this->runTime_.name(),
217  this->mesh_
218  ),
219  this->mesh_,
220  dimensioned<scalar>("C1", this->coeffDict(), 1.44)
221  ),
222  C2_
223  (
224  IOobject
225  (
226  "C2",
227  this->runTime_.name(),
228  this->mesh_
229  ),
230  this->mesh_,
231  dimensioned<scalar>("C2", this->coeffDict(), 1.92)
232  ),
233  sigmak_
234  (
235  IOobject
236  (
237  "sigmak",
238  this->runTime_.name(),
239  this->mesh_
240  ),
241  this->mesh_,
242  dimensioned<scalar>("sigmak", this->coeffDict(), 1.0)
243  ),
244  sigmaEps_
245  (
246  IOobject
247  (
248  "sigmaEps",
249  this->runTime_.name(),
250  this->mesh_
251  ),
252  this->mesh_,
253  dimensioned<scalar>("sigmaEps", this->coeffDict(), 1.3)
254  ),
255 
256  CdAv_
257  (
258  IOobject
259  (
260  "CdAv",
261  this->runTime_.name(),
262  this->mesh_
263  ),
264  this->mesh_,
266  ),
267  betap_
268  (
269  IOobject
270  (
271  "betap",
272  this->runTime_.name(),
273  this->mesh_
274  ),
275  this->mesh_,
277  ),
278  betad_
279  (
280  IOobject
281  (
282  "betad",
283  this->runTime_.name(),
284  this->mesh_
285  ),
286  this->mesh_,
288  ),
289  C4_
290  (
291  IOobject
292  (
293  "C4",
294  this->runTime_.name(),
295  this->mesh_
296  ),
297  this->mesh_,
299  ),
300  C5_
301  (
302  IOobject
303  (
304  "C5",
305  this->runTime_.name(),
306  this->mesh_
307  ),
308  this->mesh_,
310  ),
311 
312  k_
313  (
314  IOobject
315  (
316  IOobject::groupName("k", alphaRhoPhi.group()),
317  this->runTime_.name(),
318  this->mesh_,
319  IOobject::MUST_READ,
320  IOobject::AUTO_WRITE
321  ),
322  this->mesh_
323  ),
324  epsilon_
325  (
326  IOobject
327  (
328  IOobject::groupName("epsilon", alphaRhoPhi.group()),
329  this->runTime_.name(),
330  this->mesh_,
331  IOobject::MUST_READ,
332  IOobject::AUTO_WRITE
333  ),
334  this->mesh_
335  )
336 {
337  bound(k_, this->kMin_);
338  boundEpsilon();
340 }
341 
342 
343 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
344 
345 template<class BasicMomentumTransportModel>
347 {
349  {
350  return true;
351  }
352  else
353  {
354  return false;
355  }
356 }
357 
358 
359 template<class BasicMomentumTransportModel>
361 {
362  if (!this->turbulence_)
363  {
364  return;
365  }
366 
367  // Local references
368  const alphaField& alpha = this->alpha_;
369  const rhoField& rho = this->rho_;
370  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
371  const volVectorField& U = this->U_;
372  volScalarField& nut = this->nut_;
373  const Foam::fvModels& fvModels(Foam::fvModels::New(this->mesh_));
375  (
376  Foam::fvConstraints::New(this->mesh_)
377  );
378 
380 
382  (
383  fvc::div(fvc::absolute(this->phi(), U))().v()
384  );
385 
386  tmp<volTensorField> tgradU = fvc::grad(U);
388  (
389  this->GName(),
390  nut.v()*(dev(twoSymm(tgradU().v())) && tgradU().v())
391  );
392  tgradU.clear();
393 
394  // Update epsilon and G at the wall
395  epsilon_.boundaryFieldRef().updateCoeffs();
396 
398  volScalarField::Internal magU3(pow3(magU));
399 
400  // Dissipation equation
401  tmp<fvScalarMatrix> epsEqn
402  (
403  fvm::ddt(alpha, rho, epsilon_)
404  + fvm::div(alphaRhoPhi, epsilon_)
405  - fvm::laplacian(alpha*rho*DepsilonEff(), epsilon_)
406  ==
407  C1_*alpha()*rho()*G*epsilon_()/k_()
408  - fvm::SuSp(((2.0/3.0)*C1_)*alpha()*rho()*divU, epsilon_)
409  - fvm::Sp(C2_*alpha()*rho()*epsilon_()/k_(), epsilon_)
410  + epsilonSource(magU, magU3)
411  + fvModels.source(alpha, rho, epsilon_)
412  );
413 
414  epsEqn.ref().relax();
415  fvConstraints.constrain(epsEqn.ref());
416  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
417  solve(epsEqn);
418  fvConstraints.constrain(epsilon_);
419  boundEpsilon();
420 
421  // Turbulent kinetic energy equation
423  (
424  fvm::ddt(alpha, rho, k_)
425  + fvm::div(alphaRhoPhi, k_)
426  - fvm::laplacian(alpha*rho*DkEff(), k_)
427  ==
428  alpha()*rho()*G
429  - fvm::SuSp((2.0/3.0)*alpha()*rho()*divU, k_)
430  - fvm::Sp(alpha()*rho()*epsilon_()/k_(), k_)
431  + kSource(magU, magU3)
432  + fvModels.source(alpha, rho, k_)
433  );
434 
435  kEqn.ref().relax();
436  fvConstraints.constrain(kEqn.ref());
437  solve(kEqn);
439  bound(k_, this->kMin_);
440 
441  correctNut();
442 }
443 
444 
445 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
446 
447 } // End namespace RASModels
448 } // End namespace Foam
449 
450 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
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:740
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:539
static autoPtr< dictionary > New(Istream &)
Construct top-level dictionary on freestore from Istream.
Definition: dictionaryIO.C:103
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 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:253
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:197
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 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.
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
void dev(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
void twoSymm(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
const dimensionSet dimless
const dimensionSet dimLength
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
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
void pow3(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
void sqr(LagrangianPatchField< typename outerProduct< Type, Type >::type > &f, const LagrangianPatchField< Type > &f1)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
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.