SplineTwoPhaseMaterialParams.hpp
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*
4 This file is part of the Open Porous Media project (OPM).
5
6 OPM is free software: you can redistribute it and/or modify
7 it under the terms of the GNU General Public License as published by
8 the Free Software Foundation, either version 2 of the License, or
9 (at your option) any later version.
10
11 OPM is distributed in the hope that it will be useful,
12 but WITHOUT ANY WARRANTY; without even the implied warranty of
13 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 GNU General Public License for more details.
15
16 You should have received a copy of the GNU General Public License
17 along with OPM. If not, see <http://www.gnu.org/licenses/>.
18
19 Consult the COPYING file in the top-level source directory of this
20 module for the precise wording of the license and the list of
21 copyright holders.
22*/
27#ifndef OPM_SPLINE_TWO_PHASE_MATERIAL_PARAMS_HPP
28#define OPM_SPLINE_TWO_PHASE_MATERIAL_PARAMS_HPP
29
32
33#include <vector>
34#include <cassert>
35
36namespace Opm {
43template<class TraitsT>
45{
46 typedef typename TraitsT::Scalar Scalar;
47public:
49
50public:
51 typedef std::vector<Scalar> SamplePoints;
52 typedef ::Opm::Spline<Scalar> Spline;
54
55 typedef TraitsT Traits;
56
58 {
59 }
60
66 const Spline& pcnwSpline() const
67 { EnsureFinalized::check(); return pcwnSpline_; }
68
74 void setPcnwSamples(const SamplePoints& SwSamplePoints,
75 const SamplePoints& pcnwSamplePoints,
76 SplineType splineType = Spline::Monotonic)
77 {
78 assert(SwSamplePoints.size() == pcnwSamplePoints.size());
79 pcwnSpline_.setXYContainers(SwSamplePoints, pcnwSamplePoints, splineType);
80 }
81
88 const Spline& krwSpline() const
89 { EnsureFinalized::check(); return krwSpline_; }
90
97 void setKrwSamples(const SamplePoints& SwSamplePoints,
98 const SamplePoints& krwSamplePoints,
99 SplineType splineType = Spline::Monotonic)
100 {
101 assert(SwSamplePoints.size() == krwSamplePoints.size());
102 krwSpline_.setXYContainers(SwSamplePoints, krwSamplePoints, splineType);
103 }
104
111 const Spline& krnSpline() const
112 { EnsureFinalized::check(); return krnSpline_; }
113
120 void setKrnSamples(const SamplePoints& SwSamplePoints,
121 const SamplePoints& krnSamplePoints,
122 SplineType splineType = Spline::Monotonic)
123 {
124 assert(SwSamplePoints.size() == krnSamplePoints.size());
125 krnSpline_.setXYContainers(SwSamplePoints, krnSamplePoints, splineType);
126 }
127
128private:
129 Spline SwSpline_;
130 Spline pcwnSpline_;
131 Spline krwSpline_;
132 Spline krnSpline_;
133};
134} // namespace Opm
135
136#endif
Default implementation for asserting finalization of parameter objects.
Definition: EnsureFinalized.hpp:47
void finalize()
Mark the object as finalized.
Definition: EnsureFinalized.hpp:75
void check() const
Definition: EnsureFinalized.hpp:63
Class implementing cubic splines.
Definition: Spline.hpp:91
SplineType
The type of the spline to be created.
Definition: Spline.hpp:101
@ Monotonic
Definition: Spline.hpp:105
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:356
Specification of the material parameters for a two-phase material law which uses a table and spline-b...
Definition: SplineTwoPhaseMaterialParams.hpp:45
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
Spline::SplineType SplineType
Definition: SplineTwoPhaseMaterialParams.hpp:53
const Spline & krnSpline() const
Return the sampling points for the relative permeability curve of the non-wetting phase.
Definition: SplineTwoPhaseMaterialParams.hpp:111
std::vector< Scalar > SamplePoints
Definition: SplineTwoPhaseMaterialParams.hpp:51
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
SplineTwoPhaseMaterialParams()
Definition: SplineTwoPhaseMaterialParams.hpp:57
const Spline & pcnwSpline() const
Return the sampling points for the capillary pressure curve.
Definition: SplineTwoPhaseMaterialParams.hpp:66
::Opm::Spline< Scalar > Spline
Definition: SplineTwoPhaseMaterialParams.hpp:52
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
TraitsT Traits
Definition: SplineTwoPhaseMaterialParams.hpp:55
Definition: Air_Mesitylene.hpp:34