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-2014 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 
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  scalar P = 1e5;
77  scalar T = 3000.0;
78 
79  SLPtrList<thermo> EQreactions;
80 
81  EQreactions.append
82  (
83  new thermo
84  (
85  thermo(thermoData.subDict("CO2"))
86  ==
87  thermo(thermoData.subDict("CO"))
88  + 0.5*thermo(thermoData.subDict("O2"))
89  )
90  );
91 
92  EQreactions.append
93  (
94  new thermo
95  (
96  thermo(thermoData.subDict("O2"))
97  ==
98  2.0*thermo(thermoData.subDict("O"))
99  )
100  );
101 
102  EQreactions.append
103  (
104  new thermo
105  (
106  thermo(thermoData.subDict("H2O"))
107  ==
108  thermo(thermoData.subDict("H2"))
109  + 0.5*thermo(thermoData.subDict("O2"))
110  )
111  );
112 
113  EQreactions.append
114  (
115  new thermo
116  (
117  thermo(thermoData.subDict("H2O"))
118  ==
119  thermo(thermoData.subDict("H"))
120  + thermo(thermoData.subDict("OH"))
121  )
122  );
123 
124 
125  forAllConstIter(SLPtrList<thermo>, EQreactions, iter)
126  {
127  Info<< "Kc(EQreactions) = " << iter().Kc(P, T) << endl;
128  }
129 
130  Info<< nl << "End" << nl << endl;
131 
132  return 0;
133 }
134 
135 
136 // ************************************************************************* //
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
OSspecific.H
Functions used by OpenFOAM that are specific to POSIX compliant operating systems and need to be repl...
SLPtrList.H
thermo.H
Foam::absoluteEnthalpy
Thermodynamics mapping class to expose the absolute enthalpy function as the standard enthalpy functi...
Definition: absoluteEnthalpy.H:45
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)
Add newline and flush stream.
Definition: Ostream.H:251
specie.H
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
thermo
rhoThermo & thermo
Definition: setRegionFluidFields.H:3
forAllConstIter
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
Foam::IOobject::MUST_READ_IF_MODIFIED
@ MUST_READ_IF_MODIFIED
Definition: IOobject.H:109
Foam::SLPtrList
Non-intrusive singly-linked pointer list.
Definition: SLPtrList.H:47
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
janafThermo.H
argList.H
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:55
IOmanip.H
Istream and Ostream manipulators taking arguments.
main
int main(int argc, char *argv[])
Definition: postCalc.C:54
IFstream.H
perfectGas.H
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::species::thermo
Definition: thermo.H:52
setRootCase.H
dictionary.H
createTime.H
absoluteEnthalpy.H