SSG.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) 2015-2016 OpenFOAM Foundation
9  Copyright (C) 2018-2020 OpenCFD Ltd
10 -------------------------------------------------------------------------------
11 License
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 "SSG.H"
30 #include "fvOptions.H"
31 #include "wallFvPatch.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace RASModels
38 {
39 
40 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
41 
42 template<class BasicTurbulenceModel>
44 {
45  this->nut_ = this->Cmu_*sqr(k_)/epsilon_;
46  this->nut_.correctBoundaryConditions();
47  fv::options::New(this->mesh_).correct(this->nut_);
48 
49  BasicTurbulenceModel::correctNut();
50 }
51 
52 
53 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
54 
55 template<class BasicTurbulenceModel>
57 (
58  const alphaField& alpha,
59  const rhoField& rho,
60  const volVectorField& U,
61  const surfaceScalarField& alphaRhoPhi,
62  const surfaceScalarField& phi,
63  const transportModel& transport,
64  const word& propertiesName,
65  const word& type
66 )
67 :
69  (
70  type,
71  alpha,
72  rho,
73  U,
74  alphaRhoPhi,
75  phi,
76  transport,
77  propertiesName
78  ),
79 
80  Cmu_
81  (
83  (
84  "Cmu",
85  this->coeffDict_,
86  0.09
87  )
88  ),
89  C1_
90  (
92  (
93  "C1",
94  this->coeffDict_,
95  3.4
96  )
97  ),
98  C1s_
99  (
101  (
102  "C1s",
103  this->coeffDict_,
104  1.8
105  )
106  ),
107  C2_
108  (
110  (
111  "C2",
112  this->coeffDict_,
113  4.2
114  )
115  ),
116  C3_
117  (
119  (
120  "C3",
121  this->coeffDict_,
122  0.8
123  )
124  ),
125  C3s_
126  (
128  (
129  "C3s",
130  this->coeffDict_,
131  1.3
132  )
133  ),
134  C4_
135  (
137  (
138  "C4",
139  this->coeffDict_,
140  1.25
141  )
142  ),
143  C5_
144  (
146  (
147  "C5",
148  this->coeffDict_,
149  0.4
150  )
151  ),
152 
153  Ceps1_
154  (
156  (
157  "Ceps1",
158  this->coeffDict_,
159  1.44
160  )
161  ),
162  Ceps2_
163  (
165  (
166  "Ceps2",
167  this->coeffDict_,
168  1.92
169  )
170  ),
171  Cs_
172  (
174  (
175  "Cs",
176  this->coeffDict_,
177  0.25
178  )
179  ),
180  Ceps_
181  (
183  (
184  "Ceps",
185  this->coeffDict_,
186  0.15
187  )
188  ),
189 
190  k_
191  (
192  IOobject
193  (
194  "k",
195  this->runTime_.timeName(),
196  this->mesh_,
199  ),
200  0.5*tr(this->R_)
201  ),
202  epsilon_
203  (
204  IOobject
205  (
206  "epsilon",
207  this->runTime_.timeName(),
208  this->mesh_,
211  ),
212  this->mesh_
213  )
214 {
215  if (type == typeName)
216  {
217  this->printCoeffs(type);
218 
219  this->boundNormalStress(this->R_);
220  bound(epsilon_, this->epsilonMin_);
221  k_ = 0.5*tr(this->R_);
222  }
223 }
224 
225 
226 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
227 
228 template<class BasicTurbulenceModel>
230 {
232  {
233  Cmu_.readIfPresent(this->coeffDict());
234  C1_.readIfPresent(this->coeffDict());
235  C1s_.readIfPresent(this->coeffDict());
236  C2_.readIfPresent(this->coeffDict());
237  C3_.readIfPresent(this->coeffDict());
238  C3s_.readIfPresent(this->coeffDict());
239  C4_.readIfPresent(this->coeffDict());
240  C5_.readIfPresent(this->coeffDict());
241 
242  Ceps1_.readIfPresent(this->coeffDict());
243  Ceps2_.readIfPresent(this->coeffDict());
244  Cs_.readIfPresent(this->coeffDict());
245  Ceps_.readIfPresent(this->coeffDict());
246 
247  return true;
248  }
249 
250  return false;
251 }
252 
253 
254 template<class BasicTurbulenceModel>
256 {
258  (
260  (
261  "DREff",
262  (Cs_*(this->k_/this->epsilon_))*this->R_ + I*this->nu()
263  )
264  );
265 }
266 
267 
268 template<class BasicTurbulenceModel>
270 {
272  (
274  (
275  "DepsilonEff",
276  (Ceps_*(this->k_/this->epsilon_))*this->R_ + I*this->nu()
277  )
278  );
279 }
280 
281 
282 template<class BasicTurbulenceModel>
284 {
285  if (!this->turbulence_)
286  {
287  return;
288  }
289 
290  // Local references
291  const alphaField& alpha = this->alpha_;
292  const rhoField& rho = this->rho_;
293  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
294  const volVectorField& U = this->U_;
295  volSymmTensorField& R = this->R_;
296  fv::options& fvOptions(fv::options::New(this->mesh_));
297 
299 
301  const volTensorField& gradU = tgradU();
302 
303  volSymmTensorField P(-twoSymm(R & gradU));
304  volScalarField G(this->GName(), 0.5*mag(tr(P)));
305 
306  // Update epsilon and G at the wall
307  epsilon_.boundaryFieldRef().updateCoeffs();
308 
309  // Dissipation equation
310  tmp<fvScalarMatrix> epsEqn
311  (
312  fvm::ddt(alpha, rho, epsilon_)
313  + fvm::div(alphaRhoPhi, epsilon_)
314  - fvm::laplacian(alpha*rho*DepsilonEff(), epsilon_)
315  ==
316  Ceps1_*alpha*rho*G*epsilon_/k_
317  - fvm::Sp(Ceps2_*alpha*rho*epsilon_/k_, epsilon_)
318  + fvOptions(alpha, rho, epsilon_)
319  );
320 
321  epsEqn.ref().relax();
322  fvOptions.constrain(epsEqn.ref());
323  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
324  solve(epsEqn);
325  fvOptions.correct(epsilon_);
326  bound(epsilon_, this->epsilonMin_);
327 
328 
329  // Correct the trace of the tensorial production to be consistent
330  // with the near-wall generation from the wall-functions
331  const fvPatchList& patches = this->mesh_.boundary();
332 
333  forAll(patches, patchi)
334  {
335  const fvPatch& curPatch = patches[patchi];
336 
337  if (isA<wallFvPatch>(curPatch))
338  {
339  forAll(curPatch, facei)
340  {
341  label celli = curPatch.faceCells()[facei];
342  P[celli] *= min
343  (
344  G[celli]/(0.5*mag(tr(P[celli])) + SMALL),
345  1.0
346  );
347  }
348  }
349  }
350 
351  volSymmTensorField b(dev(R)/(2*k_));
352  volSymmTensorField S(symm(gradU));
353  volTensorField Omega(skew(gradU));
354 
355  // Reynolds stress equation
357  (
358  fvm::ddt(alpha, rho, R)
359  + fvm::div(alphaRhoPhi, R)
360  - fvm::laplacian(alpha*rho*DREff(), R)
361  + fvm::Sp(((C1_/2)*epsilon_ + (C1s_/2)*G)*alpha*rho/k_, R)
362  ==
363  alpha*rho*P
364  - ((1.0/3.0)*I)*(((2.0 - C1_)*epsilon_ - C1s_*G)*alpha*rho)
365  + (C2_*(alpha*rho*epsilon_))*dev(innerSqr(b))
366  + alpha*rho*k_
367  *(
368  (C3_ - C3s_*mag(b))*dev(S)
369  + C4_*dev(twoSymm(b&S))
370  + C5_*twoSymm(b&Omega)
371  )
372  + fvOptions(alpha, rho, R)
373  );
374 
375  REqn.ref().relax();
376  fvOptions.constrain(REqn.ref());
377  solve(REqn);
378  fvOptions.correct(R);
379 
380  this->boundNormalStress(R);
381 
382  k_ = 0.5*tr(R);
383  bound(k_, this->kMin_);
384 
385  correctNut();
386 
387  // Correct wall shear-stresses when applying wall-functions
388  this->correctWallShearStress(R);
389 }
390 
391 
392 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
393 
394 } // End namespace RASModels
395 } // End namespace Foam
396 
397 // ************************************************************************* //
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::symm
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:84
Foam::IOobject::AUTO_WRITE
Definition: IOobject.H:194
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:54
fvOptions.H
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::fv::optionList::correct
void correct(GeometricField< Type, fvPatchField, volMesh > &field)
Apply correction to field.
Definition: fvOptionListTemplates.C:355
Foam::skew
dimensionedTensor skew(const dimensionedTensor &dt)
Definition: dimensionedTensor.C:138
Foam::constant::universal::G
const dimensionedScalar G
Newtonian constant of gravitation.
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::RASModel
Templated abstract base class for RAS turbulence models.
Definition: RASModel.H:52
Foam::fv::options::New
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present.
Definition: fvOptions.C:103
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
Foam::bound
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:35
Foam::read
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:108
wallFvPatch.H
Foam::innerSqr
dimensionedSymmTensor innerSqr(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:62
rho
rho
Definition: readInitialConditions.H:88
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
nu
volScalarField & nu
Definition: readMechanicalProperties.H:176
R
#define R(A, B, C, D, E, F, K, M)
Foam::ReynoldsStress
Reynolds-stress turbulence model base class.
Definition: ReynoldsStress.H:53
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:227
Foam::fvm::Sp
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
Foam::constant::physicoChemical::b
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
fvOptions
fv::options & fvOptions
Definition: setRegionFluidFields.H:23
Foam::fvm::laplacian
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:48
correct
fvOptions correct(rho)
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
Foam::solve
SolverPerformance< Type > solve(faMatrix< Type > &, Istream &)
Solve returning the solution statistics given convergence tolerance.
SSG.H
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:59
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::fv::options
Finite-volume options.
Definition: fvOptions.H:55
Foam::RASModel::rhoField
BasicTurbulenceModel::rhoField rhoField
Definition: RASModel.H:98
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:42
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::RASModel::alphaField
BasicTurbulenceModel::alphaField alphaField
Definition: RASModel.H:97
Foam::RASModel::transportModel
BasicTurbulenceModel::transportModel transportModel
Definition: RASModel.H:99
U
U
Definition: pEqn.H:72
Foam::fvm::ddt
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:48
Foam::RASModels::SSG::DREff
tmp< volSymmTensorField > DREff() const
Return the effective diffusivity for R.
Definition: SSG.C:255
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::RASModels::SSG::correctNut
virtual void correctNut()
Update the eddy-viscosity.
Definition: SSG.C:43
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
patches
const polyBoundaryMesh & patches
Definition: convertProcessorPatches.H:65
Foam::RASModels::SSG::read
virtual bool read()
Read model coefficients if they have changed.
Definition: SSG.C:229
Foam::fvPatch::faceCells
virtual const labelUList & faceCells() const
Return faceCells.
Definition: fvPatch.C:113
Foam::RASModels::SSG
Speziale, Sarkar and Gatski Reynolds-stress turbulence model for incompressible and compressible flow...
Definition: SSG.H:98
Foam::tr
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
Definition: dimensionedSphericalTensor.C:51
Foam::RASModels::SSG::correct
virtual void correct()
Solve the turbulence equations and correct eddy-Viscosity and.
Definition: SSG.C:283
Foam::GeometricField< vector, fvPatchField, volMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:188
Foam::fvm::div
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:48
Foam::RASModels::SSG::DepsilonEff
tmp< volSymmTensorField > DepsilonEff() const
Return the effective diffusivity for epsilon.
Definition: SSG.C:269
Foam::twoSymm
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:95
Foam::I
static const Identity< scalar > I
Definition: Identity.H:95
Foam::IOobject::MUST_READ
Definition: IOobject.H:185
Foam::dev
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:106