eddyDissipationModelBase.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) 2016-2019 OpenCFD Ltd
9-------------------------------------------------------------------------------
10License
11
12 OpenFOAM is free software: you can redistribute it and/or modify it
13 under the terms of the GNU General Public License as published by
14 the Free Software Foundation, either version 3 of the License, or
15 (at your option) any later version.
16
17 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
18 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
19 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
20 for more details.
21
22 You should have received a copy of the GNU General Public License
23 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
24
25\*---------------------------------------------------------------------------*/
26
29
30namespace Foam
31{
32namespace combustionModels
33{
34
35// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
36
37template<class ReactionThermo, class ThermoType>
38eddyDissipationModelBase<ReactionThermo, ThermoType>::eddyDissipationModelBase
40 const word& modelType,
41 ReactionThermo& thermo,
43 const word& combustionProperties
44)
45:
47 (
48 modelType,
49 thermo,
50 turb,
51 combustionProperties
52 ),
53 CEDC_(this->coeffs().getScalar("CEDC"))
54{}
55
56
57// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
58
59template<class ReactionThermo, class ThermoType>
62{}
63
64
65// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
66
67template<class ReactionThermo, class ThermoType>
71 return
72 CEDC_*this->turbulence().epsilon()
73 / max
74 (
75 this->turbulence().k(),
76 dimensionedScalar("SMALL",sqr(dimVelocity), SMALL)
77 );
78}
79
80
81template<class ReactionThermo, class ThermoType>
85
86 if (this->active())
87 {
88 this->singleMixturePtr_->fresCorrect();
89
90 const label fuelI = this->singleMixturePtr_->fuelIndex();
91
92 const volScalarField& YFuel = this->thermo_.composition().Y()[fuelI];
93
94 const dimensionedScalar s = this->singleMixturePtr_->s();
95
96 if (this->thermo_.composition().contains("O2"))
97 {
98 const volScalarField& YO2 = this->thermo_.composition().Y("O2");
99
100 this->wFuel_ ==
101 this->rho()
102 * min(YFuel, YO2/s.value())
103 * timeScale();
104 }
105 else
106 {
108 << "You selected a combustion model that requires O2 mass"
109 << " to be present in the mixture"
110 << exit(FatalError);
111 }
112 }
113}
114
115
116template<class ReactionThermo, class ThermoType>
120 {
121 this->coeffs().readEntry("CEDC", CEDC_);
122 return true;
123 }
124
125 return false;
126}
127
128
129// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
130
131} // End namespace combustionModels
132} // End namespace Foam
133
134// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
label k
compressible::turbulenceModel & turb
Standard Eddy Dissipation Model based on the assumption that the reaction rates are controlled by the...
tmp< volScalarField > rtTurb() const
Return the reciprocal of the turbulent mixing time scale.
Base class for combustion models using singleStepReactingMixture.
Abstract base class for turbulence models (RAS, LES and laminar).
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
compressible::turbulenceModel & turbulence
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Namespace for OpenFOAM.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
const dimensionSet dimVelocity
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
error FatalError
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:60
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:51