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 | Copyright (C) 2012 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 Class
25  Foam::sensibleInternalEnergy
26 
27 Description
28  Thermodynamics mapping class to expose the sensible internal energy function
29  as the standard internal energy function e(T).
30 
31 \*---------------------------------------------------------------------------*/
32 
33 #ifndef sensibleInternalEnergy_H
34 #define sensibleInternalEnergy_H
35 
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40 
41 /*---------------------------------------------------------------------------*\
42  Class sensibleInternalEnergy Declaration
43 \*---------------------------------------------------------------------------*/
44 
45 template<class Thermo>
47 {
48 
49 public:
50 
51  // Constructors
52 
53  //- Construct
55  {}
56 
57 
58  // Member Functions
59 
60  //- Return the instantiated type name
61  static word typeName()
62  {
63  return "sensibleInternalEnergy";
64  }
65 
66  // Fundamental properties
67 
68  static word name()
69  {
70  return "e";
71  }
72 
73  //- Sensible Internal energy [J/kmol]
74  scalar he
75  (
76  const Thermo& thermo,
77  const scalar p,
78  const scalar T) const
79  {
80  return thermo.es(p, T);
81  }
82 
83  //- Heat capacity at constant volume [J/(kmol K)]
84  scalar cpv
85  (
86  const Thermo& thermo,
87  const scalar p,
88  const scalar T
89  ) const
90  {
91  return thermo.cv(p, T);
92  }
93 
94  //- cp/cv []
95  scalar cpBycpv
96  (
97  const Thermo& thermo,
98  const scalar p,
99  const scalar T
100  ) const
101  {
102  return thermo.gamma(p, T);
103  }
104 
105  //- Sensible internal energy [J/kg]
106  scalar HE
107  (
108  const Thermo& thermo,
109  const scalar p,
110  const scalar T
111  ) const
112  {
113  return thermo.Es(p, T);
114  }
115 
116  //- Temperature from sensible internal energy
117  // given an initial temperature T0
118  scalar THE
119  (
120  const Thermo& thermo,
121  const scalar e,
122  const scalar p,
123  const scalar T0
124  ) const
125  {
126  return thermo.TEs(e, p, T0);
127  }
128 };
129 
130 
131 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
132 
133 } // End namespace Foam
134 
135 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
136 
137 #endif
138 
139 // ************************************************************************* //
const double e
Elementary charge.
Definition: doubleFloat.H:78
scalar HE(const Thermo &thermo, const scalar p, const scalar T) const
Sensible internal energy [J/kg].
scalar THE(const Thermo &thermo, const scalar e, const scalar p, const scalar T0) const
Temperature from sensible internal energy.
scalar he(const Thermo &thermo, const scalar p, const scalar T) const
Sensible Internal energy [J/kmol].
scalar cpv(const Thermo &thermo, const scalar p, const scalar T) const
Heat capacity at constant volume [J/(kmol K)].
A class for handling words, derived from string.
Definition: word.H:59
scalar cpBycpv(const Thermo &thermo, const scalar p, const scalar T) const
cp/cv []
static word typeName()
Return the instantiated type name.
Thermodynamics mapping class to expose the sensible internal energy function as the standard internal...
Namespace for OpenFOAM.