singleStepCombustion.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) 2011-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 "singleStepCombustion.H"
27 #include "fvmSup.H"
28 
29 namespace Foam
30 {
31 namespace combustionModels
32 {
33 
34 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35 
36 template<class ReactionThermo, class ThermoType>
37 singleStepCombustion<ReactionThermo, ThermoType>::singleStepCombustion
38 (
39  const word& modelType,
40  ReactionThermo& thermo,
41  const compressibleTurbulenceModel& turb,
42  const word& combustionProperties
43 )
44 :
45  ThermoCombustion<ReactionThermo>(modelType, thermo, turb),
46  singleMixturePtr_(nullptr),
47  wFuel_
48  (
49  IOobject
50  (
51  this->thermo().phasePropertyName("wFuel"),
52  this->mesh().time().timeName(),
53  this->mesh(),
56  ),
57  this->mesh(),
59  ),
60  semiImplicit_(readBool(this->coeffs_.lookup("semiImplicit")))
61 {
63  {
64  singleMixturePtr_ =
66  (
67  this->thermo()
68  );
69  }
70  else
71  {
73  << "Inconsistent thermo package for " << this->type() << " model:\n"
74  << " " << this->thermo().type() << nl << nl
75  << "Please select a thermo package based on "
76  << "singleStepReactingMixture" << exit(FatalError);
77  }
78 
79  if (semiImplicit_)
80  {
81  Info<< "Combustion mode: semi-implicit" << endl;
82  }
83  else
84  {
85  Info<< "Combustion mode: explicit" << endl;
86  }
87 }
88 
89 
90 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
91 
92 template<class ReactionThermo, class ThermoType>
94 {}
95 
96 
97 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
98 
99 template<class ReactionThermo, class ThermoType>
101 (
103 ) const
104 {
105  const label specieI =
106  this->thermo().composition().species()[Y.member()];
107 
108  volScalarField wSpecie
109  (
110  wFuel_*singleMixturePtr_->specieStoichCoeffs()[specieI]
111  );
112 
113  if (semiImplicit_)
114  {
115  const label fNorm = singleMixturePtr_->specieProd()[specieI];
116  const volScalarField fres(singleMixturePtr_->fres(specieI));
117  wSpecie /= max(fNorm*(Y - fres), scalar(1e-2));
118 
119  return -fNorm*wSpecie*fres + fNorm*fvm::Sp(wSpecie, Y);
120  }
121  else
122  {
123  return wSpecie + fvm::Sp(0.0*wSpecie, Y);
124  }
125 }
126 
127 
128 template<class ReactionThermo, class ThermoType>
131 {
132  const label fuelI = singleMixturePtr_->fuelIndex();
133  volScalarField& YFuel =
134  const_cast<volScalarField&>(this->thermo().composition().Y(fuelI));
135 
136  return -singleMixturePtr_->qFuel()*(R(YFuel) & YFuel);
137 }
138 
139 
140 template<class ReactionThermo, class ThermoType>
142 {
144  {
145  return true;
146  }
147  else
148  {
149  return false;
150  }
151 }
152 
153 
154 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
155 
156 } // End namespace combustionModels
157 } // End namespace Foam
158 
159 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
bool isA(const Type &t)
Check if a dynamic_cast to typeid is possible.
Definition: typeInfo.H:134
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
Single step reacting mixture.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
bool readBool(Istream &)
Definition: boolIO.C:60
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
rhoReactionThermo & thermo
Definition: createFields.H:28
static word member(const word &name)
Return member (name without the extension)
Definition: IOobject.C:191
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
dynamicFvMesh & mesh
A class for handling words, derived from string.
Definition: word.H:59
word timeName
Definition: getTimeIndex.H:3
static const char nl
Definition: Ostream.H:265
Abstract base class for turbulence models (RAS, LES and laminar).
Thermo model wrapper for combustion models.
fileName::Type type(const fileName &, const bool followLink=true)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:481
virtual tmp< volScalarField > Qdot() const
Heat release rate [kg/m/s3].
virtual bool read()
Update properties from given dictionary.
#define R(A, B, C, D, E, F, K, M)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
PtrList< volScalarField > & Y
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
messageStream Info
const doubleScalar e
Elementary charge.
Definition: doubleScalar.H:98
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
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
Calculate the matrix for implicit and explicit sources.
virtual tmp< fvScalarMatrix > R(volScalarField &Y) const
Fuel consumption rate matrix.
Namespace for OpenFOAM.