moleFractions.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 "moleFractions.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  X_()
54 {}
55 
56 
57 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
58 
60 {}
61 
62 
63 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
64 
66 {
67  return true;
68 }
69 
70 
72 {
73  // Lookup or construct a multicomponent thermo. Lookup is likely to be
74  // appropriate if this is being used at run-time. Construction if this is
75  // being called as a one-off post-process.
76  const word thermoName =
77  IOobject::groupName(physicalProperties::typeName, phaseName_);
79  (
80  mesh_.foundObject<fluidMulticomponentThermo>(thermoName)
83  );
85  mesh_.lookupObject<fluidMulticomponentThermo>(thermoName);
86 
87  // Construct mole fraction fields corresponding to the mass fraction fields
88  const PtrList<volScalarField>& Y = thermo.Y();
89  if (X_.empty())
90  {
91  X_.setSize(Y.size());
92 
93  forAll(Y, i)
94  {
95  X_.set
96  (
97  i,
98  new volScalarField
99  (
100  IOobject
101  (
102  "X_" + Y[i].name(),
103  mesh_.time().name(),
104  mesh_
105  ),
106  mesh_,
108  )
109  );
110  }
111  }
112 
113  // Get the mixture molar mass
114  const volScalarField W(thermo.W());
115 
116  // Calculate the mole fractions
117  forAll(Y, i)
118  {
119  X_[i] = Y[i]*W/thermo.Wi(i);
120  }
121 
122  return true;
123 }
124 
125 
127 {
128  forAll(X_, i)
129  {
130  X_[i].write();
131  }
132 
133  return true;
134 }
135 
136 
137 // ************************************************************************* //
#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.
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)
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
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
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.
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 mole-fraction fields from the mass-fraction fields of the multicompon...
Definition: moleFractions.H:69
virtual ~moleFractions()
Destructor.
Definition: moleFractions.C:59
moleFractions(const word &name, const Time &t, const dictionary &dict)
Construct from Time and dictionary.
Definition: moleFractions.C:45
virtual bool execute()
Calculate the mole-fraction fields.
Definition: moleFractions.C:71
virtual bool write()
The mole-fraction fields auto-write.
virtual bool read(const dictionary &)
Read the moleFractions data.
Definition: moleFractions.C:65
A class for handling words, derived from string.
Definition: word.H:62
defineTypeNameAndDebug(adjustTimeStepToCombustion, 0)
addToRunTimeSelectionTable(functionObject, adjustTimeStepToCombustion, dictionary)
Namespace for OpenFOAM.
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
const dimensionSet dimless
scalarList W(const fluidMulticomponentThermo &thermo)
dictionary dict
PtrList< volScalarField > & Y
fluidMulticomponentThermo & thermo
Definition: createFields.H:31