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-2023 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 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
45 
47 (
48  const geometricOneField& alpha,
49  const volScalarField& rho,
50  const volVectorField& U,
51  const surfaceScalarField& alphaRhoPhi,
52  const surfaceScalarField& phi,
53  const viscosity& transport,
54  const word& modelName
55 )
56 :
57  Foam::RASModels::kEpsilon<compressible::momentumTransportModel>
58  (
60  rho,
61  U,
62  phi,
63  phi,
64  transport,
65  modelName
66  ),
67 
68  C4_
69  (
70  dimensioned<scalar>::lookupOrAddToDict
71  (
72  "C4",
73  coeffDict_,
74  0.1
75  )
76  )
77 {}
78 
79 
80 // * * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * //
81 
83 {}
84 
85 
86 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
87 
89 {
90  if (RASModel::read())
91  {
92  C4_.readIfPresent(coeffDict_);
93  return true;
94  }
95  else
96  {
97  return false;
98  }
99 }
100 
101 
103 {
104  if (!turbulence_)
105  {
106  // Re-calculate viscosity
107  nut_ = Cmu_*sqr(k_)/epsilon_;
109 
110  // Re-calculate thermal diffusivity
111  //***HGWalphat_ = mut_/Prt_;
112  // alphat_.correctBoundaryConditions();
113 
114  return;
115  }
116 
118 
119  volScalarField divU(fvc::div(phi_/fvc::interpolate(rho_)));
120 
121  if (mesh_.moving())
122  {
123  divU += fvc::div(mesh_.phi());
124  }
125 
126  tmp<volTensorField> tgradU = fvc::grad(U_);
127  volScalarField G(GName(), rho_*nut_*(tgradU() && dev(twoSymm(tgradU()))));
128  tgradU.clear();
129 
130  // Update epsilon and G at the wall
132 
133  // Add the blockage generation term so that it is included consistently
134  // in both the k and epsilon equations
135  const volScalarField& betav =
136  U_.db().lookupObject<volScalarField>("betav");
137 
138  const volScalarField& Lobs =
139  U_.db().lookupObject<volScalarField>("Lobs");
140 
141  const PDRDragModel& drag =
142  U_.db().lookupObject<PDRDragModel>("PDRDragModel");
143 
144  volScalarField GR(drag.Gk());
145 
146  volScalarField LI
147  (C4_*(Lobs + dimensionedScalar(dimLength, rootVSmall)));
148 
149  // Dissipation equation
150  tmp<fvScalarMatrix> epsEqn
151  (
152  betav*fvm::ddt(rho_, epsilon_)
153  + fvm::div(phi_, epsilon_)
155  ==
156  C1_*betav*G*epsilon_/k_
157  + 1.5*pow(Cmu_, 3.0/4.0)*GR*sqrt(k_)/LI
158  - fvm::SuSp(((2.0/3.0)*C1_)*betav*rho_*divU, epsilon_)
159  - fvm::Sp(C2_*betav*rho_*epsilon_/k_, epsilon_)
160  );
161 
162  epsEqn.ref().relax();
163 
164  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
165 
166  solve(epsEqn);
167  bound(epsilon_, epsilonMin_);
168 
169 
170  // Turbulent kinetic energy equation
171 
173  (
174  betav*fvm::ddt(rho_, k_)
175  + fvm::div(phi_, k_)
176  - fvm::laplacian(rho_*DkEff(), k_)
177  ==
178  betav*G + GR
179  - fvm::SuSp((2.0/3.0)*betav*rho_*divU, k_)
180  - fvm::Sp(betav*rho_*epsilon_/k_, k_)
181  );
182 
183  kEqn.ref().relax();
184  solve(kEqn);
185  bound(k_, kMin_);
186 
187  // Re-calculate viscosity
188  nut_ = Cmu_*sqr(k_)/epsilon_;
190 
191  // Re-calculate thermal diffusivity
192  //***HGWalphat_ = mut_/Prt_;
193  // alphat_.correctBoundaryConditions();
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.
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:207
virtual bool read()
Read model coefficients if they have changed.
Definition: RASModel.C:184
tmp< volScalarField > DkEff() const
Return the effective diffusivity for k.
Definition: kEpsilon.H:152
tmp< volScalarField > DepsilonEff() const
Return the effective diffusivity for epsilon.
Definition: kEpsilon.H:162
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:102
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:47
bool read()
Read momentumTransport dictionary.
Definition: PDRkEpsilon.C:88
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
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:194
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 > > Sp(const volScalarField::Internal &, const VolField< Type > &)
tmp< fvMatrix< Type > > SuSp(const volScalarField::Internal &, const VolField< Type > &)
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const VolField< Type > &vf, const word &name)
Definition: fvmDiv.C:48
tmp< fvMatrix< Type > > ddt(const VolField< Type > &vf)
Definition: fvmDdt.C:46
Namespace for OpenFOAM.
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
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)
bool bound(volScalarField &, const dimensionedScalar &min)
Bound the given scalar field where it is below the specified min value.
Definition: bound.C:31
SolverPerformance< Type > solve(fvMatrix< Type > &, const word &)
Solve returning the solution statistics given convergence tolerance.