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  this->nut_ = this->Cmu_*sqr(k_)/epsilon_;
44  this->nut_.correctBoundaryConditions();
45  fvConstraints::New(this->mesh_).constrain(this->nut_);
46 }
47 
48 
49 template<class BasicMomentumTransportModel>
51 {
52  return tmp<fvScalarMatrix>
53  (
54  new fvScalarMatrix
55  (
56  epsilon_,
57  dimVolume*this->rho_.dimensions()*epsilon_.dimensions()/dimTime
58  )
59  );
60 }
61 
62 
63 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
64 
65 template<class BasicMomentumTransportModel>
67 (
68  const alphaField& alpha,
69  const rhoField& rho,
70  const volVectorField& U,
71  const surfaceScalarField& alphaRhoPhi,
72  const surfaceScalarField& phi,
73  const viscosity& viscosity,
74  const word& type
75 )
76 :
77  ReynoldsStress<RASModel<BasicMomentumTransportModel>>
78  (
79  type,
80  alpha,
81  rho,
82  U,
83  alphaRhoPhi,
84  phi,
85  viscosity
86  ),
87 
88  Cmu_
89  (
90  dimensioned<scalar>::lookupOrAddToDict
91  (
92  "Cmu",
93  this->coeffDict_,
94  0.09
95  )
96  ),
97  C1_
98  (
99  dimensioned<scalar>::lookupOrAddToDict
100  (
101  "C1",
102  this->coeffDict_,
103  3.4
104  )
105  ),
106  C1s_
107  (
108  dimensioned<scalar>::lookupOrAddToDict
109  (
110  "C1s",
111  this->coeffDict_,
112  1.8
113  )
114  ),
115  C2_
116  (
117  dimensioned<scalar>::lookupOrAddToDict
118  (
119  "C2",
120  this->coeffDict_,
121  4.2
122  )
123  ),
124  C3_
125  (
126  dimensioned<scalar>::lookupOrAddToDict
127  (
128  "C3",
129  this->coeffDict_,
130  0.8
131  )
132  ),
133  C3s_
134  (
135  dimensioned<scalar>::lookupOrAddToDict
136  (
137  "C3s",
138  this->coeffDict_,
139  1.3
140  )
141  ),
142  C4_
143  (
144  dimensioned<scalar>::lookupOrAddToDict
145  (
146  "C4",
147  this->coeffDict_,
148  1.25
149  )
150  ),
151  C5_
152  (
153  dimensioned<scalar>::lookupOrAddToDict
154  (
155  "C5",
156  this->coeffDict_,
157  0.4
158  )
159  ),
160 
161  Ceps1_
162  (
163  dimensioned<scalar>::lookupOrAddToDict
164  (
165  "Ceps1",
166  this->coeffDict_,
167  1.44
168  )
169  ),
170  Ceps2_
171  (
172  dimensioned<scalar>::lookupOrAddToDict
173  (
174  "Ceps2",
175  this->coeffDict_,
176  1.92
177  )
178  ),
179  Cs_
180  (
181  dimensioned<scalar>::lookupOrAddToDict
182  (
183  "Cs",
184  this->coeffDict_,
185  0.25
186  )
187  ),
188  Ceps_
189  (
190  dimensioned<scalar>::lookupOrAddToDict
191  (
192  "Ceps",
193  this->coeffDict_,
194  0.15
195  )
196  ),
197 
198  k_
199  (
200  IOobject
201  (
202  this->groupName("k"),
203  this->runTime_.name(),
204  this->mesh_,
205  IOobject::NO_READ,
206  IOobject::AUTO_WRITE
207  ),
208  0.5*tr(this->R_)
209  ),
210  epsilon_
211  (
212  IOobject
213  (
214  this->groupName("epsilon"),
215  this->runTime_.name(),
216  this->mesh_,
217  IOobject::MUST_READ,
218  IOobject::AUTO_WRITE
219  ),
220  this->mesh_
221  )
222 {
223  if (type == typeName)
224  {
225  this->printCoeffs(type);
226 
227  this->boundNormalStress(this->R_);
228  bound(epsilon_, this->epsilonMin_);
229  k_ = 0.5*tr(this->R_);
230  }
231 }
232 
233 
234 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
235 
236 template<class BasicMomentumTransportModel>
238 {
240  {
241  Cmu_.readIfPresent(this->coeffDict());
242  C1_.readIfPresent(this->coeffDict());
243  C1s_.readIfPresent(this->coeffDict());
244  C2_.readIfPresent(this->coeffDict());
245  C3_.readIfPresent(this->coeffDict());
246  C3s_.readIfPresent(this->coeffDict());
247  C4_.readIfPresent(this->coeffDict());
248  C5_.readIfPresent(this->coeffDict());
249 
250  Ceps1_.readIfPresent(this->coeffDict());
251  Ceps2_.readIfPresent(this->coeffDict());
252  Cs_.readIfPresent(this->coeffDict());
253  Ceps_.readIfPresent(this->coeffDict());
254 
255  return true;
256  }
257  else
258  {
259  return false;
260  }
261 }
262 
263 
264 template<class BasicMomentumTransportModel>
266 {
268  (
269  "DREff",
270  (Cs_*(this->k_/this->epsilon_))*this->R_ + I*this->nu()
271  );
272 }
273 
274 
275 template<class BasicMomentumTransportModel>
277 {
279  (
280  "DepsilonEff",
281  (Ceps_*(this->k_/this->epsilon_))*this->R_ + I*this->nu()
282  );
283 }
284 
285 
286 template<class BasicMomentumTransportModel>
288 {
289  if (!this->turbulence_)
290  {
291  return;
292  }
293 
294  // Local references
295  const alphaField& alpha = this->alpha_;
296  const rhoField& rho = this->rho_;
297  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
298  const volVectorField& U = this->U_;
299  volSymmTensorField& R = this->R_;
300  const Foam::fvModels& fvModels(Foam::fvModels::New(this->mesh_));
302  (
303  Foam::fvConstraints::New(this->mesh_)
304  );
305 
307 
309  const volTensorField& gradU = tgradU();
310 
311  volSymmTensorField P(-twoSymm(R & gradU));
312  volScalarField G(this->GName(), 0.5*mag(tr(P)));
313 
314  // Update epsilon and G at the wall
315  epsilon_.boundaryFieldRef().updateCoeffs();
316 
317  // Dissipation equation
318  tmp<fvScalarMatrix> epsEqn
319  (
320  fvm::ddt(alpha, rho, epsilon_)
321  + fvm::div(alphaRhoPhi, epsilon_)
322  - fvm::laplacian(alpha*rho*DepsilonEff(), epsilon_)
323  ==
324  Ceps1_*alpha*rho*G*epsilon_/k_
325  - fvm::Sp(Ceps2_*alpha*rho*epsilon_/k_, epsilon_)
326  + epsilonSource()
327  + fvModels.source(alpha, rho, epsilon_)
328  );
329 
330  epsEqn.ref().relax();
331  fvConstraints.constrain(epsEqn.ref());
332  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
333  solve(epsEqn);
334  fvConstraints.constrain(epsilon_);
335  bound(epsilon_, this->epsilonMin_);
336 
337 
338  // Correct the trace of the tensorial production to be consistent
339  // with the near-wall generation from the wall-functions
340  const fvPatchList& patches = this->mesh_.boundary();
341 
343  {
344  const fvPatch& curPatch = patches[patchi];
345 
346  if (isA<wallFvPatch>(curPatch))
347  {
348  forAll(curPatch, facei)
349  {
350  label celli = curPatch.faceCells()[facei];
351  P[celli] *= min
352  (
353  G[celli]/(0.5*mag(tr(P[celli])) + small),
354  1.0
355  );
356  }
357  }
358  }
359 
360  volSymmTensorField b(dev(R)/(2*k_));
361  volSymmTensorField S(symm(gradU));
362  volTensorField Omega(skew(gradU));
363 
364  // Reynolds stress equation
366  (
367  fvm::ddt(alpha, rho, R)
368  + fvm::div(alphaRhoPhi, R)
369  - fvm::laplacian(alpha*rho*DREff(), R)
370  + fvm::Sp(((C1_/2)*epsilon_ + (C1s_/2)*G)*alpha*rho/k_, R)
371  ==
372  alpha*rho*P
373  - ((1.0/3.0)*I)*(((2.0 - C1_)*epsilon_ - C1s_*G)*alpha*rho)
374  + (C2_*(alpha*rho*epsilon_))*dev(innerSqr(b))
375  + alpha*rho*k_
376  *(
377  (C3_ - C3s_*mag(b))*dev(S)
378  + C4_*dev(twoSymm(b&S))
379  + C5_*twoSymm(b&Omega)
380  )
381  + this->RSource()
382  + fvModels.source(alpha, rho, R)
383  );
384 
385  REqn.ref().relax();
386  fvConstraints.constrain(REqn.ref());
387  solve(REqn);
389 
390  this->boundNormalStress(R);
391 
392  k_ = 0.5*tr(R);
393 
394  correctNut();
395 
396  // Correct wall shear-stresses when applying wall-functions
397  this->correctWallShearStress(R);
398 }
399 
400 
401 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
402 
403 } // End namespace RASModels
404 } // End namespace Foam
405 
406 // ************************************************************************* //
#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 >> &)
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:96
virtual void printCoeffs(const word &type)
Print model coefficients.
Definition: RASModel.C:32
dimensionedScalar epsilonMin_
Lower limit of epsilon.
Definition: RASModel.H:78
BasicMomentumTransportModel::rhoField rhoField
Definition: RASModel.H:97
volScalarField epsilon_
Definition: SSG.H:120
tmp< volSymmTensorField > DepsilonEff() const
Return the effective diffusivity for epsilon.
Definition: SSG.C:276
virtual tmp< fvScalarMatrix > epsilonSource() const
Source term for the epsilon equation.
Definition: SSG.C:50
volScalarField k_
Definition: SSG.H:119
virtual void correct()
Solve the turbulence equations and correct eddy-Viscosity and.
Definition: SSG.C:287
tmp< volSymmTensorField > DREff() const
Return the effective diffusivity for R.
Definition: SSG.C:265
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:67
virtual void correctNut()
Update the eddy-viscosity.
Definition: SSG.C:41
virtual bool read()
Read model coefficients if they have changed.
Definition: SSG.C:237
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:62
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
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:27
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 > > Sp(const volScalarField::Internal &, const VolField< Type > &)
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const VolField< Type > &vf, const word &name)
Definition: fvmDiv.C:48
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)
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 > &)
bool bound(volScalarField &, const dimensionedScalar &min)
Bound the given scalar field where it is below the specified min value.
Definition: bound.C:31
static scalar R(const scalar a, const scalar x)
Definition: invIncGamma.C:102
dimensionedSymmTensor innerSqr(const dimensionedSymmTensor &dt)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
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.