zoneCombustion.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) 2016-2017 OpenFOAM Foundation
9  Copyright (C) 2019-2020 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 "zoneCombustion.H"
30 
31 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
32 
33 template<class ReactionThermo>
36 (
37  const tmp<fvScalarMatrix>& tR
38 ) const
39 {
40  fvScalarMatrix& R = tR.ref();
41  scalarField& Su = R.source();
42  scalarField filteredField(Su.size(), Zero);
43 
44  forAll(zoneNames_, zonei)
45  {
46  const labelList& cells = this->mesh().cellZones()[zoneNames_[zonei]];
47 
48  forAll(cells, i)
49  {
50  filteredField[cells[i]] = Su[cells[i]];
51  }
52  }
53 
54  Su = filteredField;
55 
56  if (R.hasDiag())
57  {
58  scalarField& Sp = R.diag();
59 
60  forAll(zoneNames_, zonei)
61  {
62  const labelList& cells =
63  this->mesh().cellZones()[zoneNames_[zonei]];
64 
65  forAll(cells, i)
66  {
67  filteredField[cells[i]] = Sp[cells[i]];
68  }
69  }
70 
71  Sp = filteredField;
72  }
73 
74  return tR;
75 }
76 
77 
78 template<class ReactionThermo>
81 (
82  const tmp<volScalarField>& tS
83 ) const
84 {
85  scalarField& S = tS.ref();
86  scalarField filteredField(S.size(), Zero);
87 
88  forAll(zoneNames_, zonei)
89  {
90  const labelList& cells = this->mesh().cellZones()[zoneNames_[zonei]];
91 
92  forAll(cells, i)
93  {
94  filteredField[cells[i]] = S[cells[i]];
95  }
96  }
97 
98  S = filteredField;
99 
100  return tS;
101 }
102 
103 
104 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
105 
106 template<class ReactionThermo>
108 (
109  const word& modelType,
110  ReactionThermo& thermo,
112  const word& combustionProperties
113 )
114 :
116  (
117  modelType,
118  thermo,
119  turb,
120  combustionProperties
121  ),
122  combustionModelPtr_
123  (
125  (
126  thermo,
127  turb,
128  "zoneCombustionProperties"
129  )
130  ),
131  zoneNames_()
132 {
133  this->coeffs().readEntry("zones", zoneNames_);
134 }
135 
136 
137 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
138 
139 template<class ReactionThermo>
141 {}
142 
143 
144 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
145 
146 template<class ReactionThermo>
148 {
149  return combustionModelPtr_->thermo();
150 }
151 
152 
153 template<class ReactionThermo>
154 const ReactionThermo&
156 {
157  return combustionModelPtr_->thermo();
158 }
159 
160 
161 template<class ReactionThermo>
163 {
164  combustionModelPtr_->correct();
165 }
166 
167 
168 template<class ReactionThermo>
171 (
173 ) const
174 {
175  return filter(combustionModelPtr_->R(Y));
176 }
177 
178 
179 template<class ReactionThermo>
182 {
183  return filter(combustionModelPtr_->Qdot());
184 }
185 
186 
187 template<class ReactionThermo>
189 {
191  {
192  combustionModelPtr_->read();
193  return true;
194  }
195 
196  return false;
197 }
198 
199 
200 // ************************************************************************* //
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:67
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::CombustionModel
Combustion models for templated thermodynamics.
Definition: CombustionModel.H:55
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::combustionModels::zoneCombustion::read
virtual bool read()
Update properties from given dictionary.
Definition: zoneCombustion.C:188
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Sp
zeroField Sp
Definition: alphaSuSp.H:2
zoneCombustion.H
thermo
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
Foam::combustionModels::zoneCombustion::correct
virtual void correct()
Correct combustion rate.
Definition: zoneCombustion.C:162
Foam::combustionModels::zoneCombustion::Qdot
virtual tmp< volScalarField > Qdot() const
Heat release rate [kg/m/s3].
Definition: zoneCombustion.C:181
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::combustionModels::zoneCombustion::R
virtual tmp< fvScalarMatrix > R(volScalarField &Y) const
Fuel consumption rate matrix.
Definition: zoneCombustion.C:171
Su
zeroField Su
Definition: alphaSuSp.H:1
R
#define R(A, B, C, D, E, F, K, M)
Foam::fvScalarMatrix
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:44
Foam::dictionary::readEntry
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, bool mandatory=true) const
Definition: dictionaryTemplates.C:302
Foam::combustionModels::zoneCombustion
Zone-filtered combustion model.
Definition: zoneCombustion.H:56
Foam::combustionModels::zoneCombustion::~zoneCombustion
virtual ~zoneCombustion()
Destructor.
Definition: zoneCombustion.C:140
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Y
PtrList< volScalarField > & Y
Definition: createFieldRefs.H:7
Foam::combustionModels::zoneCombustion::thermo
virtual ReactionThermo & thermo()
Return access to the thermo package.
Definition: zoneCombustion.C:147
cells
const cellShapeList & cells
Definition: gmvOutputHeader.H:3
Foam::GeometricField< scalar, fvPatchField, volMesh >
turb
compressible::turbulenceModel & turb
Definition: setRegionFluidFields.H:10
Foam::compressibleTurbulenceModel
Abstract base class for turbulence models (RAS, LES and laminar).
Definition: compressibleTurbulenceModel.H:54