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-2025 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_, Cmu_*sqr(k_)/(nutMaxCoeff_*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_("Ceps1", coeffDict(), 1.44),
121  Ceps2_("Ceps2", coeffDict(), 1.92),
122  sigmak_("sigmak", coeffDict(), 1.0),
123  sigmaEps_("sigmaEps", coeffDict(), 1.3),
124  Cmu_("Cmu", coeffDict(), 0.09),
125  Cmu1_("Cmu1", coeffDict(), 1.25),
126  Cmu2_("Cmu2", coeffDict(), 0.9),
127  Cbeta_("Cbeta", coeffDict(), 1000.0),
128  Cbeta1_("Cbeta1", coeffDict(), 3.0),
129  Cbeta2_("Cbeta2", coeffDict(), 15.0),
130  Cbeta3_("Cbeta3", coeffDict(), -19.0),
131 
132  k_
133  (
134  IOobject
135  (
136  this->groupName("k"),
137  runTime_.name(),
138  mesh_,
141  ),
142  mesh_
143  ),
144 
145  epsilon_
146  (
147  IOobject
148  (
149  this->groupName("epsilon"),
150  runTime_.name(),
151  mesh_,
154  ),
155  mesh_
156  )
157 {
158  bound(k_, kMin_);
159  boundEpsilon();
160 }
161 
162 
163 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
164 
166 {
168  {
169  Ceps1_.readIfPresent(coeffDict());
170  Ceps2_.readIfPresent(coeffDict());
171  sigmak_.readIfPresent(coeffDict());
172  sigmaEps_.readIfPresent(coeffDict());
173  Cmu_.readIfPresent(coeffDict());
174  Cmu1_.readIfPresent(coeffDict());
175  Cmu2_.readIfPresent(coeffDict());
176  Cbeta_.readIfPresent(coeffDict());
177  Cbeta1_.readIfPresent(coeffDict());
178  Cbeta2_.readIfPresent(coeffDict());
179  Cbeta3_.readIfPresent(coeffDict());
180 
181  return true;
182  }
183  else
184  {
185  return false;
186  }
187 }
188 
189 
191 {
192  if (!turbulence_)
193  {
194  return;
195  }
196 
198 
199  tmp<volTensorField> tgradU = fvc::grad(U_);
200  const volTensorField& gradU = tgradU();
201 
203  (
204  GName(),
205  (nut_*twoSymm(gradU) - nonlinearStress_) && gradU
206  );
207 
208 
209  // Update epsilon and G at the wall
211 
212  // Dissipation equation
213  tmp<fvScalarMatrix> epsEqn
214  (
216  + fvm::div(phi_, epsilon_)
218  ==
221  );
222 
223  epsEqn.ref().relax();
224  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
225  solve(epsEqn);
226  boundEpsilon();
227 
228 
229  // Turbulent kinetic energy equation
230  tmp<fvScalarMatrix> kEqn
231  (
232  fvm::ddt(k_)
233  + fvm::div(phi_, k_)
234  - fvm::laplacian(DkEff(), k_)
235  ==
236  G
237  - fvm::Sp(epsilon_/k_, k_)
238  );
239 
240  kEqn.ref().relax();
241  solve(kEqn);
242  bound(k_, kMin_);
243 
244 
245  // Re-calculate viscosity and non-linear stress
246  correctNonlinearStress(gradU);
247 }
248 
249 
250 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
251 
252 } // End namespace RASModels
253 } // End namespace incompressible
254 } // End namespace Foam
255 
256 // ************************************************************************* //
makeMomentumTransportModelTypes(geometricOneField, geometricOneField, incompressibleMomentumTransportModel)
Bound the given scalar field where it is below the specified minimum.
void updateCoeffs()
Update the boundary condition coefficients.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
void correctBoundaryConditions()
Correct boundary field.
bool readIfPresent(const dictionary &, const unitConversion &defaultUnits=NullObjectRef< unitConversion >())
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.
void skew(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
VolField< vector > volVectorField
Definition: volFieldsFwd.H:63
addToRunTimeSelectionTable(polyPatch, mergedCyclicPolyPatch, word)
void dev(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
void twoSymm(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
SurfaceField< scalar > surfaceScalarField
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
VolField< tensor > volTensorField
Definition: volFieldsFwd.H:66
void symm(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
scalarList W(const fluidMulticomponentThermo &thermo)
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:62
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
void pow3(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
void sqr(LagrangianPatchField< typename outerProduct< Type, Type >::type > &f, const LagrangianPatchField< Type > &f1)
dimensionedSymmTensor innerSqr(const dimensionedSymmTensor &dt)
void sqrt(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
VolField< symmTensor > volSymmTensorField
Definition: volFieldsFwd.H:65
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))