hPolynomialThermo.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-2015 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 \*---------------------------------------------------------------------------*/
25 
26 #include "hPolynomialThermo.H"
27 #include "IOstreams.H"
28 
29 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
30 
31 template<class EquationOfState, int PolySize>
33 (
34  Istream& is
35 )
36 :
37  EquationOfState(is),
38  Hf_(readScalar(is)),
39  Sf_(readScalar(is)),
40  CpCoeffs_("CpCoeffs<" + Foam::name(PolySize) + '>', is),
41  hCoeffs_(),
42  sCoeffs_()
43 {
44  Hf_ *= this->W();
45  Sf_ *= this->W();
46  CpCoeffs_ *= this->W();
47 
48  hCoeffs_ = CpCoeffs_.integral();
49  sCoeffs_ = CpCoeffs_.integralMinus1();
50 
51  // Offset h poly so that it is relative to the enthalpy at Tstd
52  hCoeffs_[0] += Hf_ - hCoeffs_.value(Tstd);
53 
54  // Offset s poly so that it is relative to the entropy at Tstd
55  sCoeffs_[0] += Sf_ - sCoeffs_.value(Tstd);
56 }
57 
58 
59 template<class EquationOfState, int PolySize>
61 (
62  const dictionary& dict
63 )
64 :
65  EquationOfState(dict),
66  Hf_(readScalar(dict.subDict("thermodynamics").lookup("Hf"))),
67  Sf_(readScalar(dict.subDict("thermodynamics").lookup("Sf"))),
68  CpCoeffs_
69  (
70  dict.subDict("thermodynamics").lookup
71  (
72  "CpCoeffs<" + Foam::name(PolySize) + '>'
73  )
74  ),
75  hCoeffs_(),
76  sCoeffs_()
77 {
78  Hf_ *= this->W();
79  Sf_ *= this->W();
80  CpCoeffs_ *= this->W();
81 
82  hCoeffs_ = CpCoeffs_.integral();
83  sCoeffs_ = CpCoeffs_.integralMinus1();
84 
85  // Offset h poly so that it is relative to the enthalpy at Tstd
86  hCoeffs_[0] += Hf_ - hCoeffs_.value(Tstd);
87 
88  // Offset s poly so that it is relative to the entropy at Tstd
89  sCoeffs_[0] += Sf_ - sCoeffs_.value(Tstd);
90 }
91 
92 
93 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
94 
95 template<class EquationOfState, int PolySize>
97 (
98  Ostream& os
99 ) const
100 {
102 
103  dictionary dict("thermodynamics");
104  dict.add("Hf", Hf_/this->W());
105  dict.add("Sf", Sf_/this->W());
106  dict.add
107  (
108  word("CpCoeffs<" + Foam::name(PolySize) + '>'),
109  CpCoeffs_/this->W()
110  );
111  os << indent << dict.dictName() << dict;
112 }
113 
114 
115 // * * * * * * * * * * * * * * * Ostream Operator * * * * * * * * * * * * * //
116 
117 template<class EquationOfState, int PolySize>
118 Foam::Ostream& Foam::operator<<
119 (
120  Ostream& os,
122 )
123 {
124  os << static_cast<const EquationOfState&>(pt) << tab
125  << pt.Hf_/pt.W() << tab
126  << pt.Sf_/pt.W() << tab
127  << "CpCoeffs<" << Foam::name(PolySize) << '>' << tab
128  << pt.CpCoeffs_/pt.W();
129 
130  os.check
131  (
132  "operator<<"
133  "("
134  "Ostream&, "
135  "const hPolynomialThermo<EquationOfState, PolySize>&"
136  ")"
137  );
138 
139  return os;
140 }
141 
142 
143 // ************************************************************************* //
#define readScalar
Definition: doubleScalar.C:38
dictionary dict
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:223
static const char tab
Definition: Ostream.H:261
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
const word dictName() const
Return the local dictionary name (final part of scoped name)
Definition: dictionary.H:115
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:633
Thermodynamics package templated on the equation of state, using polynomial functions for cp...
void write(Ostream &os) const
Write to Ostream.
bool add(entry *, bool mergeEntry=false)
Add a new entry.
Definition: dictionary.C:737
Useful combination of include files which define Sin, Sout and Serr and the use of IO streams general...
A class for handling words, derived from string.
Definition: word.H:59
const dimensionedScalar Tstd
Standard temperature.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
runTime write()
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:451