Loading...
Searching...
No Matches
ShihQuadraticKE.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2019-2023 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "ShihQuadraticKE.H"
30#include "bound.H"
31#include "wallFvPatch.H"
33
34// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35
36namespace Foam
37{
40namespace RASModels
41{
42
43// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
44
47
48// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
51{
53}
54
55
57{
58 volSymmTensorField S(symm(gradU));
59 volTensorField W(skew(gradU));
60
61 volScalarField sBar((k_/epsilon_)*sqrt(2.0)*mag(S));
62 volScalarField wBar((k_/epsilon_)*sqrt(2.0)*mag(W));
63
64 volScalarField Cmu((2.0/3.0)/(Cmu1_ + sBar + Cmu2_*wBar));
65
66 nut_ = Cmu*sqr(k_)/epsilon_;
68
70 k_*sqr(k_/epsilon_)/(Cbeta_ + pow3(sBar))
71 *(
73 + Cbeta2_*twoSymm(S&W)
75 );
76}
77
78
79// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
80
82(
85 const volVectorField& U,
86 const surfaceScalarField& alphaRhoPhi,
88 const transportModel& transport,
89 const word& propertiesName,
90 const word& type
91)
92:
94 (
95 type,
96 alpha,
97 rho,
98 U,
99 alphaRhoPhi,
100 phi,
101 transport,
102 propertiesName
103 ),
104
105 Ceps1_
106 (
107 dimensioned<scalar>::getOrAddToDict
108 (
109 "Ceps1",
110 coeffDict_,
111 1.44
112 )
113 ),
114 Ceps2_
115 (
116 dimensioned<scalar>::getOrAddToDict
117 (
118 "Ceps2",
119 coeffDict_,
120 1.92
121 )
122 ),
123 sigmak_
124 (
125 dimensioned<scalar>::getOrAddToDict
126 (
127 "sigmak",
128 coeffDict_,
129 1.0
130 )
131 ),
132 sigmaEps_
133 (
134 dimensioned<scalar>::getOrAddToDict
135 (
136 "sigmaEps",
137 coeffDict_,
138 1.3
139 )
140 ),
141 Cmu1_
142 (
143 dimensioned<scalar>::getOrAddToDict
144 (
145 "Cmu1",
146 coeffDict_,
147 1.25
148 )
149 ),
150 Cmu2_
151 (
152 dimensioned<scalar>::getOrAddToDict
153 (
154 "Cmu2",
155 coeffDict_,
156 0.9
157 )
158 ),
159 Cbeta_
160 (
161 dimensioned<scalar>::getOrAddToDict
162 (
163 "Cbeta",
164 coeffDict_,
165 1000.0
166 )
167 ),
168 Cbeta1_
169 (
170 dimensioned<scalar>::getOrAddToDict
171 (
172 "Cbeta1",
173 coeffDict_,
174 3.0
175 )
176 ),
177 Cbeta2_
178 (
179 dimensioned<scalar>::getOrAddToDict
180 (
181 "Cbeta2",
182 coeffDict_,
183 15.0
184 )
185 ),
186 Cbeta3_
187 (
188 dimensioned<scalar>::getOrAddToDict
189 (
190 "Cbeta3",
191 coeffDict_,
192 -19.0
193 )
194 ),
195
196 k_
197 (
199 (
200 IOobject::groupName("k", alphaRhoPhi.group()),
201 runTime_.timeName(),
202 mesh_,
203 IOobject::MUST_READ,
204 IOobject::AUTO_WRITE
205 ),
206 mesh_
207 ),
208
209 epsilon_
210 (
212 (
213 IOobject::groupName("epsilon", alphaRhoPhi.group()),
214 runTime_.timeName(),
215 mesh_,
216 IOobject::MUST_READ,
217 IOobject::AUTO_WRITE
218 ),
219 mesh_
220 )
221{
222 bound(k_, kMin_);
223 bound(epsilon_, epsilonMin_);
224
225 if (type == typeName)
226 {
227 printCoeffs(type);
228 }
229}
230
231
232// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
233
235{
237 {
238 Ceps1_.readIfPresent(coeffDict());
239 Ceps2_.readIfPresent(coeffDict());
240 sigmak_.readIfPresent(coeffDict());
241 sigmaEps_.readIfPresent(coeffDict());
242 Cmu1_.readIfPresent(coeffDict());
243 Cmu2_.readIfPresent(coeffDict());
244 Cbeta_.readIfPresent(coeffDict());
245 Cbeta1_.readIfPresent(coeffDict());
246 Cbeta2_.readIfPresent(coeffDict());
247 Cbeta3_.readIfPresent(coeffDict());
248
249 return true;
250 }
251
252 return false;
253}
254
255
257{
258 if (!turbulence_)
259 {
260 return;
261 }
262
264
265 tmp<volTensorField> tgradU = fvc::grad(U_);
266 const volTensorField& gradU = tgradU();
267
269 (
270 GName(),
271 (nut_*twoSymm(gradU) - nonlinearStress_) && gradU
272 );
273
274
275 // Update epsilon and G at the wall
277 // Push any changed cell values to coupled neighbours
279
280 // Dissipation equation
282 (
284 + fvm::div(phi_, epsilon_)
286 ==
289 );
290
291 epsEqn.ref().relax();
292 epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
293 solve(epsEqn);
294 bound(epsilon_, epsilonMin_);
295
296
297 // Turbulent kinetic energy equation
299 (
300 fvm::ddt(k_)
301 + fvm::div(phi_, k_)
303 ==
304 G
306 );
307
308 kEqn.ref().relax();
309 solve(kEqn);
310 bound(k_, kMin_);
311
312
313 // Re-calculate viscosity and non-linear stress
315}
316
317
318// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
319
320} // End namespace RASModels
321} // End namespace incompressible
322} // End namespace Foam
323
324// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Bound the given scalar field if it has gone unbounded.
void evaluateCoupled(const UPstream::commsTypes commsType=UPstream::defaultCommsType)
Evaluate boundary conditions on coupled patches of the given type, using specified or default comms.
void updateCoeffs()
Update the boundary condition coefficients.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
void correctBoundaryConditions()
Correct boundary field.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
An abstract base class for patches that couple regions of the computational domain e....
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Generic dimensioned Type class.
bool readIfPresent(const dictionary &dict)
Update the value of dimensioned<Type> if found in the dictionary, lookup in dictionary with the name(...
virtual bool read()=0
Re-read model coefficients if they have changed.
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
Shih's quadratic algebraic Reynolds stress k-epsilon turbulence model for incompressible flows.
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)
ShihQuadraticKE(const geometricOneField &alpha, const geometricOneField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName=turbulenceModel::propertiesName, const word &type=typeName)
Construct from components.
virtual bool read()
Re-read model coefficients if they have changed.
tmp< volScalarField > DkEff() const
Return the effective diffusivity for k.
Eddy viscosity turbulence model with non-linear correction base class.
nonlinearEddyViscosity(const word &modelName, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName)
incompressible::RASModel::transportModel transportModel
A class for managing temporary objects.
Definition tmp.H:75
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition tmpI.H:235
Base-class for all transport models used by the incompressible turbulence models.
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
U
Definition pEqn.H:72
word timeName
Definition getTimeIndex.H:3
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcGrad.C:47
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition fvmDiv.C:41
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition fvmDdt.C:41
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
RASModel< turbulenceModel > RASModel
Namespace for OpenFOAM.
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
GeometricField< vector, fvPatchField, volMesh > volVectorField
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar pow3(const dimensionedScalar &ds)
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition bound.C:29
GeometricField< scalar, fvPatchField, volMesh > volScalarField
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
SymmTensor< Cmpt > devSymm(const SymmTensor< Cmpt > &st)
Return the deviatoric part of the symmetric part of a SymmTensor.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
GeometricField< tensor, fvPatchField, volMesh > volTensorField
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
dimensionedSymmTensor innerSqr(const dimensionedSymmTensor &dt)
dimensionedTensor skew(const dimensionedTensor &dt)
volScalarField & alpha
CEqn solve()