realizableKE.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-2021 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 "realizableKE.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  const volTensorField& gradU,
44  const volScalarField& S2,
45  const volScalarField& magS
46 )
47 {
48  tmp<volSymmTensorField> tS = dev(symm(gradU));
49  const volSymmTensorField& S = tS();
50 
51  const volScalarField W
52  (
53  (2*sqrt(2.0))*((S&S)&&S)
54  /(
55  magS*S2
56  + dimensionedScalar(dimensionSet(0, 0, -3, 0, 0), small)
57  )
58  );
59 
60  tS.clear();
61 
62  const volScalarField phis
63  (
64  (1.0/3.0)*acos(min(max(sqrt(6.0)*W, -scalar(1)), scalar(1)))
65  );
66  const volScalarField As(sqrt(6.0)*cos(phis));
67  const volScalarField Us(sqrt(S2/2.0 + magSqr(skew(gradU))));
68 
69  return 1.0/(A0_ + As*Us*k_/epsilon_);
70 }
71 
72 
73 template<class BasicMomentumTransportModel>
75 (
76  const volTensorField& gradU,
77  const volScalarField& S2,
78  const volScalarField& magS
79 )
80 {
81  this->nut_ = rCmu(gradU, S2, magS)*sqr(k_)/epsilon_;
82  this->nut_.correctBoundaryConditions();
83  fvConstraints::New(this->mesh_).constrain(this->nut_);
84 }
85 
86 
87 template<class BasicMomentumTransportModel>
89 {
90  const volTensorField gradU(fvc::grad(this->U_));
91  const volScalarField S2(modelName("S2"), 2*magSqr(dev(symm(gradU))));
92  const volScalarField magS(modelName("magS"), sqrt(S2));
93 
94  correctNut(gradU, S2, magS);
95 }
96 
97 
98 template<class BasicMomentumTransportModel>
100 {
101  return tmp<fvScalarMatrix>
102  (
103  new fvScalarMatrix
104  (
105  k_,
106  dimVolume*this->rho_.dimensions()*k_.dimensions()
107  /dimTime
108  )
109  );
110 }
111 
112 
113 template<class BasicMomentumTransportModel>
116 {
117  return tmp<fvScalarMatrix>
118  (
119  new fvScalarMatrix
120  (
121  epsilon_,
122  dimVolume*this->rho_.dimensions()*epsilon_.dimensions()
123  /dimTime
124  )
125  );
126 }
127 
128 
129 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
130 
131 template<class BasicMomentumTransportModel>
133 (
134  const alphaField& alpha,
135  const rhoField& rho,
136  const volVectorField& U,
137  const surfaceScalarField& alphaRhoPhi,
138  const surfaceScalarField& phi,
139  const viscosity& viscosity,
140  const word& type
141 )
142 :
144  (
145  type,
146  alpha,
147  rho,
148  U,
149  alphaRhoPhi,
150  phi,
151  viscosity
152  ),
153  A0_
154  (
156  (
157  "A0",
158  this->coeffDict_,
159  4.0
160  )
161  ),
162  C2_
163  (
165  (
166  "C2",
167  this->coeffDict_,
168  1.9
169  )
170  ),
171  sigmak_
172  (
174  (
175  "sigmak",
176  this->coeffDict_,
177  1.0
178  )
179  ),
180  sigmaEps_
181  (
183  (
184  "sigmaEps",
185  this->coeffDict_,
186  1.2
187  )
188  ),
189 
190  k_
191  (
192  IOobject
193  (
194  IOobject::groupName("k", alphaRhoPhi.group()),
195  this->runTime_.timeName(),
196  this->mesh_,
199  ),
200  this->mesh_
201  ),
202  epsilon_
203  (
204  IOobject
205  (
206  IOobject::groupName("epsilon", alphaRhoPhi.group()),
207  this->runTime_.timeName(),
208  this->mesh_,
211  ),
212  this->mesh_
213  )
214 {
215  bound(k_, this->kMin_);
216  bound(epsilon_, this->epsilonMin_);
217 
218  if (type == typeName)
219  {
220  this->printCoeffs(type);
221  }
222 }
223 
224 
225 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
226 
227 template<class BasicMomentumTransportModel>
229 {
231  {
232  A0_.readIfPresent(this->coeffDict());
233  C2_.readIfPresent(this->coeffDict());
234  sigmak_.readIfPresent(this->coeffDict());
235  sigmaEps_.readIfPresent(this->coeffDict());
236 
237  return true;
238  }
239  else
240  {
241  return false;
242  }
243 }
244 
245 
246 template<class BasicMomentumTransportModel>
248 {
249  if (!this->turbulence_)
250  {
251  return;
252  }
253 
254  // Local references
255  const alphaField& alpha = this->alpha_;
256  const rhoField& rho = this->rho_;
257  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
258  const volVectorField& U = this->U_;
259  volScalarField& nut = this->nut_;
260  const Foam::fvModels& fvModels(Foam::fvModels::New(this->mesh_));
262  (
263  Foam::fvConstraints::New(this->mesh_)
264  );
265 
267 
269  (
270  modelName("divU"),
271  fvc::div(fvc::absolute(this->phi(), U))()
272  );
273 
274  const volTensorField gradU(fvc::grad(U));
275  const volScalarField S2(modelName("S2"), 2*magSqr(dev(symm(gradU))));
276  const volScalarField magS(modelName("magS"), sqrt(S2));
277 
278  const volScalarField::Internal eta
279  (
280  modelName("eta"), magS()*k_()/epsilon_()
281  );
282  const volScalarField::Internal C1
283  (
284  modelName("C1"),
285  max(eta/(scalar(5) + eta), scalar(0.43))
286  );
287 
289  (
290  this->GName(),
291  nut*(gradU.v() && dev(twoSymm(gradU.v())))
292  );
293 
294  // Update epsilon and G at the wall
295  epsilon_.boundaryFieldRef().updateCoeffs();
296 
297  // Dissipation equation
298  tmp<fvScalarMatrix> epsEqn
299  (
300  fvm::ddt(alpha, rho, epsilon_)
301  + fvm::div(alphaRhoPhi, epsilon_)
302  - fvm::laplacian(alpha*rho*DepsilonEff(), epsilon_)
303  ==
304  C1*alpha()*rho()*magS()*epsilon_()
305  - fvm::Sp
306  (
307  C2_*alpha()*rho()*epsilon_()/(k_() + sqrt(this->nu()()*epsilon_())),
308  epsilon_
309  )
310  + epsilonSource()
311  + fvModels.source(alpha, rho, epsilon_)
312  );
313 
314  epsEqn.ref().relax();
315  fvConstraints.constrain(epsEqn.ref());
316  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
317  solve(epsEqn);
318  fvConstraints.constrain(epsilon_);
319  bound(epsilon_, this->epsilonMin_);
320 
321 
322  // Turbulent kinetic energy equation
323 
325  (
326  fvm::ddt(alpha, rho, k_)
327  + fvm::div(alphaRhoPhi, k_)
328  - fvm::laplacian(alpha*rho*DkEff(), k_)
329  ==
330  alpha()*rho()*G
331  - fvm::SuSp(2.0/3.0*alpha()*rho()*divU, k_)
332  - fvm::Sp(alpha()*rho()*epsilon_()/k_(), k_)
333  + kSource()
334  + fvModels.source(alpha, rho, k_)
335  );
336 
337  kEqn.ref().relax();
338  fvConstraints.constrain(kEqn.ref());
339  solve(kEqn);
340  fvConstraints.constrain(k_);
341  bound(k_, this->kMin_);
342 
343  correctNut(gradU, S2, magS);
344 }
345 
346 
347 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
348 
349 } // End namespace RASModels
350 } // End namespace Foam
351 
352 // ************************************************************************* //
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
dimensionedScalar acos(const dimensionedScalar &ds)
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
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
dimensionedTensor skew(const dimensionedTensor &dt)
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))
dimensionedScalar sqrt(const dimensionedScalar &ds)
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
virtual tmp< fvScalarMatrix > kSource() const
Definition: realizableKE.C:99
BasicMomentumTransportModel::alphaField alphaField
Definition: RASModel.H:96
virtual bool read()
Re-read model coefficients if they have changed.
Definition: realizableKE.C:228
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
Templated abstract base class for RAS turbulence models.
Definition: RASModel.H:53
tmp< volScalarField > rCmu(const volTensorField &gradU, const volScalarField &S2, const volScalarField &magS)
Definition: realizableKE.C:42
const dimensionSet dimTime
Dimension set for the base types.
Definition: dimensionSet.H:121
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)
dimensionedScalar cos(const dimensionedScalar &ds)
Foam::fvConstraints & fvConstraints
A class for handling words, derived from string.
Definition: word.H:59
tmp< fvMatrix< Type > > S(const Pair< tmp< volScalarField::Internal >> &, const GeometricField< Type, fvPatchField, volMesh > &)
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
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
phi
Definition: correctPhi.H:3
dimensioned< scalar > magSqr(const dimensioned< Type > &)
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
virtual tmp< fvScalarMatrix > epsilonSource() const
Definition: realizableKE.C:115
bool constrain(fvMatrix< Type > &eqn) const
Apply constraints to an equation.
Foam::fvModels & fvModels
BasicMomentumTransportModel::rhoField rhoField
Definition: RASModel.H:97
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
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
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
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: realizableKE.C:247
const scalarList W(::W(thermo))
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
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
realizableKE(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: realizableKE.C:133
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:98
scalar nut
Namespace for OpenFOAM.