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_
87  (
88  dimensioned<scalar>::lookupOrAddToDict
89  (
90  "C4",
91  coeffDict_,
92  0.1
93  )
94  )
95 {}
96 
97 
98 // * * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * //
99 
101 {}
102 
103 
104 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
105 
107 {
108  if (RASModel::read())
109  {
110  C4_.readIfPresent(coeffDict_);
111  return true;
112  }
113  else
114  {
115  return false;
116  }
117 }
118 
119 
121 {
122  if (!turbulence_)
123  {
124  // Re-calculate viscosity
125  nut_ = Cmu_*sqr(k_)/epsilon_;
127 
128  return;
129  }
130 
132 
133  volScalarField divU(fvc::div(phi_/fvc::interpolate(rho_)));
134 
135  if (mesh_.moving())
136  {
137  divU += fvc::div(mesh_.phi());
138  }
139 
140  tmp<volTensorField> tgradU = fvc::grad(U_);
141  volScalarField G(GName(), rho_*nut_*(tgradU() && dev(twoSymm(tgradU()))));
142  tgradU.clear();
143 
144  // Update epsilon and G at the wall
146 
147  // Add the blockage generation term so that it is included consistently
148  // in both the k and epsilon equations
149  const volScalarField& betav =
150  U_.db().lookupObject<volScalarField>("betav");
151 
152  const volScalarField& Lobs =
153  U_.db().lookupObject<volScalarField>("Lobs");
154 
155  const PDRDragModel& drag =
156  U_.db().lookupObject<PDRDragModel>("PDRDragModel");
157 
158  volScalarField GR(drag.Gk());
159 
160  volScalarField LI
161  (C4_*(Lobs + dimensionedScalar(dimLength, rootVSmall)));
162 
163  // Dissipation equation
164  tmp<fvScalarMatrix> epsEqn
165  (
166  betav*fvm::ddt(rho_, epsilon_)
167  + fvm::div(phi_, epsilon_)
169  ==
170  C1_*betav*G*epsilon_/k_
171  + 1.5*pow(Cmu_, 3.0/4.0)*GR*sqrt(k_)/LI
172  - fvm::SuSp(((2.0/3.0)*C1_)*betav*rho_*divU, epsilon_)
173  - fvm::Sp(C2_*betav*rho_*epsilon_/k_, epsilon_)
174  );
175 
176  epsEqn.ref().relax();
177 
178  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
179 
180  solve(epsEqn);
181  boundEpsilon();
182 
183 
184  // Turbulent kinetic energy equation
185 
187  (
188  betav*fvm::ddt(rho_, k_)
189  + fvm::div(phi_, k_)
190  - fvm::laplacian(rho_*DkEff(), k_)
191  ==
192  betav*G + GR
193  - fvm::SuSp((2.0/3.0)*betav*rho_*divU, k_)
194  - fvm::Sp(betav*rho_*epsilon_/k_, k_)
195  );
196 
197  kEqn.ref().relax();
198  solve(kEqn);
199  bound(k_, kMin_);
200 
201  correctNut();
202 }
203 
204 
205 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
206 
207 } // End namespace RASModels
208 } // End namespace compressible
209 } // End namespace Foam
210 
211 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
void updateCoeffs()
Update the boundary condition coefficients.
Generic GeometricField class.
void correctBoundaryConditions()
Correct boundary field.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
const objectRegistry & db() const
Return the local objectRegistry.
Definition: IOobject.C:312
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:195
virtual bool read()
Read model coefficients if they have changed.
Definition: RASModel.C:173
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:120
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:106
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
static autoPtr< dictionary > New(Istream &)
Construct top-level dictionary on freestore from Istream.
Definition: dictionaryIO.C:100
Generic dimensioned Type class.
bool readIfPresent(const dictionary &)
Update the value of dimensioned<Type> if found in the dictionary.
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:237
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
Abstract base class for all fluid physical properties.
Definition: viscosity.H:50
A class for handling words, derived from string.
Definition: word.H:62
Info<< "Reading strained laminar flame speed field Su\n"<< endl;volScalarField Su(IOobject("Su", runTime.name(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Reading field betav\n"<< endl;volScalarField betav(IOobject("betav", mesh.facesInstance(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE), mesh);Info<< "Reading field Lobs\n"<< endl;volScalarField Lobs(IOobject("Lobs", mesh.facesInstance(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE), mesh);Info<< "Reading field CT\n"<< endl;volSymmTensorField CT(IOobject("CT", mesh.facesInstance(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE), mesh);Info<< "Reading field Nv\n"<< endl;volScalarField Nv(IOobject("Nv", mesh.facesInstance(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE), mesh);Info<< "Reading field nsv\n"<< endl;volSymmTensorField nsv(IOobject("nsv", mesh.facesInstance(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE), mesh);IOdictionary PDRProperties(IOobject("PDRProperties", runTime.constant(), mesh, IOobject::MUST_READ_IF_MODIFIED, IOobject::NO_WRITE));autoPtr< PDRDragModel > drag
Definition: createFields.H:192
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.
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
const dimensionSet dimLength
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar sqrt(const dimensionedScalar &ds)
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
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.