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-2020,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().phasePropertyName(typeName + ":kappa"),
63  this->mesh().time().timeName(),
64  this->mesh(),
65  IOobject::NO_READ,
66  IOobject::AUTO_WRITE
67  ),
68  this->mesh(),
70  )
71 {}
72 
73 
74 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
75 
76 template<class ReactionThermo>
78 {}
79 
80 
81 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
82 
83 template<class ReactionThermo>
85 {
86  if (this->active())
87  {
88  tmp<volScalarField> tepsilon(this->turbulence().epsilon());
89  const volScalarField& epsilon = tepsilon();
90 
91  tmp<volScalarField> tmu(this->turbulence().mu());
92  const volScalarField& mu = tmu();
93 
94  tmp<volScalarField> tk(this->turbulence().k());
95  const volScalarField& k = tk();
96 
97  tmp<volScalarField> trho(this->rho());
98  const volScalarField& rho = trho();
99 
100  scalarField tauStar(epsilon.size(), Zero);
101 
102  if (version_ == EDCversions::v2016)
103  {
104  tmp<volScalarField> ttc(this->chemistryPtr_->tc());
105  const volScalarField& tc = ttc();
106 
107  forAll(tauStar, i)
108  {
109  const scalar nu = mu[i]/(rho[i] + SMALL);
110 
111  const scalar Da = clamp
112  (
113  sqrt(nu/(epsilon[i] + SMALL))/tc[i],
114  scalar(1e-10), scalar(10)
115  );
116 
117  const scalar ReT = sqr(k[i])/(nu*epsilon[i] + SMALL);
118  const scalar CtauI = min(C1_/(Da*sqrt(ReT + 1)), 2.1377);
119 
120  const scalar CgammaI =
121  clamp(C2_*sqrt(Da*(ReT + 1)), scalar(0.4082), scalar(5));
122 
123  const scalar gammaL =
124  CgammaI*pow025(nu*epsilon[i]/(sqr(k[i]) + SMALL));
125 
126  tauStar[i] = CtauI*sqrt(nu/(epsilon[i] + SMALL));
127 
128  if (gammaL >= 1)
129  {
130  kappa_[i] = 1;
131  }
132  else
133  {
134  kappa_[i] =
135  max
136  (
137  min
138  (
139  pow(gammaL, exp1_)/(1 - pow(gammaL, exp2_)),
140  1
141  ),
142  0
143  );
144  }
145  }
146 
147  // Evaluate bcs
148  kappa_.correctBoundaryConditions();
149  }
150  else
151  {
152  forAll(tauStar, i)
153  {
154  const scalar nu = mu[i]/(rho[i] + SMALL);
155  const scalar gammaL =
156  Cgamma_*pow025(nu*epsilon[i]/(sqr(k[i]) + SMALL));
157 
158  tauStar[i] = Ctau_*sqrt(nu/(epsilon[i] + SMALL));
159  if (gammaL >= 1)
160  {
161  kappa_[i] = 1;
162  }
163  else
164  {
165  kappa_[i] =
166  max
167  (
168  min
169  (
170  pow(gammaL, exp1_)/(1 - pow(gammaL, exp2_)),
171  1
172  ),
173  0
174  );
175  }
176  }
177 
178  // Evaluate bcs
179  kappa_.correctBoundaryConditions();
180  }
181 
182  Info<< "Chemistry time solved max/min : "
183  << gMax(tauStar) << " / " << gMin(tauStar) << endl;
184 
185  this->chemistryPtr_->solve(tauStar);
186  }
187 }
188 
189 
190 template<class ReactionThermo>
193 {
194  return kappa_*laminar<ReactionThermo>::R(Y);
195 }
196 
197 
198 template<class ReactionThermo>
201 {
202  tmp<volScalarField> tQdot
203  (
204  new volScalarField
205  (
206  IOobject
207  (
208  this->thermo().phasePropertyName(typeName + ":Qdot"),
209  this->mesh().time().timeName(),
210  this->mesh(),
214  ),
215  this->mesh(),
217  )
218  );
219 
220  if (this->active())
221  {
222  tQdot.ref() = kappa_*this->chemistryPtr_->Qdot();
223  }
224 
225  return tQdot;
226 }
227 
228 
229 template<class ReactionThermo>
231 {
233  {
234  version_ =
235  (
236  EDCversionNames.getOrDefault
237  (
238  "version",
239  this->coeffs(),
241  )
242  );
243  C1_ = this->coeffs().getOrDefault("C1", 0.05774);
244  C2_ = this->coeffs().getOrDefault("C2", 0.5);
245  Cgamma_ = this->coeffs().getOrDefault("Cgamma", 2.1377);
246  Ctau_ = this->coeffs().getOrDefault("Ctau", 0.4083);
247  exp1_ = this->coeffs().getOrDefault("exp1", EDCexp1[int(version_)]);
248  exp2_ = this->coeffs().getOrDefault("exp2", EDCexp2[int(version_)]);
249 
250  return true;
251  }
252 
253  return false;
254 }
255 
256 
257 // ************************************************************************* //
const EDCversions EDCdefaultVersion
virtual tmp< fvScalarMatrix > R(volScalarField &Y) const
Fuel consumption rate matrix.
Definition: EDC.C:185
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:223
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.
Ignore writing from objectRegistry::writeObject()
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...
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
Nothing to be read.
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:193
messageStream Info
Information stream (stdout output on master, null elsewhere)
virtual ~EDC()
Destructor.
Definition: EDC.C:70
virtual void correct()
Correct combustion rate.
Definition: EDC.C:77
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:172
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:167
dimensionSet clamp(const dimensionSet &a, const dimensionSet &range)
Definition: dimensionSet.C:271
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127