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  Copyright (C) 2020-2021 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "hPolynomialThermo.H"
30 #include "IOstreams.H"
31 
32 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33 
34 template<class EquationOfState, int PolySize>
36 (
37  const dictionary& dict
38 )
39 :
40  EquationOfState(dict),
41  Hf_(dict.subDict("thermodynamics").get<scalar>("Hf")),
42  Sf_(dict.subDict("thermodynamics").get<scalar>("Sf")),
43  CpCoeffs_(dict.subDict("thermodynamics").lookup(coeffsName("Cp"))),
44  Tref_(dict.subDict("thermodynamics").getOrDefault<scalar>("Tref", Tstd)),
45  Href_(dict.subDict("thermodynamics").getOrDefault<scalar>("Href", 0)),
46  Sref_(dict.subDict("thermodynamics").getOrDefault<scalar>("Sref", 0)),
47  Pref_(dict.subDict("thermodynamics").getOrDefault<scalar>("Pref", Pstd)),
48  hCoeffs_(),
49  sCoeffs_()
50 {
51  hCoeffs_ = CpCoeffs_.integral();
52  sCoeffs_ = CpCoeffs_.integralMinus1();
53  // Offset h poly so that it is relative to the enthalpy at Tref
54  hCoeffs_[0] +=
55  Href_ - hCoeffs_.value(Tref_) - EquationOfState::H(Pstd, Tref_);
56 
57  // Offset s poly so that it is relative to the entropy at Tref
58  sCoeffs_[0] +=
59  Sref_ - sCoeffs_.value(Tref_) - EquationOfState::S(Pstd, Tref_);
60 }
61 
62 
63 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
64 
65 template<class EquationOfState, int PolySize>
67 (
68  Ostream& os
69 ) const
70 {
72 
73  // Entries in dictionary format
74  {
75  os.beginBlock("thermodynamics");
76  os.writeEntry("Hf", Hf_);
77  os.writeEntry("Sf", Sf_);
78  os.writeEntry("Tref", Tref_);
79  os.writeEntry("Hsref", Href_);
80  os.writeEntry("Sref", Sref_);
81  os.writeEntry("Pref", Pref_);
82  os.writeEntry(coeffsName("Cp"), CpCoeffs_);
83  os.endBlock();
84  }
85 }
86 
87 
88 // * * * * * * * * * * * * * * * Ostream Operator * * * * * * * * * * * * * //
89 
90 template<class EquationOfState, int PolySize>
91 Foam::Ostream& Foam::operator<<
92 (
93  Ostream& os,
95 )
96 {
97  pt.write(os);
98  return os;
99 }
100 
101 
102 // ************************************************************************* //
IOstreams.H
Useful combination of include files which define Sin, Sout and Serr and the use of IO streams general...
Foam::constant::standard::Pstd
const dimensionedScalar Pstd
Standard pressure.
Definition: thermodynamicConstants.C:48
H
volScalarField H(IOobject("H", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), mesh, dimensionedScalar(dimLength, Zero))
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:67
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
os
OBJstream os(runTime.globalPath()/outputName)
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::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:36
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56