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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
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 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "hPolynomialThermo.H"
29 #include "IOstreams.H"
30 
31 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32 
33 template<class EquationOfState, int PolySize>
35 (
36  const dictionary& dict
37 )
38 :
39  EquationOfState(dict),
40  Hf_(dict.subDict("thermodynamics").get<scalar>("Hf")),
41  Sf_(dict.subDict("thermodynamics").get<scalar>("Sf")),
42  CpCoeffs_(dict.subDict("thermodynamics").lookup(coeffsName("Cp"))),
43  hCoeffs_(),
44  sCoeffs_()
45 {
46  hCoeffs_ = CpCoeffs_.integral();
47  sCoeffs_ = CpCoeffs_.integralMinus1();
48 
49  // Offset h poly so that it is relative to the enthalpy at Tstd
50  hCoeffs_[0] += Hf_ - hCoeffs_.value(Tstd);
51 
52  // Offset s poly so that it is relative to the entropy at Tstd
53  sCoeffs_[0] += Sf_ - sCoeffs_.value(Tstd);
54 }
55 
56 
57 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
58 
59 template<class EquationOfState, int PolySize>
61 (
62  Ostream& os
63 ) const
64 {
66 
67  // Entries in dictionary format
68  {
69  os.beginBlock("thermodynamics");
70  os.writeEntry("Hf", Hf_);
71  os.writeEntry("Sf", Sf_);
72  os.writeEntry(coeffsName("Cp"), CpCoeffs_);
73  os.endBlock();
74  }
75 }
76 
77 
78 // * * * * * * * * * * * * * * * Ostream Operator * * * * * * * * * * * * * //
79 
80 template<class EquationOfState, int PolySize>
81 Foam::Ostream& Foam::operator<<
82 (
83  Ostream& os,
85 )
86 {
87  pt.write(os);
88  return os;
89 }
90 
91 
92 // ************************************************************************* //
IOstreams.H
Useful combination of include files which define Sin, Sout and Serr and the use of IO streams general...
Foam::Ostream::beginBlock
virtual Ostream & beginBlock(const keyType &kw)
Write begin block group with the given name.
Definition: Ostream.C:91
Foam::constant::standard::Tstd
const dimensionedScalar Tstd
Standard temperature.
Definition: thermodynamicConstants.C:49
Foam::hPolynomialThermo::write
void write(Ostream &os) const
Write to Ostream.
Definition: hPolynomialThermo.C:61
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::Ostream::endBlock
virtual Ostream & endBlock()
Write end block group.
Definition: Ostream.C:109
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
Foam::Ostream::write
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
Foam::hPolynomialThermo
Thermodynamics package templated on the equation of state, using polynomial functions for cp,...
Definition: hPolynomialThermo.H:104
hPolynomialThermo.H
Foam::Ostream::writeEntry
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:232
Foam::vtk::write
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Component-wise write of a value (N times)
Definition: foamVtkOutputTemplates.C:35
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56