27 #ifndef OPM_SAT_CURVE_MULTIPLEXER_PARAMS_HPP
28 #define OPM_SAT_CURVE_MULTIPLEXER_PARAMS_HPP
38 #include <type_traits>
44 {
enum class SatCurveMultiplexerApproach {
45 PiecewiseLinearApproach,
57 template <
class TraitsT>
61 using Traits = TraitsT;
62 using Scalar =
typename TraitsT::Scalar;
63 enum { numPhases = 2 };
69 using LETParams =
typename LETTwoPhaseLaw::Params;
72 template <
class ParamT>
75 inline void operator () (
void* ptr )
77 delete static_cast< ParamT*
> (ptr);
81 using ParamPointerType = std::shared_ptr<void>;
95 setApproach( other.approach() );
101 setApproach( other.approach() );
105 void setApproach(SatCurveMultiplexerApproach newApproach)
107 assert(realParams_ == 0);
108 approach_ = newApproach;
110 switch (approach()) {
111 case SatCurveMultiplexerApproach::LETApproach:
112 realParams_ = ParamPointerType(
new LETParams, Deleter< LETParams > () );
115 case SatCurveMultiplexerApproach::PiecewiseLinearApproach:
116 realParams_ = ParamPointerType(
new PLParams, Deleter< PLParams > () );
121 SatCurveMultiplexerApproach approach()
const
122 {
return approach_; }
125 template <SatCurveMultiplexerApproach approachV>
126 typename std::enable_if<approachV == SatCurveMultiplexerApproach::LETApproach, LETParams>::type&
129 assert(approach() == approachV);
130 return this->
template castTo<LETParams>();
133 template <SatCurveMultiplexerApproach approachV>
134 typename std::enable_if<approachV == SatCurveMultiplexerApproach::LETApproach, const LETParams>::type&
135 getRealParams()
const
137 assert(approach() == approachV);
138 return this->
template castTo<LETParams>();
142 template <SatCurveMultiplexerApproach approachV>
143 typename std::enable_if<approachV == SatCurveMultiplexerApproach::PiecewiseLinearApproach, PLParams>::type&
146 assert(approach() == approachV);
147 return this->
template castTo<PLParams>();
150 template <SatCurveMultiplexerApproach approachV>
151 typename std::enable_if<approachV == SatCurveMultiplexerApproach::PiecewiseLinearApproach, const PLParams>::type&
152 getRealParams()
const
154 assert(approach() == approachV);
155 return this->
template castTo<PLParams>();
159 template <
class ParamT>
162 return *(
static_cast<ParamT *
> (realParams_.operator->()));
165 template <
class ParamT>
166 const ParamT& castTo()
const
168 return *(
static_cast<const ParamT *
> (realParams_.operator->()));
171 SatCurveMultiplexerApproach approach_;
172 ParamPointerType realParams_;
Default implementation for asserting finalization of parameter objects.
Specification of the material parameters for a two-phase material law which uses a table and piecewis...
Implementation of a tabulated, piecewise linear capillary pressure law.
Default implementation for asserting finalization of parameter objects.
Definition: EnsureFinalized.hpp:47
Implementation of a tabulated, piecewise linear capillary pressure law.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:50
ParamsT Params
The type of the parameter objects for this law.
Definition: PiecewiseLinearTwoPhaseMaterial.hpp:58
Specification of the material parameters for the saturation function multiplexer.
Definition: SatCurveMultiplexerParams.hpp:59
SatCurveMultiplexerParams()
The multiplexer constructor.
Definition: SatCurveMultiplexerParams.hpp:88
Implementation of the LET curve saturation functions.
Definition: TwoPhaseLETCurves.hpp:49