24#ifndef DUMUX_COMMON_CUBIC_SPLINE_HH
25#define DUMUX_COMMON_CUBIC_SPLINE_HH
31#include <dune/common/fvector.hh>
32#include <dune/common/fmatrix.hh>
33#include <dune/istl/btdmatrix.hh>
34#include <dune/istl/bvector.hh>
43template<
class Scalar =
double>
57 CubicSpline(
const std::vector<Scalar>& x,
const std::vector<Scalar> y)
60 assert (x.size() == y.size());
61 assert (x.size() >=2);
62 assert (std::is_sorted(x.begin(), x.end()));
73 void updatePoints(
const std::vector<Scalar>& x,
const std::vector<Scalar>& y)
79 numPoints_ = x.size();
81 const auto numSegments = numPoints_-1;
82 coeff_.resize(numSegments*4+2);
85 Dune::BTDMatrix<Dune::FieldMatrix<double, 1, 1>> matrix(numPoints_);
86 Dune::BlockVector<Dune::FieldVector<double, 1>> rhs(numPoints_);
87 Dune::BlockVector<Dune::FieldVector<double, 1>> d(numPoints_);
92 rhs[0] = 3.0*(y[1] - y[0]);
93 for (
int i = 1; i < numPoints_-1; ++i)
98 rhs[i] = 3.0*(y[i+1] - y[i-1]);
100 matrix[numPoints_-1][numPoints_-1] = 2.0;
101 matrix[numPoints_-1][numPoints_-2] = 1.0;
102 rhs[numPoints_-1] = 3.0*(y[numPoints_-1] - y[numPoints_-2]);
105 matrix.solve(d, rhs);
108 std::size_t offset = 0;
109 for (
int i = 0; i < numSegments; ++i, offset += 4)
111 coeff_[offset+0] = y[i];
112 coeff_[offset+1] = d[i];
113 coeff_[offset+2] = 3.0*(y[i+1]-y[i]) - 2.0*d[i] - d[i+1];
114 coeff_[offset+3] = 2.0*(y[i]-y[i+1]) + d[i] + d[i+1];
116 coeff_[offset+0] = y[numPoints_-1];
117 coeff_[offset+1] = d[numPoints_-1];
125 Scalar
eval(
const Scalar x)
const
129 else if (x > x_[numPoints_-1])
130 return coeff_[(numPoints_-1)*4] +
evalDerivative(x_[numPoints_-1])*(x - x_[numPoints_-1]);
133 const auto lookUpIndex =
std::distance(x_.begin(), std::lower_bound(x_.begin(), x_.end(), x));
134 assert(lookUpIndex != 0);
137 const auto* coeff = coeff_.data() + (lookUpIndex-1)*4;
140 const auto t = (x - x_[lookUpIndex-1])/(x_[lookUpIndex] - x_[lookUpIndex-1]);
141 return coeff[0] + t*(coeff[1] + t*coeff[2] + t*t*coeff[3]);
153 return coeff_[1]/(x_[1] - x_[0]);
154 else if (x > x_[numPoints_-1])
155 return coeff_[(numPoints_-1)*4 + 1]/(x_[numPoints_-1] - x_[numPoints_-2]);
158 const auto lookUpIndex =
std::distance(x_.begin(), std::lower_bound(x_.begin(), x_.end(), x));
159 assert(lookUpIndex != 0);
162 const auto* coeff = coeff_.data() + (lookUpIndex-1)*4;
165 const auto t = (x - x_[lookUpIndex-1])/(x_[lookUpIndex] - x_[lookUpIndex-1]);
166 const auto dtdx = 1.0/(x_[lookUpIndex] - x_[lookUpIndex-1]);
167 return dtdx*(coeff[1] + t*(2.0*coeff[2] + t*3.0*coeff[3]));
172 std::vector<Scalar> x_;
173 std::size_t numPoints_;
174 std::vector<Scalar> coeff_;
ctype distance(const Dune::FieldVector< ctype, dimWorld > &a, const Dune::FieldVector< ctype, dimWorld > &b)
Compute the shortest distance between two points.
Definition: distance.hh:138
A simple implementation of a natural cubic spline.
Definition: cubicspline.hh:45
CubicSpline(const std::vector< Scalar > &x, const std::vector< Scalar > y)
Contruct a natural cubic spline from the control points (x[i], y[i])
Definition: cubicspline.hh:57
CubicSpline()=default
Default constructor.
Scalar eval(const Scalar x) const
Evaluate the y value at a given x value.
Definition: cubicspline.hh:125
Scalar evalDerivative(const Scalar x) const
Evaluate the first derivative dy/dx at a given x value.
Definition: cubicspline.hh:150
void updatePoints(const std::vector< Scalar > &x, const std::vector< Scalar > &y)
Create a natural cubic spline from the control points (x[i], y[i])
Definition: cubicspline.hh:73