PiecewiseLinearTwoPhaseMaterialParams: use the same API as SplineTwoPhaseMaterialParams

this allows to use the two by changing the class name but without
requiring any further code changes...
This commit is contained in:
Andreas Lauser 2014-12-08 17:50:36 +01:00
parent c267cacbd6
commit 578980e16a
2 changed files with 35 additions and 41 deletions

View File

@ -46,7 +46,7 @@ namespace Opm {
template <class TraitsT, class ParamsT = PiecewiseLinearTwoPhaseMaterialParams<TraitsT> >
class PiecewiseLinearTwoPhaseMaterial : public TraitsT
{
typedef typename ParamsT::SamplePoints SamplePoints;
typedef typename ParamsT::ValueVector ValueVector;
public:
//! The traits class for this material law
@ -360,50 +360,50 @@ public:
}
private:
static Scalar eval_(const SamplePoints &xSamples,
const SamplePoints &ySamples,
static Scalar eval_(const ValueVector &xValues,
const ValueVector &yValues,
Scalar x)
{
int segIdx = findSegmentIndex_(xSamples, x);
int segIdx = findSegmentIndex_(xValues, x);
Scalar x0 = xSamples[segIdx];
Scalar x1 = xSamples[segIdx + 1];
Scalar x0 = xValues[segIdx];
Scalar x1 = xValues[segIdx + 1];
Scalar y0 = ySamples[segIdx];
Scalar y1 = ySamples[segIdx + 1];
Scalar y0 = yValues[segIdx];
Scalar y1 = yValues[segIdx + 1];
return y0 + (y1 - y0)*(x - x0)/(x1 - x0);
}
static Scalar evalDeriv_(const SamplePoints &xSamples,
const SamplePoints &ySamples,
static Scalar evalDeriv_(const ValueVector &xValues,
const ValueVector &yValues,
Scalar x)
{
int segIdx = findSegmentIndex_(xSamples, x);
int segIdx = findSegmentIndex_(xValues, x);
Scalar x0 = xSamples[segIdx];
Scalar x1 = xSamples[segIdx + 1];
Scalar x0 = xValues[segIdx];
Scalar x1 = xValues[segIdx + 1];
Scalar y0 = ySamples[segIdx];
Scalar y1 = ySamples[segIdx + 1];
Scalar y0 = yValues[segIdx];
Scalar y1 = yValues[segIdx + 1];
return (y1 - y0)/(x1 - x0);
}
static int findSegmentIndex_(const SamplePoints &xSamples, Scalar x)
static int findSegmentIndex_(const ValueVector &xValues, Scalar x)
{
int n = xSamples.size() - 1;
int n = xValues.size() - 1;
assert(n >= 1); // we need at least two sampling points!
if (xSamples[n] < x)
if (xValues[n] < x)
return n - 1;
else if (xSamples[0] > x)
else if (xValues[0] > x)
return 0;
// bisection
int lowIdx = 0, highIdx = n;
while (lowIdx + 1 < highIdx) {
int curIdx = (lowIdx + highIdx)/2;
if (xSamples[curIdx] < x)
if (xValues[curIdx] < x)
lowIdx = curIdx;
else
highIdx = curIdx;

View File

@ -38,7 +38,7 @@ class PiecewiseLinearTwoPhaseMaterialParams
typedef typename TraitsT::Scalar Scalar;
public:
typedef std::vector<Scalar> SamplePoints;
typedef std::vector<Scalar> ValueVector;
typedef TraitsT Traits;
@ -79,21 +79,15 @@ public:
/*!
* \brief Return the wetting-phase saturation values of all sampling points.
*/
const SamplePoints& SwSamples() const
const ValueVector& SwSamples() const
{ assertFinalized_(); return SwSamples_; }
/*!
* \brief Set the wetting-phase saturation values of all sampling points.
*/
void setSwSamples(const SamplePoints& samples)
{ SwSamples_ = samples; }
/*!
* \brief Return the sampling points for the capillary pressure curve.
*
* This curve is assumed to depend on the wetting phase saturation
*/
const SamplePoints& pcnwSamples() const
const ValueVector& pcnwSamples() const
{ assertFinalized_(); return pcwnSamples_; }
/*!
@ -101,8 +95,8 @@ public:
*
* This curve is assumed to depend on the wetting phase saturation
*/
void setPcnwSamples(const SamplePoints& samples)
{ pcwnSamples_ = samples; }
void setPcnwSamples(const ValueVector& SwValues, const ValueVector& values)
{ SwSamples_ = SwValues; pcwnSamples_ = values; }
/*!
* \brief Return the sampling points for the relative permeability
@ -110,7 +104,7 @@ public:
*
* This curve is assumed to depend on the wetting phase saturation
*/
const SamplePoints& krwSamples() const
const ValueVector& krwSamples() const
{ assertFinalized_(); return krwSamples_; }
/*!
@ -119,8 +113,8 @@ public:
*
* This curve is assumed to depend on the wetting phase saturation
*/
void setKrwSamples(const SamplePoints& samples)
{ krwSamples_ = samples; }
void setKrwSamples(const ValueVector& SwValues, const ValueVector& values)
{ SwSamples_ = SwValues; krwSamples_ = values; }
/*!
* \brief Return the sampling points for the relative permeability
@ -128,7 +122,7 @@ public:
*
* This curve is assumed to depend on the wetting phase saturation
*/
const SamplePoints& krnSamples() const
const ValueVector& krnSamples() const
{ assertFinalized_(); return krnSamples_; }
/*!
@ -137,8 +131,8 @@ public:
*
* This curve is assumed to depend on the wetting phase saturation
*/
void setKrnSamples(const SamplePoints& samples)
{ krnSamples_ = samples; }
void setKrnSamples(const ValueVector& SwValues, const ValueVector& values)
{ SwSamples_ = SwValues; krnSamples_ = values; }
private:
#ifndef NDEBUG
@ -151,10 +145,10 @@ private:
{ }
#endif
SamplePoints SwSamples_;
SamplePoints pcwnSamples_;
SamplePoints krwSamples_;
SamplePoints krnSamples_;
ValueVector SwSamples_;
ValueVector pcwnSamples_;
ValueVector krwSamples_;
ValueVector krnSamples_;
};
} // namespace Opm