kOmegaSSTDES.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 | Copyright (C) 2015 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 "kOmegaSSTDES.H"
27 
28 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
29 
30 namespace Foam
31 {
32 namespace LESModels
33 {
34 
35 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
36 
37 template<class BasicTurbulenceModel>
39 {
40  // Correct the turbulence viscosity
42 
43  // Correct the turbulence thermal diffusivity
44  BasicTurbulenceModel::correctNut();
45 }
46 
47 
48 template<class BasicTurbulenceModel>
50 {
51  correctNut(2*magSqr(symm(fvc::grad(this->U_))));
52 }
53 
54 
55 template<class BasicTurbulenceModel>
57 (
58  const volScalarField& magGradU,
59  const volScalarField& CDES
60 ) const
61 {
62  const volScalarField& k = this->k_;
63  const volScalarField& omega = this->omega_;
64 
65  return min(CDES*this->delta(), sqrt(k)/(this->betaStar_*omega));
66 }
67 
68 
69 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
70 
71 template<class BasicTurbulenceModel>
73 (
74  const alphaField& alpha,
75  const rhoField& rho,
76  const volVectorField& U,
77  const surfaceScalarField& alphaRhoPhi,
78  const surfaceScalarField& phi,
79  const transportModel& transport,
80  const word& propertiesName,
81  const word& type
82 )
83 :
85  (
86  type,
87  alpha,
88  rho,
89  U,
90  alphaRhoPhi,
91  phi,
92  transport,
93  propertiesName
94  ),
95 
96  kappa_
97  (
99  (
100  "kappa",
101  this->coeffDict_,
102  0.41
103  )
104  ),
105  CDESkom_
106  (
108  (
109  "CDESkom",
110  this->coeffDict_,
111  0.82
112  )
113  ),
114  CDESkeps_
115  (
117  (
118  "CDESkeps",
119  this->coeffDict_,
120  0.60
121  )
122  )
123 {
124  correctNut();
125 
126  if (type == typeName)
127  {
128  this->printCoeffs(type);
129  }
130 }
131 
132 
133 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
134 
135 template<class BasicTurbulenceModel>
137 {
139  {
140  kappa_.readIfPresent(this->coeffDict());
141  CDESkom_.readIfPresent(this->coeffDict());
142  CDESkeps_.readIfPresent(this->coeffDict());
143 
144  return true;
145  }
146  else
147  {
148  return false;
149  }
150 }
151 
152 
153 template<class BasicTurbulenceModel>
155 {
156  if (!this->turbulence_)
157  {
158  return;
159  }
160 
161  // Local references
162  const alphaField& alpha = this->alpha_;
163  const rhoField& rho = this->rho_;
164  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
165  const volVectorField& U = this->U_;
166  volScalarField& k = this->k_;
167  volScalarField& omega = this->omega_;
168  volScalarField& nut = this->nut_;
169 
171 
173 
174  tmp<volTensorField> tgradU = fvc::grad(U);
175  volScalarField magGradU(mag(tgradU()));
176  volScalarField S2(2*magSqr(symm(tgradU())));
177  volScalarField GbyNu((tgradU() && dev(twoSymm(tgradU()))));
178  volScalarField G(this->GName(), nut*GbyNu);
179  tgradU.clear();
180 
181  // Update omega and G at the wall
182  omega.boundaryField().updateCoeffs();
183 
184  volScalarField CDkOmega
185  (
186  (2*this->alphaOmega2_)*(fvc::grad(k) & fvc::grad(omega))/omega
187  );
188 
189  volScalarField F1(this->F1(CDkOmega));
190 
191  {
192  volScalarField gamma(this->gamma(F1));
193  volScalarField beta(this->beta(F1));
194 
195  // Turbulent frequency equation
196  tmp<fvScalarMatrix> omegaEqn
197  (
198  fvm::ddt(alpha, rho, omega)
199  + fvm::div(alphaRhoPhi, omega)
200  - fvm::laplacian(alpha*rho*this->DomegaEff(F1), omega)
201  ==
202  alpha*rho*gamma*GbyNu // Using unlimited GybNu
203  - fvm::SuSp((2.0/3.0)*alpha*rho*gamma*divU, omega)
204  - fvm::Sp(alpha*rho*beta*omega, omega)
205  - fvm::SuSp(alpha*rho*(F1 - scalar(1))*CDkOmega/omega, omega)
206  + this->omegaSource()
207  );
208 
209  omegaEqn().relax();
210 
211  omegaEqn().boundaryManipulate(omega.boundaryField());
212 
213  solve(omegaEqn);
214  bound(omega, this->omegaMin_);
215  }
216 
217  {
218  volScalarField CDES(this->CDES(F1));
219  volScalarField dTilda(this->dTilda(magGradU, CDES));
220 
221  // Turbulent kinetic energy equation
223  (
224  fvm::ddt(alpha, rho, k)
225  + fvm::div(alphaRhoPhi, k)
226  - fvm::laplacian(alpha*rho*this->DkEff(F1), k)
227  ==
228  min(alpha*rho*G, (this->c1_*this->betaStar_)*alpha*rho*k*omega)
229  - fvm::SuSp((2.0/3.0)*alpha*rho*divU, k)
230  - fvm::Sp(alpha*rho*sqrt(k)/dTilda, k) // modified for DES
231  + this->kSource()
232  );
233 
234  kEqn().relax();
235  solve(kEqn);
236  bound(k, this->kMin_);
237  }
238 
239  this->correctNut(S2);
240 }
241 
242 
243 template<class BasicTurbulenceModel>
245 {
246  const volScalarField& k = this->k_;
247  const volScalarField& omega = this->omega_;
248  const volVectorField& U = this->U_;
249 
250  const volScalarField CDkOmega
251  (
252  (2*this->alphaOmega2_)*(fvc::grad(k) & fvc::grad(omega))/omega
253  );
254 
255  const volScalarField F1(this->F1(CDkOmega));
256 
257  tmp<volScalarField> tLESRegion
258  (
259  new volScalarField
260  (
261  IOobject
262  (
263  "DES::LESRegion",
264  this->mesh_.time().timeName(),
265  this->mesh_,
268  ),
269  neg
270  (
271  dTilda
272  (
273  mag(fvc::grad(U)),
274  F1*CDESkom_ + (1 - F1)*CDESkeps_
275  )
276  - sqrt(k)/(this->betaStar_*omega)
277  )
278  )
279  );
280 
281  return tLESRegion;
282 }
283 
284 
285 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
286 
287 } // End namespace LESModels
288 } // End namespace Foam
289 
290 // ************************************************************************* //
beta
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
F1
#define F1(B, C, D)
Definition: SHA1.C:175
Foam::symm
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:82
Foam::GeometricField::GeometricBoundaryField::updateCoeffs
void updateCoeffs()
Update the boundary condition coefficients.
Definition: GeometricBoundaryField.C:447
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::constant::universal::G
const dimensionedScalar G
Newtonian constant of gravitation.
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:216
Foam::read
bool read(const char *, int32_t &)
Definition: int32IO.C:87
Foam::bound
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:33
Foam::LESModels::kOmegaSSTDES::correct
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: kOmegaSSTDES.C:154
Foam::fvc::div
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
Foam::GeometricField::boundaryField
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
Definition: GeometricField.C:735
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::LESModels::kOmegaSSTDES::kOmegaSSTDES
kOmegaSSTDES(const kOmegaSSTDES &)
Foam::LESModels::kOmegaSSTDES::read
virtual bool read()
Re-read model coefficients if they have changed.
Definition: kOmegaSSTDES.C:136
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
U
U
Definition: pEqn.H:46
Foam::LESModels::kEqn
One equation eddy-viscosity model.
Definition: kEqn.H:74
Foam::fvm::Sp
tmp< fvMatrix< Type > > Sp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
nut
nut
Definition: createTDFields.H:71
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
DESModel< BasicTurbulenceModel >
Foam::fvm::laplacian
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:46
delta
scalar delta
Definition: LISASMDCalcMethod2.H:8
Foam::fvm::SuSp
tmp< fvMatrix< Type > > SuSp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
Foam::LESModels::kOmegaSSTDES::dTilda
virtual tmp< volScalarField > dTilda(const volScalarField &magGradU, const volScalarField &CDES) const
Length scale.
Definition: kOmegaSSTDES.C:57
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::LESModels::kOmegaSSTDES::transportModel
BasicTurbulenceModel::transportModel transportModel
Definition: kOmegaSSTDES.H:110
rho
rho
Definition: pEqn.H:3
Foam::solve
SolverPerformance< Type > solve(fvMatrix< Type > &, const dictionary &)
Solve returning the solution statistics given convergence tolerance.
Foam::fvm::ddt
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
Foam::LESModels::kOmegaSSTDES::alphaField
BasicTurbulenceModel::alphaField alphaField
Definition: kOmegaSSTDES.H:108
Foam::LESModels::kOmegaSSTDES::correctNut
virtual void correctNut()
Definition: kOmegaSSTDES.C:49
divU
volScalarField divU(fvc::div(fvc::absolute(phi, U)))
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:142
k
label k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
Foam::eddyViscosity< LESModel< BasicTurbulenceModel > >::correct
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
Definition: eddyViscosity.C:132
Foam::type
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:588
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::fvm::div
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:46
Foam::LESModels::kOmegaSSTDES::LESRegion
virtual tmp< volScalarField > LESRegion() const
Return the LES field indicator.
Definition: kOmegaSSTDES.C:244
Foam::tmp::clear
void clear() const
If object pointer points to valid object:
Definition: tmpI.H:172
Foam::min
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::LESModels::kOmegaSSTDES::rhoField
BasicTurbulenceModel::rhoField rhoField
Definition: kOmegaSSTDES.H:109
Foam::neg
dimensionedScalar neg(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:201
Foam::twoSymm
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:93
Foam::magSqr
dimensioned< scalar > magSqr(const dimensioned< Type > &)
Foam::fvc::absolute
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:187
Foam::kOmegaSSTBase
Base class implementation of the k-omega-SST turbulence model for incompressible and compressible flo...
Definition: kOmegaSSTBase.H:115
kOmegaSSTDES.H
Foam::dev
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:104