realizableKE.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2019-2023 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "realizableKE.H"
30 #include "fvOptions.H"
31 #include "bound.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace RASModels
38 {
39 
40 // * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * * //
41 
42 template<class BasicTurbulenceModel>
44 (
45  const volTensorField& gradU,
46  const volScalarField& S2,
47  const volScalarField& magS
48 )
49 {
50  tmp<volSymmTensorField> tS = devSymm(gradU);
51  const volSymmTensorField& S = tS();
52 
54  (
55  (2*sqrt(2.0))*((S&S)&&S)
56  /(
57  magS*S2
58  + dimensionedScalar("small", dimensionSet(0, 0, -3, 0, 0), SMALL)
59  )
60  );
61 
62  tS.clear();
63 
65  (
66  (1.0/3.0)*acos(clamp(sqrt(6.0)*W, scalarMinMax(-1, 1)))
67  );
68  volScalarField As(sqrt(6.0)*cos(phis));
69  volScalarField Us(sqrt(S2/2.0 + magSqr(skew(gradU))));
70 
71  return 1.0/(A0_ + As*Us*k_/epsilon_);
72 }
73 
74 
75 template<class BasicTurbulenceModel>
77 (
78  const volTensorField& gradU,
79  const volScalarField& S2,
80  const volScalarField& magS
81 )
82 {
83  this->nut_ = rCmu(gradU, S2, magS)*sqr(k_)/epsilon_;
84  this->nut_.correctBoundaryConditions();
85  fv::options::New(this->mesh_).correct(this->nut_);
86 
87  BasicTurbulenceModel::correctNut();
88 }
89 
90 
91 template<class BasicTurbulenceModel>
93 {
94  tmp<volTensorField> tgradU = fvc::grad(this->U_);
95  volScalarField S2(2*magSqr(devSymm(tgradU())));
96  volScalarField magS(sqrt(S2));
97  correctNut(tgradU(), S2, magS);
98 }
99 
100 
101 template<class BasicTurbulenceModel>
103 {
105  (
106  k_,
107  dimVolume*this->rho_.dimensions()*k_.dimensions()/dimTime
108  );
109 }
110 
111 
112 template<class BasicTurbulenceModel>
114 {
116  (
117  epsilon_,
118  dimVolume*this->rho_.dimensions()*epsilon_.dimensions()/dimTime
119  );
120 }
121 
122 
123 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
124 
125 template<class BasicTurbulenceModel>
127 (
128  const alphaField& alpha,
129  const rhoField& rho,
130  const volVectorField& U,
131  const surfaceScalarField& alphaRhoPhi,
132  const surfaceScalarField& phi,
133  const transportModel& transport,
134  const word& propertiesName,
135  const word& type
136 )
137 :
138  eddyViscosity<RASModel<BasicTurbulenceModel>>
139  (
140  type,
141  alpha,
142  rho,
143  U,
144  alphaRhoPhi,
145  phi,
146  transport,
147  propertiesName
148  ),
149  A0_
150  (
151  dimensioned<scalar>::getOrAddToDict
152  (
153  "A0",
154  this->coeffDict_,
155  4.0
156  )
157  ),
158  C2_
159  (
160  dimensioned<scalar>::getOrAddToDict
161  (
162  "C2",
163  this->coeffDict_,
164  1.9
165  )
166  ),
167  sigmak_
168  (
169  dimensioned<scalar>::getOrAddToDict
170  (
171  "sigmak",
172  this->coeffDict_,
173  1.0
174  )
175  ),
176  sigmaEps_
177  (
178  dimensioned<scalar>::getOrAddToDict
179  (
180  "sigmaEps",
181  this->coeffDict_,
182  1.2
183  )
184  ),
185 
186  k_
187  (
188  IOobject
189  (
190  IOobject::groupName("k", U.group()),
191  this->runTime_.timeName(),
192  this->mesh_,
193  IOobject::MUST_READ,
194  IOobject::AUTO_WRITE
195  ),
196  this->mesh_
197  ),
198  epsilon_
199  (
200  IOobject
201  (
202  IOobject::groupName("epsilon", U.group()),
203  this->runTime_.timeName(),
204  this->mesh_,
205  IOobject::MUST_READ,
206  IOobject::AUTO_WRITE
207  ),
208  this->mesh_
209  )
210 {
211  bound(k_, this->kMin_);
212  bound(epsilon_, this->epsilonMin_);
213 
214  if (type == typeName)
215  {
216  this->printCoeffs(type);
217  }
218 }
219 
220 
221 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
222 
223 template<class BasicTurbulenceModel>
225 {
227  {
228  A0_.readIfPresent(this->coeffDict());
229  C2_.readIfPresent(this->coeffDict());
230  sigmak_.readIfPresent(this->coeffDict());
231  sigmaEps_.readIfPresent(this->coeffDict());
232 
233  return true;
234  }
235 
236  return false;
237 }
238 
239 
240 template<class BasicTurbulenceModel>
242 {
243  if (!this->turbulence_)
244  {
245  return;
246  }
247 
248  // Local references
249  const alphaField& alpha = this->alpha_;
250  const rhoField& rho = this->rho_;
251  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
252  const volVectorField& U = this->U_;
253  volScalarField& nut = this->nut_;
254  fv::options& fvOptions(fv::options::New(this->mesh_));
255 
256  eddyViscosity<RASModel<BasicTurbulenceModel>>::correct();
257 
259 
260  tmp<volTensorField> tgradU = fvc::grad(U);
261  volScalarField S2(2*magSqr(devSymm(tgradU())));
262  volScalarField magS(sqrt(S2));
263 
264  volScalarField eta(magS*k_/epsilon_);
265  volScalarField C1(max(eta/(scalar(5) + eta), scalar(0.43)));
266 
267  volScalarField G(this->GName(), nut*(tgradU() && devTwoSymm(tgradU())));
268 
269  // Update epsilon and G at the wall
270  epsilon_.boundaryFieldRef().updateCoeffs();
271  // Push any changed cell values to coupled neighbours
272  epsilon_.boundaryFieldRef().template evaluateCoupled<coupledFvPatch>();
273 
274  // SAF: limiting thermo->nu(). If psiThermo is used rho might be < 0
275  // temporarily when p < 0 then nu < 0 which needs limiting
276  volScalarField nuLimited
277  (
278  max
279  (
280  this->nu(),
281  dimensionedScalar(this->nu()().dimensions(), Zero)
282  )
283  );
284 
285  // Dissipation equation
286  tmp<fvScalarMatrix> epsEqn
287  (
288  fvm::ddt(alpha, rho, epsilon_)
289  + fvm::div(alphaRhoPhi, epsilon_)
290  - fvm::laplacian(alpha*rho*DepsilonEff(), epsilon_)
291  ==
292  C1*alpha*rho*magS*epsilon_
293  - fvm::Sp
294  (
295  C2_*alpha*rho*epsilon_/(k_ + sqrt(nuLimited*epsilon_)),
296  epsilon_
297  )
298  + epsilonSource()
299  + fvOptions(alpha, rho, epsilon_)
300  );
301 
302  epsEqn.ref().relax();
303  fvOptions.constrain(epsEqn.ref());
304  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
305  solve(epsEqn);
306  fvOptions.correct(epsilon_);
307  bound(epsilon_, this->epsilonMin_);
308 
309 
310  // Turbulent kinetic energy equation
311 
312  tmp<fvScalarMatrix> kEqn
313  (
314  fvm::ddt(alpha, rho, k_)
315  + fvm::div(alphaRhoPhi, k_)
316  - fvm::laplacian(alpha*rho*DkEff(), k_)
317  ==
318  alpha*rho*G
319  - fvm::SuSp(2.0/3.0*alpha*rho*divU, k_)
320  - fvm::Sp(alpha*rho*epsilon_/k_, k_)
321  + kSource()
322  + fvOptions(alpha, rho, k_)
323  );
324 
325  kEqn.ref().relax();
326  fvOptions.constrain(kEqn.ref());
327  solve(kEqn);
328  fvOptions.correct(k_);
329  bound(k_, this->kMin_);
330 
331  correctNut(tgradU(), S2, magS);
332 }
333 
334 
335 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
336 
337 } // End namespace RASModels
338 } // End namespace Foam
339 
340 // ************************************************************************* //
void correct(GeometricField< Type, PatchField, GeoMesh > &field)
Apply correction to field.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:47
dimensionedScalar acos(const dimensionedScalar &ds)
MinMax< scalar > scalarMinMax
A scalar min/max range.
Definition: MinMax.H:97
dimensionedScalar kMin_
Lower limit of k.
Definition: RASModel.H:109
dimensionedTensor skew(const dimensionedTensor &dt)
const dimensionedScalar G
Newtonian constant of gravitation.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
dimensionedSymmTensor sqr(const dimensionedVector &dv)
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:42
dimensionedScalar sqrt(const dimensionedScalar &ds)
SymmTensor< Cmpt > devSymm(const SymmTensor< Cmpt > &st)
Return the deviatoric part of the symmetric part of a SymmTensor.
Definition: SymmTensorI.H:481
Generic dimensioned Type class.
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:51
virtual tmp< fvScalarMatrix > kSource() const
Definition: realizableKE.C:95
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:76
virtual bool read()
Re-read model coefficients if they have changed.
Definition: realizableKE.C:217
fv::options & fvOptions
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:127
Templated abstract base class for RAS turbulence models.
Definition: RASModel.H:77
tmp< volScalarField > rCmu(const volTensorField &gradU, const volScalarField &S2, const volScalarField &magS)
Definition: realizableKE.C:37
constexpr const char *const group
Group name for atomic constants.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:72
word timeName
Definition: getTimeIndex.H:3
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Definition: dimensionSet.H:105
dimensionedScalar cos(const dimensionedScalar &ds)
A class for handling words, derived from Foam::string.
Definition: word.H:63
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:41
Us
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
zeroField SuSp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
virtual void printCoeffs(const word &type)
Print model coefficients.
Definition: RASModel.C:27
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
SolverPerformance< Type > solve(faMatrix< Type > &, const dictionary &solverControls)
Solve returning the solution statistics given convergence tolerance.
edgeScalarField phis(IOobject("phis", runTime.timeName(), aMesh.thisDb(), IOobject::NO_READ, IOobject::NO_WRITE), linearEdgeInterpolate(Us) &aMesh.Le())
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:41
Bound the given scalar field if it has gone unbounded.
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(psi *p+alphal *rhol0+((alphav *psiv+alphal *psil) - psi) *pSat, rhoMin);# 1 "/home/chef2/andy/OpenFOAM/release/v2406/OpenFOAM-v2406/applications/solvers/multiphase/cavitatingFoam/alphavPsi.H" 1{ alphav=clamp((rho - rholSat)/(rhovSat - rholSat), zero_one{});alphal=1.0 - alphav;Info<< "max-min alphav: "<< max(alphav).value()<< " "<< min(alphav).value()<< endl;psiModel-> correct()
Definition: pEqn.H:63
virtual tmp< fvScalarMatrix > epsilonSource() const
Definition: realizableKE.C:106
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:183
U
Definition: pEqn.H:72
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:29
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:41
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: realizableKE.C:234
Base-class for all transport models used by the incompressible turbulence models. ...
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
realizableKE(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName=turbulenceModel::propertiesName, const word &type=typeName)
Construct from components.
Definition: realizableKE.C:120
void clear() const noexcept
If object pointer points to valid object: delete object and set pointer to nullptr.
Definition: tmpI.H:289
zeroField divU
Definition: alphaSuSp.H:3
A class for managing temporary objects.
Definition: HashPtrTable.H:50
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present.
Definition: fvOptions.C:96
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:180
volScalarField & nu
dimensionedScalar epsilonMin_
Lower limit of epsilon.
Definition: RASModel.H:114
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
scalar nut
Namespace for OpenFOAM.
dimensionSet clamp(const dimensionSet &a, const dimensionSet &range)
Definition: dimensionSet.C:271
SymmTensor< Cmpt > devTwoSymm(const SymmTensor< Cmpt > &st)
Return the deviatoric part of twice the symmetric part of a SymmTensor.
Definition: SymmTensorI.H:491
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127