singleStepCombustion.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) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2019 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 "singleStepCombustion.H"
30 #include "fvmSup.H"
31 
32 namespace Foam
33 {
34 namespace combustionModels
35 {
36 
37 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
38 
39 template<class ReactionThermo, class ThermoType>
41 (
42  const word& modelType,
43  ReactionThermo& thermo,
45  const word& combustionProperties
46 )
47 :
49  singleMixturePtr_(nullptr),
50  wFuel_
51  (
52  IOobject
53  (
54  this->thermo().phasePropertyName("wFuel"),
55  this->mesh().time().timeName(),
56  this->mesh(),
59  ),
60  this->mesh(),
62  ),
63  semiImplicit_(this->coeffs_.getBool("semiImplicit"))
64 {
66  {
67  singleMixturePtr_ =
69  (
70  this->thermo()
71  );
72  }
73  else
74  {
76  << "Inconsistent thermo package for " << this->type() << " model:\n"
77  << " " << this->thermo().type() << nl << nl
78  << "Please select a thermo package based on "
79  << "singleStepReactingMixture" << exit(FatalError);
80  }
81 
82  if (semiImplicit_)
83  {
84  Info<< "Combustion mode: semi-implicit" << endl;
85  }
86  else
87  {
88  Info<< "Combustion mode: explicit" << endl;
89  }
90 }
91 
92 
93 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
94 
95 template<class ReactionThermo, class ThermoType>
97 {}
98 
99 
100 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
101 
102 template<class ReactionThermo, class ThermoType>
104 (
106 ) const
107 {
108  const label specieI =
109  this->thermo().composition().species().find(Y.member());
110 
111  volScalarField wSpecie
112  (
113  wFuel_*singleMixturePtr_->specieStoichCoeffs()[specieI]
114  );
115 
116  if (semiImplicit_)
117  {
118  const label fNorm = singleMixturePtr_->specieProd()[specieI];
119  const volScalarField fres(singleMixturePtr_->fres(specieI));
120  wSpecie /= max(fNorm*(Y - fres), scalar(1e-2));
121 
122  return -fNorm*wSpecie*fres + scalar(fNorm)*fvm::Sp(wSpecie, Y);
123  }
125  return wSpecie + fvm::Sp(0.0*wSpecie, Y);
126 }
127 
128 
129 template<class ReactionThermo, class ThermoType>
132 {
133  const label fuelI = singleMixturePtr_->fuelIndex();
134  volScalarField& YFuel =
135  const_cast<volScalarField&>(this->thermo().composition().Y(fuelI));
136 
137  return -singleMixturePtr_->qFuel()*(R(YFuel) & YFuel);
138 }
139 
140 
141 template<class ReactionThermo, class ThermoType>
143 {
145  {
146  return true;
147  }
148 
149  return false;
150 }
151 
152 
153 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
154 
155 } // End namespace combustionModels
156 } // End namespace Foam
157 
158 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
Single step reacting mixture.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
compressible::turbulenceModel & turb
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
Ignore writing from objectRegistry::writeObject()
static word member(const word &name)
Return member (name without the extension)
Definition: IOobject.C:313
psiReactionThermo & thermo
Definition: createFields.H:28
word timeName
Definition: getTimeIndex.H:3
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
dynamicFvMesh & mesh
A class for handling words, derived from Foam::string.
Definition: word.H:63
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
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).
Thermo model wrapper for combustion models.
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
Nothing to be read.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
messageStream Info
Information stream (stdout output on master, null elsewhere)
Base class for combustion models using singleStepReactingMixture.
const Type * isA(const U &obj)
Check if dynamic_cast to Type is possible.
Definition: typeInfo.H:88
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
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
Calculate the finiteVolume matrix for implicit and explicit sources.
virtual tmp< fvScalarMatrix > R(volScalarField &Y) const
Fuel consumption rate matrix.
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127