24#ifndef DUMUX_FIXED_LENGTH_SPLINE_HH
25#define DUMUX_FIXED_LENGTH_SPLINE_HH
27#include <dune/common/fvector.hh>
28#include <dune/common/fmatrix.hh>
29#include <dune/istl/btdmatrix.hh>
40template<
class ScalarT,
int nSamples>
43 FixedLengthSpline_<ScalarT, nSamples> >
47 using Scalar = ScalarT;
48 using Vector = Dune::FieldVector<Scalar, nSamples>;
49 using BlockVector = Dune::BlockVector<Dune::FieldVector<Scalar, 1> >;
50 using BTDMatrix = Dune::BTDMatrix<Dune::FieldMatrix<Scalar, 1, 1> >;
83 template <
class ScalarArrayX,
class ScalarArrayY>
85 const ScalarArrayX &x,
86 const ScalarArrayY &y,
91 for (
int i = 0; i < numberSamples; ++i) {
109 template <
class ScalarContainerX,
class ScalarContainerY>
111 const ScalarContainerY &y,
112 Scalar m0, Scalar m1)
114 assert(x.size() == y.size());
117 std::copy(x.begin(), x.end(),
xPos_.begin());
118 std::copy(y.begin(), y.end(),
yPos_.begin());
134 template <
class Po
intArray>
136 const PointArray &points,
144 for (
int i = 0; i < numberSamples; ++i) {
145 xPos_[i] = points[i][0];
146 yPos_[i] = points[i][1];
163 template <
class XYContainer>
170 typename XYContainer::const_iterator it = points.begin();
171 typename XYContainer::const_iterator endIt = points.end();
172 for (
int i = 0; it != endIt; ++i, ++it) {
196 template <
class XYContainer>
203 typename XYContainer::const_iterator it = points.begin();
204 typename XYContainer::const_iterator endIt = points.end();
205 for (
int i = 0; it != endIt; ++i, ++it) {
206 xPos_[i] = std::get<0>(*it);
207 yPos_[i] = std::get<1>(*it);
231 template <
class ScalarArrayX,
class ScalarArrayY>
233 const ScalarArrayX &x,
234 const ScalarArrayY &y)
238 for (
int i = 0; i < numberSamples; ++i) {
257 template <
class ScalarContainerX,
class ScalarContainerY>
259 const ScalarContainerY &y)
261 assert(x.size() == y.size());
262 assert(x.size() > 1);
264 setNumSamples_(x.size());
265 std::copy(x.begin(), x.end(),
xPos_.begin());
266 std::copy(y.begin(), y.end(),
yPos_.begin());
283 template <
class Po
intArray>
285 const PointArray &points)
289 for (
int i = 0; i < numberSamples; ++i) {
290 xPos_[i] = points[i][0];
291 yPos_[i] = points[i][1];
307 template <
class XYContainer>
312 typename XYContainer::const_iterator it = points.begin();
313 typename XYContainer::const_iterator endIt = points.end();
314 for (
int i = 0; it != endIt; ++ i, ++it) {
337 template <
class XYContainer>
342 typename XYContainer::const_iterator it = points.begin();
343 typename XYContainer::const_iterator endIt = points.end();
344 for (
int i = 0; it != endIt; ++i, ++it) {
345 xPos_[i] = std::get<0>(*it);
346 yPos_[i] = std::get<1>(*it);
391 Scalar
x_(
int i)
const
397 Scalar
y_(
int i)
const
The common code for all 3rd order polynomial splines.
The common code for all 3rd order polynomial splines with more than two sampling points.
Definition: fixedlengthspline_.hh:44
FixedLengthSpline_()
Definition: fixedlengthspline_.hh:53
BlockVector m_
Definition: fixedlengthspline_.hh:409
void setContainerOfTuples(const XYContainer &points, Scalar m0, Scalar m1)
Set the sampling points and the boundary slopes of a full spline using a STL-compatible container of ...
Definition: fixedlengthspline_.hh:197
void makeFullSpline_(Scalar m0, Scalar m1)
Create a full spline from the already set sampling points.
Definition: fixedlengthspline_.hh:359
void setXYArrays(int numberSamples, const ScalarArrayX &x, const ScalarArrayY &y)
Set the sampling points natural spline using C-style arrays.
Definition: fixedlengthspline_.hh:232
Vector xPos_
Definition: fixedlengthspline_.hh:407
void setContainerOfPoints(const XYContainer &points, Scalar m0, Scalar m1)
Set the sampling points and the boundary slopes of a full spline using a STL-compatible container of ...
Definition: fixedlengthspline_.hh:164
void makeNaturalSpline_()
Create a natural spline from the already set sampling points.
Definition: fixedlengthspline_.hh:376
void setXYArrays(int numberSamples, const ScalarArrayX &x, const ScalarArrayY &y, Scalar m0, Scalar m1)
Set the sampling points and the boundary slopes of a full spline using C-style arrays.
Definition: fixedlengthspline_.hh:84
void setArrayOfPoints(int numberSamples, const PointArray &points, Scalar m0, Scalar m1)
Set the sampling points and the boundary slopes of a full spline using a C-style array.
Definition: fixedlengthspline_.hh:135
Vector yPos_
Definition: fixedlengthspline_.hh:408
Scalar x_(int i) const
Returns the x coordinate of the i-th sampling point.
Definition: fixedlengthspline_.hh:391
void setXYContainers(const ScalarContainerX &x, const ScalarContainerY &y)
Set the sampling points of a natural spline using STL-compatible containers.
Definition: fixedlengthspline_.hh:258
void setContainerOfPoints(const XYContainer &points)
Set the sampling points of a natural spline using a STL-compatible container of array-like objects.
Definition: fixedlengthspline_.hh:308
Scalar y_(int i) const
Returns the y coordinate of the i-th sampling point.
Definition: fixedlengthspline_.hh:397
Scalar moment_(int i) const
Returns the moment (i.e. second derivative) of the spline at the i-th sampling point.
Definition: fixedlengthspline_.hh:404
int numSamples() const
Returns the number of sampling points.
Definition: fixedlengthspline_.hh:61
void setXYContainers(const ScalarContainerX &x, const ScalarContainerY &y, Scalar m0, Scalar m1)
Set the sampling points and the boundary slopes of a full spline using STL-compatible containers.
Definition: fixedlengthspline_.hh:110
void setArrayOfPoints(int numberSamples, const PointArray &points)
Set the sampling points of a natural spline using a C-style array.
Definition: fixedlengthspline_.hh:284
void setContainerOfTuples(const XYContainer &points)
Set the sampling points of a natural spline using a STL-compatible container of tuple-like objects.
Definition: fixedlengthspline_.hh:338
The common code for all 3rd order polynomial splines.
Definition: splinecommon_.hh:44
void makeFullSystem_(Matrix &M, Vector &d, Scalar m0, Scalar m1)
Make the linear system of equations Mx = d which results in the moments of the full spline.
Definition: splinecommon_.hh:427
void makeNaturalSystem_(Matrix &M, Vector &d)
Make the linear system of equations Mx = d which results in the moments of the natural spline....
Definition: splinecommon_.hh:452