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-2017 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& combustionProperties,
38  const word& phaseName
39 )
40 :
41  Type(modelType, mesh, combustionProperties, phaseName),
42  integrateReactionRate_
43  (
44  this->coeffs().lookupOrDefault("integrateReactionRate", true)
45  )
46 {
47  if (integrateReactionRate_)
48  {
49  Info<< " using integrated reaction rate" << endl;
50  }
51  else
52  {
53  Info<< " using instantaneous reaction rate" << endl;
54  }
55 }
56 
57 
58 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
59 
60 template<class Type>
62 {}
63 
64 
65 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
66 
67 template<class Type>
70 {
71  return this->chemistryPtr_->tc();
72 }
73 
74 
75 template<class Type>
77 {
78  if (this->active())
79  {
80  if (integrateReactionRate_)
81  {
82  if (fv::localEulerDdt::enabled(this->mesh()))
83  {
84  const scalarField& rDeltaT =
85  fv::localEulerDdt::localRDeltaT(this->mesh());
86 
87  if (this->coeffs().found("maxIntegrationTime"))
88  {
89  scalar maxIntegrationTime
90  (
91  readScalar(this->coeffs().lookup("maxIntegrationTime"))
92  );
93 
94  this->chemistryPtr_->solve
95  (
96  min(1.0/rDeltaT, maxIntegrationTime)()
97  );
98  }
99  else
100  {
101  this->chemistryPtr_->solve((1.0/rDeltaT)());
102  }
103  }
104  else
105  {
106  this->chemistryPtr_->solve(this->mesh().time().deltaTValue());
107  }
108  }
109  else
110  {
111  this->chemistryPtr_->calculate();
112  }
113  }
114 }
115 
116 
117 template<class Type>
120 {
122 
123  fvScalarMatrix& Su = tSu.ref();
124 
125  if (this->active())
126  {
127  const label specieI =
128  this->thermo().composition().species()[Y.member()];
129 
130  Su += this->chemistryPtr_->RR(specieI);
131  }
132 
133  return tSu;
134 }
135 
136 
137 template<class Type>
140 {
141  tmp<volScalarField> tQdot
142  (
143  new volScalarField
144  (
145  IOobject
146  (
147  IOobject::groupName(typeName + ":Qdot", this->phaseName_),
148  this->mesh().time().timeName(),
149  this->mesh(),
150  IOobject::NO_READ,
151  IOobject::NO_WRITE,
152  false
153  ),
154  this->mesh(),
156  )
157  );
158 
159  if (this->active())
160  {
161  tQdot.ref() = this->chemistryPtr_->Qdot();
162  }
163 
164  return tQdot;
165 }
166 
167 
168 template<class Type>
170 {
171  if (Type::read())
172  {
173  integrateReactionRate_ =
174  this->coeffs().lookupOrDefault("integrateReactionRate", true);
175  return true;
176  }
177  else
178  {
179  return false;
180  }
181 }
182 
183 
184 // ************************************************************************* //
#define readScalar
Definition: doubleScalar.C:38
virtual void correct()
Correct combustion rate.
Definition: laminar.C:76
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:253
word member() const
Return member (name without the extension)
Definition: IOobject.C:364
tmp< volScalarField > tc() const
Return the chemical time scale.
Definition: laminar.C:69
void read(Istream &, label &, const dictionary &)
In-place read with dictionary lookup.
virtual tmp< volScalarField > Qdot() const
Heat release rate [kg/m/s3].
Definition: laminar.C:139
psiReactionThermo & thermo
Definition: createFields.H:31
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
Laminar combustion model.
Definition: laminar.H:49
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
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
const dimensionSet dimEnergy
virtual ~laminar()
Destructor.
Definition: laminar.C:61
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
PtrList< volScalarField > & Y
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
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:119
Calculate the matrix for implicit and explicit sources.
virtual bool read()
Update properties from given dictionary.
Definition: laminar.C:169