FluidStateEnthalpyModules.hpp
Go to the documentation of this file.
1 // -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2 // vi: set et ts=4 sw=4 sts=4:
3 /*
4  This file is part of the Open Porous Media project (OPM).
5 
6  OPM is free software: you can redistribute it and/or modify
7  it under the terms of the GNU General Public License as published by
8  the Free Software Foundation, either version 2 of the License, or
9  (at your option) any later version.
10 
11  OPM is distributed in the hope that it will be useful,
12  but WITHOUT ANY WARRANTY; without even the implied warranty of
13  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14  GNU General Public License for more details.
15 
16  You should have received a copy of the GNU General Public License
17  along with OPM. If not, see <http://www.gnu.org/licenses/>.
18 
19  Consult the COPYING file in the top-level source directory of this
20  module for the precise wording of the license and the list of
21  copyright holders.
22 */
28 #ifndef OPM_FLUID_STATE_ENTHALPY_MODULES_HPP
29 #define OPM_FLUID_STATE_ENTHALPY_MODULES_HPP
30 
31 #include <opm/common/ErrorMacros.hpp>
32 #include <opm/common/Exceptions.hpp>
33 
35 #include <opm/common/Valgrind.hpp>
36 
37 #include <algorithm>
38 
39 namespace Opm {
44 template <class Scalar,
45  unsigned numPhases,
46  class Implementation>
48 {
49 public:
51  { Valgrind::SetUndefined(enthalpy_); }
52 
56  const Scalar& enthalpy(unsigned phaseIdx) const
57  { return enthalpy_[phaseIdx]; }
58 
62  Scalar internalEnergy(unsigned phaseIdx) const
63  { return enthalpy_[phaseIdx] - asImp_().pressure(phaseIdx)/asImp_().density(phaseIdx); }
64 
68  void setEnthalpy(unsigned phaseIdx, const Scalar& value)
69  { enthalpy_[phaseIdx] = value; }
70 
75  template <class FluidState>
76  void assign(const FluidState& fs)
77  {
78  for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
79  enthalpy_[phaseIdx] = Opm::decay<Scalar>(fs.enthalpy(phaseIdx));
80  }
81  }
82 
91  void checkDefined() const
92  {
93  Valgrind::CheckDefined(enthalpy_);
94  }
95 
96 protected:
97  const Implementation& asImp_() const
98  { return *static_cast<const Implementation*>(this); }
99 
100  Scalar enthalpy_[numPhases];
101 };
102 
108 template <class Scalar,
109  unsigned numPhases,
110  class Implementation>
112 {
113 public:
115  { }
116 
120  const Scalar& internalEnergy(unsigned /* phaseIdx */) const
121  {
122  static Scalar tmp = 0;
123  Valgrind::SetUndefined(tmp);
124  return tmp;
125  }
126 
130  const Scalar& enthalpy(unsigned /* phaseIdx */) const
131  {
132  static Scalar tmp = 0;
133  Valgrind::SetUndefined(tmp);
134  return tmp;
135  }
136 
141  template <class FluidState>
142  void assign(const FluidState& /* fs */)
143  { }
144 
153  void checkDefined() const
154  { }
155 };
156 
157 } // namespace Opm
158 
159 #endif
A traits class which provides basic mathematical functions for arbitrary scalar floating point values...
void checkDefined() const
Make sure that all attributes are defined.
Definition: FluidStateEnthalpyModules.hpp:153
const Scalar & enthalpy(unsigned) const
The specific enthalpy of a fluid phase [J/kg].
Definition: FluidStateEnthalpyModules.hpp:130
Definition: Air_Mesitylene.hpp:33
const Scalar & internalEnergy(unsigned) const
The specific internal energy of a fluid phase [J/kg].
Definition: FluidStateEnthalpyModules.hpp:120
void assign(const FluidState &fs)
Retrieve all parameters from an arbitrary fluid state.
Definition: FluidStateEnthalpyModules.hpp:76
void setEnthalpy(unsigned phaseIdx, const Scalar &value)
Set the specific enthalpy of a phase [J/kg].
Definition: FluidStateEnthalpyModules.hpp:68
Module for the modular fluid state which stores the enthalpies explicitly.
Definition: FluidStateEnthalpyModules.hpp:47
Scalar internalEnergy(unsigned phaseIdx) const
The specific internal energy of a fluid phase [J/kg].
Definition: FluidStateEnthalpyModules.hpp:62
void assign(const FluidState &)
Retrieve all parameters from an arbitrary fluid state.
Definition: FluidStateEnthalpyModules.hpp:142
const Scalar & enthalpy(unsigned phaseIdx) const
The specific enthalpy of a fluid phase [J/kg].
Definition: FluidStateEnthalpyModules.hpp:56
void checkDefined() const
Make sure that all attributes are defined.
Definition: FluidStateEnthalpyModules.hpp:91
Module for the modular fluid state which does not store the enthalpies but returns 0 instead...
Definition: FluidStateEnthalpyModules.hpp:111