All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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-2020 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 BasicMomentumTransportModel>
41 {
42  this->nut_ = this->Cmu_*sqr(k_)/epsilon_;
43  this->nut_.correctBoundaryConditions();
44  fv::options::New(this->mesh_).correct(this->nut_);
45 }
46 
47 
48 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
49 
50 template<class BasicMomentumTransportModel>
52 (
53  const alphaField& alpha,
54  const rhoField& rho,
55  const volVectorField& U,
56  const surfaceScalarField& alphaRhoPhi,
57  const surfaceScalarField& phi,
58  const transportModel& transport,
59  const word& type
60 )
61 :
63  (
64  type,
65  alpha,
66  rho,
67  U,
68  alphaRhoPhi,
69  phi,
70  transport
71  ),
72 
73  Cmu_
74  (
76  (
77  "Cmu",
78  this->coeffDict_,
79  0.09
80  )
81  ),
82  C1_
83  (
85  (
86  "C1",
87  this->coeffDict_,
88  3.4
89  )
90  ),
91  C1s_
92  (
94  (
95  "C1s",
96  this->coeffDict_,
97  1.8
98  )
99  ),
100  C2_
101  (
103  (
104  "C2",
105  this->coeffDict_,
106  4.2
107  )
108  ),
109  C3_
110  (
112  (
113  "C3",
114  this->coeffDict_,
115  0.8
116  )
117  ),
118  C3s_
119  (
121  (
122  "C3s",
123  this->coeffDict_,
124  1.3
125  )
126  ),
127  C4_
128  (
130  (
131  "C4",
132  this->coeffDict_,
133  1.25
134  )
135  ),
136  C5_
137  (
139  (
140  "C5",
141  this->coeffDict_,
142  0.4
143  )
144  ),
145 
146  Ceps1_
147  (
149  (
150  "Ceps1",
151  this->coeffDict_,
152  1.44
153  )
154  ),
155  Ceps2_
156  (
158  (
159  "Ceps2",
160  this->coeffDict_,
161  1.92
162  )
163  ),
164  Cs_
165  (
167  (
168  "Cs",
169  this->coeffDict_,
170  0.25
171  )
172  ),
173  Ceps_
174  (
176  (
177  "Ceps",
178  this->coeffDict_,
179  0.15
180  )
181  ),
182 
183  k_
184  (
185  IOobject
186  (
187  "k",
188  this->runTime_.timeName(),
189  this->mesh_,
192  ),
193  0.5*tr(this->R_)
194  ),
195  epsilon_
196  (
197  IOobject
198  (
199  "epsilon",
200  this->runTime_.timeName(),
201  this->mesh_,
204  ),
205  this->mesh_
206  )
207 {
208  if (type == typeName)
209  {
210  this->printCoeffs(type);
211 
212  this->boundNormalStress(this->R_);
213  bound(epsilon_, this->epsilonMin_);
214  k_ = 0.5*tr(this->R_);
215  }
216 }
217 
218 
219 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
220 
221 template<class BasicMomentumTransportModel>
223 {
225  {
226  Cmu_.readIfPresent(this->coeffDict());
227  C1_.readIfPresent(this->coeffDict());
228  C1s_.readIfPresent(this->coeffDict());
229  C2_.readIfPresent(this->coeffDict());
230  C3_.readIfPresent(this->coeffDict());
231  C3s_.readIfPresent(this->coeffDict());
232  C4_.readIfPresent(this->coeffDict());
233  C5_.readIfPresent(this->coeffDict());
234 
235  Ceps1_.readIfPresent(this->coeffDict());
236  Ceps2_.readIfPresent(this->coeffDict());
237  Cs_.readIfPresent(this->coeffDict());
238  Ceps_.readIfPresent(this->coeffDict());
239 
240  return true;
241  }
242  else
243  {
244  return false;
245  }
246 }
247 
248 
249 template<class BasicMomentumTransportModel>
251 {
253  (
254  "DREff",
255  (Cs_*(this->k_/this->epsilon_))*this->R_ + I*this->nu()
256  );
257 }
258 
259 
260 template<class BasicMomentumTransportModel>
262 {
264  (
265  "DepsilonEff",
266  (Ceps_*(this->k_/this->epsilon_))*this->R_ + I*this->nu()
267  );
268 }
269 
270 
271 template<class BasicMomentumTransportModel>
273 {
274  if (!this->turbulence_)
275  {
276  return;
277  }
278 
279  // Local references
280  const alphaField& alpha = this->alpha_;
281  const rhoField& rho = this->rho_;
282  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
283  const volVectorField& U = this->U_;
284  volSymmTensorField& R = this->R_;
285  fv::options& fvOptions(fv::options::New(this->mesh_));
286 
288 
289  tmp<volTensorField> tgradU(fvc::grad(U));
290  const volTensorField& gradU = tgradU();
291 
292  volSymmTensorField P(-twoSymm(R & gradU));
293  volScalarField G(this->GName(), 0.5*mag(tr(P)));
294 
295  // Update epsilon and G at the wall
296  epsilon_.boundaryFieldRef().updateCoeffs();
297 
298  // Dissipation equation
299  tmp<fvScalarMatrix> epsEqn
300  (
301  fvm::ddt(alpha, rho, epsilon_)
302  + fvm::div(alphaRhoPhi, epsilon_)
303  - fvm::laplacian(alpha*rho*DepsilonEff(), epsilon_)
304  ==
305  Ceps1_*alpha*rho*G*epsilon_/k_
306  - fvm::Sp(Ceps2_*alpha*rho*epsilon_/k_, epsilon_)
307  + fvOptions(alpha, rho, epsilon_)
308  );
309 
310  epsEqn.ref().relax();
311  fvOptions.constrain(epsEqn.ref());
312  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
313  solve(epsEqn);
314  fvOptions.correct(epsilon_);
315  bound(epsilon_, this->epsilonMin_);
316 
317 
318  // Correct the trace of the tensorial production to be consistent
319  // with the near-wall generation from the wall-functions
320  const fvPatchList& patches = this->mesh_.boundary();
321 
322  forAll(patches, patchi)
323  {
324  const fvPatch& curPatch = patches[patchi];
325 
326  if (isA<wallFvPatch>(curPatch))
327  {
328  forAll(curPatch, facei)
329  {
330  label celli = curPatch.faceCells()[facei];
331  P[celli] *= min
332  (
333  G[celli]/(0.5*mag(tr(P[celli])) + small),
334  1.0
335  );
336  }
337  }
338  }
339 
340  volSymmTensorField b(dev(R)/(2*k_));
341  volSymmTensorField S(symm(gradU));
342  volTensorField Omega(skew(gradU));
343 
344  // Reynolds stress equation
346  (
347  fvm::ddt(alpha, rho, R)
348  + fvm::div(alphaRhoPhi, R)
349  - fvm::laplacian(alpha*rho*DREff(), R)
350  + fvm::Sp(((C1_/2)*epsilon_ + (C1s_/2)*G)*alpha*rho/k_, R)
351  ==
352  alpha*rho*P
353  - ((1.0/3.0)*I)*(((2.0 - C1_)*epsilon_ - C1s_*G)*alpha*rho)
354  + (C2_*(alpha*rho*epsilon_))*dev(innerSqr(b))
355  + alpha*rho*k_
356  *(
357  (C3_ - C3s_*mag(b))*dev(S)
358  + C4_*dev(twoSymm(b&S))
359  + C5_*twoSymm(b&Omega)
360  )
361  + fvOptions(alpha, rho, R)
362  );
363 
364  REqn.ref().relax();
365  fvOptions.constrain(REqn.ref());
366  solve(REqn);
367  fvOptions.correct(R);
368 
369  this->boundNormalStress(R);
370 
371  k_ = 0.5*tr(R);
372 
373  correctNut();
374 
375  // Correct wall shear-stresses when applying wall-functions
376  this->correctWallShearStress(R);
377 }
378 
379 
380 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
381 
382 } // End namespace RASModels
383 } // End namespace Foam
384 
385 // ************************************************************************* //
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:434
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
fv::options & fvOptions
tmp< volSymmTensorField > DREff() const
Return the effective diffusivity for R.
Definition: SSG.C:250
dimensionedTensor skew(const dimensionedTensor &dt)
void correct(GeometricField< Type, fvPatchField, volMesh > &field)
Apply correction to field.
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:261
static tmp< GeometricField< symmTensor, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< symmTensor >> &)
Return a temporary field constructed from name,.
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
BasicMomentumTransportModel::alphaField alphaField
Definition: RASModel.H:88
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
phi
Definition: pEqn.H:104
virtual void correct()
Solve the turbulence equations and correct eddy-Viscosity and.
Definition: SSG.C:272
SSG(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &type=typeName)
Construct from components.
Definition: SSG.C:52
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)
BasicMomentumTransportModel::transportModel transportModel
Definition: RASModel.H:90
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
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(rho0+psi *p, rhoMin);# 1 "/home/ubuntu/OpenFOAM-8/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:68
virtual bool read()
Read model coefficients if they have changed.
Definition: SSG.C:222
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
BasicMomentumTransportModel::rhoField rhoField
Definition: RASModel.H:89
label patchi
U
Definition: pEqn.H:72
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:70
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.
dimensioned< scalar > mag(const dimensioned< Type > &)
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
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:101
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
const dimensionedScalar & G
Newtonian constant of gravitation.
Namespace for OpenFOAM.
Reynolds-stress turbulence model base class.