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
|
|
|
/*
|
|
|
|
|
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/>.
|
2016-03-14 14:11:22 +01:00
|
|
|
|
|
|
|
|
Consult the COPYING file in the top-level source directory of this
|
|
|
|
|
module for the precise wording of the license and the list of
|
|
|
|
|
copyright holders.
|
2013-12-02 16:31:53 +01:00
|
|
|
*/
|
2013-05-30 15:45:31 +02:00
|
|
|
/*!
|
|
|
|
|
* \file
|
2018-07-27 12:45:59 +02:00
|
|
|
* \copydoc Opm::SinglePhaseFluidSystem
|
2013-05-30 15:45:31 +02:00
|
|
|
*/
|
2014-05-08 16:53:53 +02:00
|
|
|
#ifndef OPM_SINGLE_PHASE_FLUIDSYSTEM_HPP
|
|
|
|
|
#define OPM_SINGLE_PHASE_FLUIDSYSTEM_HPP
|
2013-05-30 15:45:31 +02:00
|
|
|
|
2013-09-10 14:55:59 +02:00
|
|
|
#include "BaseFluidSystem.hpp"
|
|
|
|
|
#include "NullParameterCache.hpp"
|
2013-05-30 15:45:31 +02:00
|
|
|
|
2013-09-10 14:55:59 +02:00
|
|
|
#include <opm/material/fluidsystems/LiquidPhase.hpp>
|
|
|
|
|
#include <opm/material/fluidsystems/GasPhase.hpp>
|
|
|
|
|
#include <opm/material/components/SimpleH2O.hpp>
|
|
|
|
|
#include <opm/material/components/H2O.hpp>
|
|
|
|
|
#include <opm/material/components/N2.hpp>
|
|
|
|
|
#include <opm/material/components/TabulatedComponent.hpp>
|
2013-05-30 15:45:31 +02:00
|
|
|
|
2018-02-01 15:40:01 +01:00
|
|
|
#include <opm/material/common/Unused.hpp>
|
2015-09-22 11:20:55 +02:00
|
|
|
|
2013-05-30 15:45:31 +02:00
|
|
|
#include <limits>
|
|
|
|
|
#include <cassert>
|
|
|
|
|
|
|
|
|
|
namespace Opm {
|
|
|
|
|
|
|
|
|
|
/*!
|
|
|
|
|
* \ingroup Fluidsystems
|
|
|
|
|
*
|
|
|
|
|
* \brief A fluid system for single phase models.
|
|
|
|
|
*
|
|
|
|
|
* The fluid is defined as a template parameter. For existing
|
|
|
|
|
* components the Opm::LiquidPhase<Component> and
|
|
|
|
|
* Opm::GasPhase<Component> may be used.
|
|
|
|
|
*/
|
|
|
|
|
template <class Scalar, class Fluid>
|
2018-07-27 12:45:59 +02:00
|
|
|
class SinglePhaseFluidSystem
|
|
|
|
|
: public BaseFluidSystem<Scalar, SinglePhaseFluidSystem<Scalar, Fluid> >
|
2013-05-30 15:45:31 +02:00
|
|
|
{
|
2018-07-27 12:45:59 +02:00
|
|
|
typedef SinglePhaseFluidSystem<Scalar, Fluid> ThisType;
|
2013-05-30 15:45:31 +02:00
|
|
|
typedef BaseFluidSystem<Scalar, ThisType> Base;
|
|
|
|
|
|
|
|
|
|
public:
|
|
|
|
|
//! \copydoc BaseFluidSystem::ParameterCache
|
2016-04-17 11:28:04 +02:00
|
|
|
template <class Evaluation>
|
|
|
|
|
struct ParameterCache : public Opm::NullParameterCache<Evaluation>
|
|
|
|
|
{};
|
2013-05-30 15:45:31 +02:00
|
|
|
|
|
|
|
|
/****************************************
|
|
|
|
|
* Fluid phase related static parameters
|
|
|
|
|
****************************************/
|
|
|
|
|
|
|
|
|
|
//! \copydoc BaseFluidSystem::numPhases
|
2013-06-06 11:21:02 +02:00
|
|
|
static const int numPhases = 1;
|
2013-05-30 15:45:31 +02:00
|
|
|
|
|
|
|
|
//! \copydoc BaseFluidSystem::phaseName
|
2017-01-18 16:53:58 +01:00
|
|
|
static const char* phaseName(unsigned phaseIdx OPM_OPTIM_UNUSED)
|
2013-05-30 15:45:31 +02:00
|
|
|
{
|
|
|
|
|
assert(0 <= phaseIdx && phaseIdx < numPhases);
|
|
|
|
|
|
|
|
|
|
return Fluid::name();
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
//! \copydoc BaseFluidSystem::isLiquid
|
2015-09-22 11:20:55 +02:00
|
|
|
static bool isLiquid(unsigned /*phaseIdx*/)
|
2013-05-30 15:45:31 +02:00
|
|
|
{
|
|
|
|
|
//assert(0 <= phaseIdx && phaseIdx < numPhases);
|
|
|
|
|
|
|
|
|
|
return Fluid::isLiquid();
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
//! \copydoc BaseFluidSystem::isCompressible
|
2015-09-22 11:20:55 +02:00
|
|
|
static bool isCompressible(unsigned /*phaseIdx*/)
|
2013-05-30 15:45:31 +02:00
|
|
|
{
|
|
|
|
|
//assert(0 <= phaseIdx && phaseIdx < numPhases);
|
|
|
|
|
|
|
|
|
|
// let the fluid decide
|
|
|
|
|
return Fluid::isCompressible();
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
//! \copydoc BaseFluidSystem::isIdealGas
|
2015-09-22 11:20:55 +02:00
|
|
|
static bool isIdealMixture(unsigned /*phaseIdx*/)
|
2013-05-30 15:45:31 +02:00
|
|
|
{
|
|
|
|
|
//assert(0 <= phaseIdx && phaseIdx < numPhases);
|
|
|
|
|
|
|
|
|
|
// we assume immisibility
|
|
|
|
|
return true;
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
//! \copydoc BaseFluidSystem::isIdealMixture
|
2015-09-22 11:20:55 +02:00
|
|
|
static bool isIdealGas(unsigned /*phaseIdx*/)
|
2013-05-30 15:45:31 +02:00
|
|
|
{
|
|
|
|
|
//assert(0 <= phaseIdx && phaseIdx < numPhases);
|
|
|
|
|
|
|
|
|
|
// let the fluid decide
|
|
|
|
|
return Fluid::isIdealGas();
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
/****************************************
|
|
|
|
|
* Component related static parameters
|
|
|
|
|
****************************************/
|
|
|
|
|
|
|
|
|
|
//! \copydoc BaseFluidSystem::numComponents
|
2013-06-06 11:21:02 +02:00
|
|
|
static const int numComponents = 1;
|
2013-05-30 15:45:31 +02:00
|
|
|
|
|
|
|
|
//! \copydoc BaseFluidSystem::componentName
|
2017-01-18 16:53:58 +01:00
|
|
|
static const char* componentName(unsigned compIdx OPM_OPTIM_UNUSED)
|
2013-05-30 15:45:31 +02:00
|
|
|
{
|
|
|
|
|
assert(0 <= compIdx && compIdx < numComponents);
|
|
|
|
|
|
|
|
|
|
return Fluid::name();
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
//! \copydoc BaseFluidSystem::molarMass
|
2015-09-22 11:20:55 +02:00
|
|
|
static Scalar molarMass(unsigned /*compIdx*/)
|
2013-05-30 15:45:31 +02:00
|
|
|
{
|
|
|
|
|
//assert(0 <= compIdx && compIdx < numComponents);
|
|
|
|
|
|
|
|
|
|
return Fluid::molarMass();
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
/*!
|
|
|
|
|
* \brief Critical temperature of a component [K].
|
|
|
|
|
*
|
|
|
|
|
* \param compIdx The index of the component to consider
|
|
|
|
|
*/
|
2015-09-22 11:20:55 +02:00
|
|
|
static Scalar criticalTemperature(unsigned /*compIdx*/)
|
2013-05-30 15:45:31 +02:00
|
|
|
{
|
|
|
|
|
//assert(0 <= compIdx && compIdx < numComponents);
|
|
|
|
|
|
|
|
|
|
return Fluid::criticalTemperature();
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
/*!
|
|
|
|
|
* \brief Critical pressure of a component [Pa].
|
|
|
|
|
*
|
|
|
|
|
* \param compIdx The index of the component to consider
|
|
|
|
|
*/
|
2015-09-22 11:20:55 +02:00
|
|
|
static Scalar criticalPressure(unsigned /*compIdx*/)
|
2013-05-30 15:45:31 +02:00
|
|
|
{
|
|
|
|
|
//assert(0 <= compIdx && compIdx < numComponents);
|
|
|
|
|
|
|
|
|
|
return Fluid::criticalPressure();
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
/*!
|
|
|
|
|
* \brief The acentric factor of a component [].
|
|
|
|
|
*
|
|
|
|
|
* \param compIdx The index of the component to consider
|
|
|
|
|
*/
|
2015-09-22 11:20:55 +02:00
|
|
|
static Scalar acentricFactor(unsigned /*compIdx*/)
|
2013-05-30 15:45:31 +02:00
|
|
|
{
|
|
|
|
|
//assert(0 <= compIdx && compIdx < numComponents);
|
|
|
|
|
|
|
|
|
|
return Fluid::acentricFactor();
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
/****************************************
|
|
|
|
|
* thermodynamic relations
|
|
|
|
|
****************************************/
|
|
|
|
|
|
|
|
|
|
//! \copydoc BaseFluidSystem::init
|
|
|
|
|
static void init()
|
|
|
|
|
{ }
|
|
|
|
|
|
|
|
|
|
//! \copydoc BaseFluidSystem::density
|
2016-04-17 11:28:04 +02:00
|
|
|
template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCacheEval = LhsEval>
|
2017-01-18 17:04:59 +01:00
|
|
|
static LhsEval density(const FluidState& fluidState,
|
|
|
|
|
const ParameterCache<ParamCacheEval>& /*paramCache*/,
|
2015-09-22 11:20:55 +02:00
|
|
|
unsigned phaseIdx)
|
2013-05-30 15:45:31 +02:00
|
|
|
{
|
|
|
|
|
assert(0 <= phaseIdx && phaseIdx < numPhases);
|
|
|
|
|
|
2017-06-12 16:36:14 +02:00
|
|
|
const auto& T = Opm::decay<LhsEval>(fluidState.temperature(phaseIdx));
|
|
|
|
|
const auto& p = Opm::decay<LhsEval>(fluidState.pressure(phaseIdx));
|
2015-05-21 15:33:16 +02:00
|
|
|
return Fluid::density(T, p);
|
2013-05-30 15:45:31 +02:00
|
|
|
}
|
|
|
|
|
|
|
|
|
|
//! \copydoc BaseFluidSystem::viscosity
|
2016-04-17 11:28:04 +02:00
|
|
|
template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCacheEval = LhsEval>
|
2017-01-18 17:04:59 +01:00
|
|
|
static LhsEval viscosity(const FluidState& fluidState,
|
|
|
|
|
const ParameterCache<ParamCacheEval>& /*paramCache*/,
|
2015-09-22 11:20:55 +02:00
|
|
|
unsigned phaseIdx)
|
2013-05-30 15:45:31 +02:00
|
|
|
{
|
|
|
|
|
assert(0 <= phaseIdx && phaseIdx < numPhases);
|
|
|
|
|
|
2017-06-12 16:36:14 +02:00
|
|
|
const auto& T = Opm::decay<LhsEval>(fluidState.temperature(phaseIdx));
|
|
|
|
|
const auto& p = Opm::decay<LhsEval>(fluidState.pressure(phaseIdx));
|
2015-05-21 15:33:16 +02:00
|
|
|
return Fluid::viscosity(T, p);
|
2013-05-30 15:45:31 +02:00
|
|
|
}
|
|
|
|
|
|
|
|
|
|
//! \copydoc BaseFluidSystem::fugacityCoefficient
|
2016-04-17 11:28:04 +02:00
|
|
|
template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCacheEval = LhsEval>
|
2017-01-18 17:04:59 +01:00
|
|
|
static LhsEval fugacityCoefficient(const FluidState& /*fluidState*/,
|
|
|
|
|
const ParameterCache<ParamCacheEval>& /*paramCache*/,
|
2015-09-22 11:20:55 +02:00
|
|
|
unsigned phaseIdx,
|
|
|
|
|
unsigned compIdx)
|
2013-05-30 15:45:31 +02:00
|
|
|
{
|
2017-01-18 17:04:59 +01:00
|
|
|
assert(0 <= phaseIdx && phaseIdx < numPhases);
|
|
|
|
|
assert(0 <= compIdx && compIdx < numComponents);
|
2013-05-30 15:45:31 +02:00
|
|
|
|
|
|
|
|
if (phaseIdx == compIdx)
|
|
|
|
|
// TODO (?): calculate the real fugacity coefficient of
|
|
|
|
|
// the component in the fluid. Probably that's not worth
|
|
|
|
|
// the effort, since the fugacity coefficient of the other
|
|
|
|
|
// component is infinite anyway...
|
|
|
|
|
return 1.0;
|
|
|
|
|
return std::numeric_limits<Scalar>::infinity();
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
//! \copydoc BaseFluidSystem::enthalpy
|
2016-04-17 11:28:04 +02:00
|
|
|
template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCacheEval = LhsEval>
|
2017-01-18 17:04:59 +01:00
|
|
|
static LhsEval enthalpy(const FluidState& fluidState,
|
|
|
|
|
const ParameterCache<ParamCacheEval>& /*paramCache*/,
|
2015-09-22 11:20:55 +02:00
|
|
|
unsigned phaseIdx)
|
2013-05-30 15:45:31 +02:00
|
|
|
{
|
|
|
|
|
assert(0 <= phaseIdx && phaseIdx < numPhases);
|
|
|
|
|
|
2017-06-12 16:36:14 +02:00
|
|
|
const auto& T = Opm::decay<LhsEval>(fluidState.temperature(phaseIdx));
|
|
|
|
|
const auto& p = Opm::decay<LhsEval>(fluidState.pressure(phaseIdx));
|
2015-05-21 15:33:16 +02:00
|
|
|
return Fluid::enthalpy(T, p);
|
2013-05-30 15:45:31 +02:00
|
|
|
}
|
|
|
|
|
|
|
|
|
|
//! \copydoc BaseFluidSystem::thermalConductivity
|
2016-04-17 11:28:04 +02:00
|
|
|
template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCacheEval = LhsEval>
|
2017-01-18 17:04:59 +01:00
|
|
|
static LhsEval thermalConductivity(const FluidState& fluidState,
|
|
|
|
|
const ParameterCache<ParamCacheEval>& /*paramCache*/,
|
2015-09-22 11:20:55 +02:00
|
|
|
unsigned phaseIdx)
|
2013-05-30 15:45:31 +02:00
|
|
|
{
|
|
|
|
|
assert(0 <= phaseIdx && phaseIdx < numPhases);
|
|
|
|
|
|
2017-06-12 16:36:14 +02:00
|
|
|
const auto& T = Opm::decay<LhsEval>(fluidState.temperature(phaseIdx));
|
|
|
|
|
const auto& p = Opm::decay<LhsEval>(fluidState.pressure(phaseIdx));
|
2015-05-21 15:33:16 +02:00
|
|
|
return Fluid::thermalConductivity(T, p);
|
2013-05-30 15:45:31 +02:00
|
|
|
}
|
|
|
|
|
|
|
|
|
|
//! \copydoc BaseFluidSystem::heatCapacity
|
2016-04-17 11:28:04 +02:00
|
|
|
template <class FluidState, class LhsEval = typename FluidState::Scalar, class ParamCacheEval = LhsEval>
|
2017-01-18 17:04:59 +01:00
|
|
|
static LhsEval heatCapacity(const FluidState& fluidState,
|
|
|
|
|
const ParameterCache<ParamCacheEval>& /*paramCache*/,
|
2015-09-22 11:20:55 +02:00
|
|
|
unsigned phaseIdx)
|
2013-05-30 15:45:31 +02:00
|
|
|
{
|
|
|
|
|
assert(0 <= phaseIdx && phaseIdx < numPhases);
|
|
|
|
|
|
2017-06-12 16:36:14 +02:00
|
|
|
const auto& T = Opm::decay<LhsEval>(fluidState.temperature(phaseIdx));
|
|
|
|
|
const auto& p = Opm::decay<LhsEval>(fluidState.pressure(phaseIdx));
|
2015-05-21 15:33:16 +02:00
|
|
|
return Fluid::heatCapacity(T, p);
|
2013-05-30 15:45:31 +02:00
|
|
|
}
|
|
|
|
|
};
|
|
|
|
|
|
2018-07-27 12:45:59 +02:00
|
|
|
} // namespace Opm
|
2013-05-30 15:45:31 +02:00
|
|
|
|
|
|
|
|
#endif
|