basicSpecieMixture.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) 2014-2022 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 "basicSpecieMixture.H"
27 
28 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
29 
30 namespace Foam
31 {
32  defineTypeNameAndDebug(basicSpecieMixture, 0);
33 }
34 
35 
36 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
37 
39 (
40  const dictionary& thermoDict,
41  const wordList& specieNames,
42  const fvMesh& mesh,
43  const word& phaseName
44 )
45 :
46  basicMixture(thermoDict, mesh, phaseName),
47  phaseName_(phaseName),
48  species_(specieNames),
49  defaultSpecie_
50  (
51  species_.size()
52  ? thermoDict.lookupBackwardsCompatible<word>
53  (
54  {"defaultSpecie", "inertSpecie"}
55  )
56  : word("undefined")
57  ),
58  defaultSpecieIndex_(-1),
59  active_(species_.size(), true),
60  Y_(species_.size())
61 {
62  if (species_.size())
63  {
64  if (species_.found(defaultSpecie_))
65  {
66  defaultSpecieIndex_ = species_[defaultSpecie_];
67  }
68  else
69  {
70  FatalIOErrorInFunction(thermoDict)
71  << "default specie " << defaultSpecie_
72  << " not found in available species " << species_
73  << exit(FatalIOError);
74  }
75  }
76 
77  tmp<volScalarField> tYdefault;
78 
79  forAll(species_, i)
80  {
82  (
83  IOobject::groupName(species_[i], phaseName),
84  mesh.time().timeName(),
85  mesh,
87  );
88 
89  // check if field exists and can be read
90  if (header.headerOk())
91  {
92  Y_.set
93  (
94  i,
95  new volScalarField
96  (
97  IOobject
98  (
99  IOobject::groupName(species_[i], phaseName),
100  mesh.time().timeName(),
101  mesh,
104  ),
105  mesh
106  )
107  );
108  }
109  else
110  {
111  // Read Ydefault if not already read
112  if (!tYdefault.valid())
113  {
114  const word YdefaultName
115  (
116  IOobject::groupName("Ydefault", phaseName)
117  );
118 
120  (
121  YdefaultName,
122  mesh.time().timeName(),
123  mesh,
126  );
127 
129  (
130  YdefaultName,
131  mesh.time().constant(),
132  mesh,
135  );
136 
138  (
139  YdefaultName,
140  Time::timeName(0),
141  mesh,
142  IOobject::MUST_READ,
144  );
145 
146  if (timeIO.headerOk())
147  {
148  tYdefault = new volScalarField(timeIO, mesh);
149  }
150  else if (constantIO.headerOk())
151  {
152  tYdefault = new volScalarField(constantIO, mesh);
153  }
154  else
155  {
156  tYdefault = new volScalarField(time0IO, mesh);
157  }
158  }
159 
160  Y_.set
161  (
162  i,
163  new volScalarField
164  (
165  IOobject
166  (
167  IOobject::groupName(species_[i], phaseName),
168  mesh.time().timeName(),
169  mesh,
172  ),
173  tYdefault()
174  )
175  );
176  }
177  }
178 
179  correctMassFractions();
180 }
181 
182 
184 {
185  if (species_.size())
186  {
188  (
190  (
191  IOobject::groupName("Yt", phaseName_),
192  Y_[0],
193  calculatedFvPatchScalarField::typeName
194  )
195  );
196  volScalarField& Yt = tYt.ref();
197 
198  for (label i=1; i<Y_.size(); i++)
199  {
200  Yt += Y_[i];
201  }
202 
203  if (mag(min(Yt).value()) < rootVSmall)
204  {
206  << "Sum of mass fractions is zero for species " << species()
207  << exit(FatalError);
208  }
209 
210  forAll(Y_, i)
211  {
212  Y_[i] /= Yt;
213  }
214  }
215 }
216 
217 
219 {
220  if (species_.size())
221  {
223  (
225  (
226  IOobject::groupName("Yt", phaseName_),
227  Y_[0].mesh(),
229  )
230  );
231  volScalarField& Yt = tYt.ref();
232 
233  forAll(Y_, i)
234  {
235  if (solve(i))
236  {
237  Y_[i].max(scalar(0));
238  Yt += Y_[i];
239  }
240  }
241 
242  Y_[defaultSpecieIndex_] = scalar(1) - Yt;
243  Y_[defaultSpecieIndex_].max(scalar(0));
244  }
245 }
246 
247 
248 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
Templated form of IOobject providing type information for file reading and header type checking...
Definition: IOobject.H:537
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
const dimensionSet dimless
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:372
fvMesh & mesh
virtual word timeName() const
Return current time name.
Definition: Time.C:676
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
void normalise()
Normalise the mass fractions.
bool valid() const
Is this temporary object valid,.
Definition: tmpI.H:167
static word timeName(const scalar, const int precision=curPrecision_)
Return time name of given scalar time.
Definition: Time.C:666
A class for handling words, derived from string.
Definition: word.H:59
static word groupName(Name name, const word &group)
const word & constant() const
Return constant name.
Definition: TimePaths.H:123
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
Foam::basicMixture.
Definition: basicMixture.H:50
defineTypeNameAndDebug(combustionModel, 0)
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:318
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
void max(const dimensioned< Type > &)
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
SolverPerformance< Type > solve(fvMatrix< Type > &, const word &)
Solve returning the solution statistics given convergence tolerance.
dimensioned< scalar > mag(const dimensioned< Type > &)
A class for managing temporary objects.
Definition: PtrList.H:53
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:98
void correctMassFractions()
Scale the mass fractions to sum to 1.
basicSpecieMixture(const dictionary &, const wordList &specieNames, const fvMesh &, const word &)
Construct from dictionary, species names, mesh and phase name.
ITstream & lookupBackwardsCompatible(const wordList &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream, trying a list of keywords.
Definition: dictionary.C:875
Namespace for OpenFOAM.
IOerror FatalIOError