multicomponentThermo.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) 2023-2025 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 "multicomponentThermo.H"
27 
28 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
29 
30 namespace Foam
31 {
33 }
34 
35 
36 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
37 
39 {
40  if (species_.size())
41  {
43  (
45  (
46  IOobject::groupName("Yt", Y_[0].group()),
47  Y_[0],
48  calculatedFvPatchScalarField::typeName
49  )
50  );
51  volScalarField& Yt = tYt.ref();
52 
53  for (label i=1; i<Y_.size(); i++)
54  {
55  Yt += Y_[i];
56  }
57 
58  if (mag(min(Yt).value()) < rootVSmall)
59  {
61  << "Sum of mass fractions is zero for species " << species()
62  << exit(FatalError);
63  }
64 
65  forAll(Y_, i)
66  {
67  Y_[i] /= Yt;
68  }
69  }
70 }
71 
72 
73 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
74 
76 (
77  const dictionary& dict,
78  const wordList& specieNames,
79  const fvMesh& mesh,
80  const word& phaseName
81 )
82 :
83  species_(specieNames),
84  defaultSpecieName_
85  (
86  species_.size()
87  ? dict.lookupBackwardsCompatible<word>
88  (
89  {"defaultSpecie", "inertSpecie"}
90  )
91  : word("undefined")
92  ),
93  defaultSpeciei_
94  (
95  species_.size()
96  && species_.found(defaultSpecieName_)
97  ? species_[defaultSpecieName_]
98  : -1
99  ),
100  active_(species_.size(), true),
101  Y_(species_.size())
102 {
103  if (species_.size() && defaultSpeciei_ == -1)
104  {
106  << "default specie " << defaultSpecieName_
107  << " not found in available species " << species_
108  << exit(FatalIOError);
109  }
110 
111  // Read the species' mass fractions
112  forAll(species_, i)
113  {
115  (
116  IOobject::groupName(species_[i], phaseName),
117  mesh.time().name(),
118  mesh,
120  );
121 
122  if (header.headerOk())
123  {
124  // Read the mass fraction field
125  Y_.set
126  (
127  i,
128  new volScalarField
129  (
130  IOobject
131  (
132  IOobject::groupName(species_[i], phaseName),
133  mesh.time().name(),
134  mesh,
137  ),
138  mesh
139  )
140  );
141  }
142  else
143  {
144  // Read Ydefault if not already read
145  if (!Ydefault_.valid())
146  {
147  Ydefault_ = new volScalarField
148  (
149  IOobject
150  (
151  IOobject::groupName("Ydefault", phaseName),
152  mesh.time().name(),
153  mesh,
156  ),
157  mesh
158  );
159  }
160 
161  Y_.set
162  (
163  i,
164  new volScalarField
165  (
166  IOobject
167  (
168  IOobject::groupName(species_[i], phaseName),
169  mesh.time().name(),
170  mesh,
173  ),
174  Ydefault_()
175  )
176  );
177  }
178  }
179 
180  correctMassFractions();
181 }
182 
183 
184 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
185 
187 {}
188 
189 
191 {}
192 
193 
194 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
195 
196 const Foam::speciesTable&
198 {
199  return species_;
200 }
201 
202 
204 {
205  return defaultSpeciei_;
206 }
207 
208 
209 const Foam::List<bool>&
211 {
212  return active_;
213 }
214 
215 
217 {
218  if (Pstream::parRun())
219  {
221  const_cast<List<bool>&>(this->speciesActive());
222 
223  Pstream::listCombineGather(speciesActive, orEqOp<bool>());
225 
227  const_cast<PtrList<volScalarField>&>(this->Y());
228 
229  forAll(Y, speciei)
230  {
231  if (speciesActive[speciei])
232  {
233  Y[speciei].writeOpt() = IOobject::AUTO_WRITE;
234  }
235  }
236  }
237 
238  if (Ydefault_.valid())
239  {
240  Ydefault_->writeOpt() = IOobject::AUTO_WRITE;
241  }
242 }
243 
244 
247 {
248  return Y_;
249 }
250 
251 
254 {
255  return Y_;
256 }
257 
258 
260 {
261  if (species().size())
262  {
264  (
266  (
268  Y()[0].mesh(),
270  )
271  );
272  volScalarField& Yt = tYt.ref();
273 
274  forAll(Y(), i)
275  {
276  if (solveSpecie(i))
277  {
278  Y()[i].max(scalar(0));
279  Yt += Y()[i];
280  }
281  }
282 
283  Y()[defaultSpecie()] = scalar(1) - Yt;
284  Y()[defaultSpecie()].max(scalar(0));
285  }
286 
287  if (Ydefault_.valid() && Ydefault_->writeOpt() == IOobject::NO_WRITE)
288  {
289  Ydefault_.clear();
290  }
291 }
292 
293 
294 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
Generic GeometricField class.
static tmp< GeometricField< Type, GeoMesh, PrimitiveField > > New(const word &name, const Internal &, const PtrList< Patch > &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
static word groupName(Name name, const word &group)
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
static void listCombineGather(const List< commsStruct > &comms, List< T > &Value, const CombineOp &cop, const int tag, const label comm)
static void listCombineScatter(const List< commsStruct > &comms, List< T > &Value, const int tag, const label comm)
Scatter data. Reverse of combineGather.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:399
virtual const fvMesh & mesh() const =0
Return const access to the mesh.
virtual const word & phaseName() const =0
Phase name.
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
const word & name() const
Return const reference to name.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:96
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:420
A wordList with hashed indices for faster lookup by name.
void correctMassFractions()
Scale the mass fractions to sum to 1.
virtual const speciesTable & species() const
The table of species.
implementation(const dictionary &, const wordList &, const fvMesh &, const word &)
Construct from dictionary, specie names, mesh and phase name.
virtual label defaultSpecie() const
The index of the default specie.
virtual void normaliseY()
Normalise the mass fractions by clipping positive and deriving.
virtual const List< bool > & speciesActive() const
Access the specie active flags.
speciesTable species_
Table of specie names.
virtual void syncSpeciesActive() const
Synchronise the specie active flags.
virtual PtrList< volScalarField > & Y()
Access the mass-fraction fields.
PtrList< volScalarField > Y_
Species mass fractions.
Base-class for multi-component thermodynamic properties.
virtual const speciesTable & species() const =0
The table of species.
bool solveSpecie(const label speciei) const
Should the given specie be solved for? I.e., is it active and.
virtual PtrList< volScalarField > & Y()=0
Access the mass-fraction fields.
virtual label defaultSpecie() const =0
The index of the default specie.
virtual ~multicomponentThermo()
Destructor.
virtual const List< bool > & speciesActive() const =0
Access the specie active flags.
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:197
Templated form of IOobject providing type information for file reading and header type checking.
Definition: IOobject.H:530
A class for handling words, derived from string.
Definition: word.H:62
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:346
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
const char *const group
Group name for atomic constants.
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
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
const dimensionSet dimless
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:62
defineTypeNameAndDebug(combustionModel, 0)
IOerror FatalIOError
error FatalError
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dictionary dict