ShihQuadraticKE.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 "ShihQuadraticKE.H"
27 #include "bound.H"
28 #include "wallFvPatch.H"
31 
33 (
34  geometricOneField,
35  geometricOneField,
36  incompressibleMomentumTransportModel
37 )
38 
39 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
40 
41 namespace Foam
42 {
43 namespace incompressible
44 {
45 namespace RASModels
46 {
47 
48 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
49 
50 defineTypeNameAndDebug(ShihQuadraticKE, 0);
52 (
53  RASincompressibleMomentumTransportModel,
54  ShihQuadraticKE,
55  dictionary
56 );
57 
58 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
59 
61 {
63 }
64 
65 
67 {
68  volSymmTensorField S(symm(gradU));
69  volTensorField W(skew(gradU));
70 
71  volScalarField sBar((k_/epsilon_)*sqrt(2.0)*mag(S));
72  volScalarField wBar((k_/epsilon_)*sqrt(2.0)*mag(W));
73 
74  volScalarField Cmu((2.0/3.0)/(Cmu1_ + sBar + Cmu2_*wBar));
75 
76  nut_ = Cmu*sqr(k_)/epsilon_;
78 
80  k_*sqr(k_/epsilon_)/(Cbeta_ + pow3(sBar))
81  *(
83  + Cbeta2_*twoSymm(S&W)
84  + Cbeta3_*dev(symm(W&W))
85  );
86 }
87 
88 
89 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
90 
92 (
93  const geometricOneField& alpha,
94  const geometricOneField& rho,
95  const volVectorField& U,
96  const surfaceScalarField& alphaRhoPhi,
97  const surfaceScalarField& phi,
98  const viscosity& viscosity,
99  const word& type
100 )
101 :
102  nonlinearEddyViscosity<incompressible::RASModel>
103  (
104  type,
105  alpha,
106  rho,
107  U,
108  alphaRhoPhi,
109  phi,
110  viscosity
111  ),
112 
113  Ceps1_
114  (
116  (
117  "Ceps1",
118  coeffDict_,
119  1.44
120  )
121  ),
122  Ceps2_
123  (
125  (
126  "Ceps2",
127  coeffDict_,
128  1.92
129  )
130  ),
131  sigmak_
132  (
134  (
135  "sigmak",
136  coeffDict_,
137  1.0
138  )
139  ),
140  sigmaEps_
141  (
143  (
144  "sigmaEps",
145  coeffDict_,
146  1.3
147  )
148  ),
149  Cmu1_
150  (
152  (
153  "Cmu1",
154  coeffDict_,
155  1.25
156  )
157  ),
158  Cmu2_
159  (
161  (
162  "Cmu2",
163  coeffDict_,
164  0.9
165  )
166  ),
167  Cbeta_
168  (
170  (
171  "Cbeta",
172  coeffDict_,
173  1000.0
174  )
175  ),
176  Cbeta1_
177  (
179  (
180  "Cbeta1",
181  coeffDict_,
182  3.0
183  )
184  ),
185  Cbeta2_
186  (
188  (
189  "Cbeta2",
190  coeffDict_,
191  15.0
192  )
193  ),
194  Cbeta3_
195  (
197  (
198  "Cbeta3",
199  coeffDict_,
200  -19.0
201  )
202  ),
203 
204  k_
205  (
206  IOobject
207  (
208  this->groupName("k"),
209  runTime_.name(),
210  mesh_,
213  ),
214  mesh_
215  ),
216 
217  epsilon_
218  (
219  IOobject
220  (
221  this->groupName("epsilon"),
222  runTime_.name(),
223  mesh_,
226  ),
227  mesh_
228  )
229 {
230  bound(k_, kMin_);
231  bound(epsilon_, epsilonMin_);
232 
233  if (type == typeName)
234  {
235  printCoeffs(type);
236  }
237 }
238 
239 
240 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
241 
243 {
245  {
246  Ceps1_.readIfPresent(coeffDict());
247  Ceps2_.readIfPresent(coeffDict());
248  sigmak_.readIfPresent(coeffDict());
249  sigmaEps_.readIfPresent(coeffDict());
250  Cmu1_.readIfPresent(coeffDict());
251  Cmu2_.readIfPresent(coeffDict());
252  Cbeta_.readIfPresent(coeffDict());
253  Cbeta1_.readIfPresent(coeffDict());
254  Cbeta2_.readIfPresent(coeffDict());
255  Cbeta3_.readIfPresent(coeffDict());
256 
257  return true;
258  }
259  else
260  {
261  return false;
262  }
263 }
264 
265 
267 {
268  if (!turbulence_)
269  {
270  return;
271  }
272 
274 
275  tmp<volTensorField> tgradU = fvc::grad(U_);
276  const volTensorField& gradU = tgradU();
277 
279  (
280  GName(),
281  (nut_*twoSymm(gradU) - nonlinearStress_) && gradU
282  );
283 
284 
285  // Update epsilon and G at the wall
287 
288  // Dissipation equation
289  tmp<fvScalarMatrix> epsEqn
290  (
292  + fvm::div(phi_, epsilon_)
294  ==
297  );
298 
299  epsEqn.ref().relax();
300  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
301  solve(epsEqn);
302  bound(epsilon_, epsilonMin_);
303 
304 
305  // Turbulent kinetic energy equation
306  tmp<fvScalarMatrix> kEqn
307  (
308  fvm::ddt(k_)
309  + fvm::div(phi_, k_)
310  - fvm::laplacian(DkEff(), k_)
311  ==
312  G
313  - fvm::Sp(epsilon_/k_, k_)
314  );
315 
316  kEqn.ref().relax();
317  solve(kEqn);
318  bound(k_, kMin_);
319 
320 
321  // Re-calculate viscosity and non-linear stress
322  correctNonlinearStress(gradU);
323 }
324 
325 
326 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
327 
328 } // End namespace RASModels
329 } // End namespace incompressible
330 } // End namespace Foam
331 
332 // ************************************************************************* //
makeMomentumTransportModelTypes(geometricOneField, geometricOneField, incompressibleMomentumTransportModel)
Bound the given scalar field where it is below the specified minimum.
void updateCoeffs()
Update the boundary condition coefficients.
void correctBoundaryConditions()
Correct boundary field.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
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 > DkEff() const
Return the effective diffusivity for k.
ShihQuadraticKE(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.
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
tmp< volScalarField > DepsilonEff() const
Return the effective diffusivity for epsilon.
virtual void correctNonlinearStress(const volTensorField &gradU)
virtual bool read()
Read RASProperties dictionary.
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< VolField< typename outerProduct< vector, Type >::type > > grad(const SurfaceField< Type > &ssf)
Definition: fvcGrad.C:46
tmp< fvMatrix< Type > > S(const Pair< tmp< volScalarField::Internal >> &, const VolField< Type > &)
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 > > 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.
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
VolField< vector > volVectorField
Definition: volFieldsFwd.H:62
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar pow3(const dimensionedScalar &ds)
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
SurfaceField< scalar > surfaceScalarField
VolField< tensor > volTensorField
Definition: volFieldsFwd.H:65
scalarList W(const fluidMulticomponentThermo &thermo)
dimensionedScalar sqrt(const dimensionedScalar &ds)
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:61
dimensioned< scalar > mag(const dimensioned< Type > &)
defineTypeNameAndDebug(combustionModel, 0)
bool bound(volScalarField &, const dimensionedScalar &min)
Bound the given scalar field where it is below the specified min value.
Definition: bound.C:31
dimensionedSymmTensor innerSqr(const dimensionedSymmTensor &dt)
VolField< symmTensor > volSymmTensorField
Definition: volFieldsFwd.H:64
dimensionedTensor skew(const dimensionedTensor &dt)
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)
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))