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-2023 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 :
68  ReynoldsStress<RASModel<BasicTurbulenceModel>>
69  (
70  type,
71  alpha,
72  rho,
73  U,
74  alphaRhoPhi,
75  phi,
76  transport,
77  propertiesName
78  ),
79 
80  Cmu_
81  (
82  dimensioned<scalar>::getOrAddToDict
83  (
84  "Cmu",
85  this->coeffDict_,
86  0.09
87  )
88  ),
89  C1_
90  (
91  dimensioned<scalar>::getOrAddToDict
92  (
93  "C1",
94  this->coeffDict_,
95  3.4
96  )
97  ),
98  C1s_
99  (
100  dimensioned<scalar>::getOrAddToDict
101  (
102  "C1s",
103  this->coeffDict_,
104  1.8
105  )
106  ),
107  C2_
108  (
109  dimensioned<scalar>::getOrAddToDict
110  (
111  "C2",
112  this->coeffDict_,
113  4.2
114  )
115  ),
116  C3_
117  (
118  dimensioned<scalar>::getOrAddToDict
119  (
120  "C3",
121  this->coeffDict_,
122  0.8
123  )
124  ),
125  C3s_
126  (
127  dimensioned<scalar>::getOrAddToDict
128  (
129  "C3s",
130  this->coeffDict_,
131  1.3
132  )
133  ),
134  C4_
135  (
136  dimensioned<scalar>::getOrAddToDict
137  (
138  "C4",
139  this->coeffDict_,
140  1.25
141  )
142  ),
143  C5_
144  (
145  dimensioned<scalar>::getOrAddToDict
146  (
147  "C5",
148  this->coeffDict_,
149  0.4
150  )
151  ),
152 
153  Ceps1_
154  (
155  dimensioned<scalar>::getOrAddToDict
156  (
157  "Ceps1",
158  this->coeffDict_,
159  1.44
160  )
161  ),
162  Ceps2_
163  (
164  dimensioned<scalar>::getOrAddToDict
165  (
166  "Ceps2",
167  this->coeffDict_,
168  1.92
169  )
170  ),
171  Cs_
172  (
173  dimensioned<scalar>::getOrAddToDict
174  (
175  "Cs",
176  this->coeffDict_,
177  0.25
178  )
179  ),
180  Ceps_
181  (
182  dimensioned<scalar>::getOrAddToDict
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_,
197  IOobject::NO_READ,
198  IOobject::AUTO_WRITE
199  ),
200  0.5*tr(this->R_)
201  ),
202  epsilon_
203  (
204  IOobject
205  (
206  "epsilon",
207  this->runTime_.timeName(),
208  this->mesh_,
209  IOobject::MUST_READ,
210  IOobject::AUTO_WRITE
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  // Push any changed cell values to coupled neighbours
309  epsilon_.boundaryFieldRef().template evaluateCoupled<coupledFvPatch>();
310 
311  // Dissipation equation
312  tmp<fvScalarMatrix> epsEqn
313  (
314  fvm::ddt(alpha, rho, epsilon_)
315  + fvm::div(alphaRhoPhi, epsilon_)
316  - fvm::laplacian(alpha*rho*DepsilonEff(), epsilon_)
317  ==
318  Ceps1_*alpha*rho*G*epsilon_/k_
319  - fvm::Sp(Ceps2_*alpha*rho*epsilon_/k_, epsilon_)
320  + fvOptions(alpha, rho, epsilon_)
321  );
322 
323  epsEqn.ref().relax();
324  fvOptions.constrain(epsEqn.ref());
325  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
326  solve(epsEqn);
327  fvOptions.correct(epsilon_);
328  bound(epsilon_, this->epsilonMin_);
329 
330 
331  // Correct the trace of the tensorial production to be consistent
332  // with the near-wall generation from the wall-functions
333  const fvPatchList& patches = this->mesh_.boundary();
334 
335  forAll(patches, patchi)
336  {
337  const fvPatch& curPatch = patches[patchi];
338 
339  if (isA<wallFvPatch>(curPatch))
340  {
341  forAll(curPatch, facei)
342  {
343  label celli = curPatch.faceCells()[facei];
344  P[celli] *= min
345  (
346  G[celli]/(0.5*mag(tr(P[celli])) + SMALL),
347  1.0
348  );
349  }
350  }
351  }
352 
353  volSymmTensorField b(dev(R)/(2*k_));
354  volSymmTensorField S(symm(gradU));
355  volTensorField Omega(skew(gradU));
356 
357  // Reynolds stress equation
358  tmp<fvSymmTensorMatrix> REqn
359  (
360  fvm::ddt(alpha, rho, R)
361  + fvm::div(alphaRhoPhi, R)
362  - fvm::laplacian(alpha*rho*DREff(), R)
363  + fvm::Sp(((C1_/2)*epsilon_ + (C1s_/2)*G)*alpha*rho/k_, R)
364  ==
365  alpha*rho*P
366  - ((1.0/3.0)*I)*(((2.0 - C1_)*epsilon_ - C1s_*G)*alpha*rho)
367  + (C2_*(alpha*rho*epsilon_))*dev(innerSqr(b))
368  + alpha*rho*k_
369  *(
370  (C3_ - C3s_*mag(b))*dev(S)
371  + C4_*devTwoSymm(b&S)
372  + C5_*twoSymm(b&Omega)
373  )
374  + fvOptions(alpha, rho, R)
375  );
376 
377  REqn.ref().relax();
378  fvOptions.constrain(REqn.ref());
379  solve(REqn);
380  fvOptions.correct(R);
381 
382  this->boundNormalStress(R);
383 
384  k_ = 0.5*tr(R);
385  bound(k_, this->kMin_);
386 
387  correctNut();
388 
389  // Correct wall shear-stresses when applying wall-functions
390  this->correctWallShearStress(R);
391 }
392 
393 
394 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
395 
396 } // End namespace RASModels
397 } // End namespace Foam
398 
399 // ************************************************************************* //
void correct(GeometricField< Type, PatchField, GeoMesh > &field)
Apply correction to field.
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
Definition: volFieldsFwd.H:85
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:47
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
tmp< volSymmTensorField > DREff() const
Return the effective diffusivity for R.
Definition: SSG.C:248
dimensionedTensor skew(const dimensionedTensor &dt)
const dimensionedScalar G
Newtonian constant of gravitation.
GeometricField< tensor, fvPatchField, volMesh > volTensorField
Definition: volFieldsFwd.H:86
dimensionedSymmTensor sqr(const dimensionedVector &dv)
tmp< volSymmTensorField > DepsilonEff() const
Return the effective diffusivity for epsilon.
Definition: SSG.C:262
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
Generic dimensioned Type class.
Finite-volume options.
Definition: fvOptions.H:51
Speziale, Sarkar and Gatski Reynolds-stress turbulence model for incompressible and compressible flow...
Definition: SSG.H:94
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
fv::options & fvOptions
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:127
virtual void correct()
Solve the turbulence equations and correct eddy-Viscosity and.
Definition: SSG.C:276
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
Templated abstract base class for RAS turbulence models.
Definition: RASModel.H:77
word timeName
Definition: getTimeIndex.H:3
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
static const Identity< scalar > I
Definition: Identity.H:100
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
A class for handling words, derived from Foam::string.
Definition: word.H:63
virtual const labelUList & faceCells() const
Return faceCells.
Definition: fvPatch.C:107
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:41
volScalarField k_
Definition: SSG.H:134
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
virtual void printCoeffs(const word &type)
Print model coefficients.
Definition: RASModel.C:27
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
SolverPerformance< Type > solve(faMatrix< Type > &, const dictionary &solverControls)
Solve returning the solution statistics given convergence tolerance.
volScalarField epsilon_
Definition: SSG.H:135
dimensionedSymmTensor innerSqr(const dimensionedSymmTensor &dt)
void boundNormalStress(volSymmTensorField &R) const
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:41
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(psi *p+alphal *rhol0+((alphav *psiv+alphal *psil) - psi) *pSat, rhoMin);# 1 "/home/chef2/andy/OpenFOAM/release/v2312/OpenFOAM-v2312/applications/solvers/multiphase/cavitatingFoam/alphavPsi.H" 1{ alphav=clamp((rho - rholSat)/(rhovSat - rholSat), zero_one{});alphal=1.0 - alphav;Info<< "max-min alphav: "<< max(alphav).value()<< " "<< min(alphav).value()<< endl;psiModel-> correct()
Definition: pEqn.H:63
virtual bool read()
Read model coefficients if they have changed.
Definition: SSG.C:222
U
Definition: pEqn.H:72
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:29
#define R(A, B, C, D, E, F, K, M)
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:41
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
Base-class for all transport models used by the incompressible turbulence models. ...
const polyBoundaryMesh & patches
A class for managing temporary objects.
Definition: HashPtrTable.H:50
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present.
Definition: fvOptions.C:96
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
virtual void correctNut()
Update the eddy-viscosity.
Definition: SSG.C:36
volScalarField & nu
dimensionedScalar epsilonMin_
Lower limit of epsilon.
Definition: RASModel.H:114
Namespace for OpenFOAM.
Reynolds-stress turbulence model base class.
SymmTensor< Cmpt > devTwoSymm(const SymmTensor< Cmpt > &st)
Return the deviatoric part of twice the symmetric part of a SymmTensor.
Definition: SymmTensorI.H:491