kEpsilon.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) 2011-2022 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 "kEpsilon.H"
27 #include "fvModels.H"
28 #include "fvConstraints.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 {
43  this->nut_ = Cmu_*sqr(k_)/epsilon_;
44  this->nut_.correctBoundaryConditions();
45  fvConstraints::New(this->mesh_).constrain(this->nut_);
46 }
47 
48 
49 template<class BasicMomentumTransportModel>
51 {
52  return tmp<fvScalarMatrix>
53  (
54  new fvScalarMatrix
55  (
56  k_,
57  dimVolume*this->rho_.dimensions()*k_.dimensions()
58  /dimTime
59  )
60  );
61 }
62 
63 
64 template<class BasicMomentumTransportModel>
66 {
67  return tmp<fvScalarMatrix>
68  (
69  new fvScalarMatrix
70  (
71  epsilon_,
72  dimVolume*this->rho_.dimensions()*epsilon_.dimensions()
73  /dimTime
74  )
75  );
76 }
77 
78 
79 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
80 
81 template<class BasicMomentumTransportModel>
83 (
84  const alphaField& alpha,
85  const rhoField& rho,
86  const volVectorField& U,
87  const surfaceScalarField& alphaRhoPhi,
88  const surfaceScalarField& phi,
89  const viscosity& viscosity,
90  const word& type
91 )
92 :
94  (
95  type,
96  alpha,
97  rho,
98  U,
99  alphaRhoPhi,
100  phi,
101  viscosity
102  ),
103 
104  Cmu_
105  (
107  (
108  "Cmu",
109  this->coeffDict_,
110  0.09
111  )
112  ),
113  C1_
114  (
116  (
117  "C1",
118  this->coeffDict_,
119  1.44
120  )
121  ),
122  C2_
123  (
125  (
126  "C2",
127  this->coeffDict_,
128  1.92
129  )
130  ),
131  C3_
132  (
134  (
135  "C3",
136  this->coeffDict_,
137  0
138  )
139  ),
140  sigmak_
141  (
143  (
144  "sigmak",
145  this->coeffDict_,
146  1.0
147  )
148  ),
149  sigmaEps_
150  (
152  (
153  "sigmaEps",
154  this->coeffDict_,
155  1.3
156  )
157  ),
158 
159  k_
160  (
161  IOobject
162  (
163  IOobject::groupName("k", alphaRhoPhi.group()),
164  this->runTime_.timeName(),
165  this->mesh_,
168  ),
169  this->mesh_
170  ),
171  epsilon_
172  (
173  IOobject
174  (
175  IOobject::groupName("epsilon", alphaRhoPhi.group()),
176  this->runTime_.timeName(),
177  this->mesh_,
180  ),
181  this->mesh_
182  )
183 {
184  bound(k_, this->kMin_);
185  bound(epsilon_, this->epsilonMin_);
186 
187  if (type == typeName)
188  {
189  this->printCoeffs(type);
190  }
191 }
192 
193 
194 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
195 
196 template<class BasicMomentumTransportModel>
198 {
200  {
201  Cmu_.readIfPresent(this->coeffDict());
202  C1_.readIfPresent(this->coeffDict());
203  C2_.readIfPresent(this->coeffDict());
204  C3_.readIfPresent(this->coeffDict());
205  sigmak_.readIfPresent(this->coeffDict());
206  sigmaEps_.readIfPresent(this->coeffDict());
207 
208  return true;
209  }
210  else
211  {
212  return false;
213  }
214 }
215 
216 
217 template<class BasicMomentumTransportModel>
219 {
220  if (!this->turbulence_)
221  {
222  return;
223  }
224 
225  // Local references
226  const alphaField& alpha = this->alpha_;
227  const rhoField& rho = this->rho_;
228  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
229  const volVectorField& U = this->U_;
230  volScalarField& nut = this->nut_;
231  const Foam::fvModels& fvModels(Foam::fvModels::New(this->mesh_));
233  (
234  Foam::fvConstraints::New(this->mesh_)
235  );
236 
238 
240  (
241  fvc::div(fvc::absolute(this->phi(), U))()
242  );
243 
244  tmp<volTensorField> tgradU = fvc::grad(U);
246  (
247  this->GName(),
248  nut()*(dev(twoSymm(tgradU().v())) && tgradU().v())
249  );
250  tgradU.clear();
251 
252  // Update epsilon and G at the wall
253  epsilon_.boundaryFieldRef().updateCoeffs();
254 
255  // Dissipation equation
256  tmp<fvScalarMatrix> epsEqn
257  (
258  fvm::ddt(alpha, rho, epsilon_)
259  + fvm::div(alphaRhoPhi, epsilon_)
260  - fvm::laplacian(alpha*rho*DepsilonEff(), epsilon_)
261  ==
262  C1_*alpha()*rho()*G*epsilon_()/k_()
263  - fvm::SuSp(((2.0/3.0)*C1_ - C3_)*alpha()*rho()*divU, epsilon_)
264  - fvm::Sp(C2_*alpha()*rho()*epsilon_()/k_(), epsilon_)
265  + epsilonSource()
266  + fvModels.source(alpha, rho, epsilon_)
267  );
268 
269  epsEqn.ref().relax();
270  fvConstraints.constrain(epsEqn.ref());
271  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
272  solve(epsEqn);
273  fvConstraints.constrain(epsilon_);
274  bound(epsilon_, this->epsilonMin_);
275 
276  // Turbulent kinetic energy equation
278  (
279  fvm::ddt(alpha, rho, k_)
280  + fvm::div(alphaRhoPhi, k_)
281  - fvm::laplacian(alpha*rho*DkEff(), k_)
282  ==
283  alpha()*rho()*G
284  - fvm::SuSp((2.0/3.0)*alpha()*rho()*divU, k_)
285  - fvm::Sp(alpha()*rho()*epsilon_()/k_(), k_)
286  + kSource()
287  + fvModels.source(alpha, rho, k_)
288  );
289 
290  kEqn.ref().relax();
291  fvConstraints.constrain(kEqn.ref());
292  solve(kEqn);
293  fvConstraints.constrain(k_);
294  bound(k_, this->kMin_);
295 
296  correctNut();
297 }
298 
299 
300 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
301 
302 } // End namespace RASModels
303 } // End namespace Foam
304 
305 // ************************************************************************* //
virtual tmp< fvScalarMatrix > epsilonSource() const
Definition: kEpsilon.C:65
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
tmp< fvMatrix< Type > > SuSp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
void clear() const
If object pointer points to valid object:
Definition: tmpI.H:237
U
Definition: pEqn.H:72
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(rho0+psi *p, rhoMin);# 1 "/home/ubuntu/OpenFOAM-10/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
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
dimensionedSymmTensor sqr(const dimensionedVector &dv)
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
virtual tmp< fvScalarMatrix > kSource() const
Definition: kEpsilon.C:50
Generic dimensioned Type class.
const dimensionedScalar G
Newtonian constant of gravitation.
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:49
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
Templated abstract base class for RAS turbulence models.
Definition: RASModel.H:53
const dimensionSet dimTime
bool read(const char *, int32_t &)
Definition: int32IO.C:85
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
Foam::fvConstraints & fvConstraints
A class for handling words, derived from string.
Definition: word.H:59
static word groupName(Name name, const word &group)
Finite volume models.
Definition: fvModels.H:60
Info<< "Reading field p\"<< endl;volScalarField p(IOobject("p", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Reading field U\"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);pressureReference pressureReference(p, simple.dict());mesh.schemes().setFluxRequired(p.name());Info<< "Reading field pa\"<< endl;volScalarField pa(IOobject("pa", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Reading field Ua\"<< endl;volVectorField Ua(IOobject("Ua", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);# 65 "/home/ubuntu/OpenFOAM-10/applications/solvers/incompressible/adjointShapeOptimisationFoam/createFields.H" 2label paRefCell=0;scalar paRefValue=0.0;setRefCell(pa, simple.dict(), paRefCell, paRefValue);mesh.schemes().setFluxRequired(pa.name());autoPtr< viscosityModel > viscosity(viscosityModel::New(mesh))
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
phi
Definition: correctPhi.H:3
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.
static autoPtr< dictionary > New(Istream &)
Construct top-level dictionary on freestore from Istream.
Definition: dictionaryIO.C:96
tmp< volScalarField > divU
Definition: alphaSuSp.H:9
bool constrain(fvMatrix< Type > &eqn) const
Apply constraints to an equation.
Foam::fvModels & fvModels
kEpsilon(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.
Definition: kEpsilon.C:83
virtual bool read()
Re-read model coefficients if they have changed.
Definition: kEpsilon.C:197
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:202
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:33
Finite volume constraints.
Definition: fvConstraints.H:57
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
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const dimensionSet dimVolume
virtual void correctNut()
Definition: kEpsilon.C:41
SolverPerformance< Type > solve(fvMatrix< Type > &, const word &)
Solve returning the solution statistics given convergence tolerance.
A class for managing temporary objects.
Definition: PtrList.H:53
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: kEpsilon.C:218
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:98
scalar nut
Namespace for OpenFOAM.