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-2020 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  const ReactionThermo& thermo,
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 (integrateReactionRate_)
85  {
86  if (fv::localEulerDdt::enabled(this->mesh()))
87  {
88  const scalarField& rDeltaT =
89  fv::localEulerDdt::localRDeltaT(this->mesh());
90 
91  if (this->coeffs().found("maxIntegrationTime"))
92  {
93  const scalar maxIntegrationTime
94  (
95  this->coeffs().template lookup<scalar>("maxIntegrationTime")
96  );
97 
98  this->chemistryPtr_->solve
99  (
100  min(1.0/rDeltaT, maxIntegrationTime)()
101  );
102  }
103  else
104  {
105  this->chemistryPtr_->solve((1.0/rDeltaT)());
106  }
107  }
108  else
109  {
110  this->chemistryPtr_->solve(this->mesh().time().deltaTValue());
111  }
112  }
113  else
114  {
115  this->chemistryPtr_->calculate();
116  }
117 }
118 
119 
120 template<class ReactionThermo>
123 {
125  fvScalarMatrix& Su = tSu.ref();
126 
127  const label specieI = this->thermo().composition().species()[Y.member()];
128  Su += this->chemistryPtr_->RR(specieI);
129 
130  return tSu;
131 }
132 
133 
134 template<class ReactionThermo>
137 {
138  return this->chemistryPtr_->Qdot();
139 }
140 
141 
142 template<class ReactionThermo>
144 {
146  {
147  integrateReactionRate_ =
148  this->coeffs().lookupOrDefault("integrateReactionRate", true);
149  return true;
150  }
151  else
152  {
153  return false;
154  }
155 }
156 
157 
158 // ************************************************************************* //
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:251
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/s^3].
Definition: laminar.C:136
static word member(const word &name)
Return member (name without the extension)
Definition: IOobject.C:159
dynamicFvMesh & mesh
laminar(const word &modelType, const ReactionThermo &thermo, const compressibleMomentumTransportModel &turb, const word &combustionProperties)
Construct from components.
Definition: laminar.C:34
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
Chemistry model wrapper for combustion models.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
virtual ~laminar()
Destructor.
Definition: laminar.C:67
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
bool found
virtual tmp< fvScalarMatrix > R(volScalarField &Y) const
Fuel consumption rate matrix.
Definition: laminar.C:122
Abstract base class for turbulence models (RAS, LES and laminar).
Calculate the matrix for implicit and explicit sources.
virtual bool read()
Update properties from given dictionary.
Definition: laminar.C:143