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-2022 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 "fvModels.H"
28 #include "fvConstraints.H"
29 
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace LESModels
35 {
36 
37 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
38 
39 template<class BasicMomentumTransportModel>
41 {
42  this->nut_ = Ck_*sqrt(this->k())*this->delta();
43  this->nut_.correctBoundaryConditions();
44  fvConstraints::New(this->mesh_).constrain(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 viscosity& viscosity,
59  const word& type
60 )
61 :
63  (
64  type,
65  alpha,
66  rho,
67  U,
68  alphaRhoPhi,
69  phi,
70  viscosity
71  ),
72 
73  Ck_
74  (
76  (
77  "Ck",
78  this->coeffDict_,
79  0.094
80  )
81  ),
82  Cm_
83  (
85  (
86  "Cm",
87  this->coeffDict_,
88  4.13
89  )
90  ),
91  Ce_
92  (
94  (
95  "Ce",
96  this->coeffDict_,
97  1.05
98  )
99  ),
100  Cs_
101  (
103  (
104  "Cs",
105  this->coeffDict_,
106  0.25
107  )
108  )
109 {
110  if (type == typeName)
111  {
112  this->printCoeffs(type);
113  this->boundNormalStress(this->R_);
114  }
115 }
116 
117 
118 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
119 
120 template<class BasicMomentumTransportModel>
122 {
124  {
125  Ck_.readIfPresent(this->coeffDict());
126  Cm_.readIfPresent(this->coeffDict());
127  Ce_.readIfPresent(this->coeffDict());
128  Cs_.readIfPresent(this->coeffDict());
129 
130  return true;
131  }
132  else
133  {
134  return false;
135  }
136 }
137 
138 
139 template<class BasicMomentumTransportModel>
142 {
143  const volScalarField k(this->k());
144 
145  return volScalarField::New
146  (
147  IOobject::groupName("epsilon", this->alphaRhoPhi_.group()),
148  this->Ce_*k*sqrt(k)/this->delta()
149  );
150 }
151 
152 
153 template<class BasicMomentumTransportModel>
156 {
157  const volScalarField k(this->k());
158 
159  return volScalarField::New
160  (
161  IOobject::groupName("omega", this->alphaRhoPhi_.group()),
162  (this->Ce_/this->Ck_)*sqrt(k)/this->delta()
163  );
164 }
165 
166 
167 template<class BasicMomentumTransportModel>
169 {
170  if (!this->turbulence_)
171  {
172  return;
173  }
174 
175  // Local references
176  const alphaField& alpha = this->alpha_;
177  const rhoField& rho = this->rho_;
178  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
179  const volVectorField& U = this->U_;
180  volSymmTensorField& R = this->R_;
181  const Foam::fvModels& fvModels(Foam::fvModels::New(this->mesh_));
183  (
184  Foam::fvConstraints::New(this->mesh_)
185  );
186 
188 
189  tmp<volTensorField> tgradU(fvc::grad(U));
190  const volTensorField& gradU = tgradU();
191 
192  volSymmTensorField D(symm(gradU));
193 
194  volSymmTensorField P(-twoSymm(R & gradU));
195 
196  volScalarField k(this->k());
197 
199  (
200  fvm::ddt(alpha, rho, R)
201  + fvm::div(alphaRhoPhi, R)
202  - fvm::laplacian(I*this->nu() + Cs_*(k/this->epsilon())*R, R)
203  + fvm::Sp(Cm_*alpha*rho*sqrt(k)/this->delta(), R)
204  ==
205  alpha*rho*P
206  + (4.0/5.0)*alpha*rho*k*D
207  - ((2.0/3.0)*(1.0 - Cm_/this->Ce_)*I)*(alpha*rho*this->epsilon())
208  + this->RSource()
209  + fvModels.source(alpha, rho, R)
210  );
211 
212  REqn.ref().relax();
213  fvConstraints.constrain(REqn.ref());
214  REqn.ref().solve();
215  fvConstraints.constrain(R);
216  this->boundNormalStress(R);
217 
218  correctNut();
219 }
220 
221 
222 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
223 
224 } // End namespace LESModels
225 } // End namespace Foam
226 
227 // ************************************************************************* //
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:60
U
Definition: pEqn.H:72
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(rho0+psi *p, rhoMin);# 1 "/home/ubuntu/OpenFOAM-10/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
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
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 > &)
BasicMomentumTransportModel::rhoField rhoField
BasicMomentumTransportModel::alphaField alphaField
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
virtual void correct()
Correct sub-grid stress, eddy-Viscosity and related properties.
DeardorffDiffStress(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity, const word &type=typeName)
Constructor from components.
bool read(const char *, int32_t &)
Definition: int32IO.C:85
virtual tmp< volScalarField > omega() const
Return the turbulence specific dissipation rate.
static const Identity< scalar > I
Definition: Identity.H:93
Foam::fvConstraints & fvConstraints
A class for handling words, derived from string.
Definition: word.H:59
static word groupName(Name name, const word &group)
Finite volume models.
Definition: fvModels.H:60
Info<< "Reading field p\"<< endl;volScalarField p(IOobject("p", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Reading field U\"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);pressureReference pressureReference(p, simple.dict());mesh.schemes().setFluxRequired(p.name());Info<< "Reading field pa\"<< endl;volScalarField pa(IOobject("pa", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Reading field Ua\"<< endl;volVectorField Ua(IOobject("Ua", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);# 65 "/home/ubuntu/OpenFOAM-10/applications/solvers/incompressible/adjointShapeOptimisationFoam/createFields.H" 2label paRefCell=0;scalar paRefValue=0.0;setRefCell(pa, simple.dict(), paRefCell, paRefValue);mesh.schemes().setFluxRequired(pa.name());autoPtr< viscosityModel > viscosity(viscosityModel::New(mesh))
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
phi
Definition: correctPhi.H:3
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:46
static autoPtr< dictionary > New(Istream &)
Construct top-level dictionary on freestore from Istream.
Definition: dictionaryIO.C:96
bool constrain(fvMatrix< Type > &eqn) const
Apply constraints to an equation.
Foam::fvModels & fvModels
#define R(A, B, C, D, E, F, K, M)
Finite volume constraints.
Definition: fvConstraints.H:57
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
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.
A class for managing temporary objects.
Definition: PtrList.H:53
Namespace for OpenFOAM.
Reynolds-stress turbulence model base class.