My Project
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 
33 
34 #include <algorithm>
35 
36 namespace Opm {
41 template <class Scalar,
42  unsigned numPhases,
43  class Implementation>
45 {
46 public:
48  { Valgrind::SetUndefined(enthalpy_); }
49 
53  const Scalar& enthalpy(unsigned phaseIdx) const
54  { return enthalpy_[phaseIdx]; }
55 
59  Scalar internalEnergy(unsigned phaseIdx) const
60  { return enthalpy_[phaseIdx] - asImp_().pressure(phaseIdx)/asImp_().density(phaseIdx); }
61 
65  void setEnthalpy(unsigned phaseIdx, const Scalar& value)
66  { enthalpy_[phaseIdx] = value; }
67 
72  template <class FluidState>
73  void assign(const FluidState& fs)
74  {
75  for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
76  enthalpy_[phaseIdx] = decay<Scalar>(fs.enthalpy(phaseIdx));
77  }
78  }
79 
88  void checkDefined() const
89  {
90  Valgrind::CheckDefined(enthalpy_);
91  }
92 
93 protected:
94  const Implementation& asImp_() const
95  { return *static_cast<const Implementation*>(this); }
96 
97  Scalar enthalpy_[numPhases];
98 };
99 
105 template <class Scalar,
106  unsigned numPhases,
107  class Implementation>
109 {
110 public:
112  { }
113 
117  const Scalar& internalEnergy(unsigned /* phaseIdx */) const
118  {
119  static Scalar tmp = 0;
120  Valgrind::SetUndefined(tmp);
121  return tmp;
122  }
123 
127  const Scalar& enthalpy(unsigned /* phaseIdx */) const
128  {
129  static Scalar tmp = 0;
130  Valgrind::SetUndefined(tmp);
131  return tmp;
132  }
133 
138  template <class FluidState>
139  void assign(const FluidState& /* fs */)
140  { }
141 
150  void checkDefined() const
151  { }
152 };
153 
154 } // namespace Opm
155 
156 #endif
A traits class which provides basic mathematical functions for arbitrary scalar floating point values...
Some templates to wrap the valgrind client request macros.
Module for the modular fluid state which stores the enthalpies explicitly.
Definition: FluidStateEnthalpyModules.hpp:45
void checkDefined() const
Make sure that all attributes are defined.
Definition: FluidStateEnthalpyModules.hpp:88
const Scalar & enthalpy(unsigned phaseIdx) const
The specific enthalpy of a fluid phase [J/kg].
Definition: FluidStateEnthalpyModules.hpp:53
void assign(const FluidState &fs)
Retrieve all parameters from an arbitrary fluid state.
Definition: FluidStateEnthalpyModules.hpp:73
void setEnthalpy(unsigned phaseIdx, const Scalar &value)
Set the specific enthalpy of a phase [J/kg].
Definition: FluidStateEnthalpyModules.hpp:65
Scalar internalEnergy(unsigned phaseIdx) const
The specific internal energy of a fluid phase [J/kg].
Definition: FluidStateEnthalpyModules.hpp:59
Module for the modular fluid state which does not store the enthalpies but returns 0 instead.
Definition: FluidStateEnthalpyModules.hpp:109
void assign(const FluidState &)
Retrieve all parameters from an arbitrary fluid state.
Definition: FluidStateEnthalpyModules.hpp:139
void checkDefined() const
Make sure that all attributes are defined.
Definition: FluidStateEnthalpyModules.hpp:150
const Scalar & internalEnergy(unsigned) const
The specific internal energy of a fluid phase [J/kg].
Definition: FluidStateEnthalpyModules.hpp:117
const Scalar & enthalpy(unsigned) const
The specific enthalpy of a fluid phase [J/kg].
Definition: FluidStateEnthalpyModules.hpp:127