Loading...
Searching...
No Matches
hConstThermo.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 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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\*---------------------------------------------------------------------------*/
29#include "hConstThermo.H"
30#include "IOstreams.H"
31
32// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33
34template<class EquationOfState>
35Foam::hConstThermo<EquationOfState>::hConstThermo(const dictionary& dict)
36:
37 EquationOfState(dict),
38 Cp_(dict.subDict("thermodynamics").get<scalar>("Cp")),
39 Hf_(dict.subDict("thermodynamics").get<scalar>("Hf")),
40 Tref_(dict.subDict("thermodynamics").getOrDefault<scalar>("Tref", Tstd)),
41 Hsref_(dict.subDict("thermodynamics").getOrDefault<scalar>("Href", 0))
42{}
43
44
45// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
46
47template<class EquationOfState>
49{
50 EquationOfState::write(os);
51
52 // Entries in dictionary format
53 {
54 os.beginBlock("thermodynamics");
55 os.writeEntry("Cp", Cp_);
56 os.writeEntry("Hf", Hf_);
57 os.writeEntryIfDifferent<scalar>("Tref", Tstd, Tref_);
58 os.writeEntryIfDifferent<scalar>("Href", 0, Hsref_);
59 os.endBlock();
60 }
61}
62
63
64// * * * * * * * * * * * * * * * Ostream Operator * * * * * * * * * * * * * //
65
66template<class EquationOfState>
67Foam::Ostream& Foam::operator<<
68(
69 Ostream& os,
70 const hConstThermo<EquationOfState>& ct
71)
72{
73 ct.write(os);
74 return os;
75}
76
77
78// ************************************************************************* //
Useful combination of include files which define Sin, Sout and Serr and the use of IO streams general...
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Constant properties thermodynamics package templated into the EquationOfState.
void write(Ostream &os) const
Write to Ostream.
OBJstream os(runTime.globalPath()/outputName)
dictionary dict