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 {
62  epsilon_ = max(epsilon_, 0.09*sqr(k_)/(this->nutMaxCoeff_*this->nu()));
63 }
64 
65 
67 {
69 }
70 
71 
73 {
74  volSymmTensorField S(symm(gradU));
75  volTensorField W(skew(gradU));
76 
77  volScalarField sBar((k_/epsilon_)*sqrt(2.0)*mag(S));
78  volScalarField wBar((k_/epsilon_)*sqrt(2.0)*mag(W));
79 
80  volScalarField Cmu((2.0/3.0)/(Cmu1_ + sBar + Cmu2_*wBar));
81 
82  boundEpsilon();
83  nut_ = Cmu*sqr(k_)/epsilon_;
85 
87  k_*sqr(k_/epsilon_)/(Cbeta_ + pow3(sBar))
88  *(
90  + Cbeta2_*twoSymm(S&W)
91  + Cbeta3_*dev(symm(W&W))
92  );
93 }
94 
95 
96 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
97 
99 (
100  const geometricOneField& alpha,
101  const geometricOneField& rho,
102  const volVectorField& U,
103  const surfaceScalarField& alphaRhoPhi,
104  const surfaceScalarField& phi,
105  const viscosity& viscosity,
106  const word& type
107 )
108 :
109  nonlinearEddyViscosity<incompressible::RASModel>
110  (
111  type,
112  alpha,
113  rho,
114  U,
115  alphaRhoPhi,
116  phi,
117  viscosity
118  ),
119 
120  Ceps1_
121  (
123  (
124  "Ceps1",
125  coeffDict_,
126  1.44
127  )
128  ),
129  Ceps2_
130  (
132  (
133  "Ceps2",
134  coeffDict_,
135  1.92
136  )
137  ),
138  sigmak_
139  (
141  (
142  "sigmak",
143  coeffDict_,
144  1.0
145  )
146  ),
147  sigmaEps_
148  (
150  (
151  "sigmaEps",
152  coeffDict_,
153  1.3
154  )
155  ),
156  Cmu1_
157  (
159  (
160  "Cmu1",
161  coeffDict_,
162  1.25
163  )
164  ),
165  Cmu2_
166  (
168  (
169  "Cmu2",
170  coeffDict_,
171  0.9
172  )
173  ),
174  Cbeta_
175  (
177  (
178  "Cbeta",
179  coeffDict_,
180  1000.0
181  )
182  ),
183  Cbeta1_
184  (
186  (
187  "Cbeta1",
188  coeffDict_,
189  3.0
190  )
191  ),
192  Cbeta2_
193  (
195  (
196  "Cbeta2",
197  coeffDict_,
198  15.0
199  )
200  ),
201  Cbeta3_
202  (
204  (
205  "Cbeta3",
206  coeffDict_,
207  -19.0
208  )
209  ),
210 
211  k_
212  (
213  IOobject
214  (
215  this->groupName("k"),
216  runTime_.name(),
217  mesh_,
220  ),
221  mesh_
222  ),
223 
224  epsilon_
225  (
226  IOobject
227  (
228  this->groupName("epsilon"),
229  runTime_.name(),
230  mesh_,
233  ),
234  mesh_
235  )
236 {
237  bound(k_, kMin_);
238  boundEpsilon();
239 
240  if (type == typeName)
241  {
242  printCoeffs(type);
243  }
244 }
245 
246 
247 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
248 
250 {
252  {
253  Ceps1_.readIfPresent(coeffDict());
254  Ceps2_.readIfPresent(coeffDict());
255  sigmak_.readIfPresent(coeffDict());
256  sigmaEps_.readIfPresent(coeffDict());
257  Cmu1_.readIfPresent(coeffDict());
258  Cmu2_.readIfPresent(coeffDict());
259  Cbeta_.readIfPresent(coeffDict());
260  Cbeta1_.readIfPresent(coeffDict());
261  Cbeta2_.readIfPresent(coeffDict());
262  Cbeta3_.readIfPresent(coeffDict());
263 
264  return true;
265  }
266  else
267  {
268  return false;
269  }
270 }
271 
272 
274 {
275  if (!turbulence_)
276  {
277  return;
278  }
279 
281 
282  tmp<volTensorField> tgradU = fvc::grad(U_);
283  const volTensorField& gradU = tgradU();
284 
286  (
287  GName(),
288  (nut_*twoSymm(gradU) - nonlinearStress_) && gradU
289  );
290 
291 
292  // Update epsilon and G at the wall
294 
295  // Dissipation equation
296  tmp<fvScalarMatrix> epsEqn
297  (
299  + fvm::div(phi_, epsilon_)
301  ==
304  );
305 
306  epsEqn.ref().relax();
307  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
308  solve(epsEqn);
309  boundEpsilon();
310 
311 
312  // Turbulent kinetic energy equation
313  tmp<fvScalarMatrix> kEqn
314  (
315  fvm::ddt(k_)
316  + fvm::div(phi_, k_)
317  - fvm::laplacian(DkEff(), k_)
318  ==
319  G
320  - fvm::Sp(epsilon_/k_, k_)
321  );
322 
323  kEqn.ref().relax();
324  solve(kEqn);
325  bound(k_, kMin_);
326 
327 
328  // Re-calculate viscosity and non-linear stress
329  correctNonlinearStress(gradU);
330 }
331 
332 
333 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
334 
335 } // End namespace RASModels
336 } // End namespace incompressible
337 } // End namespace Foam
338 
339 // ************************************************************************* //
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.
virtual void correctNut()
Correct the eddy-viscosity nut.
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 > > div(const surfaceScalarField &flux, const VolField< Type > &vf, const word &name)
Definition: fvmDiv.C:48
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const VolField< Type > &)
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:65
addToRunTimeSelectionTable(polyPatch, mergedCyclicPolyPatch, word)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar pow3(const dimensionedScalar &ds)
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
SurfaceField< scalar > surfaceScalarField
VolField< tensor > volTensorField
Definition: volFieldsFwd.H:68
scalarList W(const fluidMulticomponentThermo &thermo)
dimensionedScalar sqrt(const dimensionedScalar &ds)
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:64
dimensioned< scalar > mag(const dimensioned< Type > &)
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
dimensionedSymmTensor innerSqr(const dimensionedSymmTensor &dt)
VolField< symmTensor > volSymmTensorField
Definition: volFieldsFwd.H:67
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.
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))