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-2025 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 
33 {"v1981", "v1996", "v2005", "v2016"};
34 
38 
39 namespace Foam
40 {
41 namespace combustionModels
42 {
45 }
46 }
47 
48 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
49 
51 (
52  const word& modelType,
55  const word& combustionProperties
56 )
57 :
58  combustionModel(modelType, thermo, turb, combustionProperties),
59  version_
60  (
62  [
63  this->coeffs().lookupOrDefault
64  (
65  "version",
67  )
68  ]
69  ),
70  C1_(this->coeffs().lookupOrDefault("C1", 0.05774)),
71  C2_(this->coeffs().lookupOrDefault("C2", 0.5)),
72  Cgamma_(this->coeffs().lookupOrDefault("Cgamma", 2.1377)),
73  Ctau_(this->coeffs().lookupOrDefault("Ctau", 0.4083)),
74  exp1_(this->coeffs().lookupOrDefault("exp1", EDCexp1[int(version_)])),
75  exp2_(this->coeffs().lookupOrDefault("exp2", EDCexp2[int(version_)])),
76  kappa_
77  (
78  IOobject
79  (
80  this->thermo().phasePropertyName(typedName("kappa")),
81  this->mesh().time().name(),
82  this->mesh(),
83  IOobject::NO_READ,
84  IOobject::AUTO_WRITE
85  ),
86  this->mesh(),
88  ),
89  outerCorrect_
90  (
91  this->coeffs().lookupOrDefault("outerCorrect", true)
92  ),
93  timeIndex_(-1),
94  chemistryPtr_(basicChemistryModel::New(thermo))
95 {}
96 
97 
98 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
99 
101 {}
102 
103 
104 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
105 
107 {
108  if (!outerCorrect_ && timeIndex_ == this->mesh().time().timeIndex())
109  {
110  return;
111  }
112 
113  tmp<volScalarField> tepsilon(this->turbulence().epsilon());
114  const volScalarField& epsilon = tepsilon();
115 
116  tmp<volScalarField> tnu(this->turbulence().nu());
117  const volScalarField& nu = tnu();
118 
119  tmp<volScalarField> tk(this->turbulence().k());
120  const volScalarField& k = tk();
121 
122  scalarField tauStar(epsilon.size(), 0);
123 
124  if (version_ == EDCversions::v2016)
125  {
126  tmp<volScalarField> ttc(chemistryPtr_->tc());
127  const volScalarField& tc = ttc();
128 
129  forAll(tauStar, i)
130  {
131  const scalar Da =
132  max(min(sqrt(nu[i]/(epsilon[i] + small))/tc[i], 10), 1e-10);
133 
134  const scalar ReT = sqr(k[i])/(nu[i]*epsilon[i] + small);
135  const scalar CtauI = min(C1_/(Da*sqrt(ReT + 1)), Ctau_);
136 
137  const scalar CgammaI =
138  max(min(C2_*sqrt(Da*(ReT + 1)), 5), Cgamma_);
139 
140  const scalar gammaL =
141  CgammaI*pow025(nu[i]*epsilon[i]/(sqr(k[i]) + small));
142 
143  tauStar[i] = CtauI*sqrt(nu[i]/(epsilon[i] + small));
144 
145  if (gammaL >= 1)
146  {
147  kappa_[i] = 1;
148  }
149  else
150  {
151  kappa_[i] =
152  max
153  (
154  min
155  (
156  pow(gammaL, exp1_)/(1 - pow(gammaL, exp2_)),
157  1
158  ),
159  0
160  );
161  }
162  }
163  }
164  else
165  {
166  forAll(tauStar, i)
167  {
168  const scalar gammaL =
169  Cgamma_*pow025(nu[i]*epsilon[i]/(sqr(k[i]) + small));
170 
171  tauStar[i] = Ctau_*sqrt(nu[i]/(epsilon[i] + small));
172  if (gammaL >= 1)
173  {
174  kappa_[i] = 1;
175  }
176  else
177  {
178  kappa_[i] =
179  max
180  (
181  min
182  (
183  pow(gammaL, exp1_)/(1 - pow(gammaL, exp2_)),
184  1
185  ),
186  0
187  );
188  }
189  }
190  }
191 
192  chemistryPtr_->solve(tauStar);
193 
194  timeIndex_ = this->mesh().time().timeIndex();
195 }
196 
197 
200 {
201  return kappa_*chemistryPtr_->RR()[speciei];
202 }
203 
204 
207 {
209  fvScalarMatrix& Su = tSu.ref();
210 
211  const label speciei = this->thermo().species()[Y.member()];
212  Su += chemistryPtr_->RR()[speciei];
213 
214  return kappa_*tSu;
215 }
216 
217 
220 {
221  return volScalarField::New
222  (
223  this->thermo().phasePropertyName(typedName("Qdot")),
224  kappa_*chemistryPtr_->Qdot()
225  );
226 }
227 
228 
230 {
231  if (combustionModel::read())
232  {
233  version_ =
234  (
236  [
237  this->coeffs().lookupOrDefault
238  (
239  "version",
241  )
242  ]
243  );
244  C1_ = this->coeffs().lookupOrDefault("C1", 0.05774);
245  C2_ = this->coeffs().lookupOrDefault("C2", 0.5);
246  Cgamma_ = this->coeffs().lookupOrDefault("Cgamma", 2.1377);
247  Ctau_ = this->coeffs().lookupOrDefault("Ctau", 0.4083);
248  exp1_ = this->coeffs().lookupOrDefault("exp1", EDCexp1[int(version_)]);
249  exp2_ = this->coeffs().lookupOrDefault("exp2", EDCexp2[int(version_)]);
250  outerCorrect_ = this->coeffs().lookupOrDefault("outerCorrect", true);
251 
252  return true;
253  }
254  else
255  {
256  return false;
257  }
258 }
259 
260 
261 // ************************************************************************* //
label k
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
static tmp< GeometricField< Type, GeoMesh, PrimitiveField > > New(const word &name, const Internal &, const PtrList< Patch > &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
static word member(const word &name)
Return member (name without the extension)
Definition: IOobject.C:146
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:55
label timeIndex() const
Return current time index.
Definition: TimeStateI.H:28
Base class for chemistry models.
Base class for combustion models.
virtual bool read()
Update properties from given dictionary.
Eddy Dissipation Concept (EDC) turbulent combustion model.
Definition: EDC.H:135
virtual void correct()
Correct combustion rate.
Definition: EDC.C:106
EDC(const word &modelType, const fluidMulticomponentThermo &type, const compressibleMomentumTransportModel &turb, const word &combustionProperties)
Construct from components.
Definition: EDC.C:51
virtual ~EDC()
Destructor.
Definition: EDC.C:100
virtual tmp< volScalarField::Internal > R(const label speciei) const
Specie consumption rate field.
Definition: EDC.C:199
virtual tmp< volScalarField > Qdot() const
Heat release rate [kg/m/s^3].
Definition: EDC.C:219
virtual bool read()
Update properties from given dictionary.
Definition: EDC.C:229
Base class for single-phase compressible turbulence models.
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
Base-class for multi-component fluid thermodynamic properties.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:420
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:197
A class for handling words, derived from string.
Definition: word.H:62
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
const scalar epsilon
defineTypeNameAndDebug(diffusion, 0)
const scalar EDCexp2[]
Definition: EDC.H:126
const scalar EDCexp1[]
Definition: EDC.H:125
const NamedEnum< EDCversions, 4 > EDCversionNames
Definition: EDC.C:33
const EDCversions EDCdefaultVersion
Definition: EDC.C:36
EDCversions
EDC model versions.
Definition: EDC.H:115
addToRunTimeSelectionTable(combustionModel, diffusion, dictionary)
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
tmp< VolField< Type > > Su(const VolField< Type > &su, const VolField< Type > &vf)
Definition: fvcSup.C:44
Namespace for OpenFOAM.
const doubleScalar e
Definition: doubleScalar.H:106
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
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
void pow025(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
const dimensionSet dimless
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
const dimensionSet dimTime
void pow(LagrangianPatchField< typename powProduct< Type, r >::type > &f, const LagrangianPatchField< Type > &f1)
word typedName(Name name)
Return the name of the object within the given type.
Definition: typeInfo.H:181
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
void sqr(LagrangianPatchField< typename outerProduct< Type, Type >::type > &f, const LagrangianPatchField< Type > &f1)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
const dimensionSet dimMass
void sqrt(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
label timeIndex
Definition: getTimeIndex.H:4
autoPtr< incompressible::momentumTransportModel > turbulence(incompressible::momentumTransportModel::New(U, phi, viscosity))
PtrList< volScalarField > & Y
fluidMulticomponentThermo & thermo
Definition: createFields.H:31