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-2021 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"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace combustionModels
36 {
37  defineTypeNameAndDebug(laminar, 0);
38  addToRunTimeSelectionTable(combustionModel, laminar, dictionary);
39 }
40 }
41 
42 
43 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44 
46 (
47  const word& modelType,
48  const fluidReactionThermo& thermo,
50  const word& combustionProperties
51 )
52 :
54  (
55  modelType,
56  thermo,
57  turb,
58  combustionProperties
59  ),
60  outerCorrect_
61  (
62  this->coeffs().lookupOrDefault("outerCorrect", false)
63  ),
64  integrateReactionRate_
65  (
66  this->coeffs().lookupOrDefault("integrateReactionRate", true)
67  ),
68  maxIntegrationTime_
69  (
70  this->coeffs().lookupOrDefault("maxIntegrationTime", vGreat)
71  ),
72  timeIndex_(-1),
73  chemistryPtr_(basicChemistryModel::New(thermo))
74 {
75  if (integrateReactionRate_)
76  {
77  Info<< " using integrated reaction rate" << endl;
78  }
79  else
80  {
81  Info<< " using instantaneous reaction rate" << endl;
82  }
83 }
84 
85 
86 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
87 
89 {}
90 
91 
92 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
93 
95 {
96  if (!outerCorrect_ && timeIndex_ == this->mesh().time().timeIndex())
97  {
98  return;
99  }
100 
101  if (integrateReactionRate_)
102  {
103  if (fv::localEulerDdt::enabled(this->mesh()))
104  {
105  const scalarField& rDeltaT =
107 
108  chemistryPtr_->solve(min(1/rDeltaT, maxIntegrationTime_)());
109  }
110  else
111  {
112  const scalar deltaT = this->mesh().time().deltaTValue();
113 
114  chemistryPtr_->solve(min(deltaT, maxIntegrationTime_));
115  }
116  }
117  else
118  {
119  chemistryPtr_->calculate();
120  }
121 
122  timeIndex_ = this->mesh().time().timeIndex();
123 }
124 
125 
128 {
130  fvScalarMatrix& Su = tSu.ref();
131 
132  const label specieI = this->thermo().composition().species()[Y.member()];
133  Su += chemistryPtr_->RR(specieI);
134 
135  return tSu;
136 }
137 
138 
140 {
141  return chemistryPtr_->Qdot();
142 }
143 
144 
146 {
147  if (combustionModel::read())
148  {
149  outerCorrect_ =
150  this->coeffs().lookupOrDefault("outerCorrect", false);
151  integrateReactionRate_ =
152  this->coeffs().lookupOrDefault("integrateReactionRate", true);
153  maxIntegrationTime_ =
154  this->coeffs().lookupOrDefault("maxIntegrationTime", vGreat);
155  return true;
156  }
157  else
158  {
159  return false;
160  }
161 }
162 
163 
164 // ************************************************************************* //
virtual void correct()
Correct combustion rate.
Definition: laminar.C:94
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
fluidReactionThermo & thermo
Definition: createFields.H:28
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Base-class for multi-component fluid thermodynamic properties.
virtual tmp< volScalarField > Qdot() const
Heat release rate [kg/m/s^3].
Definition: laminar.C:139
Macros for easy insertion into run-time selection tables.
static word member(const word &name)
Return member (name without the extension)
Definition: IOobject.C:149
static autoPtr< basicChemistryModel > New(const fluidReactionThermo &thermo)
Select based on fluid reaction thermo.
const dimensionSet dimTime
dynamicFvMesh & mesh
static bool enabled(const fvMesh &mesh)
Return true if LTS is enabled.
Definition: localEulerDdt.C:37
A class for handling words, derived from string.
Definition: word.H:59
laminar(const word &modelType, const fluidReactionThermo &thermo, const compressibleMomentumTransportModel &turb, const word &combustionProperties)
Construct from components.
Definition: laminar.C:46
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
static const volScalarField & localRDeltaT(const fvMesh &mesh)
Return the reciprocal of the local time-step.
Definition: localEulerDdt.C:46
virtual bool read()
Update properties from given dictionary.
Base class for combustion models.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
const dimensionSet dimMass
virtual ~laminar()
Destructor.
Definition: laminar.C:88
PtrList< volScalarField > & Y
label timeIndex
Definition: getTimeIndex.H:4
messageStream Info
const tmp< volScalarField::Internal > & Su
Definition: alphaSuSp.H:6
A class for managing temporary objects.
Definition: PtrList.H:53
defineTypeNameAndDebug(diffusion, 0)
virtual tmp< fvScalarMatrix > R(volScalarField &Y) const
Fuel consumption rate matrix.
Definition: laminar.C:127
Abstract base class for turbulence models (RAS, LES and laminar).
Calculate the matrix for implicit and explicit sources.
addToRunTimeSelectionTable(combustionModel, diffusion, dictionary)
Namespace for OpenFOAM.
virtual bool read()
Update properties from given dictionary.
Definition: laminar.C:145