EDC.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) 2017 OpenFOAM Foundation
9  Copyright (C) 2019-2023 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 "EDC.H"
30 
31 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32 
33 template<class ReactionThermo>
35 (
36  const word& modelType,
37  ReactionThermo& thermo,
39  const word& combustionProperties
40 )
41 :
42  laminar<ReactionThermo>(modelType, thermo, turb, combustionProperties),
43  version_
44  (
45  EDCversionNames.getOrDefault
46  (
47  "version",
48  this->coeffs(),
50  )
51  ),
52  C1_(this->coeffs().getOrDefault("C1", 0.05774)),
53  C2_(this->coeffs().getOrDefault("C2", 0.5)),
54  Cgamma_(this->coeffs().getOrDefault("Cgamma", 2.1377)),
55  Ctau_(this->coeffs().getOrDefault("Ctau", 0.4083)),
56  exp1_(this->coeffs().getOrDefault("exp1", EDCexp1[int(version_)])),
57  exp2_(this->coeffs().getOrDefault("exp2", EDCexp2[int(version_)])),
58  kappa_
59  (
60  IOobject
61  (
62  this->thermo().phaseScopedName(typeName, "kappa"),
63  this->mesh().time().timeName(),
64  this->mesh(),
65  IOobject::NO_READ,
66  IOobject::AUTO_WRITE,
67  IOobject::REGISTER
68  ),
69  this->mesh(),
71  )
72 {}
73 
74 
75 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
76 
77 template<class ReactionThermo>
79 {}
80 
81 
82 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
83 
84 template<class ReactionThermo>
86 {
87  if (this->active())
88  {
89  tmp<volScalarField> tepsilon(this->turbulence().epsilon());
90  const auto& epsilon = tepsilon();
91 
92  tmp<volScalarField> tmu(this->turbulence().mu());
93  const auto& mu = tmu();
94 
95  tmp<volScalarField> tk(this->turbulence().k());
96  const auto& k = tk();
97 
98  tmp<volScalarField> trho(this->rho());
99  const auto& rho = trho();
100 
101  scalarField tauStar(epsilon.size(), Zero);
102 
103  if (version_ == EDCversions::v2016)
104  {
105  tmp<volScalarField> ttc(this->chemistryPtr_->tc());
106  const auto& tc = ttc();
107 
108  forAll(tauStar, i)
109  {
110  const scalar nu = mu[i]/(rho[i] + SMALL);
111 
112  const scalar Da = clamp
113  (
114  sqrt(nu/(epsilon[i] + SMALL))/tc[i],
115  scalar(1e-10), scalar(10)
116  );
117 
118  const scalar ReT = sqr(k[i])/(nu*epsilon[i] + SMALL);
119  const scalar CtauI = min(C1_/(Da*sqrt(ReT + 1)), 2.1377);
120 
121  const scalar CgammaI =
122  clamp(C2_*sqrt(Da*(ReT + 1)), scalar(0.4082), scalar(5));
123 
124  const scalar gammaL =
125  CgammaI*pow025(nu*epsilon[i]/(sqr(k[i]) + SMALL));
126 
127  tauStar[i] = CtauI*sqrt(nu/(epsilon[i] + SMALL));
128 
129  if (gammaL >= 1)
130  {
131  kappa_[i] = 1;
132  }
133  else
134  {
135  kappa_[i] =
136  max
137  (
138  min
139  (
140  pow(gammaL, exp1_)/(1 - pow(gammaL, exp2_)),
141  1
142  ),
143  0
144  );
145  }
146  }
147 
148  // Evaluate bcs
149  kappa_.correctBoundaryConditions();
150  }
151  else
152  {
153  forAll(tauStar, i)
154  {
155  const scalar nu = mu[i]/(rho[i] + SMALL);
156  const scalar gammaL =
157  Cgamma_*pow025(nu*epsilon[i]/(sqr(k[i]) + SMALL));
158 
159  tauStar[i] = Ctau_*sqrt(nu/(epsilon[i] + SMALL));
160  if (gammaL >= 1)
161  {
162  kappa_[i] = 1;
163  }
164  else
165  {
166  kappa_[i] =
167  max
168  (
169  min
170  (
171  pow(gammaL, exp1_)/(1 - pow(gammaL, exp2_)),
172  1
173  ),
174  0
175  );
176  }
177  }
178 
179  // Evaluate bcs
180  kappa_.correctBoundaryConditions();
181  }
182 
183  Info<< "Chemistry time solved max/min : "
184  << gMax(tauStar) << " / " << gMin(tauStar) << endl;
185 
186  this->chemistryPtr_->solve(tauStar);
187  }
188 }
189 
190 
191 template<class ReactionThermo>
194 {
195  return kappa_*laminar<ReactionThermo>::R(Y);
196 }
197 
198 
199 template<class ReactionThermo>
202 {
203  auto tQdot = volScalarField::New
204  (
205  this->thermo().phaseScopedName(typeName, "Qdot"),
207  this->mesh(),
209  );
210 
211  if (this->active())
212  {
213  tQdot.ref() = kappa_*this->chemistryPtr_->Qdot();
214  }
215 
216  return tQdot;
217 }
218 
219 
220 template<class ReactionThermo>
222 {
224  {
225  version_ =
226  (
227  EDCversionNames.getOrDefault
228  (
229  "version",
230  this->coeffs(),
232  )
233  );
234  C1_ = this->coeffs().getOrDefault("C1", 0.05774);
235  C2_ = this->coeffs().getOrDefault("C2", 0.5);
236  Cgamma_ = this->coeffs().getOrDefault("Cgamma", 2.1377);
237  Ctau_ = this->coeffs().getOrDefault("Ctau", 0.4083);
238  exp1_ = this->coeffs().getOrDefault("exp1", EDCexp1[int(version_)]);
239  exp2_ = this->coeffs().getOrDefault("exp2", EDCexp2[int(version_)]);
240 
241  return true;
242  }
243 
244  return false;
245 }
246 
247 
248 // ************************************************************************* //
const EDCversions EDCdefaultVersion
virtual tmp< fvScalarMatrix > R(volScalarField &Y) const
Fuel consumption rate matrix.
Definition: EDC.C:186
Type gMin(const FieldField< Field, Type > &f)
const Enum< EDCversions > EDCversionNames
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
virtual bool read()
Update properties from given dictionary.
Definition: EDC.C:214
dimensionedSymmTensor sqr(const dimensionedVector &dv)
compressible::turbulenceModel & turb
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
dimensionedScalar pow025(const dimensionedScalar &ds)
tmp< volScalarField > trho
label k
Boltzmann constant.
const dimensionSet dimless
Dimensionless.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
const scalar EDCexp1[]
Definition: EDC.H:123
psiReactionThermo & thermo
Definition: createFields.H:28
word timeName
Definition: getTimeIndex.H:3
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
Laminar combustion model.
Definition: laminar.H:52
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
dynamicFvMesh & mesh
A class for handling words, derived from Foam::string.
Definition: word.H:63
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions and patch type. [Takes current timeN...
Abstract base class for turbulence models (RAS, LES and laminar).
Type gMax(const FieldField< Field, Type > &f)
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("pos", dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar("neg", dimless, -1.0));surfaceScalarField phi("phi", fvc::flux(rhoU));Info<< "Creating turbulence model\"<< endl;autoPtr< compressible::turbulenceModel > turbulence(compressible::turbulenceModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
Eddy Dissipation Concept (EDC) turbulent combustion model.
Definition: EDC.H:131
const dimensionSet dimEnergy
const scalar EDCexp2[]
Definition: EDC.H:124
const dimensionedScalar mu
Atomic mass unit.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
PtrList< volScalarField > & Y
scalar epsilon
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
virtual tmp< volScalarField > Qdot() const
Heat release rate [kg/m/s3].
Definition: EDC.C:194
messageStream Info
Information stream (stdout output on master, null elsewhere)
virtual ~EDC()
Destructor.
Definition: EDC.C:71
virtual void correct()
Correct combustion rate.
Definition: EDC.C:78
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:180
volScalarField & nu
virtual tmp< fvScalarMatrix > R(volScalarField &Y) const
Fuel consumption rate matrix.
Definition: laminar.C:117
Do not request registration (bool: false)
virtual bool read()
Update properties from given dictionary.
Definition: laminar.C:156
dimensionSet clamp(const dimensionSet &a, const dimensionSet &range)
Definition: dimensionSet.C:271
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127