laminar.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 | Copyright (C) 2013-2015 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 Type>
34 (
35  const word& modelType,
36  const fvMesh& mesh,
37  const word& phaseName
38 )
39 :
40  Type(modelType, mesh, phaseName),
41  integrateReactionRate_
42  (
43  this->coeffs().lookupOrDefault("integrateReactionRate", true)
44  )
45 {
46  if (integrateReactionRate_)
47  {
48  Info<< " using integrated reaction rate" << endl;
49  }
50  else
51  {
52  Info<< " using instantaneous reaction rate" << endl;
53  }
54 }
55 
56 
57 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
58 
59 template<class Type>
61 {}
62 
63 
64 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
65 
66 template<class Type>
69 {
70  return this->chemistryPtr_->tc();
71 }
72 
73 
74 template<class Type>
76 {
77  if (this->active())
78  {
79  if (integrateReactionRate_)
80  {
81  if (fv::localEulerDdt::enabled(this->mesh()))
82  {
83  const scalarField& rDeltaT =
84  fv::localEulerDdt::localRDeltaT(this->mesh());
85 
86  if (this->coeffs().found("maxIntegrationTime"))
87  {
88  scalar maxIntegrationTime
89  (
90  readScalar(this->coeffs().lookup("maxIntegrationTime"))
91  );
92 
93  this->chemistryPtr_->solve
94  (
95  min(1.0/rDeltaT, maxIntegrationTime)()
96  );
97  }
98  else
99  {
100  this->chemistryPtr_->solve((1.0/rDeltaT)());
101  }
102  }
103  else
104  {
105  this->chemistryPtr_->solve(this->mesh().time().deltaTValue());
106  }
107  }
108  else
109  {
110  this->chemistryPtr_->calculate();
111  }
112  }
113 }
114 
115 
116 template<class Type>
119 {
121 
122  fvScalarMatrix& Su = tSu();
123 
124  if (this->active())
125  {
126  const label specieI =
127  this->thermo().composition().species()[Y.member()];
128 
129  Su += this->chemistryPtr_->RR(specieI);
130  }
131 
132  return tSu;
133 }
134 
135 
136 template<class Type>
139 {
141  (
142  new volScalarField
143  (
144  IOobject
145  (
146  IOobject::groupName(typeName + ":dQ", this->phaseName_),
147  this->mesh().time().timeName(),
148  this->mesh(),
149  IOobject::NO_READ,
150  IOobject::NO_WRITE,
151  false
152  ),
153  this->mesh(),
154  dimensionedScalar("dQ", dimEnergy/dimTime, 0.0),
155  zeroGradientFvPatchScalarField::typeName
156  )
157  );
158 
159  if (this->active())
160  {
161  tdQ() = this->chemistryPtr_->dQ();
162  }
163 
164  return tdQ;
165 }
166 
167 
168 template<class Type>
171 {
173  (
174  new volScalarField
175  (
176  IOobject
177  (
178  IOobject::groupName(typeName + ":Sh", this->phaseName_),
179  this->mesh().time().timeName(),
180  this->mesh(),
181  IOobject::NO_READ,
182  IOobject::NO_WRITE,
183  false
184  ),
185  this->mesh(),
187  zeroGradientFvPatchScalarField::typeName
188  )
189  );
190 
191  if (this->active())
192  {
193  tSh() = this->chemistryPtr_->Sh();
194  }
195 
196  return tSh;
197 }
198 
199 
200 template<class Type>
202 {
203  if (Type::read())
204  {
205  this->coeffs().lookup("integrateReactionRate")
206  >> integrateReactionRate_;
207  return true;
208  }
209  else
210  {
211  return false;
212  }
213 }
214 
215 
216 // ************************************************************************* //
Foam::combustionModels::laminar::tc
tmp< volScalarField > tc() const
Return the chemical time scale.
Definition: laminar.C:68
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::combustionModels::laminar::R
virtual tmp< fvScalarMatrix > R(volScalarField &Y) const
Fuel consumption rate matrix.
Definition: laminar.C:118
Foam::combustionModels::laminar::correct
virtual void correct()
Correct combustion rate.
Definition: laminar.C:75
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::dimEnergy
const dimensionSet dimEnergy
Foam::read
bool read(const char *, int32_t &)
Definition: int32IO.C:87
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
thermo
rhoThermo & thermo
Definition: setRegionFluidFields.H:3
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
Foam::combustionModels::laminar::~laminar
virtual ~laminar()
Destructor.
Definition: laminar.C:60
localEulerDdtScheme.H
Foam::label
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
Foam::fvScalarMatrix
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::fvc::Su
tmp< GeometricField< Type, fvPatchField, volMesh > > Su(const GeometricField< Type, fvPatchField, volMesh > &su, const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcSup.C:44
Foam::Info
messageStream Info
Foam::combustionModels::laminar::laminar
laminar(const laminar &)
Disallow copy construct.
laminar.H
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:41
Foam::combustionModels::laminar::dQ
virtual tmp< volScalarField > dQ() const
Heat release rate calculated from fuel consumption rate matrix.
Definition: laminar.C:138
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam::dimMass
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
Foam::combustionModels::laminar::Sh
virtual tmp< volScalarField > Sh() const
Return source for enthalpy equation [kg/m/s3].
Definition: laminar.C:170
fvmSup.H
Calculate the matrix for implicit and explicit sources.
readScalar
#define readScalar
Definition: doubleScalar.C:38
found
bool found
Definition: TABSMDCalcMethod2.H:32
Foam::combustionModels::laminar::read
virtual bool read()
Update properties from given dictionary.
Definition: laminar.C:201
timeName
word timeName
Definition: getTimeIndex.H:3
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:68
Foam::dimVolume
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::min
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
Y
PtrList< volScalarField > & Y
Definition: createFields.H:36
lookup
stressControl lookup("compactNormalStress") >> compactNormalStress