sensibleInternalEnergy.H
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) 2012-2017 OpenFOAM Foundation
9 Copyright (C) 2018 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
27Class
28 Foam::sensibleInternalEnergy
29
30Group
31 grpSpecieThermo
32
33Description
34 Thermodynamics mapping class to expose the sensible internal energy
35 functions.
36
37\*---------------------------------------------------------------------------*/
38
39#ifndef sensibleInternalEnergy_H
40#define sensibleInternalEnergy_H
41
42// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43
44namespace Foam
45{
46
47/*---------------------------------------------------------------------------*\
48 Class sensibleInternalEnergy Declaration
49\*---------------------------------------------------------------------------*/
50
51template<class Thermo>
53{
54
55public:
56
57 //- Constructor
59 {}
60
61
62 // Member Functions
63
64 //- Return the instantiated type name
65 static word typeName()
66 {
67 return "sensibleInternalEnergy";
68 }
69
70
71 // Fundamental properties
73 static word energyName()
74 {
75 return "e";
76 }
77
78 //- Heat capacity at constant volume [J/(kg K)]
79 scalar Cpv
80 (
81 const Thermo& thermo,
82 const scalar p,
83 const scalar T
84 ) const
85 {
86 #ifdef __clang__
87 // Using volatile to prevent compiler optimisations leading to
88 // a sigfpe
89 volatile const scalar cv = thermo.Cv(p, T);
90 return cv;
91 #else
92 return thermo.Cv(p, T);
93 #endif
94 }
95
96 //- Ratio of specific heats Cp/Cv []
97 scalar CpByCpv
98 (
99 const Thermo& thermo,
100 const scalar p,
101 const scalar T
102 ) const
103 {
104 #ifdef __clang__
105 // Using volatile to prevent compiler optimisations leading to
106 // a sigfpe
107 volatile const scalar gamma = thermo.gamma(p, T);
108 return gamma;
109 #else
110 return thermo.gamma(p, T);
111 #endif
112 }
113
114 //- Sensible internal energy [J/kg]
115 scalar HE
116 (
117 const Thermo& thermo,
118 const scalar p,
119 const scalar T
120 ) const
121 {
122 #ifdef __clang__
123 // Using volatile to prevent compiler optimisations leading to
124 // a sigfpe
125 volatile const scalar es = thermo.Es(p, T);
126 return es;
127 #else
128 return thermo.Es(p, T);
129 #endif
130 }
131
132 //- Temperature from sensible internal energy given an initial
133 //- temperature T0 [K]
134 scalar THE
135 (
136 const Thermo& thermo,
137 const scalar e,
138 const scalar p,
139 const scalar T0
140 ) const
141 {
142 #ifdef __clang__
143 // Using volatile to prevent compiler optimisations leading to
144 // a sigfpe
145 volatile const scalar tes = thermo.TEs(e, p, T0);
146 return tes;
147 #else
148 return thermo.TEs(e, p, T0);
149 #endif
150 }
151};
152
153
154// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
155
156} // End namespace Foam
157
158// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
159
160#endif
161
162// ************************************************************************* //
Thermodynamics mapping class to expose the sensible internal energy functions.
scalar HE(const Thermo &thermo, const scalar p, const scalar T) const
Sensible internal energy [J/kg].
static word typeName()
Return the instantiated type name.
scalar Cpv(const Thermo &thermo, const scalar p, const scalar T) const
Heat capacity at constant volume [J/(kg K)].
scalar THE(const Thermo &thermo, const scalar e, const scalar p, const scalar T0) const
scalar CpByCpv(const Thermo &thermo, const scalar p, const scalar T) const
Ratio of specific heats Cp/Cv [].
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
A class for handling words, derived from Foam::string.
Definition: word.H:68
volScalarField & p
const volScalarField & T
const scalar gamma
Definition: EEqn.H:9
Namespace for OpenFOAM.
scalar T0
Definition: createFields.H:22
volScalarField & e
Definition: createFields.H:11