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