27#ifndef OPM_SPLINE_TWO_PHASE_MATERIAL_PARAMS_HPP
28#define OPM_SPLINE_TWO_PHASE_MATERIAL_PARAMS_HPP
43template<
class TraitsT>
46 typedef typename TraitsT::Scalar Scalar;
51 typedef std::vector<Scalar> SamplePoints;
52 typedef ::Opm::Spline<Scalar> Spline;
55 typedef TraitsT Traits;
57 SplineTwoPhaseMaterialParams()
67 { EnsureFinalized::check();
return pcwnSpline_; }
75 const SamplePoints& pcnwSamplePoints,
76 SplineType splineType = Spline::Monotonic)
78 assert(SwSamplePoints.size() == pcnwSamplePoints.size());
79 pcwnSpline_.setXYContainers(SwSamplePoints, pcnwSamplePoints, splineType);
89 { EnsureFinalized::check();
return krwSpline_; }
98 const SamplePoints& krwSamplePoints,
99 SplineType splineType = Spline::Monotonic)
101 assert(SwSamplePoints.size() == krwSamplePoints.size());
102 krwSpline_.setXYContainers(SwSamplePoints, krwSamplePoints, splineType);
112 { EnsureFinalized::check();
return krnSpline_; }
121 const SamplePoints& krnSamplePoints,
122 SplineType splineType = Spline::Monotonic)
124 assert(SwSamplePoints.size() == krnSamplePoints.size());
125 krnSpline_.setXYContainers(SwSamplePoints, krnSamplePoints, splineType);
Default implementation for asserting finalization of parameter objects.
Class implementing cubic splines.
OPM_HOST_DEVICE EnsureFinalized()
The default constructor.
Definition EnsureFinalized.hpp:58
OPM_HOST_DEVICE void finalize()
Mark the object as finalized.
Definition EnsureFinalized.hpp:82
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
const Spline & krwSpline() const
Return the sampling points for the relative permeability curve of the wetting phase.
Definition SplineTwoPhaseMaterialParams.hpp:88
const Spline & krnSpline() const
Return the sampling points for the relative permeability curve of the non-wetting phase.
Definition SplineTwoPhaseMaterialParams.hpp:111
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
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
Class implementing cubic splines.
Definition Spline.hpp:91
SplineType
The type of the spline to be created.
Definition Spline.hpp:101
This class implements a small container which holds the transmissibility mulitpliers for all the face...
Definition Exceptions.hpp:30