All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
zoneCombustion.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) 2016-2023 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 "zoneCombustion.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace combustionModels
34 {
37 }
38 }
39 
40 
41 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
42 
44 Foam::combustionModels::zoneCombustion::filter
45 (
46  const tmp<fvScalarMatrix>& tR
47 ) const
48 {
49  fvScalarMatrix& R = tR.ref();
50  scalarField& Su = R.source();
51  scalarField filteredField(Su.size(), 0);
52 
53  forAll(zoneNames_, zonei)
54  {
55  const labelList& cells = this->mesh().cellZones()[zoneNames_[zonei]];
56 
57  forAll(cells, i)
58  {
59  filteredField[cells[i]] = Su[cells[i]];
60  }
61  }
62 
63  Su = filteredField;
64 
65  if (R.hasDiag())
66  {
67  scalarField& Sp = R.diag();
68 
69  forAll(zoneNames_, zonei)
70  {
71  const labelList& cells =
72  this->mesh().cellZones()[zoneNames_[zonei]];
73 
74  forAll(cells, i)
75  {
76  filteredField[cells[i]] = Sp[cells[i]];
77  }
78  }
79 
80  Sp = filteredField;
81  }
82 
83  return tR;
84 }
85 
86 
87 template<class GeoField>
89 Foam::combustionModels::zoneCombustion::filter
90 (
91  const tmp<GeoField>& tS
92 ) const
93 {
94  scalarField& S = tS.ref();
95  scalarField filteredField(S.size(), 0);
96 
97  forAll(zoneNames_, zonei)
98  {
99  const labelList& cells = this->mesh().cellZones()[zoneNames_[zonei]];
100 
101  forAll(cells, i)
102  {
103  filteredField[cells[i]] = S[cells[i]];
104  }
105  }
106 
107  S = filteredField;
108 
109  return tS;
110 }
111 
112 
113 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
114 
116 (
117  const word& modelType,
120  const word& combustionProperties
121 )
122 :
124  (
125  modelType,
126  thermo,
127  turb,
128  combustionProperties
129  ),
130  combustionModelPtr_
131  (
133  (
134  thermo,
135  turb,
136  "zoneCombustionProperties"
137  )
138  ),
139  zoneNames_(this->coeffs().lookup("zones"))
140 {}
141 
142 
143 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
144 
146 {}
147 
148 
149 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
150 
152 {
153  combustionModelPtr_->correct();
154 }
155 
156 
159 {
160  return filter(combustionModelPtr_->R(speciei));
161 }
162 
163 
166 {
167  return filter(combustionModelPtr_->R(Y));
168 }
169 
170 
173 {
174  return filter(combustionModelPtr_->Qdot());
175 }
176 
177 
179 {
180  if (combustionModel::read())
181  {
182  combustionModelPtr_->read();
183  return true;
184  }
185  else
186  {
187  return false;
188  }
189 }
190 
191 
192 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
Base class for combustion models.
const fvMesh & mesh() const
Return const access to the mesh database.
virtual bool read()
Update properties from given dictionary.
Zone-filtered combustion model.
virtual void correct()
Correct combustion rate.
zoneCombustion(const word &modelType, const fluidMulticomponentThermo &thermo, const compressibleMomentumTransportModel &turb, const word &combustionProperties)
Construct from components.
virtual tmp< volScalarField::Internal > R(const label speciei) const
Specie consumption rate field.
virtual tmp< volScalarField > Qdot() const
Heat release rate [kg/m/s^3].
virtual bool read()
Update properties from given dictionary.
Base class for single-phase compressible turbulence models.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Base-class for multi-component fluid thermodynamic properties.
const cellZoneList & cellZones() const
Return cell zones.
Definition: polyMesh.H:452
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
A class for handling words, derived from string.
Definition: word.H:62
const cellShapeList & cells
defineTypeNameAndDebug(diffusion, 0)
addToRunTimeSelectionTable(combustionModel, diffusion, dictionary)
tmp< VolField< Type > > Su(const VolField< Type > &su, const VolField< Type > &vf)
Definition: fvcSup.C:44
tmp< VolField< Type > > Sp(const volScalarField &sp, const VolField< Type > &vf)
Definition: fvcSup.C:67
tmp< fvMatrix< Type > > S(const Pair< tmp< volScalarField::Internal >> &, const VolField< Type > &)
Namespace for OpenFOAM.
List< label > labelList
A List of labels.
Definition: labelList.H:56
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
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
PtrList< volScalarField > & Y
fluidMulticomponentThermo & thermo
Definition: createFields.H:31