2015-06-18 13:47:07 +02:00
|
|
|
// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
|
|
|
|
|
// vi: set et ts=4 sw=4 sts=4:
|
2013-12-02 16:31:53 +01:00
|
|
|
/*
|
2013-12-02 16:41:45 +01:00
|
|
|
Copyright (C) 2011-2013 by Andreas Lauser
|
2013-12-02 16:31:53 +01:00
|
|
|
|
|
|
|
|
This file is part of the Open Porous Media project (OPM).
|
|
|
|
|
|
|
|
|
|
OPM is free software: you can redistribute it and/or modify
|
|
|
|
|
it under the terms of the GNU General Public License as published by
|
|
|
|
|
the Free Software Foundation, either version 2 of the License, or
|
|
|
|
|
(at your option) any later version.
|
|
|
|
|
|
|
|
|
|
OPM is distributed in the hope that it will be useful,
|
|
|
|
|
but WITHOUT ANY WARRANTY; without even the implied warranty of
|
|
|
|
|
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
|
|
|
|
|
GNU General Public License for more details.
|
|
|
|
|
|
|
|
|
|
You should have received a copy of the GNU General Public License
|
|
|
|
|
along with OPM. If not, see <http://www.gnu.org/licenses/>.
|
|
|
|
|
*/
|
2013-05-30 15:45:31 +02:00
|
|
|
/*!
|
|
|
|
|
* \file
|
|
|
|
|
*
|
|
|
|
|
* \brief Modules for the ModularFluidState which represent enthalpy.
|
|
|
|
|
*/
|
2013-11-13 18:45:52 +01:00
|
|
|
#ifndef OPM_FLUID_STATE_ENTHALPY_MODULES_HPP
|
|
|
|
|
#define OPM_FLUID_STATE_ENTHALPY_MODULES_HPP
|
2013-05-30 15:45:31 +02:00
|
|
|
|
2015-04-27 16:34:36 +02:00
|
|
|
#include <opm/material/common/ErrorMacros.hpp>
|
|
|
|
|
#include <opm/material/common/Exceptions.hpp>
|
2013-05-30 15:45:31 +02:00
|
|
|
|
2015-05-21 15:33:14 +02:00
|
|
|
#include <opm/material/common/MathToolbox.hpp>
|
|
|
|
|
#include <opm/material/common/Valgrind.hpp>
|
|
|
|
|
|
2013-05-30 15:45:31 +02:00
|
|
|
#include <algorithm>
|
|
|
|
|
|
|
|
|
|
namespace Opm {
|
|
|
|
|
/*!
|
|
|
|
|
* \brief Module for the modular fluid state which stores the
|
2013-12-03 11:56:45 +01:00
|
|
|
* enthalpies explicitly.
|
2013-05-30 15:45:31 +02:00
|
|
|
*/
|
|
|
|
|
template <class Scalar,
|
2015-07-28 17:24:27 +02:00
|
|
|
int numPhases,
|
2013-05-30 15:45:31 +02:00
|
|
|
class Implementation>
|
|
|
|
|
class FluidStateExplicitEnthalpyModule
|
|
|
|
|
{
|
|
|
|
|
public:
|
|
|
|
|
FluidStateExplicitEnthalpyModule()
|
|
|
|
|
{ Valgrind::SetUndefined(enthalpy_); }
|
|
|
|
|
|
|
|
|
|
/*!
|
|
|
|
|
* \brief The specific enthalpy of a fluid phase [J/kg]
|
|
|
|
|
*/
|
2015-05-21 15:33:14 +02:00
|
|
|
const Scalar& enthalpy(int phaseIdx) const
|
2013-05-30 15:45:31 +02:00
|
|
|
{ return enthalpy_[phaseIdx]; }
|
|
|
|
|
|
|
|
|
|
/*!
|
|
|
|
|
* \brief The specific internal energy of a fluid phase [J/kg]
|
|
|
|
|
*/
|
|
|
|
|
Scalar internalEnergy(int phaseIdx) const
|
|
|
|
|
{ return enthalpy_[phaseIdx] - asImp_().pressure(phaseIdx)/asImp_().density(phaseIdx); }
|
|
|
|
|
|
|
|
|
|
/*!
|
|
|
|
|
* \brief Set the specific enthalpy of a phase [J/kg]
|
|
|
|
|
*/
|
2015-05-21 15:33:14 +02:00
|
|
|
void setEnthalpy(int phaseIdx, const Scalar& value)
|
2013-05-30 15:45:31 +02:00
|
|
|
{ enthalpy_[phaseIdx] = value; }
|
|
|
|
|
|
|
|
|
|
/*!
|
|
|
|
|
* \brief Retrieve all parameters from an arbitrary fluid
|
|
|
|
|
* state.
|
|
|
|
|
*/
|
|
|
|
|
template <class FluidState>
|
|
|
|
|
void assign(const FluidState& fs)
|
|
|
|
|
{
|
2015-05-21 15:33:14 +02:00
|
|
|
typedef typename FluidState::Scalar FsScalar;
|
|
|
|
|
typedef Opm::MathToolbox<FsScalar> FsToolbox;
|
2013-05-30 15:45:31 +02:00
|
|
|
for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
|
2015-05-21 15:33:14 +02:00
|
|
|
enthalpy_[phaseIdx] = FsToolbox::template toLhs<Scalar>(fs.enthalpy(phaseIdx));
|
2013-05-30 15:45:31 +02:00
|
|
|
}
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
/*!
|
|
|
|
|
* \brief Make sure that all attributes are defined.
|
|
|
|
|
*
|
|
|
|
|
* This method does not do anything if the program is not run
|
|
|
|
|
* under valgrind. If it is, then valgrind will print an error
|
|
|
|
|
* message if some attributes of the object have not been properly
|
|
|
|
|
* defined.
|
|
|
|
|
*/
|
|
|
|
|
void checkDefined() const
|
|
|
|
|
{
|
|
|
|
|
Valgrind::CheckDefined(enthalpy_);
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
protected:
|
|
|
|
|
const Implementation &asImp_() const
|
|
|
|
|
{ return *static_cast<const Implementation*>(this); }
|
|
|
|
|
|
|
|
|
|
Scalar enthalpy_[numPhases];
|
|
|
|
|
};
|
|
|
|
|
|
|
|
|
|
/*!
|
|
|
|
|
* \brief Module for the modular fluid state which does not store the
|
2015-05-21 15:33:14 +02:00
|
|
|
* enthalpies but returns 0 instead.
|
|
|
|
|
*
|
|
|
|
|
* Also, the returned values are marked as undefined in Valgrind... */
|
2013-05-30 15:45:31 +02:00
|
|
|
template <class Scalar,
|
2015-07-28 17:24:27 +02:00
|
|
|
int numPhases,
|
2013-05-30 15:45:31 +02:00
|
|
|
class Implementation>
|
|
|
|
|
class FluidStateNullEnthalpyModule
|
|
|
|
|
{
|
|
|
|
|
public:
|
|
|
|
|
FluidStateNullEnthalpyModule()
|
|
|
|
|
{ }
|
|
|
|
|
|
|
|
|
|
/*!
|
|
|
|
|
* \brief The specific internal energy of a fluid phase [J/kg]
|
|
|
|
|
*/
|
2015-05-21 15:33:14 +02:00
|
|
|
const Scalar& internalEnergy(int phaseIdx) const
|
|
|
|
|
{
|
|
|
|
|
static Scalar tmp = 0;
|
|
|
|
|
Valgrind::SetUndefined(tmp);
|
|
|
|
|
return tmp;
|
|
|
|
|
}
|
2013-05-30 15:45:31 +02:00
|
|
|
|
|
|
|
|
/*!
|
|
|
|
|
* \brief The specific enthalpy of a fluid phase [J/kg]
|
|
|
|
|
*/
|
2015-05-21 15:33:14 +02:00
|
|
|
const Scalar& enthalpy(int phaseIdx) const
|
|
|
|
|
{
|
|
|
|
|
static Scalar tmp = 0;
|
|
|
|
|
Valgrind::SetUndefined(tmp);
|
|
|
|
|
return tmp;
|
|
|
|
|
}
|
2013-05-30 15:45:31 +02:00
|
|
|
|
|
|
|
|
/*!
|
|
|
|
|
* \brief Retrieve all parameters from an arbitrary fluid
|
|
|
|
|
* state.
|
|
|
|
|
*/
|
|
|
|
|
template <class FluidState>
|
|
|
|
|
void assign(const FluidState& fs)
|
|
|
|
|
{ }
|
|
|
|
|
|
|
|
|
|
/*!
|
|
|
|
|
* \brief Make sure that all attributes are defined.
|
|
|
|
|
*
|
|
|
|
|
* This method does not do anything if the program is not run
|
|
|
|
|
* under valgrind. If it is, then valgrind will print an error
|
|
|
|
|
* message if some attributes of the object have not been properly
|
|
|
|
|
* defined.
|
|
|
|
|
*/
|
|
|
|
|
void checkDefined() const
|
|
|
|
|
{ }
|
|
|
|
|
};
|
|
|
|
|
|
2013-11-04 14:15:53 +01:00
|
|
|
} // namespace Opm
|
2013-05-30 15:45:31 +02:00
|
|
|
|
|
|
|
|
#endif
|