Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
10License
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
26Application
27 equilibriumCO
28
29Group
30 grpThermophysicalUtilities
31
32Description
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 "PtrDynList.H"
51#include "IOdictionary.H"
52
53using namespace Foam;
54
56 thermo;
57
58// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
59
60int 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 (
78 (
79 "thermoData",
80 runTime.constant(),
81 runTime,
85 )
86 );
87
88
89 const scalar P = 1e5;
90 const scalar T = 3000.0;
91
92 // Oxidant (mole-based)
93 thermo O2(thermoData.subDict("O2")); O2 *= O2.W();
94 thermo N2(thermoData.subDict("N2")); N2 *= N2.W();
95
96 // Intermediates (mole-based)
97 thermo H2(thermoData.subDict("H2")); H2 *= H2.W();
98 thermo OH(thermoData.subDict("OH")); OH *= OH.W();
99 thermo H(thermoData.subDict("H")); H *= H.W();
100 thermo O(thermoData.subDict("O")); O *= O.W();
101
102 // Products (mole-based)
103 thermo CO2(thermoData.subDict("CO2")); CO2 *= CO2.W();
104 thermo H2O(thermoData.subDict("H2O")); H2O *= H2O.W();
105 thermo CO(thermoData.subDict("CO")); CO *= CO.W();
106
107 PtrDynList<thermo> EQreactions(4);
108
109 EQreactions.emplace_back((CO2 == CO + 0.5*O2));
110 EQreactions.emplace_back((O2 == 2*O));
111 EQreactions.emplace_back((H2O == H2 + 0.5*O2));
112 EQreactions.emplace_back((H2O == H + OH));
113
114
115 for (const thermo& react : EQreactions)
116 {
117 Info<< "Kc(EQreactions) = " << react.Kc(P, T) << endl;
118 }
119
120 Info<< nl << "End" << nl << endl;
121
122 return 0;
123}
124
125
126// ************************************************************************* //
Istream and Ostream manipulators taking arguments.
Functions used by OpenFOAM that are specific to POSIX compliant operating systems and need to be repl...
water
Definition H2O.H:58
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
@ NO_REGISTER
Do not request registration (bool: false).
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
Liquid N2.
Definition N2.H:59
A dynamically resizable PtrList with allocation management.
Definition PtrDynList.H:58
Thermodynamics mapping class to expose the absolute enthalpy functions.
static void noFunctionObjects(bool addWithOption=false)
Remove '-noFunctionObjects' option and ignore any occurrences.
Definition argList.C:562
static void noParallel()
Remove the parallel options.
Definition argList.C:599
static void addNote(const string &note)
Add extra notes for the usage information.
Definition argList.C:477
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
scalar W() const
Molecular weight [kg/kmol].
const volScalarField & T
engineTime & runTime
volScalarField H(IOobject("H", runTime.timeName(), mesh.thisDb(), IOobject::NO_READ, IOobject::AUTO_WRITE), mesh, dimensionedScalar(dimLength, Zero))
Namespace for OpenFOAM.
messageStream Info
Information stream (stdout output on master, null elsewhere).
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50