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-2016 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 "fvOptions.H"
28 #include "wallFvPatch.H"
29 
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace RASModels
35 {
36 
37 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
38 
39 template<class BasicTurbulenceModel>
41 {
42  this->nut_ = this->Cmu_*sqr(k_)/epsilon_;
43  this->nut_.correctBoundaryConditions();
44  fv::options::New(this->mesh_).correct(this->nut_);
45 
46  BasicTurbulenceModel::correctNut();
47 }
48 
49 
50 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
51 
52 template<class BasicTurbulenceModel>
54 (
55  const alphaField& alpha,
56  const rhoField& rho,
57  const volVectorField& U,
58  const surfaceScalarField& alphaRhoPhi,
59  const surfaceScalarField& phi,
60  const transportModel& transport,
61  const word& propertiesName,
62  const word& type
63 )
64 :
66  (
67  type,
68  alpha,
69  rho,
70  U,
71  alphaRhoPhi,
72  phi,
73  transport,
74  propertiesName
75  ),
76 
77  Cmu_
78  (
80  (
81  "Cmu",
82  this->coeffDict_,
83  0.09
84  )
85  ),
86  C1_
87  (
89  (
90  "C1",
91  this->coeffDict_,
92  3.4
93  )
94  ),
95  C1s_
96  (
98  (
99  "C1s",
100  this->coeffDict_,
101  1.8
102  )
103  ),
104  C2_
105  (
107  (
108  "C2",
109  this->coeffDict_,
110  4.2
111  )
112  ),
113  C3_
114  (
116  (
117  "C3",
118  this->coeffDict_,
119  0.8
120  )
121  ),
122  C3s_
123  (
125  (
126  "C3s",
127  this->coeffDict_,
128  1.3
129  )
130  ),
131  C4_
132  (
134  (
135  "C4",
136  this->coeffDict_,
137  1.25
138  )
139  ),
140  C5_
141  (
143  (
144  "C5",
145  this->coeffDict_,
146  0.4
147  )
148  ),
149 
150  Ceps1_
151  (
153  (
154  "Ceps1",
155  this->coeffDict_,
156  1.44
157  )
158  ),
159  Ceps2_
160  (
162  (
163  "Ceps2",
164  this->coeffDict_,
165  1.92
166  )
167  ),
168  Cs_
169  (
171  (
172  "Cs",
173  this->coeffDict_,
174  0.25
175  )
176  ),
177  Ceps_
178  (
180  (
181  "Ceps",
182  this->coeffDict_,
183  0.15
184  )
185  ),
186 
187  k_
188  (
189  IOobject
190  (
191  "k",
192  this->runTime_.timeName(),
193  this->mesh_,
196  ),
197  0.5*tr(this->R_)
198  ),
199  epsilon_
200  (
201  IOobject
202  (
203  "epsilon",
204  this->runTime_.timeName(),
205  this->mesh_,
208  ),
209  this->mesh_
210  )
211 {
212  if (type == typeName)
213  {
214  this->printCoeffs(type);
215 
216  this->boundNormalStress(this->R_);
217  bound(epsilon_, this->epsilonMin_);
218  k_ = 0.5*tr(this->R_);
219  }
220 }
221 
222 
223 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
224 
225 template<class BasicTurbulenceModel>
227 {
229  {
230  Cmu_.readIfPresent(this->coeffDict());
231  C1_.readIfPresent(this->coeffDict());
232  C1s_.readIfPresent(this->coeffDict());
233  C2_.readIfPresent(this->coeffDict());
234  C3_.readIfPresent(this->coeffDict());
235  C3s_.readIfPresent(this->coeffDict());
236  C4_.readIfPresent(this->coeffDict());
237  C5_.readIfPresent(this->coeffDict());
238 
239  Ceps1_.readIfPresent(this->coeffDict());
240  Ceps2_.readIfPresent(this->coeffDict());
241  Cs_.readIfPresent(this->coeffDict());
242  Ceps_.readIfPresent(this->coeffDict());
243 
244  return true;
245  }
246  else
247  {
248  return false;
249  }
250 }
251 
252 
253 template<class BasicTurbulenceModel>
255 {
257  (
259  (
260  "DREff",
261  (Cs_*(this->k_/this->epsilon_))*this->R_ + I*this->nu()
262  )
263  );
264 }
265 
266 
267 template<class BasicTurbulenceModel>
269 {
271  (
273  (
274  "DepsilonEff",
275  (Ceps_*(this->k_/this->epsilon_))*this->R_ + I*this->nu()
276  )
277  );
278 }
279 
280 
281 template<class BasicTurbulenceModel>
283 {
284  if (!this->turbulence_)
285  {
286  return;
287  }
288 
289  // Local references
290  const alphaField& alpha = this->alpha_;
291  const rhoField& rho = this->rho_;
292  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
293  const volVectorField& U = this->U_;
294  volSymmTensorField& R = this->R_;
295  fv::options& fvOptions(fv::options::New(this->mesh_));
296 
298 
299  tmp<volTensorField> tgradU(fvc::grad(U));
300  const volTensorField& gradU = tgradU();
301 
302  volSymmTensorField P(-twoSymm(R & gradU));
303  volScalarField G(this->GName(), 0.5*mag(tr(P)));
304 
305  // Update epsilon and G at the wall
306  epsilon_.boundaryFieldRef().updateCoeffs();
307 
308  // Dissipation equation
309  tmp<fvScalarMatrix> epsEqn
310  (
311  fvm::ddt(alpha, rho, epsilon_)
312  + fvm::div(alphaRhoPhi, epsilon_)
313  - fvm::laplacian(alpha*rho*DepsilonEff(), epsilon_)
314  ==
315  Ceps1_*alpha*rho*G*epsilon_/k_
316  - fvm::Sp(Ceps2_*alpha*rho*epsilon_/k_, epsilon_)
317  + fvOptions(alpha, rho, epsilon_)
318  );
319 
320  epsEqn.ref().relax();
321  fvOptions.constrain(epsEqn.ref());
322  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
323  solve(epsEqn);
324  fvOptions.correct(epsilon_);
325  bound(epsilon_, this->epsilonMin_);
326 
327 
328  // Correct the trace of the tensorial production to be consistent
329  // with the near-wall generation from the wall-functions
330  const fvPatchList& patches = this->mesh_.boundary();
331 
332  forAll(patches, patchi)
333  {
334  const fvPatch& curPatch = patches[patchi];
335 
336  if (isA<wallFvPatch>(curPatch))
337  {
338  forAll(curPatch, facei)
339  {
340  label celli = curPatch.faceCells()[facei];
341  P[celli] *= min
342  (
343  G[celli]/(0.5*mag(tr(P[celli])) + SMALL),
344  1.0
345  );
346  }
347  }
348  }
349 
350  volSymmTensorField b(dev(R)/(2*k_));
351  volSymmTensorField S(symm(gradU));
352  volTensorField Omega(skew(gradU));
353 
354  // Reynolds stress equation
356  (
357  fvm::ddt(alpha, rho, R)
358  + fvm::div(alphaRhoPhi, R)
359  - fvm::laplacian(alpha*rho*DREff(), R)
360  + fvm::Sp(((C1_/2)*epsilon_ + (C1s_/2)*G)*alpha*rho/k_, R)
361  ==
362  alpha*rho*P
363  - ((1.0/3.0)*I)*(((2.0 - C1_)*epsilon_ - C1s_*G)*alpha*rho)
364  + (C2_*(alpha*rho*epsilon_))*dev(innerSqr(b))
365  + alpha*rho*k_
366  *(
367  (C3_ - C3s_*mag(b))*dev(S)
368  + C4_*dev(twoSymm(b&S))
369  + C5_*twoSymm(b&Omega)
370  )
371  + fvOptions(alpha, rho, R)
372  );
373 
374  REqn.ref().relax();
375  fvOptions.constrain(REqn.ref());
376  solve(REqn);
377  fvOptions.correct(R);
378 
379  this->boundNormalStress(R);
380 
381  k_ = 0.5*tr(R);
382 
383  correctNut();
384 
385  // Correct wall shear-stresses when applying wall-functions
386  this->correctWallShearStress(R);
387 }
388 
389 
390 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
391 
392 } // End namespace RASModels
393 } // End namespace Foam
394 
395 // ************************************************************************* //
surfaceScalarField & phi
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
Definition: volFieldsFwd.H:58
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
U
Definition: pEqn.H:83
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
tmp< volSymmTensorField > DREff() const
Return the effective diffusivity for R.
Definition: SSG.C:254
dimensionedTensor skew(const dimensionedTensor &dt)
void correct(GeometricField< Type, fvPatchField, volMesh > &field)
Apply correction to field.
const dimensionedScalar G
Newtonian constant of gravitation.
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
dimensionedSymmTensor sqr(const dimensionedVector &dv)
tmp< volSymmTensorField > DepsilonEff() const
Return the effective diffusivity for epsilon.
Definition: SSG.C:268
BasicTurbulenceModel::rhoField rhoField
Definition: RASModel.H:100
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
Generic dimensioned Type class.
patches[0]
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
Finite-volume options.
Definition: fvOptions.H:52
Speziale, Sarkar and Gatski Reynolds-stress turbulence model for incompressible and compressible flow...
Definition: SSG.H:95
BasicTurbulenceModel::transportModel transportModel
Definition: RASModel.H:101
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
fv::options & fvOptions
virtual void correct()
Solve the turbulence equations and correct eddy-Viscosity and.
Definition: SSG.C:282
Templated abstract base class for RAS turbulence models.
Definition: RASModel.H:49
bool read(const char *, int32_t &)
Definition: int32IO.C:85
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
static const Identity< scalar > I
Definition: Identity.H:93
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
A class for handling words, derived from string.
Definition: word.H:59
SolverPerformance< Type > solve(fvMatrix< Type > &, const dictionary &)
Solve returning the solution statistics given convergence tolerance.
virtual const labelUList & faceCells() const
Return faceCells.
Definition: fvPatch.C:93
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
void constrain(fvMatrix< Type > &eqn)
Apply constraints to equation.
dimensionedSymmTensor innerSqr(const dimensionedSymmTensor &dt)
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:46
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/ubuntu/OpenFOAM-5.0/applications/solvers/multiphase/cavitatingFoam/alphavPsi.H" 1{ alphav=max(min((rho - rholSat)/(rhovSat - rholSat), scalar(1)), scalar(0));alphal=1.0 - alphav;Info<< "max-min alphav: "<< max(alphav).value()<< " "<< min(alphav).value()<< endl;psiModel-> correct()
Definition: pEqn.H:72
virtual bool read()
Read model coefficients if they have changed.
Definition: SSG.C:226
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
fileName::Type type(const fileName &, const bool followLink=true)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:485
label patchi
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:33
#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:46
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:63
dimensioned< scalar > mag(const dimensioned< Type > &)
A class for managing temporary objects.
Definition: PtrList.H:53
static options & New(const fvMesh &mesh)
Construct fvOptions and register to datbase if not present.
Definition: fvOptions.C:103
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
virtual void correctNut()
Update the eddy-viscosity.
Definition: SSG.C:40
volScalarField & nu
Namespace for OpenFOAM.
BasicTurbulenceModel::alphaField alphaField
Definition: RASModel.H:99
Reynolds-stress turbulence model base class.