moleFractions.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2016-2020 OpenFOAM Foundation
9  Copyright (C) 2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "moleFractions.H"
30 #include "basicThermo.H"
31 
32 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
33 
34 template<class ThermoType>
36 {
37  const auto& thermo =
38  mesh_.lookupObject<ThermoType>
39  (
40  IOobject::groupName(basicThermo::dictName, phaseName_)
41  );
42 
43  const PtrList<volScalarField>& Y = thermo.composition().Y();
44 
45  const volScalarField W(thermo.W());
46 
47  forAll(Y, i)
48  {
49  const dimensionedScalar Wi
50  (
52  thermo.composition().W(i)
53  );
54 
55  X_[i] = W*Y[i]/Wi;
56  }
57 }
58 
59 
60 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
61 
62 template<class ThermoType>
64 (
65  const word& name,
66  const Time& runTime,
67  const dictionary& dict
68 )
69 :
70  fvMeshFunctionObject(name, runTime, dict),
71  phaseName_(dict.getOrDefault<word>("phase", word::null))
72 {
73  const word dictName
74  (
76  );
77 
78  if (mesh_.foundObject<ThermoType>(dictName))
79  {
80  const auto& thermo = mesh_.lookupObject<ThermoType>(dictName);
81 
82  const PtrList<volScalarField>& Y = thermo.composition().Y();
83 
84  X_.setSize(Y.size());
85 
86  forAll(Y, i)
87  {
88  X_.set
89  (
90  i,
91  new volScalarField
92  (
93  IOobject
94  (
95  "X_" + Y[i].name(),
96  mesh_.time().timeName(),
97  mesh_,
100  ),
101  mesh_,
103  )
104  );
105  }
106 
107  calcMoleFractions();
108  }
109  else
110  {
112  << "Cannot find thermodynamics model of type "
113  << ThermoType::typeName;
114 
115  if (!phaseName_.empty())
116  {
117  FatalError
118  << " for phase " << phaseName_;
119  }
120 
122  }
123 }
124 
125 
126 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
127 
128 template<class ThermoType>
130 (
131  const dictionary& dict
132 )
133 {
135  {
136  phaseName_ = dict.getOrDefault<word>("phase", word::null);
137 
138  return true;
139  }
141  return false;
142 }
143 
144 
145 
146 template<class ThermoType>
148 {
149  calcMoleFractions();
150 
151  return true;
152 }
153 
154 
155 // ************************************************************************* //
word dictName() const
The local dictionary name (final part of scoped name)
Definition: dictionaryI.H:53
dictionary dict
const Type & lookupObject(const word &name, const bool recursive=false) const
Lookup and return const reference to the object of the given Type. Fatal if not found or the wrong ty...
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
const word dictName("faMeshDefinition")
engineTime & runTime
const dimensionSet dimless
Dimensionless.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
const word & name() const noexcept
Return the name of this functionObject.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
virtual bool read(const dictionary &dict)
Read the moleFractions data.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
psiReactionThermo & thermo
Definition: createFields.H:28
moleFractions(const word &name, const Time &t, const dictionary &dict)
Construct from Time and dictionary.
Definition: moleFractions.C:57
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
A class for handling words, derived from Foam::string.
Definition: word.H:63
static const word null
An empty word.
Definition: word.H:84
Calculates mole-fraction fields from the mass-fraction fields of the psi/rhoReactionThermo and caches...
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition: Time.C:714
const dimensionSet dimMoles(0, 0, 0, 0, 1, 0, 0)
Definition: dimensionSets.H:53
virtual bool execute()
Calculate the mole-fraction fields.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
PtrList< volScalarField > & Y
Nothing to be read.
Automatically write from objectRegistry::writeObject()
virtual bool read(const dictionary &dict)
Read optional controls.
bool foundObject(const word &name, const bool recursive=false) const
Is the named Type found?
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
const fvMesh & mesh_
Reference to the fvMesh.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127