Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
11License
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
33template<class ReactionThermo>
34Foam::combustionModels::EDC<ReactionThermo>::EDC
35(
36 const word& modelType,
37 ReactionThermo& thermo,
39 const word& combustionProperties
40)
41:
42 laminar<ReactionThermo>(modelType, thermo, turb, combustionProperties),
43 version_
44 (
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 (
61 (
62 this->thermo().phaseScopedName(typeName, "kappa"),
63 this->mesh().time().timeName(),
64 this->mesh(),
68 ),
69 this->mesh(),
71 )
72{}
73
74
75// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
76
77template<class ReactionThermo>
79{}
80
81
82// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
83
84template<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
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 auto limits = gMinMax(tauStar);
184 Info<< "Chemistry time solved min/max : "
185 << limits.min() << ", " << limits.max() << endl;
186
187 this->chemistryPtr_->solve(tauStar);
188 }
189}
190
191
192template<class ReactionThermo>
196 return kappa_*laminar<ReactionThermo>::R(Y);
197}
198
199
200template<class ReactionThermo>
203{
204 auto tQdot = volScalarField::New
205 (
206 this->thermo().phaseScopedName(typeName, "Qdot"),
208 this->mesh(),
210 );
211
212 if (this->active())
213 {
214 tQdot.ref() = kappa_*this->chemistryPtr_->Qdot();
216
217 return tQdot;
218}
219
220
221template<class ReactionThermo>
223{
225 {
226 version_ =
227 (
228 EDCversionNames.getOrDefault
229 (
230 "version",
231 this->coeffs(),
233 )
234 );
235 C1_ = this->coeffs().getOrDefault("C1", 0.05774);
236 C2_ = this->coeffs().getOrDefault("C2", 0.5);
237 Cgamma_ = this->coeffs().getOrDefault("Cgamma", 2.1377);
238 Ctau_ = this->coeffs().getOrDefault("Ctau", 0.4083);
239 exp1_ = this->coeffs().getOrDefault("exp1", EDCexp1[int(version_)]);
240 exp2_ = this->coeffs().getOrDefault("exp2", EDCexp2[int(version_)]);
241
242 return true;
243 }
244
245 return false;
246}
247
248
249// ************************************************************************* //
label k
compressible::turbulenceModel & turb
autoPtr< BasicChemistryModel< ReactionThermo > > chemistryPtr_
Pointer to chemistry model.
virtual ReactionThermo & thermo()
Return access to the thermo package.
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())
@ NO_REGISTER
Do not request registration (bool: false).
@ REGISTER
Request registration (bool: true).
@ NO_READ
Nothing to be read.
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const Time & time() const noexcept
Return Time associated with the objectRegistry.
Definition IOobject.C:456
const dictionary & coeffs() const
Return const dictionary of the model.
const volScalarField & rho() const
Return const access to rho.
Switch active() const noexcept
Is combustion active?
const fvMesh & mesh() const
Return const access to the mesh database.
virtual void correct()
Correct combustion rate.
Definition EDC.C:78
virtual tmp< fvScalarMatrix > R(volScalarField &Y) const
Fuel consumption rate matrix.
Definition EDC.C:187
virtual ~EDC()
Destructor.
Definition EDC.C:71
virtual tmp< volScalarField > Qdot() const
Heat release rate [kg/m/s3].
Definition EDC.C:195
virtual bool read()
Update properties from given dictionary.
Definition EDC.C:215
tmp< volScalarField > tc() const
Return the chemical time scale.
Definition laminar.C:71
virtual tmp< fvScalarMatrix > R(volScalarField &Y) const
Fuel consumption rate matrix.
Definition laminar.C:117
virtual bool read()
Update properties from given dictionary.
Definition laminar.C:156
Abstract base class for turbulence models (RAS, LES and laminar).
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
A class for managing temporary objects.
Definition tmp.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
PtrList< volScalarField > & Y
auto limits
Definition setRDeltaT.H:186
dynamicFvMesh & mesh
scalar epsilon
word timeName
Definition getTimeIndex.H:3
const Enum< EDCversions > EDCversionNames
const scalar EDCexp2[]
Definition EDC.H:124
const scalar EDCexp1[]
Definition EDC.H:123
const EDCversions EDCdefaultVersion
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
const dimensionSet dimless
Dimensionless.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimEnergy
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
dimensionSet clamp(const dimensionSet &a, const dimensionSet &range)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensionedScalar sqrt(const dimensionedScalar &ds)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
MinMax< Type > gMinMax(const FieldField< Field, Type > &f)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
const dimensionSet dimVolume(pow3(dimLength))
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dimensionedScalar pow025(const dimensionedScalar &ds)
volScalarField & nu
tmp< volScalarField > trho
psiReactionThermo & thermo
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299