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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2013-2017 OpenFOAM Foundation
9  Copyright (C) 2019-2023 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "laminar.H"
30 #include "fvmSup.H"
31 #include "localEulerDdtScheme.H"
32 
33 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
34 
35 template<class ReactionThermo>
37 (
38  const word& modelType,
39  ReactionThermo& thermo,
41  const word& combustionProperties
42 )
43 :
44  ChemistryCombustion<ReactionThermo>
45  (
46  modelType,
47  thermo,
48  turb,
49  combustionProperties
50  ),
51  integrateReactionRate_
52  (
53  this->coeffs().getOrDefault("integrateReactionRate", true)
54  )
55 {
56  if (integrateReactionRate_)
57  {
58  Info<< " using integrated reaction rate" << endl;
59  }
60  else
61  {
62  Info<< " using instantaneous reaction rate" << endl;
63  }
64 }
65 
66 
67 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
68 
69 template<class ReactionThermo>
71 {}
72 
73 
74 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
75 
76 template<class ReactionThermo>
79 {
80  return this->chemistryPtr_->tc();
81 }
82 
83 
84 template<class ReactionThermo>
86 {
87  if (this->active())
88  {
89  if (integrateReactionRate_)
90  {
91  if (fv::localEulerDdt::enabled(this->mesh()))
92  {
93  const scalarField& rDeltaT =
95 
96  scalar maxTime;
97  if (this->coeffs().readIfPresent("maxIntegrationTime", maxTime))
98  {
99  this->chemistryPtr_->solve
100  (
101  min(1.0/rDeltaT, maxTime)()
102  );
103  }
104  else
105  {
106  this->chemistryPtr_->solve((1.0/rDeltaT)());
107  }
108  }
109  else
110  {
111  this->chemistryPtr_->solve(this->mesh().time().deltaTValue());
112  }
113  }
114  else
115  {
116  this->chemistryPtr_->calculate();
117  }
118  }
119 }
120 
121 
122 template<class ReactionThermo>
125 {
127  auto& Su = tSu.ref();
128 
129  if (this->active())
130  {
131  const label specieI =
132  this->thermo().composition().species().find(Y.member());
133 
134  Su += this->chemistryPtr_->RR(specieI);
135  }
137  return tSu;
138 }
139 
140 
141 template<class ReactionThermo>
144 {
145  auto tQdot = volScalarField::New
146  (
147  this->thermo().phaseScopedName(typeName, "Qdot"),
149  this->mesh(),
151  );
152 
153  if (this->active())
154  {
155  tQdot.ref() = this->chemistryPtr_->Qdot();
156  }
157 
158  return tQdot;
159 }
160 
161 
162 template<class ReactionThermo>
164 {
166  {
167  integrateReactionRate_ =
168  this->coeffs().getOrDefault("integrateReactionRate", true);
169  return true;
170  }
171 
172  return false;
173 }
174 
175 
176 // ************************************************************************* //
virtual void correct()
Correct combustion rate.
Definition: laminar.C:78
zeroField Su
Definition: alphaSuSp.H:1
compressible::turbulenceModel & turb
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
tmp< volScalarField > tc() const
Return the chemical time scale.
Definition: laminar.C:71
virtual tmp< volScalarField > Qdot() const
Heat release rate [kg/m/s3].
Definition: laminar.C:136
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:72
psiReactionThermo & thermo
Definition: createFields.H:28
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
Laminar combustion model.
Definition: laminar.H:52
dynamicFvMesh & mesh
static bool enabled(const fvMesh &mesh)
Return true if LTS is enabled.
Definition: localEulerDdt.C:32
A class for handling words, derived from Foam::string.
Definition: word.H:63
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
static const volScalarField & localRDeltaT(const fvMesh &mesh)
Return the reciprocal of the local time-step.
Definition: localEulerDdt.C:41
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
Chemistry model wrapper for combustion models.
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions and patch type. [Takes current timeN...
Abstract base class for turbulence models (RAS, LES and laminar).
const dimensionSet dimEnergy
virtual ~laminar()
Destructor.
Definition: laminar.C:63
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
PtrList< volScalarField > & Y
meshDefDict readIfPresent("polyMeshPatches", polyPatchNames)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
messageStream Info
Information stream (stdout output on master, null elsewhere)
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
A class for managing temporary objects.
Definition: HashPtrTable.H:50
virtual tmp< fvScalarMatrix > R(volScalarField &Y) const
Fuel consumption rate matrix.
Definition: laminar.C:117
Do not request registration (bool: false)
Calculate the finiteVolume matrix for implicit and explicit sources.
virtual bool read()
Update properties from given dictionary.
Definition: laminar.C:156
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127