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