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-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 "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 BasicMomentumTransportModel>
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 
46 
47 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48 
49 template<class BasicMomentumTransportModel>
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& type
59 )
60 :
62  (
63  type,
64  alpha,
65  rho,
66  U,
67  alphaRhoPhi,
68  phi,
69  transport
70  ),
71 
72  Ck_
73  (
75  (
76  "Ck",
77  this->coeffDict_,
78  0.094
79  )
80  ),
81  Cm_
82  (
84  (
85  "Cm",
86  this->coeffDict_,
87  4.13
88  )
89  ),
90  Ce_
91  (
93  (
94  "Ce",
95  this->coeffDict_,
96  1.05
97  )
98  ),
99  Cs_
100  (
102  (
103  "Cs",
104  this->coeffDict_,
105  0.25
106  )
107  )
108 {
109  if (type == typeName)
110  {
111  this->printCoeffs(type);
112  this->boundNormalStress(this->R_);
113  }
114 }
115 
116 
117 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
118 
119 template<class BasicMomentumTransportModel>
121 {
123  {
124  Ck_.readIfPresent(this->coeffDict());
125  Cm_.readIfPresent(this->coeffDict());
126  Ce_.readIfPresent(this->coeffDict());
127  Cs_.readIfPresent(this->coeffDict());
128 
129  return true;
130  }
131  else
132  {
133  return false;
134  }
135 }
136 
137 
138 template<class BasicMomentumTransportModel>
141 {
142  volScalarField k(this->k());
143 
144  return volScalarField::New
145  (
146  IOobject::groupName("epsilon", this->alphaRhoPhi_.group()),
147  this->Ce_*k*sqrt(k)/this->delta()
148  );
149 }
150 
151 
152 template<class BasicMomentumTransportModel>
154 {
155  if (!this->turbulence_)
156  {
157  return;
158  }
159 
160  // Local references
161  const alphaField& alpha = this->alpha_;
162  const rhoField& rho = this->rho_;
163  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
164  const volVectorField& U = this->U_;
165  volSymmTensorField& R = this->R_;
166  fv::options& fvOptions(fv::options::New(this->mesh_));
167 
169 
170  tmp<volTensorField> tgradU(fvc::grad(U));
171  const volTensorField& gradU = tgradU();
172 
173  volSymmTensorField D(symm(gradU));
174 
175  volSymmTensorField P(-twoSymm(R & gradU));
176 
177  volScalarField k(this->k());
178 
180  (
181  fvm::ddt(alpha, rho, R)
182  + fvm::div(alphaRhoPhi, R)
183  - fvm::laplacian(I*this->nu() + Cs_*(k/this->epsilon())*R, R)
184  + fvm::Sp(Cm_*alpha*rho*sqrt(k)/this->delta(), R)
185  ==
186  alpha*rho*P
187  + (4.0/5.0)*alpha*rho*k*D
188  - ((2.0/3.0)*(1.0 - Cm_/this->Ce_)*I)*(alpha*rho*this->epsilon())
189  + fvOptions(alpha, rho, R)
190  );
191 
192  REqn.ref().relax();
193  fvOptions.constrain(REqn.ref());
194  REqn.ref().solve();
195  fvOptions.correct(R);
196  this->boundNormalStress(R);
197 
198  correctNut();
199 }
200 
201 
202 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
203 
204 } // End namespace LESModels
205 } // End namespace Foam
206 
207 // ************************************************************************* //
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
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
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
dimensionedScalar sqrt(const dimensionedScalar &ds)
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
BasicMomentumTransportModel::rhoField rhoField
BasicMomentumTransportModel::alphaField alphaField
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
phi
Definition: pEqn.H:104
DeardorffDiffStress(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &type=typeName)
Constructor from components.
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(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
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
BasicMomentumTransportModel::transportModel transportModel
virtual void correctNut()
Update the eddy-viscosity.
scalar epsilon
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
virtual bool read()
Read model coefficients if they have changed.
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
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
Namespace for OpenFOAM.
Reynolds-stress turbulence model base class.