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 | Copyright (C) 2011-2017 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 Application
25  equilibriumCO
26 
27 Description
28  Calculates the equilibrium level of carbon monoxide.
29 
30 \*---------------------------------------------------------------------------*/
31 
32 #include "argList.H"
33 #include "Time.H"
34 #include "dictionary.H"
35 #include "IFstream.H"
36 #include "OSspecific.H"
37 #include "IOmanip.H"
38 
39 #include "specie.H"
40 #include "perfectGas.H"
41 #include "thermo.H"
42 #include "janafThermo.H"
43 #include "absoluteEnthalpy.H"
44 
45 #include "SLPtrList.H"
46 #include "IOdictionary.H"
47 
48 using namespace Foam;
49 
51  thermo;
52 
53 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
54 
55 int main(int argc, char *argv[])
56 {
57  #include "setRootCase.H"
58  #include "createTime.H"
59 
60  Info<< nl << "Reading thermodynamic data IOdictionary" << endl;
61 
62  IOdictionary thermoData
63  (
64  IOobject
65  (
66  "thermoData",
67  runTime.constant(),
68  runTime,
71  false
72  )
73  );
74 
75 
76 
77  const scalar P = 1e5;
78  const scalar T = 3000.0;
79 
80  // Oxidant (mole-based)
81  thermo O2(thermoData.subDict("O2")); O2 *= O2.W();
82  thermo N2(thermoData.subDict("N2")); N2 *= N2.W();
83 
84  // Intermediates (mole-based)
85  thermo H2(thermoData.subDict("H2")); H2 *= H2.W();
86  thermo OH(thermoData.subDict("OH")); OH *= OH.W();
87  thermo H(thermoData.subDict("H")); H *= H.W();
88  thermo O(thermoData.subDict("O")); O *= O.W();
89 
90  // Products (mole-based)
91  thermo CO2(thermoData.subDict("CO2")); CO2 *= CO2.W();
92  thermo H2O(thermoData.subDict("H2O")); H2O *= H2O.W();
93  thermo CO(thermoData.subDict("CO")); CO *= CO.W();
94 
95  SLPtrList<thermo> EQreactions;
96 
97  EQreactions.append
98  (
99  new thermo(CO2 == CO + 0.5*O2)
100  );
101 
102  EQreactions.append
103  (
104  new thermo(O2 == 2*O)
105  );
106 
107  EQreactions.append
108  (
109  new thermo(H2O == H2 + 0.5*O2)
110  );
111 
112  EQreactions.append
113  (
114  new thermo(H2O == H + OH)
115  );
116 
117 
118  forAllConstIter(SLPtrList<thermo>, EQreactions, iter)
119  {
120  Info<< "Kc(EQreactions) = " << iter().Kc(P, T) << endl;
121  }
122 
123  Info<< nl << "end" << endl;
124 
125  return 0;
126 }
127 
128 
129 // ************************************************************************* //
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
Thermodynamics mapping class to expose the absolute enthalpy functions.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:52
psiReactionThermo & thermo
Definition: createFields.H:31
Functions used by OpenFOAM that are specific to POSIX compliant operating systems and need to be repl...
Template class for non-intrusive linked PtrLists.
Definition: LPtrList.H:47
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
Istream and Ostream manipulators taking arguments.
static const char nl
Definition: Ostream.H:262
void append(const T * &a)
Add at tail of list.
Definition: LList.H:169
Non-intrusive singly-linked pointer list.
messageStream Info
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
Namespace for OpenFOAM.