equilibriumCO.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) 2011-2017 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 Application
27  equilibriumCO
28 
29 Group
30  grpThermophysicalUtilities
31 
32 Description
33  Calculate the equilibrium level of carbon monoxide.
34 
35 \*---------------------------------------------------------------------------*/
36 
37 #include "argList.H"
38 #include "Time.H"
39 #include "dictionary.H"
40 #include "IFstream.H"
41 #include "OSspecific.H"
42 #include "IOmanip.H"
43 
44 #include "specie.H"
45 #include "perfectGas.H"
46 #include "thermo.H"
47 #include "janafThermo.H"
48 #include "absoluteEnthalpy.H"
49 
50 #include "SLPtrList.H"
51 #include "IOdictionary.H"
52 
53 using namespace Foam;
54 
56  thermo;
57 
58 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
59 
60 int main(int argc, char *argv[])
61 {
63  (
64  "Calculate the equilibrium level of carbon monoxide."
65  );
66 
68  argList::noFunctionObjects(); // Never use function objects
69 
70  #include "setRootCase.H"
71  #include "createTime.H"
72 
73  Info<< nl << "Reading thermodynamic data IOdictionary" << endl;
74 
75  IOdictionary thermoData
76  (
77  IOobject
78  (
79  "thermoData",
80  runTime.constant(),
81  runTime,
84  false
85  )
86  );
87 
88 
89 
90  const scalar P = 1e5;
91  const scalar T = 3000.0;
92 
93  // Oxidant (mole-based)
94  thermo O2(thermoData.subDict("O2")); O2 *= O2.W();
95  thermo N2(thermoData.subDict("N2")); N2 *= N2.W();
96 
97  // Intermediates (mole-based)
98  thermo H2(thermoData.subDict("H2")); H2 *= H2.W();
99  thermo OH(thermoData.subDict("OH")); OH *= OH.W();
100  thermo H(thermoData.subDict("H")); H *= H.W();
101  thermo O(thermoData.subDict("O")); O *= O.W();
102 
103  // Products (mole-based)
104  thermo CO2(thermoData.subDict("CO2")); CO2 *= CO2.W();
105  thermo H2O(thermoData.subDict("H2O")); H2O *= H2O.W();
106  thermo CO(thermoData.subDict("CO")); CO *= CO.W();
107 
108  SLPtrList<thermo> EQreactions;
109 
110  EQreactions.append
111  (
112  new thermo(CO2 == CO + 0.5*O2)
113  );
114 
115  EQreactions.append
116  (
117  new thermo(O2 == 2*O)
118  );
119 
120  EQreactions.append
121  (
122  new thermo(H2O == H2 + 0.5*O2)
123  );
124 
125  EQreactions.append
126  (
127  new thermo(H2O == H + OH)
128  );
129 
130 
131  for (const thermo& react : EQreactions)
132  {
133  Info<< "Kc(EQreactions) = " << react.Kc(P, T) << endl;
134  }
135 
136  Info<< nl << "End" << nl << endl;
137 
138  return 0;
139 }
140 
141 
142 // ************************************************************************* //
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:191
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:50
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:165
OSspecific.H
Functions used by OpenFOAM that are specific to POSIX compliant operating systems and need to be repl...
SLPtrList.H
Non-intrusive singly-linked pointer list.
Foam::N2
Liquid N2.
Definition: N2.H:56
thermo
psiReactionThermo & thermo
Definition: createFields.H:28
H
volScalarField H(IOobject("H", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), mesh, dimensionedScalar(dimLength, Zero))
Foam::absoluteEnthalpy
Thermodynamics mapping class to expose the absolute enthalpy functions.
Definition: absoluteEnthalpy.H:43
Foam::argList::addNote
static void addNote(const string &note)
Definition: argList.C:405
thermo
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
Foam::endl
Ostream & endl(Ostream &os)
Definition: Ostream.H:381
specie.H
Foam::thermophysicalProperties::W
scalar W() const
Definition: thermophysicalPropertiesI.H:29
Foam::LPtrList
Template class for non-intrusive linked PtrLists.
Definition: LPtrList.H:46
Foam::argList::noFunctionObjects
static void noFunctionObjects(bool addWithOption=false)
Definition: argList.C:466
Foam::Info
messageStream Info
janafThermo.H
Foam::H2O
water
Definition: H2O.H:55
argList.H
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:51
IOmanip.H
Istream and Ostream manipulators taking arguments.
IFstream.H
perfectGas.H
Foam
Definition: atmBoundaryLayer.C:26
Foam::species::thermo
Definition: thermo.H:51
IOdictionary.H
Time.H
setRootCase.H
Foam::nl
constexpr char nl
Definition: Ostream.H:424
dictionary.H
createTime.H
Foam::IOobject::MUST_READ_IF_MODIFIED
@ MUST_READ_IF_MODIFIED
Definition: IOobject.H:182
Foam::argList::noParallel
static void noParallel()
Definition: argList.C:503
Foam::TimePaths::constant
const word & constant() const
Definition: TimePathsI.H:89
absoluteEnthalpy.H