laminar.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) 2013-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 "laminar.H"
27 #include "fvmSup.H"
28 #include "localEulerDdtScheme.H"
29 
30 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
31 
32 template<class ReactionThermo>
34 (
35  const word& modelType,
36  ReactionThermo& thermo,
37  const compressibleTurbulenceModel& turb,
38  const word& combustionProperties
39 )
40 :
42  (
43  modelType,
44  thermo,
45  turb,
46  combustionProperties
47  ),
48  integrateReactionRate_
49  (
50  this->coeffs().lookupOrDefault("integrateReactionRate", true)
51  )
52 {
53  if (integrateReactionRate_)
54  {
55  Info<< " using integrated reaction rate" << endl;
56  }
57  else
58  {
59  Info<< " using instantaneous reaction rate" << endl;
60  }
61 }
62 
63 
64 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
65 
66 template<class ReactionThermo>
68 {}
69 
70 
71 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
72 
73 template<class ReactionThermo>
76 {
77  return this->chemistryPtr_->tc();
78 }
79 
80 
81 template<class ReactionThermo>
83 {
84  if (this->active())
85  {
86  if (integrateReactionRate_)
87  {
88  if (fv::localEulerDdt::enabled(this->mesh()))
89  {
90  const scalarField& rDeltaT =
91  fv::localEulerDdt::localRDeltaT(this->mesh());
92 
93  if (this->coeffs().found("maxIntegrationTime"))
94  {
95  scalar maxIntegrationTime
96  (
97  readScalar(this->coeffs().lookup("maxIntegrationTime"))
98  );
99 
100  this->chemistryPtr_->solve
101  (
102  min(1.0/rDeltaT, maxIntegrationTime)()
103  );
104  }
105  else
106  {
107  this->chemistryPtr_->solve((1.0/rDeltaT)());
108  }
109  }
110  else
111  {
112  this->chemistryPtr_->solve(this->mesh().time().deltaTValue());
113  }
114  }
115  else
116  {
117  this->chemistryPtr_->calculate();
118  }
119  }
120 }
121 
122 
123 template<class ReactionThermo>
126 {
128 
129  fvScalarMatrix& Su = tSu.ref();
130 
131  if (this->active())
132  {
133  const label specieI =
134  this->thermo().composition().species()[Y.member()];
135 
136  Su += this->chemistryPtr_->RR(specieI);
137  }
138 
139  return tSu;
140 }
141 
142 
143 template<class ReactionThermo>
146 {
147  tmp<volScalarField> tQdot
148  (
149  new volScalarField
150  (
151  IOobject
152  (
153  this->thermo().phasePropertyName(typeName + ":Qdot"),
154  this->mesh().time().timeName(),
155  this->mesh(),
156  IOobject::NO_READ,
157  IOobject::NO_WRITE,
158  false
159  ),
160  this->mesh(),
162  )
163  );
164 
165  if (this->active())
166  {
167  tQdot.ref() = this->chemistryPtr_->Qdot();
168  }
169 
170  return tQdot;
171 }
172 
173 
174 template<class ReactionThermo>
176 {
178  {
179  integrateReactionRate_ =
180  this->coeffs().lookupOrDefault("integrateReactionRate", true);
181  return true;
182  }
183  else
184  {
185  return false;
186  }
187 }
188 
189 
190 // ************************************************************************* //
#define readScalar
Definition: doubleScalar.C:38
virtual void correct()
Correct combustion rate.
Definition: laminar.C:82
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
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
zeroField Su
Definition: alphaSuSp.H:1
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
tmp< volScalarField > tc() const
Return the chemical time scale.
Definition: laminar.C:75
rhoReactionThermo & thermo
Definition: createFields.H:28
virtual tmp< volScalarField > Qdot() const
Heat release rate [kg/m/s3].
Definition: laminar.C:145
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
Laminar combustion model.
Definition: laminar.H:51
stressControl lookup("compactNormalStress") >> compactNormalStress
dynamicFvMesh & mesh
A class for handling words, derived from string.
Definition: word.H:59
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:72
word timeName
Definition: getTimeIndex.H:3
Chemistry model wrapper for combustion models.
Abstract base class for turbulence models (RAS, LES and laminar).
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
const dimensionSet dimEnergy
virtual ~laminar()
Destructor.
Definition: laminar.C:67
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 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
bool found
virtual tmp< fvScalarMatrix > R(volScalarField &Y) const
Fuel consumption rate matrix.
Definition: laminar.C:125
Calculate the matrix for implicit and explicit sources.
virtual bool read()
Update properties from given dictionary.
Definition: laminar.C:175