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  (
75  IOobject::groupName(basicThermo::dictName, phaseName_)
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_,
98  IOobject::NO_READ,
99  IOobject::AUTO_WRITE
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 
130 template<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 
148 template<class ThermoType>
150 {
151  calcMoleFractions();
152 
153  return true;
154 }
155 
156 
157 // ************************************************************************* //
runTime
engineTime & runTime
Definition: createEngineTime.H:13
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
basicThermo.H
Foam::Time
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:73
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
dictName
const word dictName("faMeshDefinition")
thermo
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
Foam::dimMoles
const dimensionSet dimMoles(0, 0, 0, 0, 1, 0, 0)
Definition: dimensionSets.H:55
Foam::moleFractions::moleFractions
moleFractions(const word &name, const Time &t, const dictionary &dict)
Construct from Time and dictionary.
Definition: moleFractions.C:64
Foam::moleFractions::execute
virtual bool execute()
Calculate the mole-fraction fields.
Definition: moleFractions.C:149
Foam::moleFractions
Calculates mole-fraction fields from the mass-fraction fields of the psi/rhoReactionThermo and caches...
Definition: moleFractions.H:155
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::blockMeshTools::read
void read(Istream &, label &val, const dictionary &)
In-place read with dictionary lookup.
Definition: blockMeshTools.C:57
Foam::moleFractions::read
virtual bool read(const dictionary &dict)
Read the moleFractions data.
Definition: moleFractions.C:132
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:59
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
moleFractions.H
Foam::dimMass
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:51
Y
PtrList< volScalarField > & Y
Definition: createFieldRefs.H:7
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::dimless
const dimensionSet dimless
Dimensionless.
Definition: dimensionSets.C:189