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-------------------------------------------------------------------------------
11License
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
34template<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 (
51 dimMass/dimMoles,
52 thermo.composition().W(i)
53 );
54
55 X_[i] = W*Y[i]/Wi;
56 }
57}
58
59
60// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
61
62template<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,
92 (
94 (
95 "X_" + Y[i].name(),
97 mesh_,
100 ),
101 mesh_,
103 )
104 );
105 }
106
107 calcMoleFractions();
108 }
109 else
110 {
111 if (phaseName_ != word::null)
112 {
114 << "Cannot find thermodynamics model of type "
115 << ThermoType::typeName
116 << " for phase " << phaseName_
117 << exit(FatalError);
118 }
119
121 << "Cannot find thermodynamics model of type "
122 << ThermoType::typeName
123 << exit(FatalError);
124 }
125}
126
127
128// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
129
130template<class ThermoType>
132(
133 const dictionary& dict
134)
135{
137 {
138 phaseName_ = dict.getOrDefault<word>("phase", word::null);
139
140 return true;
141 }
142
143 return false;
144}
145
146
147
148template<class ThermoType>
150{
151 calcMoleFractions();
152
153 return true;
154}
155
156
157// ************************************************************************* //
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
virtual bool read()
Re-read model coefficients if they have changed.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:80
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
static const word dictName
Definition: basicThermo.H:256
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
const word & name() const noexcept
Return the name of this functionObject.
const fvMesh & mesh_
Reference to the fvMesh.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:290
Calculates mole-fraction fields from the mass-fraction fields of the psi/rhoReactionThermo and caches...
virtual bool execute()
Calculate the mole-fraction fields.
bool foundObject(const word &name, const bool recursive=false) const
Is the named Type found?
const Type & lookupObject(const word &name, const bool recursive=false) const
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
A class for handling words, derived from Foam::string.
Definition: word.H:68
static const word null
An empty word.
Definition: word.H:80
PtrList< volScalarField > & Y
engineTime & runTime
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
const word dictName("faMeshDefinition")
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimless
Dimensionless.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
error FatalError
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333