PDRkEpsilon.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2019-2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "PDRkEpsilon.H"
30 #include "PDRDragModel.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace compressible
38 {
39 namespace RASModels
40 {
41 
42 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
43 
44 defineTypeNameAndDebug(PDRkEpsilon, 0);
45 addToRunTimeSelectionTable(RASModel, PDRkEpsilon, dictionary);
46 
47 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48 
50 (
51  const geometricOneField& alpha,
52  const volScalarField& rho,
53  const volVectorField& U,
54  const surfaceScalarField& alphaRhoPhi,
55  const surfaceScalarField& phi,
56  const fluidThermo& thermophysicalModel,
57  const word& turbulenceModelName,
58  const word& modelName
59 )
60 :
61  Foam::RASModels::kEpsilon<EddyDiffusivity<compressible::turbulenceModel>>
62  (
63  geometricOneField(),
64  rho,
65  U,
66  phi,
67  phi,
68  thermophysicalModel,
69  turbulenceModelName,
70  modelName
71  ),
72 
73  C4_
74  (
75  dimensioned<scalar>::getOrAddToDict
76  (
77  "C4",
78  coeffDict_,
79  0.1
80  )
81  )
82 {}
83 
84 
85 // * * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * //
86 
87 PDRkEpsilon::~PDRkEpsilon()
88 {}
89 
90 
91 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
92 
93 bool PDRkEpsilon::read()
94 {
95  if (RASModel::read())
96  {
97  C4_.readIfPresent(coeffDict_);
98  return true;
99  }
100 
101  return false;
102 }
103 
104 
106 {
107  if (!turbulence_)
108  {
109  // Re-calculate viscosity
110  nut_ = Cmu_*sqr(k_)/epsilon_;
111  nut_.correctBoundaryConditions();
112 
113  // Re-calculate thermal diffusivity
114  //***HGWalphat_ = mut_/Prt_;
115  //alphat_.correctBoundaryConditions();
116 
117  return;
118  }
119 
121 
123 
124  if (mesh_.moving())
125  {
126  divU += fvc::div(mesh_.phi());
127  }
128 
129  tmp<volTensorField> tgradU = fvc::grad(U_);
130  volScalarField G(GName(), rho_*nut_*(tgradU() && dev(twoSymm(tgradU()))));
131  tgradU.clear();
132 
133  // Update epsilon and G at the wall
134  epsilon_.boundaryFieldRef().updateCoeffs();
135 
136  // Add the blockage generation term so that it is included consistently
137  // in both the k and epsilon equations
138  const volScalarField& betav =
139  U_.db().lookupObject<volScalarField>("betav");
140 
141  const volScalarField& Lobs =
142  U_.db().lookupObject<volScalarField>("Lobs");
143 
144  const PDRDragModel& drag =
145  U_.db().lookupObject<PDRDragModel>("PDRDragModel");
146 
147  volScalarField GR(drag.Gk());
148 
149  volScalarField LI
150  (C4_*(Lobs + dimensionedScalar("minLength", dimLength, VSMALL)));
151 
152  // Dissipation equation
153  tmp<fvScalarMatrix> epsEqn
154  (
155  betav*fvm::ddt(rho_, epsilon_)
156  + fvm::div(phi_, epsilon_)
157  - fvm::laplacian(rho_*DepsilonEff(), epsilon_)
158  ==
159  C1_*betav*G*epsilon_/k_
160  + 1.5*pow(Cmu_, 3.0/4.0)*GR*sqrt(k_)/LI
161  - fvm::SuSp(((2.0/3.0)*C1_)*betav*rho_*divU, epsilon_)
162  - fvm::Sp(C2_*betav*rho_*epsilon_/k_, epsilon_)
163  );
164 
165  epsEqn.ref().relax();
166 
167  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
168 
169  solve(epsEqn);
170  bound(epsilon_, epsilonMin_);
171 
172 
173  // Turbulent kinetic energy equation
174 
175  tmp<fvScalarMatrix> kEqn
176  (
177  betav*fvm::ddt(rho_, k_)
178  + fvm::div(phi_, k_)
179  - fvm::laplacian(rho_*DkEff(), k_)
180  ==
181  betav*G + GR
182  - fvm::SuSp((2.0/3.0)*betav*rho_*divU, k_)
183  - fvm::Sp(betav*rho_*epsilon_/k_, k_)
184  );
185 
186  kEqn.ref().relax();
187  solve(kEqn);
188  bound(k_, kMin_);
189 
190  // Re-calculate viscosity
191  nut_ = Cmu_*sqr(k_)/epsilon_;
192  nut_.correctBoundaryConditions();
193 
194  // Re-calculate thermal diffusivity
195  //***HGWalphat_ = mut_/Prt_;
196  //alphat_.correctBoundaryConditions();
197 }
198 
199 
200 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
201 
202 } // End namespace RASModels
203 } // End namespace compressible
204 } // End namespace Foam
205 
206 // ************************************************************************* //
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
drag
Info<< "Reading strained laminar flame speed field Su\n"<< endl;volScalarField Su(IOobject("Su", runTime.timeName(), 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:167
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:47
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
Foam::constant::universal::G
const dimensionedScalar G
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Definition: readThermalProperties.H:212
Foam::bound
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Definition: bound.C:28
Foam::compressible::RASModels::PDRkEpsilon::PDRkEpsilon
PDRkEpsilon(const geometricOneField &alpha, const volScalarField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const fluidThermo &thermophysicalModel, const word &turbulenceModelName=turbulenceModel::typeName, const word &modelName=typeName)
Foam::fvc::div
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:42
Foam::blockMeshTools::read
void read(Istream &, label &val, const dictionary &)
Definition: blockMeshTools.C:50
Foam::fvm::SuSp
tmp< fvMatrix< Type > > SuSp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
Foam::fvm::Sp
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
Foam::compressible::turbulenceModel
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
Definition: turbulentFluidThermoModel.H:59
Foam::fvm::laplacian
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:41
rho
rho
Definition: readInitialConditions.H:88
correct
fvOptions correct(rho)
Foam::solve
SolverPerformance< Type > solve(faMatrix< Type > &, Istream &)
divU
zeroField divU
Definition: alphaSuSp.H:3
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:36
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
PDRkEpsilon.H
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:53
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:68
compressible
bool compressible
Definition: pEqn.H:2
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Definition: atmBoundaryLayer.C:26
Foam::fvc::interpolate
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
betav
const volScalarField & betav
Definition: setRegionSolidFields.H:35
Foam::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:58
Foam::RASModel::correct
virtual void correct()
Definition: RASModel.C:223
U
U
Definition: pEqn.H:72
Foam::fvm::ddt
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:41
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:44
Foam::compressible::RASModel
RASModel< EddyDiffusivity< turbulenceModel > > RASModel
Definition: turbulentFluidThermoModel.H:62
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:50
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:137
PDRDragModel.H
Foam::fvm::div
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:41
Foam::twoSymm
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:88
Foam::compressible::defineTypeNameAndDebug
defineTypeNameAndDebug(alphatPhaseChangeWallFunctionFvPatchScalarField, 0)
Foam::RASModel::read
virtual bool read()
Definition: RASModel.C:164
Foam::dev
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:99