massFractions.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 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 "massFractions.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace functionObjects
35 {
38 }
39 }
40 
41 
42 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
43 
45 (
46  const word& name,
47  const Time& runTime,
48  const dictionary& dict
49 )
50 :
51  fvMeshFunctionObject(name, runTime, dict),
52  phaseName_(dict.lookupOrDefault<word>("phase", word::null)),
53  Y_()
54 {}
55 
56 
57 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
58 
60 {}
61 
62 
63 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
64 
66 {
67  return true;
68 }
69 
70 
72 {
73  // A multicomponent thermo should not exist for this context, as the
74  // following would create a different, inconsistent definition of the
75  // composition of an existing thermo model
76  const word thermoName =
77  IOobject::groupName(physicalProperties::typeName, phaseName_);
78  if (mesh_.foundObject<fluidMulticomponentThermo>(thermoName))
79  {
81  << "Cannot create mass fractions. Mass fractions already exist "
82  << "within the multicomponent thermo model, \"" << thermoName
83  << "\"." << exit(FatalError);
84  }
85 
86  // Read Ydefault if it exists
88  (
89  "Ydefault",
90  time_.name(),
91  mesh_,
94  false
95  );
96  const bool YdefaultIoExists = YdefaultIo.headerOk();
97  const volScalarField Ydefault
98  (
99  YdefaultIo,
100  mesh_,
102  );
103 
104  // Back up Ydefault if it exists
105  if (YdefaultIoExists)
106  {
107  fileHandler().cp
108  (
109  YdefaultIo.filePath(),
110  YdefaultIo.path(false)/"molesToMassFractions:" + YdefaultIo.name()
111  );
112  }
113 
114  // Write Ydefault out again, but limited to a value of small. This prevents
115  // errors in construction of the thermo if no non-default fractions exist
116  {
117  volScalarField YdefaultLim(Ydefault);
118  YdefaultLim.max(small);
119  YdefaultLim.write();
120  }
121 
122  // Construct a multicomponent thermo
125  fluidMulticomponentThermo& thermo = thermoPtr();
126  const PtrList<volScalarField>& Y = thermo.Y();
127 
128  // Restore the original Ydefault if it exists, and create a new Ydefault if
129  // it does not
130  if (YdefaultIoExists)
131  {
132  fileHandler().mv
133  (
134  YdefaultIo.path(false)/"molesToMassFractions:" + YdefaultIo.name(),
135  YdefaultIo.filePath()
136  );
137  }
138  else
139  {
140  Ydefault.write();
141  }
142 
143  // One-mole constant for conversions
144  static const dimensionedScalar oneMole(dimMoles, 1);
145 
146  // Construct lists of specie molar mass, fields of specie mass, and a field
147  // of total mass
149  PtrList<volScalarField> m(Y.size());
150  volScalarField mTotal
151  (
152  IOobject("mTotal", time_.name(), mesh_),
153  mesh_,
155  );
156  bool fromMoleFractions = false, fromMoles = false;
157  forAll(Y, i)
158  {
159  W.set(i, new dimensionedScalar(thermo.Wi(i)));
160 
162  (
163  Y[i].name(),
164  time_.name(),
165  mesh_,
167  );
169  (
170  "X_" + Y[i].name(),
171  time_.name(),
172  mesh_,
174  );
176  (
177  "n_" + Y[i].name(),
178  time_.name(),
179  mesh_,
181  );
182 
183  // Check consistency of input
184  if (YIo.headerOk())
185  {
187  << "Mass fraction field " << YIo.name()
188  << " already exists on disk" << exit(FatalError);
189  }
190  fromMoleFractions = fromMoleFractions || XIo.headerOk();
191  fromMoles = fromMoles || nIo.headerOk();
192  if (fromMoleFractions && fromMoles)
193  {
195  << "Mole fraction fields and moles fields "
196  << " both found on disk"
197  << exit(FatalError);
198  }
199 
200  // Sum the contributions
201  if (XIo.headerOk())
202  {
203  m.set(i, W[i]*oneMole*volScalarField(XIo, mesh_));
204  }
205  if (nIo.headerOk())
206  {
207  m.set(i, W[i]*volScalarField(nIo, mesh_));
208  }
209  if (XIo.headerOk() || nIo.headerOk())
210  {
211  mTotal += m[i];
212  }
213  }
214 
215  // Steal the thermo's mass fraction fields and delete the thermo
216  Y_.transfer(thermo.Y());
217  thermoPtr.clear();
218 
219  // Divide the specie masses by the total mass to get the mass fractions
220  forAll(Y_, i)
221  {
222  if (m.set(i))
223  {
224  Y_[i] == m[i]/mTotal;
225  }
226  else
227  {
228  Y_.set(i, nullptr);
229  }
230  }
231 
232  return true;
233 }
234 
235 
237 {
238  forAll(Y_, i)
239  {
240  if (Y_.set(i))
241  {
242  Y_[i].write();
243  }
244  }
245 
246  return true;
247 }
248 
249 
250 // ************************************************************************* //
#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.
void max(const dimensioned< Type > &)
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
const word & name() const
Return name.
Definition: IOobject.H:310
static word groupName(Name name, const word &group)
fileName path(const bool global) const
Return complete path including the processor sub-directory.
Definition: IOobject.C:380
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
bool set(const label) const
Is element set.
Definition: PtrListI.H:62
void transfer(PtrList< T > &)
Transfer the contents of the argument PtrList into this PtrList.
Definition: PtrList.C:213
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:76
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
void clear()
Delete object (if the pointer is valid) and set pointer to.
Definition: autoPtrI.H:126
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
virtual bool cp(const fileName &src, const fileName &dst, const bool followLink=true) const =0
Copy, recursively if necessary, the source to the destination.
virtual bool mv(const fileName &src, const fileName &dst, const bool followLink=false) const =0
Rename src to dst.
Base-class for multi-component fluid thermodynamic properties.
static autoPtr< fluidMulticomponentThermo > New(const fvMesh &, const word &phaseName=word::null)
Standard selection based on fvMesh.
Abstract base-class for Time/database functionObjects.
Specialisation of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
This function object calculates mass-fraction fields from mole-fraction or moles fields present on di...
Definition: massFractions.H:74
virtual ~massFractions()
Destructor.
Definition: massFractions.C:59
massFractions(const word &name, const Time &t, const dictionary &dict)
Construct from Time and dictionary.
Definition: massFractions.C:45
virtual bool execute()
Calculate the mass-fraction fields.
Definition: massFractions.C:71
virtual bool write()
The mass-fraction fields auto-write.
virtual bool read(const dictionary &)
Read the massFractions data.
Definition: massFractions.C:65
virtual bool write(const bool write=true) const
Write using setting from DB.
Templated form of IOobject providing type information for file reading and header type checking.
Definition: IOobject.H:531
bool headerOk()
Read header (uses typeGlobalFile to find file) and check.
fileName filePath() const
Return the path for the file for this Type.
Definition: IOobject.H:563
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
defineTypeNameAndDebug(adjustTimeStepToCombustion, 0)
addToRunTimeSelectionTable(functionObject, adjustTimeStepToCombustion, dictionary)
Namespace for OpenFOAM.
const fileOperation & fileHandler()
Get current file handler.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
const dimensionSet dimless
scalarList W(const fluidMulticomponentThermo &thermo)
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:64
const dimensionSet dimMoles
const dimensionSet dimMass
error FatalError
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dictionary dict
PtrList< volScalarField > & Y
fluidMulticomponentThermo & thermo
Definition: createFields.H:31