3.3.0
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
splinecommon_.hh
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program 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 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program 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 this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_SPLINE_COMMON__HH
25#define DUMUX_SPLINE_COMMON__HH
26
27#include <algorithm>
28#include <iostream>
29#include <cassert>
30
31#include <dune/common/exceptions.hh>
32#include <dune/common/float_cmp.hh>
33
34#include "math.hh"
35
36namespace Dumux {
37
42template<class ScalarT, class ImplementationT>
44{
45 using Scalar = ScalarT;
46 using Implementation = ImplementationT;
47
48 Implementation &asImp_()
49 { return *static_cast<Implementation*>(this); }
50 const Implementation &asImp_() const
51 { return *static_cast<const Implementation*>(this); }
52
53public:
57 bool applies(Scalar x) const
58 {
59 return x_(0) <= x && x <= x_(numSamples_() - 1);
60 }
61
65 Scalar xMin() const
66 { return x_(0); }
67
71 Scalar xMax() const
72 { return x_(numSamples_() - 1); }
73
91 void printCSV(Scalar xi0, Scalar xi1, int k) const
92 {
93 using std::max;
94 using std::min;
95 Scalar x0 = min(xi0, xi1);
96 Scalar x1 = max(xi0, xi1);
97 const int n = numSamples_() - 1;
98 for (int i = 0; i <= k; ++i) {
99 double x = i*(x1 - x0)/k + x0;
100 double x_p1 = x + (x1 - x0)/k;
101 double y;
102 double dy_dx;
103 double mono = 1;
104 if (!applies(x)) {
105 if (x < 0) {
106 dy_dx = evalDerivative(x_(0));
107 y = (x - x_(0))*dy_dx + y_(0);
108 mono = (dy_dx>0)?1:-1;
109 }
110 else if (x > x_(n)) {
111 dy_dx = evalDerivative(x_(n));
112 y = (x - x_(n))*dy_dx + y_(n);
113 mono = (dy_dx>0)?1:-1;
114 }
115 else {
116 std::cerr << "ooops: " << x << "\n";
117 exit(1);
118 }
119 }
120 else {
121 y = eval(x);
122 dy_dx = evalDerivative(x);
123 mono = monotonic(max<Scalar>(x_(0), x), min<Scalar>(x_(n), x_p1));
124 }
125
126 std::cout << x << " " << y << " " << dy_dx << " " << mono << "\n";
127 }
128 }
129
141 Scalar eval(Scalar x, bool extrapolate=false) const
142 {
143 assert(extrapolate || applies(x));
144
145 // handle extrapolation
146 if (extrapolate) {
147 if (x < xMin()) {
148 Scalar m = evalDerivative_(xMin(), /*segmentIdx=*/0);
149 Scalar y0 = y_(0);
150 return y0 + m*(x - xMin());
151 }
152 else if (x > xMax()) {
153 Scalar m = evalDerivative_(xMax(), /*segmentIdx=*/numSamples_()-2);
154 Scalar y0 = y_(numSamples_()-1);
155 return y0 + m*(x - xMax());
156 }
157 }
158
159 return eval_(x, segmentIdx_(x));
160 }
161
175 Scalar evalDerivative(Scalar x, bool extrapolate=false) const
176 {
177 assert(extrapolate || applies(x));
178 if (extrapolate) {
179 if (Dune::FloatCmp::le(x, xMin()))
180 return evalDerivative_(xMin(), 0);
181 else if (Dune::FloatCmp::ge(x, xMax()))
182 return evalDerivative_(xMax(), numSamples_() - 2);
183 }
184
185 return evalDerivative_(x, segmentIdx_(x));
186 }
187
194 Scalar intersect(Scalar a, Scalar b, Scalar c, Scalar d) const
195 {
196 return intersectIntervall(xMin(), xMax(), a, b, c, d);
197 }
198
205 Scalar intersectInterval(Scalar x0, Scalar x1,
206 Scalar a, Scalar b, Scalar c, Scalar d) const
207 {
208 assert(applies(x0) && applies(x1));
209
210 Scalar tmpSol[3];
211 int nSol = 0;
212 int iFirst = segmentIdx_(x0);
213 int iLast = segmentIdx_(x1);
214 for (int i = iFirst; i <= iLast; ++i)
215 {
216 nSol += intersectSegment_(tmpSol, i, a, b, c, d, x0, x1);
217
218 if (nSol > 1) {
219 DUNE_THROW(Dune::MathError,
220 "Spline has more than one intersection"); //<<a<<"x^3 + "<<b<"x^2 + "<<c<"x + "<<d);
221 }
222 }
223
224 if (nSol != 1)
225 DUNE_THROW(Dune::MathError,
226 "Spline has no intersection"); //<<a<"x^3 + " <<b<"x^2 + "<<c<"x + "<<d<<"!");
227
228 return tmpSol[0];
229 }
230
239 int monotonic(Scalar x0, Scalar x1) const
240 {
241 assert(applies(x0));
242 assert(applies(x1));
243 assert(Dune::FloatCmp::ne<Scalar>(x0, x1));
244
245 // make sure that x0 is smaller than x1
246 if (x0 > x1)
247 {
248 using std::swap;
249 swap(x0, x1);
250 }
251
252 // corner case where the whole spline is a constant
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)))
256 {
257 // actually the is monotonically increasing as well as
258 // monotonously decreasing
259 return 2;
260 }
261
262
263 int i = segmentIdx_(x0);
264 if (x_(i) <= x0 && x1 <= x_(i+1)) {
265 // the interval to check is completely included in a
266 // single segment
267 return monotonic_(i, x0, x1);
268 }
269
270 // make sure that the segments which are completly in the
271 // interval [x0, x1] all exhibit the same monotonicity.
272 int iEnd = segmentIdx_(x1);
273 int r = monotonic_(i, x0, x_(1));
274 for (; i < iEnd - 1; ++i)
275 if (r != monotonic_(i, x_(i), x_(i + 1)))
276 return 0;
277
278 // check for the last segment
279 if (x_(iEnd) < x1 && r != monotonic_(iEnd, x_(iEnd), x1))
280 { return 0; }
281
282 return r;
283 }
284
289 int monotonic() const
290 { return monotonic(xMin(), xMax()); }
291
292protected:
293 // this is an internal class, so everything is protected!
294 SplineCommon_() = default;
295
302 template <class DestVector, class SourceVector>
303 void assignSamplingPoints_(DestVector &destX,
304 DestVector &destY,
305 const SourceVector &srcX,
306 const SourceVector &srcY,
307 int numSamples)
308 {
309 assert(numSamples >= 2);
310
311 // copy sample points, make sure that the first x value is
312 // smaller than the last one
313 for (int i = 0; i < numSamples; ++i) {
314 int idx = i;
315 if (srcX[0] > srcX[numSamples - 1])
316 idx = numSamples - i - 1;
317 destX[i] = srcX[idx];
318 destY[i] = srcY[idx];
319 }
320 }
321
322 template <class DestVector, class ListIterator>
323 void assignFromArrayList_(DestVector &destX,
324 DestVector &destY,
325 const ListIterator &srcBegin,
326 const ListIterator &srcEnd,
327 int numSamples)
328 {
329 assert(numSamples >= 2);
330
331 // find out wether the x values are in reverse order
332 ListIterator it = srcBegin;
333 ++it;
334 bool reverse = false;
335 if ((*srcBegin)[0] > (*it)[0])
336 reverse = true;
337 --it;
338
339 // loop over all sampling points
340 for (int i = 0; it != srcEnd; ++i, ++it) {
341 int idx = i;
342 if (reverse)
343 idx = numSamples - i - 1;
344 destX[i] = (*it)[0];
345 destY[i] = (*it)[1];
346 }
347 }
348
356 template <class DestVector, class ListIterator>
357 void assignFromTupleList_(DestVector &destX,
358 DestVector &destY,
359 ListIterator srcBegin,
360 ListIterator srcEnd,
361 int numSamples)
362 {
363 assert(numSamples >= 2);
364
365 // copy sample points, make sure that the first x value is
366 // smaller than the last one
367
368 // find out wether the x values are in reverse order
369 ListIterator it = srcBegin;
370 ++it;
371 bool reverse = false;
372 if (std::get<0>(*srcBegin) > std::get<0>(*it))
373 reverse = true;
374 --it;
375
376 // loop over all sampling points
377 for (int i = 0; it != srcEnd; ++i, ++it) {
378 int idx = i;
379 if (reverse)
380 idx = numSamples - i - 1;
381 destX[i] = std::get<0>(*it);
382 destY[i] = std::get<1>(*it);
383 }
384 }
385
386
394 template <class Vector, class Matrix>
395 void makePeriodicSystem_(Matrix &M, Vector &d)
396 {
397 // a periodic spline only makes sense if the first and the
398 // last sampling point have the same y value
399 assert(y_(0) == y_(numSamples_() - 1));
400
401 makeNaturalSystem(M, d);
402
403 int n = numSamples_() - 1;
404
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;
409
410 Scalar d_n =
411 6 / (h_(n) + h_(1))
412 *
413 ( (y_(1) - y_(n))/h_(1) - (y_(n) - y_(n-1))/h_(n));
414
415 // insert lambda_n and mu_1
416 M[1][n] = mu_1;
417 M[n][1] = lambda_n;
418 M[n][n-1] = mu_n;
419 d[n] = d_n;
420 }
421
426 template <class Vector, class Matrix>
427 void makeFullSystem_(Matrix &M, Vector &d, Scalar m0, Scalar m1)
428 {
429 makeNaturalSystem_(M, d);
430
431 int n = numSamples_() - 1;
432 // first row
433 M[0][1] = 1;
434 d[0] = 6/h_(1) * ( (y_(1) - y_(0))/h_(1) - m0);
435
436 // last row
437 M[n][n - 1] = 1;
438
439 // right hand side
440 d[n] =
441 6/h_(n)
442 *
443 (m1 - (y_(n) - y_(n - 1))/h_(n));
444 }
445
451 template <class Vector, class Matrix>
452 void makeNaturalSystem_(Matrix &M, Vector &d)
453 {
454 M = 0;
455 d = 0;
456
457 const int n = numSamples_() - 1;
458
459 // second to next to last rows
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;
463 const Scalar d_i =
464 6 / (h_(i) + h_(i+1))
465 *
466 ( (y_(i+1) - y_(i))/h_(i+1) - (y_(i) - y_(i-1))/h_(i));
467
468 M[i][i-1] = mu_i;
469 M[i][i] = 2;
470 M[i][i + 1] = lambda_i;
471 d[i] = d_i;
472 }
473
474 // first row
475 M[0][0] = 2;
476
477 // last row
478 M[n][n] = 2;
479 // right hand side
480 d[0] = 0.0;
481 d[n] = 0.0;
482 }
483
484 // evaluate the spline at a given the position and given the
485 // segment index
486 Scalar eval_(Scalar x, int i) const
487 {
488 // See: J. Stoer: "Numerische Mathematik 1", 9th edition,
489 // Springer, 2005, p. 109
490 Scalar h_i1 = h_(i + 1);
491 Scalar x_i = x - x_(i);
492 Scalar x_i1 = x_(i+1) - x;
493
494 Scalar A_i =
495 (y_(i+1) - y_(i))/h_i1
496 -
497 h_i1/6*(moment_(i+1) - moment_(i));
498 Scalar B_i = y_(i) - moment_(i)* (h_i1*h_i1) / 6;
499
500 return
501 moment_(i)* x_i1*x_i1*x_i1 / (6 * h_i1)
502 +
503 moment_(i + 1)* x_i*x_i*x_i / (6 * h_i1)
504 +
505 A_i*x_i
506 +
507 B_i;
508 }
509
510 // evaluate the derivative of a spline given the actual position
511 // and the segment index
512 Scalar evalDerivative_(Scalar x, int i) const
513 {
514 // See: J. Stoer: "Numerische Mathematik 1", 9th edition,
515 // Springer, 2005, p. 109
516 Scalar h_i1 = h_(i + 1);
517 Scalar x_i = x - x_(i);
518 Scalar x_i1 = x_(i+1) - x;
519
520 Scalar A_i =
521 (y_(i+1) - y_(i))/h_i1
522 -
523 h_i1/6*(moment_(i+1) - moment_(i));
524
525 return
526 -moment_(i) * x_i1*x_i1 / (2 * h_i1)
527 +
528 moment_(i + 1) * x_i*x_i / (2 * h_i1)
529 +
530 A_i;
531 }
532
533
534 // returns the monotonicality of an interval of a spline segment
535 //
536 // The return value have the following meaning:
537 //
538 // 1: spline is monotonously increasing in the specified interval
539 // 0: spline is not monotonic in the specified interval
540 // -1: spline is monotonously decreasing in the specified interval
541 int monotonic_(int i, Scalar x0, Scalar x1) const
542 {
543 // shift the interval so that it is consistent with the
544 // definitions by Stoer
545 x0 = x0 - x_(i);
546 x1 = x1 - x_(i);
547
548 Scalar a = a_(i);
549 Scalar b = b_(i);
550 Scalar c = c_(i);
551
552 Scalar disc = 4*b*b - 12*a*c;
553 if (disc < 0) {
554 // discriminant is smaller than 0, i.e. the segment does
555 // not exhibit any extrema.
556 return (x0*(x0*3*a + 2*b) + c > 0) ? 1 : -1;
557 }
558 using std::sqrt;
559 disc = sqrt(disc);
560 Scalar xE1 = (-2*b + disc)/(6*a);
561 Scalar xE2 = (-2*b - disc)/(6*a);
562
563 if (Dune::FloatCmp::eq<Scalar>(disc, 0)) {
564 // saddle point -> no extrema
565 if (Dune::FloatCmp::eq<Scalar>(xE1, x0))
566 // make sure that we're not picking the saddle point
567 // to determine whether we're monotonically increasing
568 // or decreasing
569 x0 = x1;
570 return sign(x0*(x0*3*a + 2*b) + c);
571 }
572 if ((x0 < xE1 && xE1 < x1) ||
573 (x0 < xE2 && xE2 < x1))
574 {
575 // there is an extremum in the range (x0, x1)
576 return 0;
577 }
578 // no extremum in range (x0, x1)
579 x0 = (x0 + x1)/2; // pick point in the middle of the interval
580 // to avoid extrema on the boundaries
581 return sign(x0*(x0*3*a + 2*b) + c);
582 }
583
588 int intersectSegment_(Scalar *sol,
589 int segIdx,
590 Scalar a, Scalar b, Scalar c, Scalar d,
591 Scalar x0 = -1e100, Scalar x1 = 1e100) const
592 {
593 int n = invertCubicPolynomial(sol,
594 a_(segIdx) - a,
595 b_(segIdx) - b,
596 c_(segIdx) - c,
597 d_(segIdx) - d);
598 using std::max;
599 x0 = max(x_(segIdx), x0);
600 x1 = max(x_(segIdx+1), x1);
601
602 // filter the intersections outside of the specified intervall
603 int k = 0;
604 for (int j = 0; j < n; ++j) {
605 sol[j] += x_(segIdx); // add the offset of the intervall. For details see Stoer
606 if (x0 <= sol[j] && sol[j] <= x1) {
607 sol[k] = sol[j];
608 ++k;
609 }
610 }
611 return k;
612 }
613
614 // find the segment index for a given x coordinate
615 int segmentIdx_(Scalar x) const
616 {
617 // bisection
618 int iLow = 0;
619 int iHigh = numSamples_() - 1;
620
621 while (iLow + 1 < iHigh) {
622 int i = (iLow + iHigh) / 2;
623 if (x_(i) > x)
624 iHigh = i;
625 else
626 iLow = i;
627 }
628 return iLow;
629 }
630
634 Scalar h_(int i) const
635 {
636 assert(x_(i) > x_(i-1)); // the sampling points must be given
637 // in ascending order
638 return x_(i) - x_(i - 1);
639 }
640
644 Scalar x_(int i) const
645 { return asImp_().x_(i); }
646
650 Scalar y_(int i) const
651 { return asImp_().y_(i); }
652
657 Scalar moment_(int i) const
658 { return asImp_().moment_(i); }
659
660 // returns the coefficient in front of the x^0 term. In Stoer this
661 // is delta.
662 Scalar a_(int i) const
663 { return (moment_(i+1) - moment_(i))/(6*h_(i+1)); }
664
665 // returns the coefficient in front of the x^2 term In Stoer this
666 // is gamma.
667 Scalar b_(int i) const
668 { return moment_(i)/2; }
669
670 // returns the coefficient in front of the x^1 term. In Stoer this
671 // is beta.
672 Scalar c_(int i) const
673 { return (y_(i+1) - y_(i))/h_(i + 1) - h_(i+1)/6*(2*moment_(i) + moment_(i+1)); }
674
675 // returns the coefficient in front of the x^0 term. In Stoer this
676 // is alpha.
677 Scalar d_(int i) const
678 { return y_(i); }
679
683 int numSamples_() const
684 { return asImp_().numSamples(); }
685};
686
687} // end namespace Dumux
688
689#endif
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
Definition: adapt.hh:29
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