All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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-2018 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"
27 #include "basicThermo.H"
28 
29 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
30 
31 template<class ThermoType>
33 {
34  const ThermoType& thermo =
35  mesh_.lookupObject<ThermoType>(basicThermo::dictName);
36 
37  const PtrList<volScalarField>& Y = thermo.composition().Y();
38 
39  const volScalarField W(thermo.W());
40 
41  forAll(Y, i)
42  {
43  const dimensionedScalar Wi
44  (
45  "W",
47  thermo.composition().W(i)
48  );
49 
50  X_[i] = W*Y[i]/Wi;
51  }
52 }
53 
54 
55 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
56 
57 template<class ThermoType>
59 (
60  const word& name,
61  const Time& runTime,
62  const dictionary& dict
63 )
64 :
65  fvMeshFunctionObject(name, runTime, dict)
66 {
67  if (mesh_.foundObject<ThermoType>(basicThermo::dictName))
68  {
69  const ThermoType& thermo =
70  mesh_.lookupObject<ThermoType>(basicThermo::dictName);
71 
72  const PtrList<volScalarField>& Y = thermo.composition().Y();
73 
74  X_.setSize(Y.size());
75 
76  forAll(Y, i)
77  {
78  X_.set
79  (
80  i,
81  new volScalarField
82  (
83  IOobject
84  (
85  "X_" + Y[i].name(),
86  mesh_.time().timeName(),
87  mesh_,
88  IOobject::NO_READ,
89  IOobject::AUTO_WRITE
90  ),
91  mesh_,
92  dimensionedScalar("X", dimless, 0)
93  )
94  );
95  }
96 
97  calculateMoleFractions();
98  }
99  else
100  {
102  << "Cannot find thermodynamics model of type "
103  << ThermoType::typeName
104  << exit(FatalError);
105  }
106 }
107 
108 
109 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
110 
111 template<class ThermoType>
113 {}
114 
115 
116 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
117 
118 template<class ThermoType>
120 (
121  const dictionary& dict
122 )
123 {
124  return true;
125 }
126 
127 
128 template<class ThermoType>
130 {
131  calculateMoleFractions();
132  return true;
133 }
134 
135 
136 template<class ThermoType>
138 {
139  return true;
140 }
141 
142 
143 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
virtual ~moleFractions()
Destructor.
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
const ObjectType & lookupObject(const word &fieldName) const
Lookup object from the objectRegistry.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
A class for handling words, derived from string.
Definition: word.H:59
virtual bool write()
The mole-fraction fields auto-write.
const word dictName("particleTrackDict")
This function object calculates mole-fraction fields from the mass-fraction fields of the psi/rhoReac...
Definition: moleFractions.H:74
void setSize(const label)
Reset size of PtrList. If extending the PtrList, new entries are.
Definition: PtrList.C:131
const dimensionSet dimMoles(0, 0, 0, 0, 1, 0, 0)
Definition: dimensionSets.H:53
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
virtual bool execute()
Calculate the mole-fraction fields.
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:63
const scalarList W(::W(thermo))
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
virtual bool read(const dictionary &)
Read the moleFractions data.