equilibriumCO.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 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 "IOmanip.H"
37 
38 #include "specie.H"
39 #include "perfectGas.H"
40 #include "thermo.H"
41 #include "janafThermo.H"
42 #include "absoluteEnthalpy.H"
43 
44 #include "SLPtrList.H"
45 #include "IOdictionary.H"
46 
47 using namespace Foam;
48 
50  thermo;
51 
52 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
53 
54 int main(int argc, char *argv[])
55 {
56  #include "setRootCase.H"
57  #include "createTime.H"
58 
59  Info<< nl << "Reading thermodynamic data IOdictionary" << endl;
60 
61  IOdictionary thermoData
62  (
63  IOobject
64  (
65  "thermoData",
66  runTime.constant(),
67  runTime,
70  false
71  )
72  );
73 
74 
75 
76  const scalar P = 1e5;
77  const scalar T = 3000.0;
78 
79  // Oxidant (mole-based)
80  thermo O2(thermoData.subDict("O2")); O2 *= O2.W();
81  thermo N2(thermoData.subDict("N2")); N2 *= N2.W();
82 
83  // Intermediates (mole-based)
84  thermo H2(thermoData.subDict("H2")); H2 *= H2.W();
85  thermo OH(thermoData.subDict("OH")); OH *= OH.W();
86  thermo H(thermoData.subDict("H")); H *= H.W();
87  thermo O(thermoData.subDict("O")); O *= O.W();
88 
89  // Products (mole-based)
90  thermo CO2(thermoData.subDict("CO2")); CO2 *= CO2.W();
91  thermo H2O(thermoData.subDict("H2O")); H2O *= H2O.W();
92  thermo CO(thermoData.subDict("CO")); CO *= CO.W();
93 
94  SLPtrList<thermo> EQreactions;
95 
96  EQreactions.append
97  (
98  new thermo(CO2 == CO + 0.5*O2)
99  );
100 
101  EQreactions.append
102  (
103  new thermo(O2 == 2*O)
104  );
105 
106  EQreactions.append
107  (
108  new thermo(H2O == H2 + 0.5*O2)
109  );
110 
111  EQreactions.append
112  (
113  new thermo(H2O == H + OH)
114  );
115 
116 
117  forAllConstIter(SLPtrList<thermo>, EQreactions, iter)
118  {
119  Info<< "Kc(EQreactions) = " << iter().Kc(P, T) << endl;
120  }
121 
122  Info<< nl << "end" << endl;
123 
124  return 0;
125 }
126 
127 
128 // ************************************************************************* //
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
Thermodynamics mapping class to expose the absolute enthalpy functions.
rhoReactionThermo & thermo
Definition: createFields.H:28
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:52
Template class for non-intrusive linked PtrLists.
Definition: LPtrList.H:47
const word & constant() const
Return constant name.
Definition: TimePaths.H:124
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
Istream and Ostream manipulators taking arguments.
static const char nl
Definition: Ostream.H:265
void append(const T * &a)
Add at tail of list.
Definition: LList.H:169
Non-intrusive singly-linked pointer list.
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
Definition: thermo.H:52
messageStream Info
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
Namespace for OpenFOAM.