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-2018 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"
27 
28 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
29 
30 template<class ReactionThermo>
32 (
33  const word& modelType,
34  ReactionThermo& thermo,
35  const compressibleTurbulenceModel& turb,
36  const word& combustionProperties
37 )
38 :
39  laminar<ReactionThermo>(modelType, thermo, turb, combustionProperties),
40  version_
41  (
43  [
44  this->coeffs().lookupOrDefault
45  (
46  "version",
47  word(EDCversionNames[EDCdefaultVersion])
48  )
49  ]
50  ),
51  C1_(this->coeffs().lookupOrDefault("C1", 0.05774)),
52  C2_(this->coeffs().lookupOrDefault("C2", 0.5)),
53  Cgamma_(this->coeffs().lookupOrDefault("Cgamma", 2.1377)),
54  Ctau_(this->coeffs().lookupOrDefault("Ctau", 0.4083)),
55  exp1_(this->coeffs().lookupOrDefault("exp1", EDCexp1[int(version_)])),
56  exp2_(this->coeffs().lookupOrDefault("exp2", EDCexp2[int(version_)])),
57  kappa_
58  (
59  IOobject
60  (
61  this->thermo().phasePropertyName(typeName + ":kappa"),
62  this->mesh().time().timeName(),
63  this->mesh(),
64  IOobject::NO_READ,
65  IOobject::AUTO_WRITE
66  ),
67  this->mesh(),
68  dimensionedScalar("kappa", dimless, 0)
69  )
70 {}
71 
72 
73 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
74 
75 template<class ReactionThermo>
77 {}
78 
79 
80 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
81 
82 template<class ReactionThermo>
84 {
85  if (this->active())
86  {
87  tmp<volScalarField> tepsilon(this->turbulence().epsilon());
88  const volScalarField& epsilon = tepsilon();
89 
90  tmp<volScalarField> tmu(this->turbulence().mu());
91  const volScalarField& mu = tmu();
92 
93  tmp<volScalarField> tk(this->turbulence().k());
94  const volScalarField& k = tk();
95 
96  tmp<volScalarField> trho(this->rho());
97  const volScalarField& rho = trho();
98 
99  scalarField tauStar(epsilon.size(), 0);
100 
101  if (version_ == EDCversions::v2016)
102  {
103  tmp<volScalarField> ttc(this->chemistryPtr_->tc());
104  const volScalarField& tc = ttc();
105 
106  forAll(tauStar, i)
107  {
108  const scalar nu = mu[i]/(rho[i] + small);
109 
110  const scalar Da =
111  max(min(sqrt(nu/(epsilon[i] + small))/tc[i], 10), 1e-10);
112 
113  const scalar ReT = sqr(k[i])/(nu*epsilon[i] + small);
114  const scalar CtauI = min(C1_/(Da*sqrt(ReT + 1)), 2.1377);
115 
116  const scalar CgammaI =
117  max(min(C2_*sqrt(Da*(ReT + 1)), 5), 0.4082);
118 
119  const scalar gammaL =
120  CgammaI*pow025(nu*epsilon[i]/(sqr(k[i]) + small));
121 
122  tauStar[i] = CtauI*sqrt(nu/(epsilon[i] + small));
123 
124  if (gammaL >= 1)
125  {
126  kappa_[i] = 1;
127  }
128  else
129  {
130  kappa_[i] =
131  max
132  (
133  min
134  (
135  pow(gammaL, exp1_)/(1 - pow(gammaL, exp2_)),
136  1
137  ),
138  0
139  );
140  }
141  }
142  }
143  else
144  {
145  forAll(tauStar, i)
146  {
147  const scalar nu = mu[i]/(rho[i] + small);
148  const scalar gammaL =
149  Cgamma_*pow025(nu*epsilon[i]/(sqr(k[i]) + small));
150 
151  tauStar[i] = Ctau_*sqrt(nu/(epsilon[i] + small));
152  if (gammaL >= 1)
153  {
154  kappa_[i] = 1;
155  }
156  else
157  {
158  kappa_[i] =
159  max
160  (
161  min
162  (
163  pow(gammaL, exp1_)/(1 - pow(gammaL, exp2_)),
164  1
165  ),
166  0
167  );
168  }
169  }
170  }
171 
172  this->chemistryPtr_->solve(tauStar);
173  }
174 }
175 
176 
177 template<class ReactionThermo>
180 {
181  return kappa_*laminar<ReactionThermo>::R(Y);
182 }
183 
184 
185 template<class ReactionThermo>
188 {
189  tmp<volScalarField> tQdot
190  (
191  new volScalarField
192  (
193  IOobject
194  (
195  this->thermo().phasePropertyName(typeName + ":Qdot"),
196  this->mesh().time().timeName(),
197  this->mesh(),
198  IOobject::NO_READ,
199  IOobject::NO_WRITE,
200  false
201  ),
202  this->mesh(),
204  )
205  );
206 
207  if (this->active())
208  {
209  tQdot.ref() = kappa_*this->chemistryPtr_->Qdot();
210  }
211 
212  return tQdot;
213 }
214 
215 
216 template<class ReactionThermo>
218 {
220  {
221  version_ =
222  (
224  [
225  this->coeffs().lookupOrDefault
226  (
227  "version",
228  word(EDCversionNames[EDCdefaultVersion])
229  )
230  ]
231  );
232  C1_ = this->coeffs().lookupOrDefault("C1", 0.05774);
233  C2_ = this->coeffs().lookupOrDefault("C2", 0.5);
234  Cgamma_ = this->coeffs().lookupOrDefault("Cgamma", 2.1377);
235  Ctau_ = this->coeffs().lookupOrDefault("Ctau", 0.4083);
236  exp1_ = this->coeffs().lookupOrDefault("exp1", EDCexp1[int(version_)]);
237  exp2_ = this->coeffs().lookupOrDefault("exp2", EDCexp2[int(version_)]);
238 
239  return true;
240  }
241  else
242  {
243  return false;
244  }
245 }
246 
247 
248 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
virtual tmp< fvScalarMatrix > R(volScalarField &Y) const
Fuel consumption rate matrix.
Definition: EDC.C:179
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
virtual bool read()
Update properties from given dictionary.
Definition: EDC.C:217
dimensionedSymmTensor sqr(const dimensionedVector &dv)
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:163
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensionedScalar pow025(const dimensionedScalar &ds)
tmp< volScalarField > trho
label k
Boltzmann constant.
rhoReactionThermo & thermo
Definition: createFields.H:28
const scalar EDCexp1[]
Definition: EDC.H:123
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
Laminar combustion model.
Definition: laminar.H:51
const dimensionedScalar e
Elementary charge.
Definition: doubleScalar.H:98
dynamicFvMesh & mesh
A class for handling words, derived from string.
Definition: word.H:59
word timeName
Definition: getTimeIndex.H:3
Abstract base class for turbulence models (RAS, LES and laminar).
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
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
const dimensionSet dimEnergy
const NamedEnum< EDCversions, 4 > EDCversionNames
Definition: EDCs.C:48
const scalar EDCexp2[]
Definition: EDC.H:124
const dimensionedScalar mu
Atomic mass unit.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
#define R(A, B, C, D, E, F, K, M)
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:187
virtual ~EDC()
Destructor.
Definition: EDC.C:76
virtual void correct()
Correct combustion rate.
Definition: EDC.C:83
A class for managing temporary objects.
Definition: PtrList.H:53
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
volScalarField & nu