kOmega2006.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) 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 "kOmega2006.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 template<class BasicMomentumTransportModel>
40 (
41  const volTensorField& gradU
42 )
43 {
44  this->nut_ = k_/max(omega_, Clim_*sqrt(2/betaStar_)*mag(dev(symm(gradU))));
45  this->nut_.correctBoundaryConditions();
46  fvConstraints::New(this->mesh_).constrain(this->nut_);
47 }
48 
49 
50 template<class BasicMomentumTransportModel>
51 tmp<volScalarField::Internal> kOmega2006<BasicMomentumTransportModel>::beta
52 (
53  const volTensorField& gradU
54 ) const
55 {
56  const volSymmTensorField::Internal S(symm(gradU()));
57  const volSymmTensorField::Internal Shat(S - 0.5*tr(S)*I);
58  const volTensorField::Internal Omega(skew(gradU.v()));
59 
60  const volScalarField::Internal ChiOmega
61  (
62  modelName("ChiOmega"),
63  mag((Omega & Omega) && Shat)/pow3(betaStar_*omega_.v())
64  );
65 
66  const volScalarField::Internal fBeta
67  (
68  modelName("fBeta"),
69  (1 + 85*ChiOmega)/(1 + 100*ChiOmega)
70  );
71 
72  return beta0_*fBeta;
73 }
74 
75 
76 template<class BasicMomentumTransportModel>
77 tmp<volScalarField::Internal>
78 kOmega2006<BasicMomentumTransportModel>::CDkOmega() const
79 {
80  return max
81  (
82  sigmaDo_*(fvc::grad(k_)().v() & fvc::grad(omega_)().v())/omega_(),
84  );
85 }
86 
87 
88 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
89 
90 template<class BasicMomentumTransportModel>
92 {
93  correctNut(fvc::grad(this->U_));
94 }
95 
96 
97 template<class BasicMomentumTransportModel>
99 {
100  return tmp<fvScalarMatrix>
101  (
102  new fvScalarMatrix
103  (
104  k_,
105  dimVolume*this->rho_.dimensions()*k_.dimensions()
106  /dimTime
107  )
108  );
109 }
110 
111 
112 template<class BasicMomentumTransportModel>
114 {
115  return tmp<fvScalarMatrix>
116  (
117  new fvScalarMatrix
118  (
119  omega_,
120  dimVolume*this->rho_.dimensions()*omega_.dimensions()/dimTime
121  )
122  );
123 }
124 
125 
126 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
127 
128 template<class BasicMomentumTransportModel>
130 (
131  const alphaField& alpha,
132  const rhoField& rho,
133  const volVectorField& U,
134  const surfaceScalarField& alphaRhoPhi,
135  const surfaceScalarField& phi,
136  const viscosity& viscosity,
137  const word& type
138 )
139 :
141  (
142  type,
143  alpha,
144  rho,
145  U,
146  alphaRhoPhi,
147  phi,
148  viscosity
149  ),
150 
151  betaStar_
152  (
154  (
155  "betaStar",
156  this->coeffDict_,
157  0.09
158  )
159  ),
160  beta0_
161  (
163  (
164  "beta0",
165  this->coeffDict_,
166  0.0708
167  )
168  ),
169  gamma_
170  (
172  (
173  "gamma",
174  this->coeffDict_,
175  0.52
176  )
177  ),
178  Clim_
179  (
181  (
182  "Clim",
183  this->coeffDict_,
184  0.875
185  )
186  ),
187  sigmaDo_
188  (
190  (
191  "sigmaDo",
192  this->coeffDict_,
193  0.125
194  )
195  ),
196  alphaK_
197  (
199  (
200  "alphaK",
201  this->coeffDict_,
202  0.6
203  )
204  ),
205  alphaOmega_
206  (
208  (
209  "alphaOmega",
210  this->coeffDict_,
211  0.5
212  )
213  ),
214 
215  k_
216  (
217  IOobject
218  (
219  IOobject::groupName("k", alphaRhoPhi.group()),
220  this->runTime_.timeName(),
221  this->mesh_,
224  ),
225  this->mesh_
226  ),
227  omega_
228  (
229  IOobject
230  (
231  IOobject::groupName("omega", alphaRhoPhi.group()),
232  this->runTime_.timeName(),
233  this->mesh_,
236  ),
237  this->mesh_
238  )
239 {
240  bound(k_, this->kMin_);
241  bound(omega_, this->omegaMin_);
242 
243  if (type == typeName)
244  {
245  this->printCoeffs(type);
246  }
247 }
248 
249 
250 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
251 
252 template<class BasicMomentumTransportModel>
254 {
256  {
257  betaStar_.readIfPresent(this->coeffDict());
258  beta0_.readIfPresent(this->coeffDict());
259  gamma_.readIfPresent(this->coeffDict());
260  sigmaDo_.readIfPresent(this->coeffDict());
261  alphaK_.readIfPresent(this->coeffDict());
262  alphaOmega_.readIfPresent(this->coeffDict());
263 
264  return true;
265  }
266  else
267  {
268  return false;
269  }
270 }
271 
272 
273 template<class BasicMomentumTransportModel>
275 {
276  if (!this->turbulence_)
277  {
278  return;
279  }
280 
281  // Local references
282  const alphaField& alpha = this->alpha_;
283  const rhoField& rho = this->rho_;
284  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
285  const volVectorField& U = this->U_;
286  volScalarField& nut = this->nut_;
287  const Foam::fvModels& fvModels(Foam::fvModels::New(this->mesh_));
289  (
290  Foam::fvConstraints::New(this->mesh_)
291  );
292 
294 
296  (
297  modelName("divU"),
298  fvc::div(fvc::absolute(this->phi(), U))().v()
299  );
300 
301  const volTensorField gradU(fvc::grad(U));
303  (
304  this->GName(),
305  nut.v()*(dev(twoSymm(gradU.v())) && gradU.v())
306  );
307 
308  // Update omega and G at the wall
309  omega_.boundaryFieldRef().updateCoeffs();
310 
311  // Turbulence specific dissipation rate equation
312  tmp<fvScalarMatrix> omegaEqn
313  (
314  fvm::ddt(alpha, rho, omega_)
315  + fvm::div(alphaRhoPhi, omega_)
316  - fvm::laplacian(alpha*rho*DomegaEff(), omega_)
317  ==
318  gamma_*alpha()*rho()*G*omega_()/k_()
319  - fvm::SuSp(((2.0/3.0)*gamma_)*alpha()*rho()*divU, omega_)
320  - fvm::Sp(beta(gradU)*alpha()*rho()*omega_(), omega_)
321  + alpha()*rho()*CDkOmega()
322  + omegaSource()
323  + fvModels.source(alpha, rho, omega_)
324  );
325 
326  omegaEqn.ref().relax();
327  fvConstraints.constrain(omegaEqn.ref());
328  omegaEqn.ref().boundaryManipulate(omega_.boundaryFieldRef());
329  solve(omegaEqn);
330  fvConstraints.constrain(omega_);
331  bound(omega_, this->omegaMin_);
332 
333 
334  // Turbulent kinetic energy equation
336  (
337  fvm::ddt(alpha, rho, k_)
338  + fvm::div(alphaRhoPhi, k_)
339  - fvm::laplacian(alpha*rho*DkEff(), k_)
340  ==
341  alpha()*rho()*G
342  - fvm::SuSp((2.0/3.0)*alpha()*rho()*divU, k_)
343  - fvm::Sp(betaStar_*alpha()*rho()*omega_(), k_)
344  + kSource()
345  + fvModels.source(alpha, rho, k_)
346  );
347 
348  kEqn.ref().relax();
349  fvConstraints.constrain(kEqn.ref());
350  solve(kEqn);
351  fvConstraints.constrain(k_);
352  bound(k_, this->kMin_);
353 
354  correctNut(gradU);
355 }
356 
357 
358 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
359 
360 } // End namespace RASModels
361 } // End namespace Foam
362 
363 // ************************************************************************* //
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
virtual tmp< fvScalarMatrix > kSource() const
Definition: kOmega2006.C:98
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 > &)
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
virtual bool read()
Read RASProperties dictionary.
Definition: kOmega2006.C:253
dimensionedTensor skew(const dimensionedTensor &dt)
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
GeometricField< tensor, fvPatchField, volMesh > volTensorField
Definition: volFieldsFwd.H:63
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)
const dimensionSet dimless
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
DimensionedField< symmTensor, volMesh > Internal
Type of the internal field from which this GeometricField is derived.
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
const Internal & v() const
Return a const-reference to the dimensioned internal field.
Definition: volFieldsI.H:31
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
static const Identity< scalar > I
Definition: Identity.H:93
Foam::fvConstraints & fvConstraints
BasicMomentumTransportModel::alphaField alphaField
Definition: kOmega2006.H:122
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
phi
Definition: correctPhi.H:3
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
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
virtual void correctNut()
Definition: kOmega2006.C:91
tmp< volScalarField > divU
Definition: alphaSuSp.H:9
bool constrain(fvMatrix< Type > &eqn) const
Apply constraints to an equation.
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: kOmega2006.C:274
Foam::fvModels & fvModels
dimensionedScalar pow3(const dimensionedScalar &ds)
virtual tmp< fvScalarMatrix > omegaSource() const
Definition: kOmega2006.C:113
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
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...
BasicMomentumTransportModel::rhoField rhoField
Definition: kOmega2006.H:123
const dimensionSet dimVolume
SolverPerformance< Type > solve(fvMatrix< Type > &, const word &)
Solve returning the solution statistics given convergence tolerance.
dimensioned< scalar > mag(const dimensioned< Type > &)
A class for managing temporary objects.
Definition: PtrList.H:53
kOmega2006(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: kOmega2006.C:130
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:98
scalar nut
Namespace for OpenFOAM.