Loading...
Searching...
No Matches
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,2023 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
35namespace Foam
36{
37namespace compressible
38{
39namespace RASModels
40{
41
42// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
43
46
47// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48
50(
51 const geometricOneField& alpha,
52 const volScalarField& rho,
53 const volVectorField& U,
54 const surfaceScalarField& alphaRhoPhi,
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
87PDRkEpsilon::~PDRkEpsilon()
88{}
89
90
91// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
92
93bool PDRkEpsilon::read()
94{
95 if (RASModel::read())
96 {
97 C4_.readIfPresent(coeffDict_);
98 return true;
99 }
100
101 return false;
102}
103
104
105void PDRkEpsilon::correct()
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
120 RASModel::correct();
121
122 volScalarField divU(fvc::div(phi_/fvc::interpolate(rho_)));
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 // Push new cell values to
136 // coupled neighbours. Note that we want to avoid the re-updateCoeffs
137 // of the wallFunctions so make sure to bypass the evaluate on
138 // those patches and only do the coupled ones.
139 epsilon_.boundaryFieldRef().evaluateCoupled<coupledFvPatch>();
140
141 // Add the blockage generation term so that it is included consistently
142 // in both the k and epsilon equations
143 const volScalarField& betav =
144 U_.db().lookupObject<volScalarField>("betav");
145
146 const volScalarField& Lobs =
147 U_.db().lookupObject<volScalarField>("Lobs");
148
149 const PDRDragModel& drag =
150 U_.db().lookupObject<PDRDragModel>("PDRDragModel");
151
152 volScalarField GR(drag.Gk());
153
155 (C4_*(Lobs + dimensionedScalar("minLength", dimLength, VSMALL)));
156
157 // Dissipation equation
158 tmp<fvScalarMatrix> epsEqn
159 (
160 betav*fvm::ddt(rho_, epsilon_)
161 + fvm::div(phi_, epsilon_)
162 - fvm::laplacian(rho_*DepsilonEff(), epsilon_)
163 ==
164 C1_*betav*G*epsilon_/k_
165 + 1.5*pow(Cmu_, 3.0/4.0)*GR*sqrt(k_)/LI
166 - fvm::SuSp(((2.0/3.0)*C1_)*betav*rho_*divU, epsilon_)
167 - fvm::Sp(C2_*betav*rho_*epsilon_/k_, epsilon_)
168 );
169
170 epsEqn.ref().relax();
171
172 epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
173
174 solve(epsEqn);
175 bound(epsilon_, epsilonMin_);
176
177
178 // Turbulent kinetic energy equation
179
180 tmp<fvScalarMatrix> kEqn
181 (
182 betav*fvm::ddt(rho_, k_)
183 + fvm::div(phi_, k_)
184 - fvm::laplacian(rho_*DkEff(), k_)
185 ==
186 betav*G + GR
187 - fvm::SuSp((2.0/3.0)*betav*rho_*divU, k_)
188 - fvm::Sp(betav*rho_*epsilon_/k_, k_)
189 );
190
191 kEqn.ref().relax();
192 solve(kEqn);
193 bound(k_, kMin_);
194
195 // Re-calculate viscosity
196 nut_ = Cmu_*sqr(k_)/epsilon_;
197 nut_.correctBoundaryConditions();
198
199 // Re-calculate thermal diffusivity
200 //***HGWalphat_ = mut_/Prt_;
201 //alphat_.correctBoundaryConditions();
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.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
const objectRegistry & db() const noexcept
Return the local objectRegistry.
Definition IOobject.C:450
Standard k-epsilon turbulence model with additional source terms corresponding to PDR basic drag mode...
Definition PDRkEpsilon.H:85
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)
Construct from components.
const Type & lookupObject(const word &name, const bool recursive=false) const
Lookup and return const reference to the object of the given Type. Fatal if not found or the wrong ty...
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
U
Definition pEqn.H:72
bool compressible
Definition pEqn.H:2
zeroField divU
Definition alphaSuSp.H:3
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
RASModel< EddyDiffusivity< turbulenceModel > > RASModel
const dimensionedScalar G
Newtonian constant of gravitation.
Namespace for OpenFOAM.
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
GeometricField< vector, fvPatchField, volMesh > volVectorField
dimensionedSymmTensor sqr(const dimensionedVector &dv)
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition bound.C:29
GeometricField< scalar, fvPatchField, volMesh > volScalarField
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
volScalarField & alpha
const volScalarField & betav
CEqn solve()
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