EDC.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) 2017-2021 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 "EDC.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 template<>
33  {"v1981", "v1996", "v2005", "v2016"};
34 
37 
41 
42 namespace Foam
43 {
44 namespace combustionModels
45 {
48 }
49 }
50 
51 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
52 
54 (
55  const word& modelType,
56  const fluidReactionThermo& thermo,
58  const word& combustionProperties
59 )
60 :
61  combustionModel(modelType, thermo, turb, combustionProperties),
62  version_
63  (
65  [
66  this->coeffs().lookupOrDefault
67  (
68  "version",
70  )
71  ]
72  ),
73  C1_(this->coeffs().lookupOrDefault("C1", 0.05774)),
74  C2_(this->coeffs().lookupOrDefault("C2", 0.5)),
75  Cgamma_(this->coeffs().lookupOrDefault("Cgamma", 2.1377)),
76  Ctau_(this->coeffs().lookupOrDefault("Ctau", 0.4083)),
77  exp1_(this->coeffs().lookupOrDefault("exp1", EDCexp1[int(version_)])),
78  exp2_(this->coeffs().lookupOrDefault("exp2", EDCexp2[int(version_)])),
79  kappa_
80  (
81  IOobject
82  (
83  this->thermo().phasePropertyName(typeName + ":kappa"),
84  this->mesh().time().timeName(),
85  this->mesh(),
88  ),
89  this->mesh(),
91  ),
92  outerCorrect_
93  (
94  this->coeffs().lookupOrDefault("outerCorrect", true)
95  ),
96  timeIndex_(-1),
97  chemistryPtr_(basicChemistryModel::New(thermo))
98 {}
99 
100 
101 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
102 
104 {}
105 
106 
107 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
108 
110 {
111  if (!outerCorrect_ && timeIndex_ == this->mesh().time().timeIndex())
112  {
113  return;
114  }
115 
116  tmp<volScalarField> tepsilon(this->turbulence().epsilon());
117  const volScalarField& epsilon = tepsilon();
118 
119  tmp<volScalarField> tmu(this->turbulence().mu());
120  const volScalarField& mu = tmu();
121 
122  tmp<volScalarField> tk(this->turbulence().k());
123  const volScalarField& k = tk();
124 
125  tmp<volScalarField> trho(this->rho());
126  const volScalarField& rho = trho();
127 
128  scalarField tauStar(epsilon.size(), 0);
129 
130  if (version_ == EDCversions::v2016)
131  {
132  tmp<volScalarField> ttc(chemistryPtr_->tc());
133  const volScalarField& tc = ttc();
134 
135  forAll(tauStar, i)
136  {
137  const scalar nu = mu[i]/(rho[i] + small);
138 
139  const scalar Da =
140  max(min(sqrt(nu/(epsilon[i] + small))/tc[i], 10), 1e-10);
141 
142  const scalar ReT = sqr(k[i])/(nu*epsilon[i] + small);
143  const scalar CtauI = min(C1_/(Da*sqrt(ReT + 1)), 2.1377);
144 
145  const scalar CgammaI =
146  max(min(C2_*sqrt(Da*(ReT + 1)), 5), 0.4082);
147 
148  const scalar gammaL =
149  CgammaI*pow025(nu*epsilon[i]/(sqr(k[i]) + small));
150 
151  tauStar[i] = CtauI*sqrt(nu/(epsilon[i] + small));
152 
153  if (gammaL >= 1)
154  {
155  kappa_[i] = 1;
156  }
157  else
158  {
159  kappa_[i] =
160  max
161  (
162  min
163  (
164  pow(gammaL, exp1_)/(1 - pow(gammaL, exp2_)),
165  1
166  ),
167  0
168  );
169  }
170  }
171  }
172  else
173  {
174  forAll(tauStar, i)
175  {
176  const scalar nu = mu[i]/(rho[i] + small);
177  const scalar gammaL =
178  Cgamma_*pow025(nu*epsilon[i]/(sqr(k[i]) + small));
179 
180  tauStar[i] = Ctau_*sqrt(nu/(epsilon[i] + small));
181  if (gammaL >= 1)
182  {
183  kappa_[i] = 1;
184  }
185  else
186  {
187  kappa_[i] =
188  max
189  (
190  min
191  (
192  pow(gammaL, exp1_)/(1 - pow(gammaL, exp2_)),
193  1
194  ),
195  0
196  );
197  }
198  }
199  }
200 
201  chemistryPtr_->solve(tauStar);
202 
203  timeIndex_ = this->mesh().time().timeIndex();
204 }
205 
206 
209 {
211  fvScalarMatrix& Su = tSu.ref();
212 
213  const label specieI = this->thermo().composition().species()[Y.member()];
214  Su += chemistryPtr_->RR(specieI);
215 
216  return kappa_*tSu;
217 }
218 
219 
222 {
223  return volScalarField::New
224  (
225  this->thermo().phasePropertyName(typeName + ":Qdot"),
226  kappa_*chemistryPtr_->Qdot()
227  );
228 }
229 
230 
232 {
233  if (combustionModel::read())
234  {
235  version_ =
236  (
238  [
239  this->coeffs().lookupOrDefault
240  (
241  "version",
243  )
244  ]
245  );
246  C1_ = this->coeffs().lookupOrDefault("C1", 0.05774);
247  C2_ = this->coeffs().lookupOrDefault("C2", 0.5);
248  Cgamma_ = this->coeffs().lookupOrDefault("Cgamma", 2.1377);
249  Ctau_ = this->coeffs().lookupOrDefault("Ctau", 0.4083);
250  exp1_ = this->coeffs().lookupOrDefault("exp1", EDCexp1[int(version_)]);
251  exp2_ = this->coeffs().lookupOrDefault("exp2", EDCexp2[int(version_)]);
252  outerCorrect_ = this->coeffs().lookupOrDefault("outerCorrect", true);
253 
254  return true;
255  }
256  else
257  {
258  return false;
259  }
260 }
261 
262 
263 // ************************************************************************* //
const EDCversions EDCdefaultVersion
Definition: EDC.C:39
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
fluidReactionThermo & thermo
Definition: createFields.H:28
EDCversions
EDC model versions.
Definition: EDC.H:114
virtual tmp< fvScalarMatrix > R(volScalarField &Y) const
Fuel consumption rate matrix.
Definition: EDC.C:208
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
virtual bool read()
Update properties from given dictionary.
Definition: EDC.C:231
dimensionedSymmTensor sqr(const dimensionedVector &dv)
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensionedScalar pow025(const dimensionedScalar &ds)
tmp< volScalarField > trho
Base-class for multi-component fluid thermodynamic properties.
const dimensionSet dimless
label k
Boltzmann constant.
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:51
Macros for easy insertion into run-time selection tables.
static word member(const word &name)
Return member (name without the extension)
Definition: IOobject.C:149
const scalar EDCexp1[]
Definition: EDC.H:125
static autoPtr< basicChemistryModel > New(const fluidReactionThermo &thermo)
Select based on fluid reaction thermo.
const dimensionSet dimTime
dynamicFvMesh & mesh
EDC(const word &modelType, const fluidReactionThermo &type, const compressibleMomentumTransportModel &turb, const word &combustionProperties)
Construct from components.
Definition: EDC.C:54
A class for handling words, derived from string.
Definition: word.H:59
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:72
word timeName
Definition: getTimeIndex.H:3
virtual bool read()
Update properties from given dictionary.
Info<< "Reading field U\"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);volScalarField rho(IOobject("rho", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), thermo.rho());volVectorField rhoU(IOobject("rhoU", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *U);volScalarField rhoE(IOobject("rhoE", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *(e+0.5 *magSqr(U)));surfaceScalarField pos(IOobject("pos", runTime.timeName(), mesh), mesh, dimensionedScalar(dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar(dimless, -1.0));surfaceScalarField phi("phi", fvc::flux(rhoU));Info<< "Creating turbulence model\"<< endl;autoPtr< compressible::momentumTransportModel > turbulence(compressible::momentumTransportModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
const dimensionedScalar mu
Atomic mass unit.
Base class for combustion models.
Eddy Dissipation Concept (EDC) turbulent combustion model.
Definition: EDC.H:132
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
const dimensionSet dimMass
const NamedEnum< EDCversions, 4 > EDCversionNames
Definition: EDC.C:36
const scalar EDCexp2[]
Definition: EDC.H:126
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
PtrList< volScalarField > & Y
scalar epsilon
label timeIndex
Definition: getTimeIndex.H:4
virtual tmp< volScalarField > Qdot() const
Heat release rate [kg/m/s^3].
Definition: EDC.C:221
virtual ~EDC()
Destructor.
Definition: EDC.C:103
virtual void correct()
Correct combustion rate.
Definition: EDC.C:109
const doubleScalar e
Elementary charge.
Definition: doubleScalar.H:105
const tmp< volScalarField::Internal > & Su
Definition: alphaSuSp.H:6
A class for managing temporary objects.
Definition: PtrList.H:53
defineTypeNameAndDebug(diffusion, 0)
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
Abstract base class for turbulence models (RAS, LES and laminar).
addToRunTimeSelectionTable(combustionModel, diffusion, dictionary)
Namespace for OpenFOAM.