qZeta.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-2023 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 "qZeta.H"
27 #include "fvcMagSqrGradGrad.H"
28 #include "bound.H"
30 
32 (
33  geometricOneField,
34  geometricOneField,
35  incompressibleMomentumTransportModel
36 )
37 
38 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
39 
40 namespace Foam
41 {
42 namespace incompressible
43 {
44 namespace RASModels
45 {
46 
47 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
48 
49 defineTypeNameAndDebug(qZeta, 0);
51 (
52  RASincompressibleMomentumTransportModel,
53  qZeta,
54  dictionary
55 );
56 
57 // * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * * //
58 
59 tmp<volScalarField> qZeta::fMu() const
60 {
61  const volScalarField Rt(q_*k_/(2.0*nu()*zeta_));
62 
63  if (anisotropic_)
64  {
65  return exp((-scalar(2.5) + Rt/20.0)/pow3(scalar(1) + Rt/130.0));
66  }
67  else
68  {
69  return
70  exp(-6.0/sqr(scalar(1) + Rt/50.0))
71  *(scalar(1) + 3.0*exp(-Rt/10.0));
72  }
73 }
74 
75 
76 tmp<volScalarField> qZeta::f2() const
77 {
78  tmp<volScalarField> Rt = q_*k_/(2.0*nu()*zeta_);
79  return scalar(1) - 0.3*exp(-sqr(Rt));
80 }
81 
82 
83 void qZeta::correctNut()
84 {
85  nut_ = Cmu_*fMu()*sqr(k_)/epsilon_;
87 }
88 
89 
90 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
91 
93 (
94  const geometricOneField& alpha,
95  const geometricOneField& rho,
96  const volVectorField& U,
97  const surfaceScalarField& alphaRhoPhi,
98  const surfaceScalarField& phi,
99  const viscosity& viscosity,
100  const word& type
101 )
102 :
103  eddyViscosity<incompressible::RASModel>
104  (
105  type,
106  alpha,
107  rho,
108  U,
109  alphaRhoPhi,
110  phi,
111  viscosity
112  ),
113 
114  Cmu_
115  (
117  (
118  "Cmu",
119  coeffDict_,
120  0.09
121  )
122  ),
123  C1_
124  (
126  (
127  "C1",
128  coeffDict_,
129  1.44
130  )
131  ),
132  C2_
133  (
135  (
136  "C2",
137  coeffDict_,
138  1.92
139  )
140  ),
141  sigmaZeta_
142  (
144  (
145  "sigmaZeta",
146  coeffDict_,
147  1.3
148  )
149  ),
150  anisotropic_
151  (
153  (
154  "anisotropic",
155  coeffDict_,
156  false
157  )
158  ),
159 
160  qMin_("qMin", sqrt(kMin_)),
161  zetaMin_("zetaMin", epsilonMin_/(2*qMin_)),
162 
163  k_
164  (
165  IOobject
166  (
167  this->groupName("k"),
168  runTime_.name(),
169  mesh_,
172  ),
173  mesh_
174  ),
175 
176  epsilon_
177  (
178  IOobject
179  (
180  this->groupName("epsilon"),
181  runTime_.name(),
182  mesh_,
185  ),
186  mesh_
187  ),
188 
189  q_
190  (
191  IOobject
192  (
193  this->groupName("q"),
194  runTime_.name(),
195  mesh_,
198  ),
199  sqrt(max(k_, kMin_)),
200  k_.boundaryField().types()
201  ),
202 
203  zeta_
204  (
205  IOobject
206  (
207  this->groupName("zeta"),
208  runTime_.name(),
209  mesh_,
212  ),
213  max(epsilon_, epsilonMin_)/(2.0*q_),
214  epsilon_.boundaryField().types()
215  )
216 {
217  bound(zeta_, zetaMin_);
218 
219  if (type == typeName)
220  {
221  printCoeffs(type);
222  }
223 }
224 
225 
226 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
227 
228 bool qZeta::read()
229 {
231  {
232  Cmu_.readIfPresent(coeffDict());
233  C1_.readIfPresent(coeffDict());
234  C2_.readIfPresent(coeffDict());
235  sigmaZeta_.readIfPresent(coeffDict());
236  anisotropic_.readIfPresent("anisotropic", coeffDict());
237 
238  qMin_.readIfPresent(*this);
239  zetaMin_.readIfPresent(*this);
240 
241  return true;
242  }
243  else
244  {
245  return false;
246  }
247 }
248 
249 
250 void qZeta::correct()
251 {
252  if (!turbulence_)
253  {
254  return;
255  }
256 
258 
259  volScalarField G(GName(), nut_/(2.0*q_)*2*magSqr(symm(fvc::grad(U_))));
260  const volScalarField E(nu()*nut_/q_*fvc::magSqrGradGrad(U_));
261 
262  // Zeta equation
263  tmp<fvScalarMatrix> zetaEqn
264  (
265  fvm::ddt(zeta_)
266  + fvm::div(phi_, zeta_)
268  ==
269  (2.0*C1_ - 1)*G*zeta_/q_
270  - fvm::SuSp((2.0*C2_*f2() - dimensionedScalar(1.0))*zeta_/q_, zeta_)
271  + E
272  );
273 
274  zetaEqn.ref().relax();
275  solve(zetaEqn);
276  bound(zeta_, zetaMin_);
277 
278 
279  // q equation
280  tmp<fvScalarMatrix> qEqn
281  (
282  fvm::ddt(q_)
283  + fvm::div(phi_, q_)
284  - fvm::laplacian(DqEff(), q_)
285  ==
286  G - fvm::Sp(zeta_/q_, q_)
287  );
288 
289  qEqn.ref().relax();
290  solve(qEqn);
291  bound(q_, qMin_);
292 
293 
294  // Re-calculate k and epsilon
295  k_ = sqr(q_);
297 
298  epsilon_ = 2*q_*zeta_;
300 
301  correctNut();
302 }
303 
304 
305 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
306 
307 } // End namespace RASModels
308 } // End namespace incompressible
309 } // End namespace Foam
310 
311 // ************************************************************************* //
Bound the given scalar field where it is below the specified minimum.
void correctBoundaryConditions()
Correct boundary field.
bool readIfPresent(const word &, const dictionary &)
Update the value of the Switch if it is found in the dictionary.
Definition: Switch.C:135
static Switch lookupOrAddToDict(const word &, dictionary &, const Switch &defaultValue=false)
Construct from dictionary, supplying default value so that if the.
Definition: Switch.C:111
static dimensioned< scalar > lookupOrAddToDict(const word &, dictionary &, const dimensionSet &dims=dimless, const scalar &defaultValue=pTraits< scalar >::zero)
Construct from dictionary, with default value.
bool readIfPresent(const dictionary &)
Update the value of dimensioned<Type> if found in the dictionary.
virtual bool read()=0
Re-read model coefficients if they have changed.
Definition: eddyViscosity.C:74
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
tmp< volScalarField > f2() const
qZeta(const geometricOneField &alpha, const geometricOneField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity, const word &type=typeName)
Construct from components.
tmp< volScalarField > DqEff() const
Return the effective diffusivity for q.
Definition: qZeta.H:164
dimensionedScalar sigmaZeta_
Definition: qZeta.H:84
dimensionedScalar qMin_
Lower limit of q.
Definition: qZeta.H:88
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
virtual bool read()
Read RASProperties dictionary.
dimensionedScalar zetaMin_
Lower limit of zeta.
Definition: qZeta.H:91
tmp< volScalarField > fMu() const
tmp< volScalarField > DzetaEff() const
Return the effective diffusivity for epsilon.
Definition: qZeta.H:174
Calculate the magnitude of the square of the gradient of the gradient of the given volField.
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
const dimensionedScalar G
Newtonian constant of gravitation.
tmp< volScalarField > magSqrGradGrad(const VolField< Type > &vf)
tmp< VolField< typename outerProduct< vector, Type >::type > > grad(const SurfaceField< Type > &ssf)
Definition: fvcGrad.C:46
tmp< fvMatrix< Type > > laplacian(const VolField< Type > &vf, const word &name)
Definition: fvmLaplacian.C:47
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const VolField< Type > &)
tmp< fvMatrix< Type > > SuSp(const volScalarField::Internal &, const VolField< Type > &)
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const VolField< Type > &vf, const word &name)
Definition: fvmDiv.C:48
tmp< fvMatrix< Type > > ddt(const VolField< Type > &vf)
Definition: fvmDdt.C:46
Namespace for OpenFOAM.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
dimensionedScalar exp(const dimensionedScalar &ds)
VolField< vector > volVectorField
Definition: volFieldsFwd.H:62
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar pow3(const dimensionedScalar &ds)
SurfaceField< scalar > surfaceScalarField
dimensionedScalar sqrt(const dimensionedScalar &ds)
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:61
defineTypeNameAndDebug(combustionModel, 0)
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
dimensioned< scalar > magSqr(const dimensioned< Type > &)
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.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
makeMomentumTransportModelTypes(geometricOneField, geometricOneField, incompressibleMomentumTransportModel)
Definition: qZeta.C:32
Info<< "Reading field U\n"<< endl;volVectorField U(IOobject("U", runTime.name(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Creating face flux\n"<< endl;surfaceScalarField phi(IOobject("phi", runTime.name(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), mesh, dimensionedScalar(mesh.Sf().dimensions() *U.dimensions(), 0));autoPtr< viscosityModel > viscosity(viscosityModel::New(mesh))