SSG.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) 2015-2023 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 "fvModels.H"
28 #include "fvConstraints.H"
29 #include "wallFvPatch.H"
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace RASModels
36 {
37 
38 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
39 
40 template<class BasicMomentumTransportModel>
42 {
43  tmp<volScalarField> tCmuk2(Cmu_*sqr(k_));
44  epsilon_ = max(epsilon_, tCmuk2()/(this->nutMaxCoeff_*this->nu()));
45  return tCmuk2;
46 }
47 
48 
49 template<class BasicMomentumTransportModel>
51 {
52  this->nut_ = boundEpsilon()/epsilon_;
53  this->nut_.correctBoundaryConditions();
54  fvConstraints::New(this->mesh_).constrain(this->nut_);
55 }
56 
57 
58 template<class BasicMomentumTransportModel>
60 {
61  return tmp<fvScalarMatrix>
62  (
63  new fvScalarMatrix
64  (
65  epsilon_,
66  dimVolume*this->rho_.dimensions()*epsilon_.dimensions()/dimTime
67  )
68  );
69 }
70 
71 
72 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
73 
74 template<class BasicMomentumTransportModel>
76 (
77  const alphaField& alpha,
78  const rhoField& rho,
79  const volVectorField& U,
80  const surfaceScalarField& alphaRhoPhi,
81  const surfaceScalarField& phi,
82  const viscosity& viscosity,
83  const word& type
84 )
85 :
86  ReynoldsStress<RASModel<BasicMomentumTransportModel>>
87  (
88  type,
89  alpha,
90  rho,
91  U,
92  alphaRhoPhi,
93  phi,
94  viscosity
95  ),
96 
97  Cmu_
98  (
99  dimensioned<scalar>::lookupOrAddToDict
100  (
101  "Cmu",
102  this->coeffDict_,
103  0.09
104  )
105  ),
106  C1_
107  (
108  dimensioned<scalar>::lookupOrAddToDict
109  (
110  "C1",
111  this->coeffDict_,
112  3.4
113  )
114  ),
115  C1s_
116  (
117  dimensioned<scalar>::lookupOrAddToDict
118  (
119  "C1s",
120  this->coeffDict_,
121  1.8
122  )
123  ),
124  C2_
125  (
126  dimensioned<scalar>::lookupOrAddToDict
127  (
128  "C2",
129  this->coeffDict_,
130  4.2
131  )
132  ),
133  C3_
134  (
135  dimensioned<scalar>::lookupOrAddToDict
136  (
137  "C3",
138  this->coeffDict_,
139  0.8
140  )
141  ),
142  C3s_
143  (
144  dimensioned<scalar>::lookupOrAddToDict
145  (
146  "C3s",
147  this->coeffDict_,
148  1.3
149  )
150  ),
151  C4_
152  (
153  dimensioned<scalar>::lookupOrAddToDict
154  (
155  "C4",
156  this->coeffDict_,
157  1.25
158  )
159  ),
160  C5_
161  (
162  dimensioned<scalar>::lookupOrAddToDict
163  (
164  "C5",
165  this->coeffDict_,
166  0.4
167  )
168  ),
169 
170  Ceps1_
171  (
172  dimensioned<scalar>::lookupOrAddToDict
173  (
174  "Ceps1",
175  this->coeffDict_,
176  1.44
177  )
178  ),
179  Ceps2_
180  (
181  dimensioned<scalar>::lookupOrAddToDict
182  (
183  "Ceps2",
184  this->coeffDict_,
185  1.92
186  )
187  ),
188  Cs_
189  (
190  dimensioned<scalar>::lookupOrAddToDict
191  (
192  "Cs",
193  this->coeffDict_,
194  0.25
195  )
196  ),
197  Ceps_
198  (
199  dimensioned<scalar>::lookupOrAddToDict
200  (
201  "Ceps",
202  this->coeffDict_,
203  0.15
204  )
205  ),
206 
207  k_
208  (
209  IOobject
210  (
211  this->groupName("k"),
212  this->runTime_.name(),
213  this->mesh_,
214  IOobject::NO_READ,
215  IOobject::AUTO_WRITE
216  ),
217  0.5*tr(this->R_)
218  ),
219  epsilon_
220  (
221  IOobject
222  (
223  this->groupName("epsilon"),
224  this->runTime_.name(),
225  this->mesh_,
226  IOobject::MUST_READ,
227  IOobject::AUTO_WRITE
228  ),
229  this->mesh_
230  )
231 {
232  if (type == typeName)
233  {
234  this->printCoeffs(type);
235 
236  this->boundNormalStress(this->R_);
237  boundEpsilon();
238  k_ = 0.5*tr(this->R_);
239  }
240 }
241 
242 
243 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
244 
245 template<class BasicMomentumTransportModel>
247 {
249  {
250  Cmu_.readIfPresent(this->coeffDict());
251  C1_.readIfPresent(this->coeffDict());
252  C1s_.readIfPresent(this->coeffDict());
253  C2_.readIfPresent(this->coeffDict());
254  C3_.readIfPresent(this->coeffDict());
255  C3s_.readIfPresent(this->coeffDict());
256  C4_.readIfPresent(this->coeffDict());
257  C5_.readIfPresent(this->coeffDict());
258 
259  Ceps1_.readIfPresent(this->coeffDict());
260  Ceps2_.readIfPresent(this->coeffDict());
261  Cs_.readIfPresent(this->coeffDict());
262  Ceps_.readIfPresent(this->coeffDict());
263 
264  return true;
265  }
266  else
267  {
268  return false;
269  }
270 }
271 
272 
273 template<class BasicMomentumTransportModel>
275 {
277  (
278  "DREff",
279  (Cs_*(this->k_/this->epsilon_))*this->R_ + I*this->nu()
280  );
281 }
282 
283 
284 template<class BasicMomentumTransportModel>
286 {
288  (
289  "DepsilonEff",
290  (Ceps_*(this->k_/this->epsilon_))*this->R_ + I*this->nu()
291  );
292 }
293 
294 
295 template<class BasicMomentumTransportModel>
297 {
298  if (!this->turbulence_)
299  {
300  return;
301  }
302 
303  // Local references
304  const alphaField& alpha = this->alpha_;
305  const rhoField& rho = this->rho_;
306  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
307  const volVectorField& U = this->U_;
308  volSymmTensorField& R = this->R_;
309  const Foam::fvModels& fvModels(Foam::fvModels::New(this->mesh_));
311  (
312  Foam::fvConstraints::New(this->mesh_)
313  );
314 
316 
318  const volTensorField& gradU = tgradU();
319 
320  volSymmTensorField P(-twoSymm(R & gradU));
321  volScalarField G(this->GName(), 0.5*mag(tr(P)));
322 
323  // Update epsilon and G at the wall
324  epsilon_.boundaryFieldRef().updateCoeffs();
325 
326  // Dissipation equation
327  tmp<fvScalarMatrix> epsEqn
328  (
329  fvm::ddt(alpha, rho, epsilon_)
330  + fvm::div(alphaRhoPhi, epsilon_)
331  - fvm::laplacian(alpha*rho*DepsilonEff(), epsilon_)
332  ==
333  Ceps1_*alpha*rho*G*epsilon_/k_
334  - fvm::Sp(Ceps2_*alpha*rho*epsilon_/k_, epsilon_)
335  + epsilonSource()
336  + fvModels.source(alpha, rho, epsilon_)
337  );
338 
339  epsEqn.ref().relax();
340  fvConstraints.constrain(epsEqn.ref());
341  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
342  solve(epsEqn);
343  fvConstraints.constrain(epsilon_);
344  boundEpsilon();
345 
346 
347  // Correct the trace of the tensorial production to be consistent
348  // with the near-wall generation from the wall-functions
349  const fvPatchList& patches = this->mesh_.boundary();
350 
352  {
353  const fvPatch& curPatch = patches[patchi];
354 
355  if (isA<wallFvPatch>(curPatch))
356  {
357  forAll(curPatch, facei)
358  {
359  label celli = curPatch.faceCells()[facei];
360  P[celli] *= min
361  (
362  G[celli]/(0.5*mag(tr(P[celli])) + small),
363  1.0
364  );
365  }
366  }
367  }
368 
369  volSymmTensorField b(dev(R)/(2*k_));
370  volSymmTensorField S(symm(gradU));
371  volTensorField Omega(skew(gradU));
372 
373  // Reynolds stress equation
375  (
376  fvm::ddt(alpha, rho, R)
377  + fvm::div(alphaRhoPhi, R)
378  - fvm::laplacian(alpha*rho*DREff(), R)
379  + fvm::Sp(((C1_/2)*epsilon_ + (C1s_/2)*G)*alpha*rho/k_, R)
380  ==
381  alpha*rho*P
382  - ((1.0/3.0)*I)*(((2.0 - C1_)*epsilon_ - C1s_*G)*alpha*rho)
383  + (C2_*(alpha*rho*epsilon_))*dev(innerSqr(b))
384  + alpha*rho*k_
385  *(
386  (C3_ - C3s_*mag(b))*dev(S)
387  + C4_*dev(twoSymm(b&S))
388  + C5_*twoSymm(b&Omega)
389  )
390  + this->RSource()
391  + fvModels.source(alpha, rho, R)
392  );
393 
394  REqn.ref().relax();
395  fvConstraints.constrain(REqn.ref());
396  solve(REqn);
398 
399  this->boundNormalStress(R);
400 
401  k_ = 0.5*tr(R);
402 
403  correctNut();
404 
405  // Correct wall shear-stresses when applying wall-functions
406  this->correctWallShearStress(R);
407 }
408 
409 
410 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
411 
412 } // End namespace RASModels
413 } // End namespace Foam
414 
415 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Generic GeometricField class.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
Templated abstract base class for RAS turbulence models.
Definition: RASModel.H:56
BasicMomentumTransportModel::alphaField alphaField
Definition: RASModel.H:93
virtual void printCoeffs(const word &type)
Print model coefficients.
Definition: RASModel.C:32
BasicMomentumTransportModel::rhoField rhoField
Definition: RASModel.H:94
tmp< volSymmTensorField > DepsilonEff() const
Return the effective diffusivity for epsilon.
Definition: SSG.C:285
virtual tmp< fvScalarMatrix > epsilonSource() const
Source term for the epsilon equation.
Definition: SSG.C:59
volScalarField k_
Definition: SSG.H:119
virtual void correct()
Solve the turbulence equations and correct eddy-Viscosity and.
Definition: SSG.C:296
tmp< volScalarField > boundEpsilon()
Bound epsilon and return Cmu*sqr(k) for nut.
Definition: SSG.C:41
tmp< volSymmTensorField > DREff() const
Return the effective diffusivity for R.
Definition: SSG.C:274
SSG(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity, const word &type=typeName)
Construct from components.
Definition: SSG.C:76
virtual void correctNut()
Correct the eddy-viscosity nut.
Definition: SSG.C:50
virtual bool read()
Read model coefficients if they have changed.
Definition: SSG.C:246
Reynolds-stress turbulence model base class.
static autoPtr< dictionary > New(Istream &)
Construct top-level dictionary on freestore from Istream.
Definition: dictionaryIO.C:100
Generic dimensioned Type class.
Finite volume constraints.
Definition: fvConstraints.H:67
bool constrain(fvMatrix< Type > &eqn) const
Apply constraints to an equation.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
Finite volume models.
Definition: fvModels.H:65
tmp< fvMatrix< Type > > source(const VolField< Type > &field) const
Return source for an equation.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
virtual const labelUList & faceCells() const
Return faceCells.
Definition: fvPatch.C:99
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
Abstract base class for all fluid physical properties.
Definition: viscosity.H:50
A class for handling words, derived from string.
Definition: word.H:62
Foam::fvConstraints & fvConstraints(Foam::fvConstraints::New(mesh))
Foam::fvModels & fvModels(Foam::fvModels::New(mesh))
label patchi
const fvPatchList & patches
volScalarField & b
Definition: createFields.H:25
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
void correct(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiCorr, const SpType &Sp, const SuType &Su)
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.
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
bool read(const char *, int32_t &)
Definition: int32IO.C:85
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
dimensionedSymmTensor sqr(const dimensionedVector &dv)
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
static const Identity< scalar > I
Definition: Identity.H:93
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
const dimensionSet dimTime
const dimensionSet dimVolume
dimensioned< scalar > mag(const dimensioned< Type > &)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
static scalar R(const scalar a, const scalar x)
Definition: invIncGamma.C:102
dimensionedSymmTensor innerSqr(const dimensionedSymmTensor &dt)
dimensionedTensor skew(const dimensionedTensor &dt)
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.