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 | 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 "SSG.H"
27 #include "wallFvPatch.H"
28 
29 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace RASModels
34 {
35 
36 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
37 
38 template<class BasicTurbulenceModel>
40 {
41  this->nut_ = this->Cmu_*sqr(k_)/epsilon_;
42  this->nut_.correctBoundaryConditions();
43 
44  BasicTurbulenceModel::correctNut();
45 }
46 
47 
48 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
49 
50 template<class BasicTurbulenceModel>
52 (
53  const alphaField& alpha,
54  const rhoField& rho,
55  const volVectorField& U,
56  const surfaceScalarField& alphaRhoPhi,
57  const surfaceScalarField& phi,
58  const transportModel& transport,
59  const word& propertiesName,
60  const word& type
61 )
62 :
64  (
65  type,
66  alpha,
67  rho,
68  U,
69  alphaRhoPhi,
70  phi,
71  transport,
72  propertiesName
73  ),
74 
75  Cmu_
76  (
78  (
79  "Cmu",
80  this->coeffDict_,
81  0.09
82  )
83  ),
84  C1_
85  (
87  (
88  "C1",
89  this->coeffDict_,
90  3.4
91  )
92  ),
93  C1s_
94  (
96  (
97  "C1s",
98  this->coeffDict_,
99  1.8
100  )
101  ),
102  C2_
103  (
105  (
106  "C2",
107  this->coeffDict_,
108  4.2
109  )
110  ),
111  C3_
112  (
114  (
115  "C3",
116  this->coeffDict_,
117  0.8
118  )
119  ),
120  C3s_
121  (
123  (
124  "C3s",
125  this->coeffDict_,
126  1.3
127  )
128  ),
129  C4_
130  (
132  (
133  "C4",
134  this->coeffDict_,
135  1.25
136  )
137  ),
138  C5_
139  (
141  (
142  "C5",
143  this->coeffDict_,
144  0.4
145  )
146  ),
147 
148  Ceps1_
149  (
151  (
152  "Ceps1",
153  this->coeffDict_,
154  1.44
155  )
156  ),
157  Ceps2_
158  (
160  (
161  "Ceps2",
162  this->coeffDict_,
163  1.92
164  )
165  ),
166  Cs_
167  (
169  (
170  "Cs",
171  this->coeffDict_,
172  0.25
173  )
174  ),
175  Ceps_
176  (
178  (
179  "Ceps",
180  this->coeffDict_,
181  0.15
182  )
183  ),
184 
185  k_
186  (
187  IOobject
188  (
189  "k",
190  this->runTime_.timeName(),
191  this->mesh_,
194  ),
195  0.5*tr(this->R_)
196  ),
197  epsilon_
198  (
199  IOobject
200  (
201  "epsilon",
202  this->runTime_.timeName(),
203  this->mesh_,
206  ),
207  this->mesh_
208  )
209 {
210  if (type == typeName)
211  {
212  this->printCoeffs(type);
213 
214  this->boundNormalStress(this->R_);
215  bound(epsilon_, this->epsilonMin_);
216  k_ = 0.5*tr(this->R_);
217  }
218 }
219 
220 
221 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
222 
223 template<class BasicTurbulenceModel>
225 {
227  {
228  Cmu_.readIfPresent(this->coeffDict());
229  C1_.readIfPresent(this->coeffDict());
230  C1s_.readIfPresent(this->coeffDict());
231  C2_.readIfPresent(this->coeffDict());
232  C3_.readIfPresent(this->coeffDict());
233  C3s_.readIfPresent(this->coeffDict());
234  C4_.readIfPresent(this->coeffDict());
235  C5_.readIfPresent(this->coeffDict());
236 
237  Ceps1_.readIfPresent(this->coeffDict());
238  Ceps2_.readIfPresent(this->coeffDict());
239  Cs_.readIfPresent(this->coeffDict());
240  Ceps_.readIfPresent(this->coeffDict());
241 
242  return true;
243  }
244  else
245  {
246  return false;
247  }
248 }
249 
250 
251 template<class BasicTurbulenceModel>
253 {
255  (
257  (
258  "DREff",
259  (Cs_*(this->k_/this->epsilon_))*this->R_ + I*this->nu()
260  )
261  );
262 }
263 
264 
265 template<class BasicTurbulenceModel>
267 {
269  (
271  (
272  "DepsilonEff",
273  (Ceps_*(this->k_/this->epsilon_))*this->R_ + I*this->nu()
274  )
275  );
276 }
277 
278 
279 template<class BasicTurbulenceModel>
281 {
282  if (!this->turbulence_)
283  {
284  return;
285  }
286 
287  // Local references
288  const alphaField& alpha = this->alpha_;
289  const rhoField& rho = this->rho_;
290  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
291  const volVectorField& U = this->U_;
292  volSymmTensorField& R = this->R_;
293 
295 
297  const volTensorField& gradU = tgradU();
298 
299  volSymmTensorField P(-twoSymm(R & gradU));
300  volScalarField G(this->GName(), 0.5*mag(tr(P)));
301 
302  // Update epsilon and G at the wall
303  epsilon_.boundaryField().updateCoeffs();
304 
305  // Dissipation equation
306  tmp<fvScalarMatrix> epsEqn
307  (
308  fvm::ddt(alpha, rho, epsilon_)
309  + fvm::div(alphaRhoPhi, epsilon_)
310  - fvm::laplacian(alpha*rho*DepsilonEff(), epsilon_)
311  ==
312  Ceps1_*alpha*rho*G*epsilon_/k_
313  - fvm::Sp(Ceps2_*alpha*rho*epsilon_/k_, epsilon_)
314  );
315 
316  epsEqn().relax();
317 
318  epsEqn().boundaryManipulate(epsilon_.boundaryField());
319 
320  solve(epsEqn);
321  bound(epsilon_, this->epsilonMin_);
322 
323 
324  // Correct the trace of the tensorial production to be consistent
325  // with the near-wall generation from the wall-functions
326  const fvPatchList& patches = this->mesh_.boundary();
327 
329  {
330  const fvPatch& curPatch = patches[patchi];
331 
332  if (isA<wallFvPatch>(curPatch))
333  {
334  forAll(curPatch, facei)
335  {
336  label faceCelli = curPatch.faceCells()[facei];
337  P[faceCelli] *= min
338  (
339  G[faceCelli]/(0.5*mag(tr(P[faceCelli])) + SMALL),
340  1.0
341  );
342  }
343  }
344  }
345 
346  volSymmTensorField b(dev(R)/(2*k_));
347  volSymmTensorField S(symm(gradU));
348  volTensorField Omega(skew(gradU));
349 
350  // Reynolds stress equation
352  (
353  fvm::ddt(alpha, rho, R)
354  + fvm::div(alphaRhoPhi, R)
355  - fvm::laplacian(alpha*rho*DREff(), R)
356  + fvm::Sp(((C1_/2)*epsilon_ + (C1s_/2)*G)*alpha*rho/k_, R)
357  ==
358  alpha*rho*P
359  - ((1.0/3.0)*I)*(((2.0 - C1_)*epsilon_ - C1s_*G)*alpha*rho)
360  + (C2_*(alpha*rho*epsilon_))*dev(innerSqr(b))
361  + alpha*rho*k_
362  *(
363  (C3_ - C3s_*mag(b))*dev(S)
364  + C4_*dev(twoSymm(b&S))
365  + C5_*twoSymm(b&Omega)
366  )
367  );
368 
369  REqn().relax();
370  solve(REqn);
371 
372  this->boundNormalStress(R);
373 
374  k_ = 0.5*tr(R);
375 
376  correctNut();
377 
378  // Correct wall shear-stresses when applying wall-functions
379  this->correctWallShearStress(R);
380 }
381 
382 
383 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
384 
385 } // End namespace RASModels
386 } // End namespace Foam
387 
388 // ************************************************************************* //
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::symm
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:82
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::skew
dimensionedTensor skew(const dimensionedTensor &dt)
Definition: dimensionedTensor.C:135
Foam::IOobject::AUTO_WRITE
@ AUTO_WRITE
Definition: IOobject.H:117
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::constant::universal::G
const dimensionedScalar G
Newtonian constant of gravitation.
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::RASModel
Templated abstract base class for RAS turbulence models.
Definition: RASModel.H:49
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
wallFvPatch.H
Foam::IOobject::MUST_READ
@ MUST_READ
Definition: IOobject.H:108
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::innerSqr
dimensionedSymmTensor innerSqr(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:60
U
U
Definition: pEqn.H:46
Foam::fvm::Sp
tmp< fvMatrix< Type > > Sp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
nu
volScalarField & nu
Definition: readMechanicalProperties.H:179
R
#define R(A, B, C, D, E, F, K, M)
Foam::ReynoldsStress
Reynolds-stress turbulence model base class.
Definition: ReynoldsStress.H:50
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Foam::constant::physicoChemical::b
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:28
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
Foam::fvm::laplacian
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:46
correct
fvOptions correct(rho)
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
SSG.H
Foam::PtrList
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:61
Foam::I
static const sphericalTensor I(1)
Foam::RASModel::rhoField
BasicTurbulenceModel::rhoField rhoField
Definition: RASModel.H:100
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::RASModel::alphaField
BasicTurbulenceModel::alphaField alphaField
Definition: RASModel.H:99
rho
rho
Definition: pEqn.H:3
Foam::RASModel::transportModel
BasicTurbulenceModel::transportModel transportModel
Definition: RASModel.H:101
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::RASModels::SSG::DREff
tmp< volSymmTensorField > DREff() const
Return the effective diffusivity for R.
Definition: SSG.C:252
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:49
Foam::RASModels::SSG::correctNut
virtual void correctNut()
Update the eddy-viscosity.
Definition: SSG.C:39
patchi
label patchi
Definition: getPatchFieldScalar.H:1
Foam::RASModels::SSG::read
virtual bool read()
Read model coefficients if they have changed.
Definition: SSG.C:224
Foam::fvPatch::faceCells
virtual const labelUList & faceCells() const
Return faceCells.
Definition: fvPatch.C:93
patches
patches[0]
Definition: createSingleCellMesh.H:36
Foam::tr
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
Definition: dimensionedSphericalTensor.C:49
Foam::RASModels::SSG::correct
virtual void correct()
Solve the turbulence equations and correct eddy-Viscosity and.
Definition: SSG.C:280
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::min
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::RASModels::SSG::DepsilonEff
tmp< volSymmTensorField > DepsilonEff() const
Return the effective diffusivity for epsilon.
Definition: SSG.C:266
Foam::twoSymm
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:93
Foam::RASModels::SSG::SSG
SSG(const SSG &)
Foam::dev
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:104