24#ifndef DUMUX_SPLINE_COMMON__HH
25#define DUMUX_SPLINE_COMMON__HH
31#include <dune/common/exceptions.hh>
32#include <dune/common/float_cmp.hh>
42template<
class ScalarT,
class ImplementationT>
45 using Scalar = ScalarT;
46 using Implementation = ImplementationT;
48 Implementation &asImp_()
49 {
return *
static_cast<Implementation*
>(
this); }
50 const Implementation &asImp_()
const
51 {
return *
static_cast<const Implementation*
>(
this); }
91 void printCSV(Scalar xi0, Scalar xi1,
int k)
const
95 Scalar x0 = min(xi0, xi1);
96 Scalar x1 = max(xi0, xi1);
98 for (
int i = 0; i <= k; ++i) {
99 double x = i*(x1 - x0)/k + x0;
100 double x_p1 = x + (x1 - x0)/k;
107 y = (x -
x_(0))*dy_dx +
y_(0);
108 mono = (dy_dx>0)?1:-1;
110 else if (x >
x_(n)) {
112 y = (x -
x_(n))*dy_dx +
y_(n);
113 mono = (dy_dx>0)?1:-1;
116 std::cerr <<
"ooops: " << x <<
"\n";
123 mono =
monotonic(max<Scalar>(
x_(0), x), min<Scalar>(
x_(n), x_p1));
126 std::cout << x <<
" " << y <<
" " << dy_dx <<
" " << mono <<
"\n";
141 Scalar
eval(Scalar x,
bool extrapolate=
false)
const
143 assert(extrapolate ||
applies(x));
150 return y0 + m*(x -
xMin());
152 else if (x >
xMax()) {
155 return y0 + m*(x -
xMax());
177 assert(extrapolate ||
applies(x));
179 if (Dune::FloatCmp::le(x,
xMin()))
181 else if (Dune::FloatCmp::ge(x,
xMax()))
194 Scalar
intersect(Scalar a, Scalar b, Scalar c, Scalar d)
const
196 return intersectIntervall(
xMin(),
xMax(), a, b, c, d);
206 Scalar a, Scalar b, Scalar c, Scalar d)
const
214 for (
int i = iFirst; i <= iLast; ++i)
219 DUNE_THROW(Dune::MathError,
220 "Spline has more than one intersection");
225 DUNE_THROW(Dune::MathError,
226 "Spline has no intersection");
243 assert(Dune::FloatCmp::ne<Scalar>(x0, x1));
253 if (Dune::FloatCmp::eq<Scalar>(
moment_(0), 0) &&
254 Dune::FloatCmp::eq<Scalar>(
moment_(1), 0) &&
255 Dune::FloatCmp::eq<Scalar>(
y_(0),
y_(1)))
264 if (
x_(i) <= x0 && x1 <=
x_(i+1)) {
274 for (; i < iEnd - 1; ++i)
302 template <
class DestVector,
class SourceVector>
305 const SourceVector &srcX,
306 const SourceVector &srcY,
309 assert(numSamples >= 2);
313 for (
int i = 0; i < numSamples; ++i) {
315 if (srcX[0] > srcX[numSamples - 1])
316 idx = numSamples - i - 1;
317 destX[i] = srcX[idx];
318 destY[i] = srcY[idx];
322 template <
class DestVector,
class ListIterator>
325 const ListIterator &srcBegin,
326 const ListIterator &srcEnd,
329 assert(numSamples >= 2);
332 ListIterator it = srcBegin;
334 bool reverse =
false;
335 if ((*srcBegin)[0] > (*it)[0])
340 for (
int i = 0; it != srcEnd; ++i, ++it) {
343 idx = numSamples - i - 1;
356 template <
class DestVector,
class ListIterator>
359 ListIterator srcBegin,
363 assert(numSamples >= 2);
369 ListIterator it = srcBegin;
371 bool reverse =
false;
372 if (std::get<0>(*srcBegin) > std::get<0>(*it))
377 for (
int i = 0; it != srcEnd; ++i, ++it) {
380 idx = numSamples - i - 1;
381 destX[i] = std::get<0>(*it);
382 destY[i] = std::get<1>(*it);
394 template <
class Vector,
class Matrix>
401 makeNaturalSystem(M, d);
405 Scalar lambda_n =
h_(1)/(
h_(n) +
h_(1));
406 Scalar lambda_1 = M[1][2];
407 Scalar mu_1 = 1 - lambda_1;
408 Scalar mu_n = 1 - lambda_n;
426 template <
class Vector,
class Matrix>
434 d[0] = 6/
h_(1) * ( (
y_(1) -
y_(0))/
h_(1) - m0);
443 (m1 - (
y_(n) -
y_(n - 1))/
h_(n));
451 template <
class Vector,
class Matrix>
460 for (
int i = 1; i < n; ++i) {
461 const Scalar lambda_i =
h_(i + 1) / (
h_(i) +
h_(i+1));
462 const Scalar mu_i = 1 - lambda_i;
464 6 / (
h_(i) +
h_(i+1))
470 M[i][i + 1] = lambda_i;
490 Scalar h_i1 =
h_(i + 1);
491 Scalar x_i = x -
x_(i);
492 Scalar x_i1 =
x_(i+1) - x;
495 (
y_(i+1) -
y_(i))/h_i1
498 Scalar B_i =
y_(i) -
moment_(i)* (h_i1*h_i1) / 6;
501 moment_(i)* x_i1*x_i1*x_i1 / (6 * h_i1)
503 moment_(i + 1)* x_i*x_i*x_i / (6 * h_i1)
516 Scalar h_i1 =
h_(i + 1);
517 Scalar x_i = x -
x_(i);
518 Scalar x_i1 =
x_(i+1) - x;
521 (
y_(i+1) -
y_(i))/h_i1
526 -
moment_(i) * x_i1*x_i1 / (2 * h_i1)
528 moment_(i + 1) * x_i*x_i / (2 * h_i1)
552 Scalar
disc = 4*b*b - 12*a*c;
556 return (x0*(x0*3*a + 2*b) + c > 0) ? 1 : -1;
560 Scalar xE1 = (-2*b +
disc)/(6*a);
561 Scalar xE2 = (-2*b -
disc)/(6*a);
563 if (Dune::FloatCmp::eq<Scalar>(
disc, 0)) {
565 if (Dune::FloatCmp::eq<Scalar>(xE1, x0))
570 return sign(x0*(x0*3*a + 2*b) + c);
572 if ((x0 < xE1 && xE1 < x1) ||
573 (x0 < xE2 && xE2 < x1))
581 return sign(x0*(x0*3*a + 2*b) + c);
590 Scalar a, Scalar b, Scalar c, Scalar d,
591 Scalar x0 = -1e100, Scalar x1 = 1e100)
const
599 x0 = max(
x_(segIdx), x0);
600 x1 = max(
x_(segIdx+1), x1);
604 for (
int j = 0; j < n; ++j) {
605 sol[j] +=
x_(segIdx);
606 if (x0 <= sol[j] && sol[j] <= x1) {
621 while (iLow + 1 < iHigh) {
622 int i = (iLow + iHigh) / 2;
634 Scalar
h_(
int i)
const
636 assert(
x_(i) >
x_(i-1));
638 return x_(i) -
x_(i - 1);
644 Scalar
x_(
int i)
const
645 {
return asImp_().x_(i); }
650 Scalar
y_(
int i)
const
651 {
return asImp_().y_(i); }
658 {
return asImp_().moment_(i); }
662 Scalar
a_(
int i)
const
667 Scalar
b_(
int i)
const
672 Scalar
c_(
int i)
const
677 Scalar
d_(
int i)
const
684 {
return asImp_().numSamples(); }
Define some often used mathematical functions.
int invertCubicPolynomial(SolContainer *sol, Scalar a, Scalar b, Scalar c, Scalar d)
Invert a cubic polynomial analytically.
Definition: math.hh:242
constexpr int sign(const ValueType &value) noexcept
Sign or signum function.
Definition: math.hh:627
The common code for all 3rd order polynomial splines.
Definition: splinecommon_.hh:44
Scalar evalDerivative_(Scalar x, int i) const
Definition: splinecommon_.hh:512
Scalar moment_(int i) const
Returns the moment (i.e. second derivative) of the spline at the i-th sampling point.
Definition: splinecommon_.hh:657
Scalar evalDerivative(Scalar x, bool extrapolate=false) const
Evaluate the spline's derivative at a given position.
Definition: splinecommon_.hh:175
Scalar y_(int i) const
Returns the y coordinate of the i-th sampling point.
Definition: splinecommon_.hh:650
Scalar intersectInterval(Scalar x0, Scalar x1, Scalar a, Scalar b, Scalar c, Scalar d) const
Find the intersections of the spline with a cubic polynomial in a sub-intervall of the spline,...
Definition: splinecommon_.hh:205
int intersectSegment_(Scalar *sol, int segIdx, Scalar a, Scalar b, Scalar c, Scalar d, Scalar x0=-1e100, Scalar x1=1e100) const
Find all the intersections of a segment of the spline with a cubic polynomial within a specified inte...
Definition: splinecommon_.hh:588
Scalar eval(Scalar x, bool extrapolate=false) const
Evaluate the spline at a given position.
Definition: splinecommon_.hh:141
int segmentIdx_(Scalar x) const
Definition: splinecommon_.hh:615
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 assignFromArrayList_(DestVector &destX, DestVector &destY, const ListIterator &srcBegin, const ListIterator &srcEnd, int numSamples)
Definition: splinecommon_.hh:323
void assignSamplingPoints_(DestVector &destX, DestVector &destY, const SourceVector &srcX, const SourceVector &srcY, int numSamples)
Set the sampling point vectors.
Definition: splinecommon_.hh:303
Scalar a_(int i) const
Definition: splinecommon_.hh:662
int monotonic_(int i, Scalar x0, Scalar x1) const
Definition: splinecommon_.hh:541
void assignFromTupleList_(DestVector &destX, DestVector &destY, ListIterator srcBegin, ListIterator srcEnd, int numSamples)
Set the sampling points.
Definition: splinecommon_.hh:357
Scalar xMax() const
Return the x value of the rightmost sampling point.
Definition: splinecommon_.hh:71
int numSamples_() const
Returns the number of sampling points.
Definition: splinecommon_.hh:683
Scalar xMin() const
Return the x value of the leftmost sampling point.
Definition: splinecommon_.hh:65
Scalar c_(int i) const
Definition: splinecommon_.hh:672
Scalar b_(int i) const
Definition: splinecommon_.hh:667
int monotonic() const
Same as monotonic(x0, x1), but with the entire range of the spline as interval.
Definition: splinecommon_.hh:289
int monotonic(Scalar x0, Scalar x1) const
Returns 1 if the spline is monotonically increasing, -1 if the spline is mononously decreasing and 0 ...
Definition: splinecommon_.hh:239
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
void makePeriodicSystem_(Matrix &M, Vector &d)
Make the linear system of equations Mx = d which results in the moments of the periodic spline.
Definition: splinecommon_.hh:395
bool applies(Scalar x) const
Return true if the given x is in range [x1, xn].
Definition: splinecommon_.hh:57
Scalar intersect(Scalar a, Scalar b, Scalar c, Scalar d) const
Find the intersections of the spline with a cubic polynomial in the whole intervall,...
Definition: splinecommon_.hh:194
Scalar h_(int i) const
Returns x[i] - x[i - 1].
Definition: splinecommon_.hh:634
void printCSV(Scalar xi0, Scalar xi1, int k) const
Prints k tuples of the format (x, y, dx/dy, isMonotonic) to stdout.
Definition: splinecommon_.hh:91
Scalar x_(int i) const
Returns the y coordinate of the i-th sampling point.
Definition: splinecommon_.hh:644
Scalar eval_(Scalar x, int i) const
Definition: splinecommon_.hh:486
Scalar d_(int i) const
Definition: splinecommon_.hh:677