PDRkEpsilon.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-2024 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 "PDRkEpsilon.H"
27 #include "PDRDragModel.H"
29 
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace compressible
35 {
36 namespace RASModels
37 {
38 
39 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40 
43 
44 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
45 
46 tmp<volScalarField> PDRkEpsilon::boundEpsilon()
47 {
48  tmp<volScalarField> tCmuk2(Cmu_*sqr(k_));
49  epsilon_ = max(epsilon_, tCmuk2()/(nutMaxCoeff_*nu()));
50  return tCmuk2;
51 }
52 
53 
54 void PDRkEpsilon::correctNut()
55 {
56  nut_ = Cmu_*sqr(k_)/epsilon_;
58  fvConstraints::New(mesh_).constrain(nut_);
59 }
60 
61 
62 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
63 
65 (
66  const geometricOneField& alpha,
67  const volScalarField& rho,
68  const volVectorField& U,
69  const surfaceScalarField& alphaRhoPhi,
70  const surfaceScalarField& phi,
71  const viscosity& transport,
72  const word& modelName
73 )
74 :
75  Foam::RASModels::kEpsilon<compressible::momentumTransportModel>
76  (
78  rho,
79  U,
80  phi,
81  phi,
82  transport,
83  modelName
84  ),
85 
86  C4_("C4", coeffDict(), 0.1)
87 {}
88 
89 
90 // * * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * //
91 
93 {}
94 
95 
96 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
97 
99 {
100  if (RASModel::read())
101  {
102  C4_.readIfPresent(coeffDict());
103  return true;
104  }
105  else
106  {
107  return false;
108  }
109 }
110 
111 
113 {
114  if (!turbulence_)
115  {
116  // Re-calculate viscosity
117  nut_ = Cmu_*sqr(k_)/epsilon_;
119 
120  return;
121  }
122 
124 
125  volScalarField divU(fvc::div(phi_/fvc::interpolate(rho_)));
126 
127  if (mesh_.moving())
128  {
129  divU += fvc::div(mesh_.phi());
130  }
131 
132  tmp<volTensorField> tgradU = fvc::grad(U_);
133  volScalarField G(GName(), rho_*nut_*(tgradU() && dev(twoSymm(tgradU()))));
134  tgradU.clear();
135 
136  // Update epsilon and G at the wall
138 
139  // Add the blockage generation term so that it is included consistently
140  // in both the k and epsilon equations
141  const volScalarField& betav =
142  U_.db().lookupObject<volScalarField>("betav");
143 
144  const volScalarField& Lobs =
145  U_.db().lookupObject<volScalarField>("Lobs");
146 
147  const PDRDragModel& drag =
148  U_.db().lookupObject<PDRDragModel>("PDRDragModel");
149 
150  volScalarField GR(drag.Gk());
151 
152  volScalarField LI
153  (C4_*(Lobs + dimensionedScalar(dimLength, rootVSmall)));
154 
155  // Dissipation equation
156  tmp<fvScalarMatrix> epsEqn
157  (
158  betav*fvm::ddt(rho_, epsilon_)
159  + fvm::div(phi_, epsilon_)
161  ==
162  C1_*betav*G*epsilon_/k_
163  + 1.5*pow(Cmu_, 3.0/4.0)*GR*sqrt(k_)/LI
164  - fvm::SuSp(((2.0/3.0)*C1_)*betav*rho_*divU, epsilon_)
165  - fvm::Sp(C2_*betav*rho_*epsilon_/k_, epsilon_)
166  );
167 
168  epsEqn.ref().relax();
169 
170  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
171 
172  solve(epsEqn);
173  boundEpsilon();
174 
175 
176  // Turbulent kinetic energy equation
177 
179  (
180  betav*fvm::ddt(rho_, k_)
181  + fvm::div(phi_, k_)
182  - fvm::laplacian(rho_*DkEff(), k_)
183  ==
184  betav*G + GR
185  - fvm::SuSp((2.0/3.0)*betav*rho_*divU, k_)
186  - fvm::Sp(betav*rho_*epsilon_/k_, k_)
187  );
188 
189  kEqn.ref().relax();
190  solve(kEqn);
191  bound(k_, kMin_);
192 
193  correctNut();
194 }
195 
196 
197 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
198 
199 } // End namespace RASModels
200 } // End namespace compressible
201 } // End namespace Foam
202 
203 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
void updateCoeffs()
Update the boundary condition coefficients.
Generic GeometricField class.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
void correctBoundaryConditions()
Correct boundary field.
const objectRegistry & db() const
Return the local objectRegistry.
Definition: IOobject.C:309
Base-class for sub-grid obstacle drag models. The available drag model is at basic....
Definition: PDRDragModel.H:56
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: RASModel.C:179
virtual bool read()
Read model coefficients if they have changed.
Definition: RASModel.C:161
tmp< volScalarField > DkEff() const
Return the effective diffusivity for k.
Definition: kEpsilon.H:160
tmp< volScalarField > DepsilonEff() const
Return the effective diffusivity for epsilon.
Definition: kEpsilon.H:170
Base class for single-phase compressible turbulence models.
Standard k-epsilon turbulence model with additional source terms corresponding to PDR basic drag mode...
Definition: PDRkEpsilon.H:81
void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: PDRkEpsilon.C:112
PDRkEpsilon(const geometricOneField &alpha, const volScalarField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &transport, const word &modelName=typeName)
Construct from components.
Definition: PDRkEpsilon.C:65
bool read()
Read momentumTransport dictionary.
Definition: PDRkEpsilon.C:98
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
static autoPtr< dictionary > New(Istream &)
Construct top-level dictionary on freestore from Istream.
Definition: dictionaryIO.C:103
bool readIfPresent(const dictionary &, const unitConversion &defaultUnits=NullObjectRef< unitConversion >())
Update the value of dimensioned<Type> if found in the dictionary.
Base class for Lagrangian drag models.
volScalarField nut_
Definition: eddyViscosity.H:60
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type and name.
A class for managing temporary objects.
Definition: tmp.H:55
void clear() const
If object pointer points to valid object:
Definition: tmpI.H:253
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:197
Abstract base class for all fluid physical properties.
Definition: viscosity.H:50
A class for handling words, derived from string.
Definition: word.H:62
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
defineTypeNameAndDebug(PDRkEpsilon, 0)
addToRunTimeSelectionTable(RASModel, PDRkEpsilon, dictionary)
RASModel< momentumTransportModel > RASModel
Typedefs for turbulence, RAS and LES models for compressible flow based on the standard laminar trans...
const dimensionedScalar G
Newtonian constant of gravitation.
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
tmp< VolField< typename outerProduct< vector, Type >::type > > grad(const SurfaceField< Type > &ssf)
Definition: fvcGrad.C:46
tmp< VolField< Type > > div(const SurfaceField< Type > &ssf)
Definition: fvcDiv.C:47
tmp< fvMatrix< Type > > laplacian(const VolField< Type > &vf, const word &name)
Definition: fvmLaplacian.C:47
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const VolField< Type > &vf, const word &name)
Definition: fvmDiv.C:48
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const VolField< Type > &)
tmp< fvMatrix< Type > > SuSp(const volScalarField::Internal &, const VolField< Type > &)
tmp< fvMatrix< Type > > ddt(const VolField< Type > &vf)
Definition: fvmDdt.C:46
Namespace for OpenFOAM.
void dev(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
void twoSymm(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
const dimensionSet dimLength
void pow(LagrangianPatchField< typename powProduct< Type, r >::type > &f, const LagrangianPatchField< Type > &f1)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
bool bound(volScalarField &, const dimensionedScalar &min)
Bound the given scalar field where it is below the specified min value.
Definition: bound.C:31
void sqr(LagrangianPatchField< typename outerProduct< Type, Type >::type > &f, const LagrangianPatchField< Type > &f1)
void sqrt(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
SolverPerformance< Type > solve(fvMatrix< Type > &, const word &)
Solve returning the solution statistics given convergence tolerance.