kEqn.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-2021 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 "kEqn.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(k_)*this->delta();
43  this->nut_.correctBoundaryConditions();
44  fvConstraints::New(this->mesh_).constrain(this->nut_);
45 }
46 
47 
48 template<class BasicMomentumTransportModel>
50 {
51  return tmp<fvScalarMatrix>
52  (
53  new fvScalarMatrix
54  (
55  k_,
56  dimVolume*this->rho_.dimensions()*k_.dimensions()
57  /dimTime
58  )
59  );
60 }
61 
62 
63 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
64 
65 template<class BasicMomentumTransportModel>
67 (
68  const alphaField& alpha,
69  const rhoField& rho,
70  const volVectorField& U,
71  const surfaceScalarField& alphaRhoPhi,
72  const surfaceScalarField& phi,
73  const viscosity& viscosity,
74  const word& type
75 )
76 :
78  (
79  type,
80  alpha,
81  rho,
82  U,
83  alphaRhoPhi,
84  phi,
85  viscosity
86  ),
87 
88  k_
89  (
90  IOobject
91  (
92  IOobject::groupName("k", this->alphaRhoPhi_.group()),
93  this->runTime_.timeName(),
94  this->mesh_,
97  ),
98  this->mesh_
99  ),
100 
101  Ck_
102  (
104  (
105  "Ck",
106  this->coeffDict_,
107  0.094
108  )
109  )
110 {
111  bound(k_, this->kMin_);
112 
113  if (type == typeName)
114  {
115  this->printCoeffs(type);
116  }
117 }
118 
119 
120 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
121 
122 template<class BasicMomentumTransportModel>
124 {
126  {
127  Ck_.readIfPresent(this->coeffDict());
128 
129  return true;
130  }
131  else
132  {
133  return false;
134  }
135 }
136 
137 
138 template<class BasicMomentumTransportModel>
140 {
141  return volScalarField::New
142  (
143  IOobject::groupName("epsilon", this->alphaRhoPhi_.group()),
144  this->Ce_*k()*sqrt(k())/this->delta()
145  );
146 }
147 
148 
149 template<class BasicMomentumTransportModel>
151 {
152  if (!this->turbulence_)
153  {
154  return;
155  }
156 
157  // Local references
158  const alphaField& alpha = this->alpha_;
159  const rhoField& rho = this->rho_;
160  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
161  const volVectorField& U = this->U_;
162  volScalarField& nut = this->nut_;
163  const Foam::fvModels& fvModels(Foam::fvModels::New(this->mesh_));
165  (
166  Foam::fvConstraints::New(this->mesh_)
167  );
168 
170 
172 
173  tmp<volTensorField> tgradU(fvc::grad(U));
174  volScalarField G(this->GName(), nut*(tgradU() && dev(twoSymm(tgradU()))));
175  tgradU.clear();
176 
178  (
179  fvm::ddt(alpha, rho, k_)
180  + fvm::div(alphaRhoPhi, k_)
181  - fvm::laplacian(alpha*rho*DkEff(), k_)
182  ==
183  alpha*rho*G
184  - fvm::SuSp((2.0/3.0)*alpha*rho*divU, k_)
185  - fvm::Sp(this->Ce_*alpha*rho*sqrt(k_)/this->delta(), k_)
186  + kSource()
187  + fvModels.source(alpha, rho, k_)
188  );
189 
190  kEqn.ref().relax();
191  fvConstraints.constrain(kEqn.ref());
192  solve(kEqn);
193  fvConstraints.constrain(k_);
194  bound(k_, this->kMin_);
195 
196  correctNut();
197 }
198 
199 
200 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
201 
202 } // End namespace LESModels
203 } // End namespace Foam
204 
205 // ************************************************************************* //
Eddy viscosity LES SGS model base class.
virtual tmp< volScalarField > epsilon() const
Return sub-grid disipation rate.
Definition: kEqn.C:139
scalar delta
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
tmp< fvMatrix< Type > > SuSp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
void clear() const
If object pointer points to valid object:
Definition: tmpI.H:237
U
Definition: pEqn.H:72
kEqn(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.
Definition: kEqn.C:67
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
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.
const dimensionedScalar G
Newtonian constant of gravitation.
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
One equation eddy-viscosity model.
Definition: kEqn.H:71
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
const dimensionSet dimTime
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
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
virtual void correct()
Correct eddy-Viscosity and related properties.
Definition: kEqn.C:150
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
virtual tmp< fvScalarMatrix > kSource() const
Definition: kEqn.C:49
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
virtual bool read()
Read model coefficients if they have changed.
Definition: kEqn.C:123
BasicMomentumTransportModel::rhoField rhoField
Definition: kEqn.H:98
tmp< volScalarField > divU
Definition: alphaSuSp.H:9
bool constrain(fvMatrix< Type > &eqn) const
Apply constraints to an equation.
Foam::fvModels & fvModels
virtual void correctNut()
Definition: kEqn.C:40
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:202
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:33
Finite volume constraints.
Definition: fvConstraints.H:57
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:46
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
const dimensionSet dimVolume
SolverPerformance< Type > solve(fvMatrix< Type > &, const word &)
Solve returning the solution statistics given convergence tolerance.
A class for managing temporary objects.
Definition: PtrList.H:53
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:98
BasicMomentumTransportModel::alphaField alphaField
Definition: kEqn.H:97
scalar nut
Namespace for OpenFOAM.