implement a material law for twophase ECL simulations

the basic idea is that implements the threephase API, but only
calculates the quantities for the selected fluid phases.
This commit is contained in:
Andreas Lauser 2015-08-21 13:39:52 +02:00
parent 0db938f199
commit 9bfc26a2d8
6 changed files with 614 additions and 44 deletions

View File

@ -30,6 +30,7 @@
#ifndef OPM_ECL_MATERIAL_LAW_MANAGER_HPP
#define OPM_ECL_MATERIAL_LAW_MANAGER_HPP
#include <opm/material/fluidmatrixinteractions/EclTwoPhaseMaterialParams.hpp>
#include <opm/material/fluidmatrixinteractions/PiecewiseLinearTwoPhaseMaterial.hpp>
#include <opm/material/fluidmatrixinteractions/EclEpsTwoPhaseLaw.hpp>
#include <opm/material/fluidmatrixinteractions/EclHysteresisTwoPhaseLaw.hpp>
@ -48,6 +49,7 @@
#include <algorithm>
namespace Opm {
/*!
@ -286,11 +288,37 @@ private:
void readGlobalThreePhaseOptions_(Opm::DeckConstPtr deck)
{
threePhaseApproach_ = Opm::EclDefaultApproach;
if (deck->hasKeyword("STONE") || deck->hasKeyword("STONE2"))
threePhaseApproach_ = Opm::EclStone2Approach;
else if (deck->hasKeyword("STONE1"))
threePhaseApproach_ = Opm::EclStone1Approach;
bool gasEnabled = deck->hasKeyword("GAS");
bool oilEnabled = deck->hasKeyword("OIL");
bool waterEnabled = deck->hasKeyword("WATER");
int numEnabled =
(gasEnabled?1:0)
+ (oilEnabled?1:0)
+ (waterEnabled?1:0);
if (numEnabled < 2)
OPM_THROW(std::runtime_error,
"At least two fluid phases must be enabled. (Is: " << numEnabled << ")");
if (numEnabled == 2) {
threePhaseApproach_ = Opm::EclTwoPhaseApproach;
if (!gasEnabled)
twoPhaseApproach_ = Opm::EclTwoPhaseOilWater;
else if (!oilEnabled)
twoPhaseApproach_ = Opm::EclTwoPhaseGasWater;
else if (!waterEnabled)
twoPhaseApproach_ = Opm::EclTwoPhaseGasOil;
}
else {
assert(numEnabled == 3);
threePhaseApproach_ = Opm::EclDefaultApproach;
if (deck->hasKeyword("STONE") || deck->hasKeyword("STONE2"))
threePhaseApproach_ = Opm::EclStone2Approach;
else if (deck->hasKeyword("STONE1"))
threePhaseApproach_ = Opm::EclStone1Approach;
}
}
void initNonElemSpecific_(DeckConstPtr deck, EclipseStateConstPtr eclState)
@ -786,6 +814,15 @@ private:
realParams.finalize();
break;
}
case EclTwoPhaseApproach: {
auto& realParams = materialParams.template getRealParams<Opm::EclTwoPhaseApproach>();
realParams.setGasOilParams(gasOilParams);
realParams.setOilWaterParams(oilWaterParams);
realParams.setApproach(twoPhaseApproach_);
realParams.finalize();
break;
}
}
}
@ -807,6 +844,11 @@ private:
auto& realParams = materialParams.template getRealParams<Opm::EclDefaultApproach>();
return realParams.oilWaterParams().drainageParams().scaledPoints();
}
case EclTwoPhaseApproach: {
auto& realParams = materialParams.template getRealParams<Opm::EclTwoPhaseApproach>();
return realParams.oilWaterParams().drainageParams().scaledPoints();
}
}
}
@ -817,7 +859,11 @@ private:
std::vector<Opm::EclEpsScalingPointsInfo<Scalar>> unscaledEpsInfo_;
OilWaterScalingInfoVector oilWaterScaledEpsInfoDrainage_;
EclMultiplexerApproach threePhaseApproach_;
Opm::EclMultiplexerApproach threePhaseApproach_;
// this attribute only makes sense for twophase simulations!
enum EclTwoPhaseApproach twoPhaseApproach_;
std::vector<std::shared_ptr<MaterialLawParams> > materialLawParams_;
std::vector<int> compressedToCartesianElemIdx_;

View File

@ -26,10 +26,13 @@
#define OPM_ECL_MULTIPLEXER_MATERIAL_HPP
#include "EclMultiplexerMaterialParams.hpp"
#include "EclDefaultMaterial.hpp"
#include "EclStone1Material.hpp"
#include "EclStone2Material.hpp"
#include "EclTwoPhaseMaterial.hpp"
#include <opm/material/common/Valgrind.hpp>
#include <opm/material/common/MathToolbox.hpp>
#include <opm/material/common/Exceptions.hpp>
#include <opm/material/common/ErrorMacros.hpp>
@ -58,6 +61,7 @@ public:
typedef Opm::EclStone1Material<TraitsT, GasOilMaterialLaw, OilWaterMaterialLaw> Stone1Material;
typedef Opm::EclStone2Material<TraitsT, GasOilMaterialLaw, OilWaterMaterialLaw> Stone2Material;
typedef Opm::EclDefaultMaterial<TraitsT, GasOilMaterialLaw, OilWaterMaterialLaw> DefaultMaterial;
typedef Opm::EclTwoPhaseMaterial<TraitsT, GasOilMaterialLaw, OilWaterMaterialLaw> TwoPhaseMaterial;
// some safety checks
static_assert(TraitsT::numPhases == 3,
@ -144,6 +148,12 @@ public:
params.template getRealParams<EclDefaultApproach>(),
fluidState);
break;
case EclTwoPhaseApproach:
TwoPhaseMaterial::capillaryPressures(values,
params.template getRealParams<EclTwoPhaseApproach>(),
fluidState);
break;
}
}
@ -258,6 +268,12 @@ public:
params.template getRealParams<EclDefaultApproach>(),
fluidState);
break;
case EclTwoPhaseApproach:
TwoPhaseMaterial::relativePermeabilities(values,
params.template getRealParams<EclTwoPhaseApproach>(),
fluidState);
break;
}
}
@ -317,6 +333,11 @@ public:
DefaultMaterial::updateHysteresis(params.template getRealParams<EclDefaultApproach>(),
fluidState);
break;
case EclTwoPhaseApproach:
TwoPhaseMaterial::updateHysteresis(params.template getRealParams<EclTwoPhaseApproach>(),
fluidState);
break;
}
}
};

View File

@ -28,6 +28,7 @@
#include "EclStone1Material.hpp"
#include "EclStone2Material.hpp"
#include "EclDefaultMaterial.hpp"
#include "EclTwoPhaseMaterial.hpp"
#include <type_traits>
#include <cassert>
@ -38,7 +39,8 @@ namespace Opm {
enum EclMultiplexerApproach {
EclDefaultApproach,
EclStone1Approach,
EclStone2Approach
EclStone2Approach,
EclTwoPhaseApproach
};
/*!
@ -57,16 +59,14 @@ class EclMultiplexerMaterialParams : public Traits
typedef Opm::EclStone1Material<Traits, GasOilMaterialLawT, OilWaterMaterialLawT> Stone1Material;
typedef Opm::EclStone2Material<Traits, GasOilMaterialLawT, OilWaterMaterialLawT> Stone2Material;
typedef Opm::EclDefaultMaterial<Traits, GasOilMaterialLawT, OilWaterMaterialLawT> DefaultMaterial;
typedef Opm::EclTwoPhaseMaterial<Traits, GasOilMaterialLawT, OilWaterMaterialLawT> TwoPhaseMaterial;
typedef typename Stone1Material::Params Stone1Params;
typedef typename Stone2Material::Params Stone2Params;
typedef typename DefaultMaterial::Params DefaultParams;
typedef typename TwoPhaseMaterial::Params TwoPhaseParams;
public:
typedef typename GasOilMaterialLawT::Params GasOilParams;
typedef typename OilWaterMaterialLawT::Params OilWaterParams;
/*!
* \brief The multiplexer constructor.
*/
@ -83,15 +83,19 @@ public:
{
switch (approach()) {
case EclStone1Approach:
delete static_cast<Stone1Material*>(realParams_);
delete static_cast<Stone1Params*>(realParams_);
break;
case EclStone2Approach:
delete static_cast<Stone2Material*>(realParams_);
delete static_cast<Stone2Params*>(realParams_);
break;
case EclDefaultApproach:
delete static_cast<DefaultMaterial*>(realParams_);
delete static_cast<DefaultParams*>(realParams_);
break;
case EclTwoPhaseApproach:
delete static_cast<TwoPhaseParams*>(realParams_);
break;
}
}
@ -106,34 +110,6 @@ public:
#endif
}
OilWaterParams& oilWaterParams()
{
switch (approach()) {
case EclStone1Approach:
return getRealParams<EclStone1Approach>().oilWaterParams();
case EclStone2Approach:
return getRealParams<EclStone2Approach>().oilWaterParams();
case EclDefaultApproach:
return getRealParams<EclDefaultApproach>().oilWaterParams();
}
}
GasOilParams& gasOilParams()
{
switch (approach()) {
case EclStone1Approach:
return getRealParams<EclStone1Approach>().gasOilParams();
case EclStone2Approach:
return getRealParams<EclStone2Approach>().gasOilParams();
case EclDefaultApproach:
return getRealParams<EclDefaultApproach>().gasOilParams();
}
}
void setApproach(EclMultiplexerApproach newApproach)
{
assert(realParams_ == 0);
@ -151,6 +127,10 @@ public:
case EclDefaultApproach:
realParams_ = new DefaultParams;
break;
case EclTwoPhaseApproach:
realParams_ = new TwoPhaseParams;
break;
}
}
@ -191,7 +171,7 @@ public:
return *static_cast<const Stone2Params*>(realParams_);
}
// get the parameter object for the Default case
// get the parameter object for the default case
template <EclMultiplexerApproach approachV>
typename std::enable_if<approachV == EclDefaultApproach, DefaultParams>::type&
getRealParams()
@ -208,6 +188,23 @@ public:
return *static_cast<const DefaultParams*>(realParams_);
}
// get the parameter object for the twophase case
template <EclMultiplexerApproach approachV>
typename std::enable_if<approachV == EclTwoPhaseApproach, TwoPhaseParams>::type&
getRealParams()
{
assert(approach() == approachV);
return *static_cast<TwoPhaseParams*>(realParams_);
}
template <EclMultiplexerApproach approachV>
typename std::enable_if<approachV == EclTwoPhaseApproach, const TwoPhaseParams>::type&
getRealParams() const
{
assert(approach() == approachV);
return *static_cast<const TwoPhaseParams*>(realParams_);
}
private:
#ifndef NDEBUG
void assertFinalized_() const

View File

@ -0,0 +1,361 @@
// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
// vi: set et ts=4 sw=4 sts=4:
/*
Copyright (C) 2015 by Andreas Lauser
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/>.
*/
/*!
* \file
* \copydoc Opm::EclTwoPhaseMaterial
*/
#ifndef OPM_ECL_TWO_PHASE_MATERIAL_HPP
#define OPM_ECL_TWO_PHASE_MATERIAL_HPP
#include "EclTwoPhaseMaterialParams.hpp"
#include <opm/material/common/Valgrind.hpp>
#include <opm/material/common/MathToolbox.hpp>
#include <opm/material/common/Exceptions.hpp>
#include <opm/material/common/ErrorMacros.hpp>
#include <algorithm>
namespace Opm {
/*!
* \ingroup FluidMatrixInteractions
*
* \brief Implements a multiplexer class that provides ECL saturation functions for
* twophase simulations.
*
* The basic idea is that all inputs and outputs are still done on three phases, but only
* the quanties for active phases are calculated.
*/
template <class TraitsT,
class GasOilMaterialLawT,
class OilWaterMaterialLawT,
class ParamsT = EclTwoPhaseMaterialParams<TraitsT,
typename GasOilMaterialLawT::Params,
typename OilWaterMaterialLawT::Params> >
class EclTwoPhaseMaterial : public TraitsT
{
public:
typedef GasOilMaterialLawT GasOilMaterialLaw;
typedef OilWaterMaterialLawT OilWaterMaterialLaw;
// some safety checks
static_assert(TraitsT::numPhases == 3,
"The number of phases considered by this capillary pressure "
"law is always three!");
static_assert(GasOilMaterialLaw::numPhases == 2,
"The number of phases considered by the gas-oil capillary "
"pressure law must be two!");
static_assert(OilWaterMaterialLaw::numPhases == 2,
"The number of phases considered by the oil-water capillary "
"pressure law must be two!");
static_assert(std::is_same<typename GasOilMaterialLaw::Scalar,
typename OilWaterMaterialLaw::Scalar>::value,
"The two two-phase capillary pressure laws must use the same "
"type of floating point values.");
typedef TraitsT Traits;
typedef ParamsT Params;
typedef typename Traits::Scalar Scalar;
static const int numPhases = 3;
static const int waterPhaseIdx = Traits::wettingPhaseIdx;
static const int oilPhaseIdx = Traits::nonWettingPhaseIdx;
static const int gasPhaseIdx = Traits::gasPhaseIdx;
//! Specify whether this material law implements the two-phase
//! convenience API
static const bool implementsTwoPhaseApi = false;
//! Specify whether this material law implements the two-phase
//! convenience API which only depends on the phase saturations
static const bool implementsTwoPhaseSatApi = false;
//! Specify whether the quantities defined by this material law
//! are saturation dependent
static const bool isSaturationDependent = true;
//! Specify whether the quantities defined by this material law
//! are dependent on the absolute pressure
static const bool isPressureDependent = false;
//! Specify whether the quantities defined by this material law
//! are temperature dependent
static const bool isTemperatureDependent = false;
//! Specify whether the quantities defined by this material law
//! are dependent on the phase composition
static const bool isCompositionDependent = false;
/*!
* \brief Implements the multiplexer three phase capillary pressure law
* used by the ECLipse simulator.
*
* This material law is valid for three fluid phases and only
* depends on the saturations.
*
* The required two-phase relations are supplied by means of template
* arguments and can be an arbitrary other material laws.
*
* \param values Container for the return values
* \param params Parameters
* \param state The fluid state
*/
template <class ContainerT, class FluidState>
static void capillaryPressures(ContainerT &values,
const Params &params,
const FluidState &fluidState)
{
typedef typename std::remove_reference<decltype(values[0])>::type Evaluation;
typedef MathToolbox<typename FluidState::Scalar> FsToolbox;
switch (params.approach()) {
case EclTwoPhaseGasOil: {
const Evaluation& So =
FsToolbox::template toLhs<Evaluation>(fluidState.saturation(oilPhaseIdx));
values[oilPhaseIdx] = 0.0;
values[gasPhaseIdx] = GasOilMaterialLaw::twoPhaseSatPcnw(params.gasOilParams(), So);
break;
}
case EclTwoPhaseOilWater: {
const Evaluation& Sw =
FsToolbox::template toLhs<Evaluation>(fluidState.saturation(waterPhaseIdx));
values[waterPhaseIdx] = 0.0;
values[oilPhaseIdx] = OilWaterMaterialLaw::twoPhaseSatPcnw(params.oilWaterParams(), Sw);
break;
}
case EclTwoPhaseGasWater: {
const Evaluation& Sw =
FsToolbox::template toLhs<Evaluation>(fluidState.saturation(waterPhaseIdx));
values[waterPhaseIdx] = 0.0;
values[gasPhaseIdx] =
OilWaterMaterialLaw::twoPhaseSatPcnw(params.oilWaterParams(), Sw)
+ GasOilMaterialLaw::twoPhaseSatPcnw(params.gasOilParams(), 0.0);
break;
}
}
}
/*!
* \brief Capillary pressure between the gas and the non-wetting
* liquid (i.e., oil) phase.
*
* This is defined as
* \f[
* p_{c,gn} = p_g - p_n
* \f]
*/
template <class FluidState, class Evaluation = typename FluidState::Scalar>
static Evaluation pcgn(const Params &params,
const FluidState &fs)
{
OPM_THROW(std::logic_error, "Not implemented: pcgn()");
}
/*!
* \brief Capillary pressure between the non-wetting liquid (i.e.,
* oil) and the wetting liquid (i.e., water) phase.
*
* This is defined as
* \f[
* p_{c,nw} = p_n - p_w
* \f]
*/
template <class FluidState, class Evaluation = typename FluidState::Scalar>
static Evaluation pcnw(const Params &params,
const FluidState &fs)
{
OPM_THROW(std::logic_error, "Not implemented: pcnw()");
}
/*!
* \brief The inverse of the capillary pressure
*/
template <class ContainerT, class FluidState>
static void saturations(ContainerT &values,
const Params &params,
const FluidState &fs)
{
OPM_THROW(std::logic_error, "Not implemented: saturations()");
}
/*!
* \brief The saturation of the gas phase.
*/
template <class FluidState, class Evaluation = typename FluidState::Scalar>
static Evaluation Sg(const Params &params,
const FluidState &fluidState)
{
OPM_THROW(std::logic_error, "Not implemented: Sg()");
}
/*!
* \brief The saturation of the non-wetting (i.e., oil) phase.
*/
template <class FluidState, class Evaluation = typename FluidState::Scalar>
static Evaluation Sn(const Params &params,
const FluidState &fluidState)
{
OPM_THROW(std::logic_error, "Not implemented: Sn()");
}
/*!
* \brief The saturation of the wetting (i.e., water) phase.
*/
template <class FluidState, class Evaluation = typename FluidState::Scalar>
static Evaluation Sw(const Params &params,
const FluidState &fluidState)
{
OPM_THROW(std::logic_error, "Not implemented: Sw()");
}
/*!
* \brief The relative permeability of all phases.
*
* The relative permeability of the water phase it uses the same
* value as the relative permeability for water in the water-oil
* law with \f$S_o = 1 - S_w\f$. The gas relative permebility is
* taken from the gas-oil material law, but with \f$S_o = 1 -
* S_g\f$. The relative permeability of the oil phase is
* calculated using the relative permeabilities of the oil phase
* in the two two-phase systems.
*
* A more detailed description can be found in the "Three phase
* oil relative permeability models" section of the ECLipse
* technical description.
*/
template <class ContainerT, class FluidState>
static void relativePermeabilities(ContainerT &values,
const Params &params,
const FluidState &fluidState)
{
typedef typename std::remove_reference<decltype(values[0])>::type Evaluation;
typedef MathToolbox<typename FluidState::Scalar> FsToolbox;
switch (params.approach()) {
case EclTwoPhaseGasOil: {
const Evaluation& So =
FsToolbox::template toLhs<Evaluation>(fluidState.saturation(oilPhaseIdx));
values[oilPhaseIdx] = GasOilMaterialLaw::twoPhaseSatKrw(params.gasOilParams(), So);
values[gasPhaseIdx] = GasOilMaterialLaw::twoPhaseSatKrn(params.gasOilParams(), So);
break;
}
case EclTwoPhaseOilWater: {
const Evaluation& Sw =
FsToolbox::template toLhs<Evaluation>(fluidState.saturation(waterPhaseIdx));
values[waterPhaseIdx] = OilWaterMaterialLaw::twoPhaseSatKrw(params.oilWaterParams(), Sw);
values[oilPhaseIdx] = OilWaterMaterialLaw::twoPhaseSatKrn(params.oilWaterParams(), Sw);
break;
}
case EclTwoPhaseGasWater: {
const Evaluation& Sw =
FsToolbox::template toLhs<Evaluation>(fluidState.saturation(waterPhaseIdx));
values[waterPhaseIdx] = OilWaterMaterialLaw::twoPhaseSatKrw(params.oilWaterParams(), Sw);
values[gasPhaseIdx] = GasOilMaterialLaw::twoPhaseSatKrn(params.gasOilParams(), Sw);
break;
}
}
}
/*!
* \brief The relative permeability of the gas phase.
*/
template <class FluidState, class Evaluation = typename FluidState::Scalar>
static Evaluation krg(const Params &params,
const FluidState &fluidState)
{
OPM_THROW(std::logic_error, "Not implemented: krg()");
}
/*!
* \brief The relative permeability of the wetting phase.
*/
template <class FluidState, class Evaluation = typename FluidState::Scalar>
static Evaluation krw(const Params &params,
const FluidState &fluidState)
{
OPM_THROW(std::logic_error, "Not implemented: krw()");
}
/*!
* \brief The relative permeability of the non-wetting (i.e., oil) phase.
*/
template <class FluidState, class Evaluation = typename FluidState::Scalar>
static Evaluation krn(const Params &params,
const FluidState &fluidState)
{
OPM_THROW(std::logic_error, "Not implemented: krn()");
}
/*!
* \brief Update the hysteresis parameters after a time step.
*
* This assumes that the nested two-phase material laws are parameters for
* EclHysteresisLaw. If they are not, calling this methid will cause a compiler
* error. (But not calling it will still work.)
*/
template <class FluidState>
static void updateHysteresis(Params &params, const FluidState &fluidState)
{
typedef MathToolbox<typename FluidState::Scalar> FsToolbox;
switch (params.approach()) {
case EclTwoPhaseGasOil: {
Scalar So = FsToolbox::value(fluidState.saturation(oilPhaseIdx));
params.oilWaterParams().update(/*pcSw=*/0.0, /*krwSw=*/0.0, /*krnSw=*/0.0);
params.gasOilParams().update(/*pcSw=*/So, /*krwSw=*/So, /*krnSw=*/So);
break;
}
case EclTwoPhaseOilWater: {
Scalar Sw = FsToolbox::value(fluidState.saturation(waterPhaseIdx));
params.oilWaterParams().update(/*pcSw=*/Sw, /*krwSw=*/Sw, /*krnSw=*/Sw);
params.gasOilParams().update(/*pcSw=*/0.0, /*krwSw=*/0.0, /*krnSw=*/0.0);
break;
}
case EclTwoPhaseGasWater: {
Scalar Sw = FsToolbox::value(fluidState.saturation(waterPhaseIdx));
params.oilWaterParams().update(/*pcSw=*/Sw, /*krwSw=*/Sw, /*krnSw=*/0);
params.gasOilParams().update(/*pcSw=*/1.0, /*krwSw=*/0.0, /*krnSw=*/Sw);
break;
}
}
}
};
} // namespace Opm
#endif

View File

@ -0,0 +1,135 @@
// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
// vi: set et ts=4 sw=4 sts=4:
/*
Copyright (C) 2013 by Andreas Lauser
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/>.
*/
/*!
* \file
* \copydoc Opm::EclTwoPhaseMaterialParams
*/
#ifndef OPM_ECL_TWO_PHASE_MATERIAL_PARAMS_HPP
#define OPM_ECL_TWO_PHASE_MATERIAL_PARAMS_HPP
#include <type_traits>
#include <cassert>
#include <memory>
namespace Opm {
enum EclTwoPhaseApproach {
EclTwoPhaseGasOil,
EclTwoPhaseOilWater,
EclTwoPhaseGasWater
};
/*!
* \brief Implementation for the parameters required by the material law for two-phase
* simulations.
*
* Essentially, this class just stores the two parameter objects for
* the twophase capillary pressure laws.
*/
template<class Traits, class GasOilParamsT, class OilWaterParamsT>
class EclTwoPhaseMaterialParams
{
typedef typename Traits::Scalar Scalar;
enum { numPhases = 3 };
public:
typedef GasOilParamsT GasOilParams;
typedef OilWaterParamsT OilWaterParams;
/*!
* \brief The default constructor.
*/
EclTwoPhaseMaterialParams()
{
#ifndef NDEBUG
finalized_ = false;
#endif
}
/*!
* \brief Finish the initialization of the parameter object.
*/
void finalize()
{
#ifndef NDEBUG
finalized_ = true;
#endif
}
void setApproach(EclTwoPhaseApproach newApproach)
{ approach_ = newApproach; }
EclTwoPhaseApproach approach() const
{ return approach_; }
/*!
* \brief The parameter object for the gas-oil twophase law.
*/
const GasOilParams& gasOilParams() const
{ assertFinalized_(); return *gasOilParams_; }
/*!
* \brief The parameter object for the gas-oil twophase law.
*/
GasOilParams& gasOilParams()
{ assertFinalized_(); return *gasOilParams_; }
/*!
* \brief Set the parameter object for the gas-oil twophase law.
*/
void setGasOilParams(std::shared_ptr<GasOilParams> val)
{ gasOilParams_ = val; }
/*!
* \brief The parameter object for the oil-water twophase law.
*/
const OilWaterParams& oilWaterParams() const
{ assertFinalized_(); return *oilWaterParams_; }
/*!
* \brief The parameter object for the oil-water twophase law.
*/
OilWaterParams& oilWaterParams()
{ assertFinalized_(); return *oilWaterParams_; }
/*!
* \brief Set the parameter object for the oil-water twophase law.
*/
void setOilWaterParams(std::shared_ptr<OilWaterParams> val)
{ oilWaterParams_ = val; }
private:
#ifndef NDEBUG
void assertFinalized_() const
{ assert(finalized_); }
bool finalized_;
#else
void assertFinalized_() const
{ }
#endif
EclTwoPhaseApproach approach_;
std::shared_ptr<GasOilParams> gasOilParams_;
std::shared_ptr<OilWaterParams> oilWaterParams_;
};
} // namespace Opm
#endif

View File

@ -47,6 +47,7 @@
#include <opm/material/fluidmatrixinteractions/EclDefaultMaterial.hpp>
#include <opm/material/fluidmatrixinteractions/EclStone1Material.hpp>
#include <opm/material/fluidmatrixinteractions/EclStone2Material.hpp>
#include <opm/material/fluidmatrixinteractions/EclTwoPhaseMaterial.hpp>
#include <opm/material/fluidmatrixinteractions/EclMultiplexerMaterial.hpp>
// include the helper classes to construct traits
@ -344,6 +345,15 @@ int main(int argc, char **argv)
testThreePhaseApi<MaterialLaw, ThreePhaseFluidState>();
//testThreePhaseSatApi<MaterialLaw, ThreePhaseFluidState>();
}
{
typedef Opm::BrooksCorey<TwoPhaseTraits> TwoPhaseMaterial;
typedef Opm::EclTwoPhaseMaterial<ThreePhaseTraits,
/*GasOilMaterial=*/TwoPhaseMaterial,
/*OilWaterMaterial=*/TwoPhaseMaterial> MaterialLaw;
testGenericApi<MaterialLaw, ThreePhaseFluidState>();
testThreePhaseApi<MaterialLaw, ThreePhaseFluidState>();
//testThreePhaseSatApi<MaterialLaw, ThreePhaseFluidState>();
}
{
typedef Opm::BrooksCorey<TwoPhaseTraits> TwoPhaseMaterial;
typedef Opm::EclMultiplexerMaterial<ThreePhaseTraits,