DeardorffDiffStress.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) 2011-2018 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 "DeardorffDiffStress.H"
27 #include "fvOptions.H"
28 
29 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace LESModels
34 {
35 
36 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
37 
38 template<class BasicTurbulenceModel>
40 {
41  this->nut_ = Ck_*sqrt(this->k())*this->delta();
42  this->nut_.correctBoundaryConditions();
43  fv::options::New(this->mesh_).correct(this->nut_);
44 
45  BasicTurbulenceModel::correctNut();
46 }
47 
48 
49 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
50 
51 template<class BasicTurbulenceModel>
53 (
54  const alphaField& alpha,
55  const rhoField& rho,
56  const volVectorField& U,
57  const surfaceScalarField& alphaRhoPhi,
58  const surfaceScalarField& phi,
59  const transportModel& transport,
60  const word& propertiesName,
61  const word& type
62 )
63 :
65  (
66  type,
67  alpha,
68  rho,
69  U,
70  alphaRhoPhi,
71  phi,
72  transport,
73  propertiesName
74  ),
75 
76  Ck_
77  (
79  (
80  "Ck",
81  this->coeffDict_,
82  0.094
83  )
84  ),
85  Cm_
86  (
88  (
89  "Cm",
90  this->coeffDict_,
91  4.13
92  )
93  ),
94  Ce_
95  (
97  (
98  "Ce",
99  this->coeffDict_,
100  1.05
101  )
102  ),
103  Cs_
104  (
106  (
107  "Cs",
108  this->coeffDict_,
109  0.25
110  )
111  )
112 {
113  if (type == typeName)
114  {
115  this->printCoeffs(type);
116  this->boundNormalStress(this->R_);
117  }
118 }
119 
120 
121 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
122 
123 template<class BasicTurbulenceModel>
125 {
127  {
128  Ck_.readIfPresent(this->coeffDict());
129  Cm_.readIfPresent(this->coeffDict());
130  Ce_.readIfPresent(this->coeffDict());
131  Cs_.readIfPresent(this->coeffDict());
132 
133  return true;
134  }
135  else
136  {
137  return false;
138  }
139 }
140 
141 
142 template<class BasicTurbulenceModel>
144 {
145  volScalarField k(this->k());
146 
147  return tmp<volScalarField>
148  (
149  new volScalarField
150  (
151  IOobject
152  (
153  IOobject::groupName("epsilon", this->alphaRhoPhi_.group()),
154  this->runTime_.timeName(),
155  this->mesh_,
158  ),
159  this->Ce_*k*sqrt(k)/this->delta()
160  )
161  );
162 }
163 
164 
165 template<class BasicTurbulenceModel>
167 {
168  if (!this->turbulence_)
169  {
170  return;
171  }
172 
173  // Local references
174  const alphaField& alpha = this->alpha_;
175  const rhoField& rho = this->rho_;
176  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
177  const volVectorField& U = this->U_;
178  volSymmTensorField& R = this->R_;
179  fv::options& fvOptions(fv::options::New(this->mesh_));
180 
182 
183  tmp<volTensorField> tgradU(fvc::grad(U));
184  const volTensorField& gradU = tgradU();
185 
186  volSymmTensorField D(symm(gradU));
187 
188  volSymmTensorField P(-twoSymm(R & gradU));
189 
190  volScalarField k(this->k());
191 
193  (
194  fvm::ddt(alpha, rho, R)
195  + fvm::div(alphaRhoPhi, R)
196  - fvm::laplacian(I*this->nu() + Cs_*(k/this->epsilon())*R, R)
197  + fvm::Sp(Cm_*alpha*rho*sqrt(k)/this->delta(), R)
198  ==
199  alpha*rho*P
200  + (4.0/5.0)*alpha*rho*k*D
201  - ((2.0/3.0)*(1.0 - Cm_/this->Ce_)*I)*(alpha*rho*this->epsilon())
202  + fvOptions(alpha, rho, R)
203  );
204 
205  REqn.ref().relax();
206  fvOptions.constrain(REqn.ref());
207  REqn.ref().solve();
208  fvOptions.correct(R);
209  this->boundNormalStress(R);
210 
211  correctNut();
212 }
213 
214 
215 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
216 
217 } // End namespace LESModels
218 } // End namespace Foam
219 
220 // ************************************************************************* //
scalar delta
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
Templated abstract base class for LES SGS models.
Definition: LESModel.H:56
fv::options & fvOptions
surfaceScalarField & phi
Differential SGS Stress Equation Model for incompressible and compressible flows. ...
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
BasicTurbulenceModel::rhoField rhoField
dimensionedScalar sqrt(const dimensionedScalar &ds)
BasicTurbulenceModel::alphaField alphaField
label k
Boltzmann constant.
Generic dimensioned Type class.
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
Finite-volume options.
Definition: fvOptions.H:52
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
BasicTurbulenceModel::transportModel transportModel
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
virtual void correct()
Correct sub-grid stress, eddy-Viscosity and related properties.
bool read(const char *, int32_t &)
Definition: int32IO.C:85
static const Identity< scalar > I
Definition: Identity.H:93
A class for handling words, derived from string.
Definition: word.H:59
static word groupName(Name name, const word &group)
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
void constrain(fvMatrix< Type > &eqn)
Apply constraints to equation.
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-6/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
fileName::Type type(const fileName &, const bool followLink=true)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:481
U
Definition: pEqn.H:72
#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
virtual void correctNut()
Update the eddy-viscosity.
scalar epsilon
virtual bool read()
Read model coefficients if they have changed.
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
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
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
volScalarField & nu
Namespace for OpenFOAM.
Reynolds-stress turbulence model base class.