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-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 "zoneCombustion.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace combustionModels
34 {
35  defineTypeNameAndDebug(zoneCombustion, 0);
36  addToRunTimeSelectionTable(combustionModel, zoneCombustion, dictionary);
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 
88 Foam::combustionModels::zoneCombustion::filter
89 (
90  const tmp<volScalarField>& tS
91 ) const
92 {
93  scalarField& S = tS.ref();
94  scalarField filteredField(S.size(), 0);
95 
96  forAll(zoneNames_, zonei)
97  {
98  const labelList& cells = this->mesh().cellZones()[zoneNames_[zonei]];
99 
100  forAll(cells, i)
101  {
102  filteredField[cells[i]] = S[cells[i]];
103  }
104  }
105 
106  S = filteredField;
107 
108  return tS;
109 }
110 
111 
112 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
113 
115 (
116  const word& modelType,
117  const fluidReactionThermo& thermo,
119  const word& combustionProperties
120 )
121 :
123  (
124  modelType,
125  thermo,
126  turb,
127  combustionProperties
128  ),
129  combustionModelPtr_
130  (
132  (
133  thermo,
134  turb,
135  "zoneCombustionProperties"
136  )
137  ),
138  zoneNames_(this->coeffs().lookup("zones"))
139 {}
140 
141 
142 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
143 
145 {}
146 
147 
148 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
149 
151 {
152  combustionModelPtr_->correct();
153 }
154 
155 
158 {
159  return filter(combustionModelPtr_->R(Y));
160 }
161 
162 
165 {
166  return filter(combustionModelPtr_->Qdot());
167 }
168 
169 
171 {
172  if (combustionModel::read())
173  {
174  combustionModelPtr_->read();
175  return true;
176  }
177  else
178  {
179  return false;
180  }
181 }
182 
183 
184 // ************************************************************************* //
virtual tmp< fvScalarMatrix > R(volScalarField &Y) const
Fuel consumption rate matrix.
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
virtual bool read()
Update properties from given dictionary.
Base-class for multi-component fluid thermodynamic properties.
fvMesh & mesh
Macros for easy insertion into run-time selection tables.
virtual tmp< volScalarField > Qdot() const
Heat release rate [kg/m/s^3].
stressControl lookup("compactNormalStress") >> compactNormalStress
static autoPtr< combustionModel > New(const fluidReactionThermo &thermo, const compressibleMomentumTransportModel &turb, const word &combustionProperties=combustionPropertiesName)
Select using thermo and turbulence.
A class for handling words, derived from string.
Definition: word.H:59
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
List< label > labelList
A List of labels.
Definition: labelList.H:56
virtual bool read()
Update properties from given dictionary.
Base class for combustion models.
PtrList< volScalarField > & Y
A class for managing temporary objects.
Definition: PtrList.H:53
zoneCombustion(const word &modelType, const fluidReactionThermo &thermo, const compressibleMomentumTransportModel &turb, const word &combustionProperties)
Construct from components.
defineTypeNameAndDebug(diffusion, 0)
virtual void correct()
Correct combustion rate.
Base class for single-phase compressible turbulence models.
addToRunTimeSelectionTable(combustionModel, diffusion, dictionary)
Namespace for OpenFOAM.