27 #ifndef OPM_SPLINE_TWO_PHASE_MATERIAL_PARAMS_HPP 28 #define OPM_SPLINE_TWO_PHASE_MATERIAL_PARAMS_HPP 43 template<
class TraitsT>
46 typedef typename TraitsT::Scalar Scalar;
51 typedef std::vector<Scalar> SamplePoints;
55 typedef TraitsT Traits;
67 { EnsureFinalized::check();
return pcwnSpline_; }
75 const SamplePoints& pcnwSamplePoints,
78 assert(SwSamplePoints.size() == pcnwSamplePoints.size());
79 pcwnSpline_.
setXYContainers(SwSamplePoints, pcnwSamplePoints, splineType);
89 { EnsureFinalized::check();
return krwSpline_; }
98 const SamplePoints& krwSamplePoints,
101 assert(SwSamplePoints.size() == krwSamplePoints.size());
102 krwSpline_.
setXYContainers(SwSamplePoints, krwSamplePoints, splineType);
112 { EnsureFinalized::check();
return krnSpline_; }
121 const SamplePoints& krnSamplePoints,
124 assert(SwSamplePoints.size() == krnSamplePoints.size());
125 krnSpline_.
setXYContainers(SwSamplePoints, krnSamplePoints, splineType);
SplineType
The type of the spline to be created.
Definition: Spline.hpp:103
const Spline & krnSpline() const
Return the sampling points for the relative permeability curve of the non-wetting phase...
Definition: SplineTwoPhaseMaterialParams.hpp:111
void setXYContainers(const ScalarContainerX &x, const ScalarContainerY &y, Scalar m0, Scalar m1, bool sortInputs=true)
Set the sampling points and the boundary slopes of a full spline using STL-compatible containers...
Definition: Spline.hpp:358
Definition: Air_Mesitylene.hpp:33
Class implementing cubic splines.
void setPcnwSamples(const SamplePoints &SwSamplePoints, const SamplePoints &pcnwSamplePoints, SplineType splineType=Spline::Monotonic)
Set the sampling points for the capillary pressure curve.
Definition: SplineTwoPhaseMaterialParams.hpp:74
Specification of the material parameters for a two-phase material law which uses a table and spline-b...
Definition: SplineTwoPhaseMaterialParams.hpp:44
const Spline & krwSpline() const
Return the sampling points for the relative permeability curve of the wetting phase.
Definition: SplineTwoPhaseMaterialParams.hpp:88
Default implementation for asserting finalization of parameter objects.
void finalize()
Mark the object as finalized.
Definition: EnsureFinalized.hpp:77
void setKrwSamples(const SamplePoints &SwSamplePoints, const SamplePoints &krwSamplePoints, SplineType splineType=Spline::Monotonic)
Set the sampling points for the relative permeability curve of the wetting phase. ...
Definition: SplineTwoPhaseMaterialParams.hpp:97
void setKrnSamples(const SamplePoints &SwSamplePoints, const SamplePoints &krnSamplePoints, SplineType splineType=Spline::Monotonic)
Set the sampling points for the relative permeability curve of the non-wetting phase.
Definition: SplineTwoPhaseMaterialParams.hpp:120
const Spline & pcnwSpline() const
Return the sampling points for the capillary pressure curve.
Definition: SplineTwoPhaseMaterialParams.hpp:66
Default implementation for asserting finalization of parameter objects.
Definition: EnsureFinalized.hpp:46
Class implementing cubic splines.
Definition: Spline.hpp:92