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-2024 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_("Cmu", this->coeffDict(), 0.09),
98  C1_("C1", this->coeffDict(), 3.4),
99  C1s_("C1s", this->coeffDict(), 1.8),
100  C2_("C2", this->coeffDict(), 4.2),
101  C3_("C3", this->coeffDict(), 0.8),
102  C3s_("C3s", this->coeffDict(), 1.3),
103  C4_("C4", this->coeffDict(), 1.25),
104  C5_("C5", this->coeffDict(), 0.4),
105 
106  Ceps1_("Ceps1", this->coeffDict(), 1.44),
107  Ceps2_("Ceps2", this->coeffDict(), 1.92),
108  Cs_("Cs", this->coeffDict(), 0.25),
109  Ceps_("Ceps", this->coeffDict(), 0.15),
110 
111  k_
112  (
113  IOobject
114  (
115  this->groupName("k"),
116  this->runTime_.name(),
117  this->mesh_,
118  IOobject::NO_READ,
119  IOobject::AUTO_WRITE
120  ),
121  0.5*tr(this->R_)
122  ),
123  epsilon_
124  (
125  IOobject
126  (
127  this->groupName("epsilon"),
128  this->runTime_.name(),
129  this->mesh_,
130  IOobject::MUST_READ,
131  IOobject::AUTO_WRITE
132  ),
133  this->mesh_
134  )
135 {
136  if (type == typeName)
137  {
138  this->boundNormalStress(this->R_);
139  boundEpsilon();
140  k_ = 0.5*tr(this->R_);
141  }
142 }
143 
144 
145 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
146 
147 template<class BasicMomentumTransportModel>
149 {
151  {
152  Cmu_.readIfPresent(this->coeffDict());
153  C1_.readIfPresent(this->coeffDict());
154  C1s_.readIfPresent(this->coeffDict());
155  C2_.readIfPresent(this->coeffDict());
156  C3_.readIfPresent(this->coeffDict());
157  C3s_.readIfPresent(this->coeffDict());
158  C4_.readIfPresent(this->coeffDict());
159  C5_.readIfPresent(this->coeffDict());
160 
161  Ceps1_.readIfPresent(this->coeffDict());
162  Ceps2_.readIfPresent(this->coeffDict());
163  Cs_.readIfPresent(this->coeffDict());
164  Ceps_.readIfPresent(this->coeffDict());
165 
166  return true;
167  }
168  else
169  {
170  return false;
171  }
172 }
173 
174 
175 template<class BasicMomentumTransportModel>
177 {
179  (
180  "DREff",
181  (Cs_*(this->k_/this->epsilon_))*this->R_ + I*this->nu()
182  );
183 }
184 
185 
186 template<class BasicMomentumTransportModel>
188 {
190  (
191  "DepsilonEff",
192  (Ceps_*(this->k_/this->epsilon_))*this->R_ + I*this->nu()
193  );
194 }
195 
196 
197 template<class BasicMomentumTransportModel>
199 {
200  if (!this->turbulence_)
201  {
202  return;
203  }
204 
205  // Local references
206  const alphaField& alpha = this->alpha_;
207  const rhoField& rho = this->rho_;
208  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
209  const volVectorField& U = this->U_;
210  volSymmTensorField& R = this->R_;
211  const Foam::fvModels& fvModels(Foam::fvModels::New(this->mesh_));
213  (
214  Foam::fvConstraints::New(this->mesh_)
215  );
216 
218 
220  const volTensorField& gradU = tgradU();
221 
222  volSymmTensorField P(-twoSymm(R & gradU));
223  volScalarField G(this->GName(), 0.5*mag(tr(P)));
224 
225  // Update epsilon and G at the wall
226  epsilon_.boundaryFieldRef().updateCoeffs();
227 
228  // Dissipation equation
229  tmp<fvScalarMatrix> epsEqn
230  (
231  fvm::ddt(alpha, rho, epsilon_)
232  + fvm::div(alphaRhoPhi, epsilon_)
233  - fvm::laplacian(alpha*rho*DepsilonEff(), epsilon_)
234  ==
235  Ceps1_*alpha*rho*G*epsilon_/k_
236  - fvm::Sp(Ceps2_*alpha*rho*epsilon_/k_, epsilon_)
237  + epsilonSource()
238  + fvModels.source(alpha, rho, epsilon_)
239  );
240 
241  epsEqn.ref().relax();
242  fvConstraints.constrain(epsEqn.ref());
243  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
244  solve(epsEqn);
245  fvConstraints.constrain(epsilon_);
246  boundEpsilon();
247 
248 
249  // Correct the trace of the tensorial production to be consistent
250  // with the near-wall generation from the wall-functions
251  const fvPatchList& patches = this->mesh_.boundary();
252 
254  {
255  const fvPatch& curPatch = patches[patchi];
256 
257  if (isA<wallFvPatch>(curPatch))
258  {
259  forAll(curPatch, facei)
260  {
261  label celli = curPatch.faceCells()[facei];
262  P[celli] *= min
263  (
264  G[celli]/(0.5*mag(tr(P[celli])) + small),
265  1.0
266  );
267  }
268  }
269  }
270 
271  volSymmTensorField b(dev(R)/(2*k_));
272  volSymmTensorField S(symm(gradU));
273  volTensorField Omega(skew(gradU));
274 
275  // Reynolds stress equation
277  (
278  fvm::ddt(alpha, rho, R)
279  + fvm::div(alphaRhoPhi, R)
280  - fvm::laplacian(alpha*rho*DREff(), R)
281  + fvm::Sp(((C1_/2)*epsilon_ + (C1s_/2)*G)*alpha*rho/k_, R)
282  ==
283  alpha*rho*P
284  - ((1.0/3.0)*I)*(((2.0 - C1_)*epsilon_ - C1s_*G)*alpha*rho)
285  + (C2_*(alpha*rho*epsilon_))*dev(innerSqr(b))
286  + alpha*rho*k_
287  *(
288  (C3_ - C3s_*mag(b))*dev(S)
289  + C4_*dev(twoSymm(b&S))
290  + C5_*twoSymm(b&Omega)
291  )
292  + this->RSource()
293  + fvModels.source(alpha, rho, R)
294  );
295 
296  REqn.ref().relax();
297  fvConstraints.constrain(REqn.ref());
298  solve(REqn);
300 
301  this->boundNormalStress(R);
302 
303  k_ = 0.5*tr(R);
304 
305  correctNut();
306 
307  // Correct wall shear-stresses when applying wall-functions
308  this->correctWallShearStress(R);
309 }
310 
311 
312 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
313 
314 } // End namespace RASModels
315 } // End namespace Foam
316 
317 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
Generic GeometricField class.
static tmp< GeometricField< Type, GeoMesh, PrimitiveField > > New(const word &name, const Internal &, const PtrList< Patch > &, 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:87
BasicMomentumTransportModel::rhoField rhoField
Definition: RASModel.H:88
tmp< volSymmTensorField > DepsilonEff() const
Return the effective diffusivity for epsilon.
Definition: SSG.C:187
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:198
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:176
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:148
Reynolds-stress turbulence model base class.
static autoPtr< dictionary > New(Istream &)
Construct top-level dictionary on freestore from Istream.
Definition: dictionaryIO.C:103
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:197
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 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.
void skew(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
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
void dev(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
void twoSymm(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
static const Identity< scalar > I
Definition: Identity.H:93
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
const dimensionSet dimTime
void symm(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
const dimensionSet dimVolume
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
static scalar R(const scalar a, const scalar x)
Definition: invIncGamma.C:102
void sqr(LagrangianPatchField< typename outerProduct< Type, Type >::type > &f, const LagrangianPatchField< Type > &f1)
void tr(LagrangianPatchField< scalar > &f, const LagrangianPatchField< tensor > &f1)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
dimensionedSymmTensor innerSqr(const dimensionedSymmTensor &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.