opm-common
Loading...
Searching...
No Matches
Spline.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_HPP
28#define OPM_SPLINE_HPP
29
31
34
35#include <iosfwd>
36#include <vector>
37
38namespace Opm
39{
89template<class Scalar>
90class Spline
91{
92 typedef TridiagonalMatrix<Scalar> Matrix;
93 typedef std::vector<Scalar> Vector;
94
95public:
102 Full,
103 Natural,
104 Periodic,
105 Monotonic
106 };
107
114 { }
115
126 Spline(Scalar x0, Scalar x1,
127 Scalar y0, Scalar y1,
128 Scalar m0, Scalar m1)
129 { set(x0, x1, y0, y1, m0, m1); }
130
140 template <class ScalarArrayX, class ScalarArrayY>
141 Spline(size_t nSamples,
142 const ScalarArrayX& x,
143 const ScalarArrayY& y,
144 SplineType splineType = Natural,
145 bool sortInputs = true)
146 { this->setXYArrays(nSamples, x, y, splineType, sortInputs); }
147
156 template <class PointArray>
157 Spline(size_t nSamples,
158 const PointArray& points,
159 SplineType splineType = Natural,
160 bool sortInputs = true)
161 { this->setArrayOfPoints(nSamples, points, splineType, sortInputs); }
162
171 template <class ScalarContainer>
172 Spline(const ScalarContainer& x,
173 const ScalarContainer& y,
174 SplineType splineType = Natural,
175 bool sortInputs = true)
176 { this->setXYContainers(x, y, splineType, sortInputs); }
177
185 template <class PointContainer>
186 explicit Spline(const PointContainer& points,
187 SplineType splineType = Natural,
188 bool sortInputs = true)
189 { this->setContainerOfPoints(points, splineType, sortInputs); }
190
201 template <class ScalarArray>
202 Spline(size_t nSamples,
203 const ScalarArray& x,
204 const ScalarArray& y,
205 Scalar m0,
206 Scalar m1,
207 bool sortInputs = true)
208 { this->setXYArrays(nSamples, x, y, m0, m1, sortInputs); }
209
219 template <class PointArray>
220 Spline(size_t nSamples,
221 const PointArray& points,
222 Scalar m0,
223 Scalar m1,
224 bool sortInputs = true)
225 { this->setArrayOfPoints(nSamples, points, m0, m1, sortInputs); }
226
236 template <class ScalarContainerX, class ScalarContainerY>
237 Spline(const ScalarContainerX& x,
238 const ScalarContainerY& y,
239 Scalar m0,
240 Scalar m1,
241 bool sortInputs = true)
242 { this->setXYContainers(x, y, m0, m1, sortInputs); }
243
252 template <class PointContainer>
253 Spline(const PointContainer& points,
254 Scalar m0,
255 Scalar m1,
256 bool sortInputs = true)
257 { this->setContainerOfPoints(points, m0, m1, sortInputs); }
258
270 void set(Scalar x0, Scalar x1,
271 Scalar y0, Scalar y1,
272 Scalar m0, Scalar m1)
273 {
274 slopeVec_.resize(2);
275 xPos_.resize(2);
276 yPos_.resize(2);
277
278 if (x0 > x1) {
279 xPos_[0] = x1;
280 xPos_[1] = x0;
281 yPos_[0] = y1;
282 yPos_[1] = y0;
283 }
284 else {
285 xPos_[0] = x0;
286 xPos_[1] = x1;
287 yPos_[0] = y0;
288 yPos_[1] = y1;
289 }
290
291 slopeVec_[0] = m0;
292 slopeVec_[1] = m1;
293
294 Matrix M(numSamples());
295 Vector d(numSamples());
296 Vector moments(numSamples());
297 this->makeFullSystem_(M, d, m0, m1);
298
299 // solve for the moments
300 M.solve(moments, d);
301
302 this->setSlopesFromMoments_(slopeVec_, moments);
303 }
304
305
309 // Full splines //
313
325 template <class ScalarArrayX, class ScalarArrayY>
326 void setXYArrays(size_t nSamples,
327 const ScalarArrayX& x,
328 const ScalarArrayY& y,
329 Scalar m0, Scalar m1,
330 bool sortInputs = true)
331 {
332 assert(nSamples > 1);
333
334 setNumSamples_(nSamples);
335 for (size_t i = 0; i < nSamples; ++i) {
336 xPos_[i] = x[i];
337 yPos_[i] = y[i];
338 }
339
340 if (sortInputs)
341 sortInput_();
342 else if (xPos_[0] > xPos_[numSamples() - 1])
344
345 makeFullSpline_(m0, m1);
346 }
347
359 template <class ScalarContainerX, class ScalarContainerY>
360 void setXYContainers(const ScalarContainerX& x,
361 const ScalarContainerY& y,
362 Scalar m0, Scalar m1,
363 bool sortInputs = true)
364 {
365 assert(x.size() == y.size());
366 assert(x.size() > 1);
367
368 setNumSamples_(x.size());
369
370 std::ranges::copy(x, xPos_.begin());
371 std::ranges::copy(y, yPos_.begin());
372
373 if (sortInputs)
374 sortInput_();
375 else if (xPos_[0] > xPos_[numSamples() - 1])
377
378 makeFullSpline_(m0, m1);
379 }
380
393 template <class PointArray>
394 void setArrayOfPoints(size_t nSamples,
395 const PointArray& points,
396 Scalar m0,
397 Scalar m1,
398 bool sortInputs = true)
399 {
400 // a spline with no or just one sampling points? what an
401 // incredible bad idea!
402 assert(nSamples > 1);
403
404 setNumSamples_(nSamples);
405 for (size_t i = 0; i < nSamples; ++i) {
406 xPos_[i] = points[i][0];
407 yPos_[i] = points[i][1];
408 }
409
410 if (sortInputs)
411 sortInput_();
412 else if (xPos_[0] > xPos_[numSamples() - 1])
414
415 makeFullSpline_(m0, m1);
416 }
417
430 template <class XYContainer>
431 void setContainerOfPoints(const XYContainer& points,
432 Scalar m0,
433 Scalar m1,
434 bool sortInputs = true)
435 {
436 // a spline with no or just one sampling points? what an
437 // incredible bad idea!
438 assert(points.size() > 1);
439
440 setNumSamples_(points.size());
441 typename XYContainer::const_iterator it = points.begin();
442 typename XYContainer::const_iterator endIt = points.end();
443 for (size_t i = 0; it != endIt; ++i, ++it) {
444 xPos_[i] = (*it)[0];
445 yPos_[i] = (*it)[1];
446 }
447
448 if (sortInputs)
449 sortInput_();
450 else if (xPos_[0] > xPos_[numSamples() - 1])
452
453 // make a full spline
454 makeFullSpline_(m0, m1);
455 }
456
472 template <class XYContainer>
473 void setContainerOfTuples(const XYContainer& points,
474 Scalar m0,
475 Scalar m1,
476 bool sortInputs = true)
477 {
478 // resize internal arrays
479 setNumSamples_(points.size());
480 typename XYContainer::const_iterator it = points.begin();
481 typename XYContainer::const_iterator endIt = points.end();
482 for (unsigned i = 0; it != endIt; ++i, ++it) {
483 xPos_[i] = std::get<0>(*it);
484 yPos_[i] = std::get<1>(*it);
485 }
486
487 if (sortInputs)
488 sortInput_();
489 else if (xPos_[0] > xPos_[numSamples() - 1])
491
492 // make a full spline
493 makeFullSpline_(m0, m1);
494 }
495
499 // Natural/Periodic splines //
503
513 template <class ScalarArrayX, class ScalarArrayY>
514 void setXYArrays(size_t nSamples,
515 const ScalarArrayX& x,
516 const ScalarArrayY& y,
517 SplineType splineType = Natural,
518 bool sortInputs = true)
519 {
520 assert(nSamples > 1);
521
522 setNumSamples_(nSamples);
523 for (size_t i = 0; i < nSamples; ++i) {
524 xPos_[i] = x[i];
525 yPos_[i] = y[i];
526 }
527
528 if (sortInputs)
529 sortInput_();
530 else if (xPos_[0] > xPos_[numSamples() - 1])
532
533 if (splineType == Periodic)
535 else if (splineType == Natural)
537 else if (splineType == Monotonic)
538 this->makeMonotonicSpline_(slopeVec_);
539 else
540 throw std::runtime_error("Spline type "+std::to_string(int(splineType))+" not supported at this place");
541 }
542
554 template <class ScalarContainerX, class ScalarContainerY>
555 void setXYContainers(const ScalarContainerX& x,
556 const ScalarContainerY& y,
557 SplineType splineType = Natural,
558 bool sortInputs = true)
559 {
560 assert(x.size() == y.size());
561 assert(x.size() > 1);
562
563 setNumSamples_(x.size());
564 std::ranges::copy(x, xPos_.begin());
565 std::ranges::copy(y, yPos_.begin());
566
567 if (sortInputs)
568 sortInput_();
569 else if (xPos_[0] > xPos_[numSamples() - 1])
571
572 if (splineType == Periodic)
574 else if (splineType == Natural)
576 else if (splineType == Monotonic)
577 this->makeMonotonicSpline_(slopeVec_);
578 else
579 throw std::runtime_error("Spline type "+std::to_string(int(splineType))+" not supported at this place");
580 }
581
594 template <class PointArray>
595 void setArrayOfPoints(size_t nSamples,
596 const PointArray& points,
597 SplineType splineType = Natural,
598 bool sortInputs = true)
599 {
600 // a spline with no or just one sampling points? what an
601 // incredible bad idea!
602 assert(nSamples > 1);
603
604 setNumSamples_(nSamples);
605 for (size_t i = 0; i < nSamples; ++i) {
606 xPos_[i] = points[i][0];
607 yPos_[i] = points[i][1];
608 }
609
610 if (sortInputs)
611 sortInput_();
612 else if (xPos_[0] > xPos_[numSamples() - 1])
614
615 if (splineType == Periodic)
617 else if (splineType == Natural)
619 else if (splineType == Monotonic)
620 this->makeMonotonicSpline_(slopeVec_);
621 else
622 throw std::runtime_error("Spline type "+std::to_string(int(splineType))+" not supported at this place");
623 }
624
636 template <class XYContainer>
637 void setContainerOfPoints(const XYContainer& points,
638 SplineType splineType = Natural,
639 bool sortInputs = true)
640 {
641 // a spline with no or just one sampling points? what an
642 // incredible bad idea!
643 assert(points.size() > 1);
644
645 setNumSamples_(points.size());
646 typename XYContainer::const_iterator it = points.begin();
647 typename XYContainer::const_iterator endIt = points.end();
648 for (size_t i = 0; it != endIt; ++ i, ++it) {
649 xPos_[i] = (*it)[0];
650 yPos_[i] = (*it)[1];
651 }
652
653 if (sortInputs)
654 sortInput_();
655 else if (xPos_[0] > xPos_[numSamples() - 1])
657
658 if (splineType == Periodic)
660 else if (splineType == Natural)
662 else if (splineType == Monotonic)
663 this->makeMonotonicSpline_(slopeVec_);
664 else
665 throw std::runtime_error("Spline type "+std::to_string(int(splineType))+" not supported at this place");
666 }
667
682 template <class XYContainer>
683 void setContainerOfTuples(const XYContainer& points,
684 SplineType splineType = Natural,
685 bool sortInputs = true)
686 {
687 // resize internal arrays
688 setNumSamples_(points.size());
689 typename XYContainer::const_iterator it = points.begin();
690 typename XYContainer::const_iterator endIt = points.end();
691 for (unsigned i = 0; it != endIt; ++i, ++it) {
692 xPos_[i] = std::get<0>(*it);
693 yPos_[i] = std::get<1>(*it);
694 }
695
696 if (sortInputs)
697 sortInput_();
698 else if (xPos_[0] > xPos_[numSamples() - 1])
700
701 if (splineType == Periodic)
703 else if (splineType == Natural)
705 else if (splineType == Monotonic)
706 this->makeMonotonicSpline_(slopeVec_);
707 else
708 throw std::runtime_error("Spline type "+std::to_string(int(splineType))+" not supported at this place");
709 }
710
714 template <class Evaluation>
715 bool applies(const Evaluation& x) const
716 { return x_(0) <= x && x <= x_(numSamples() - 1); }
717
721 size_t numSamples() const
722 { return xPos_.size(); }
723
727 Scalar xAt(size_t sampleIdx) const
728 { return x_(sampleIdx); }
729
733 Scalar valueAt(size_t sampleIdx) const
734 { return y_(sampleIdx); }
735
753 void printCSV(Scalar xi0, Scalar xi1, size_t k, std::ostream& os) const;
754
766 template <class Evaluation>
767 Evaluation eval(const Evaluation& x, bool extrapolate = false) const
768 {
769 if (!extrapolate && !applies(x))
770 throw NumericalProblem("Tried to evaluate a spline outside of its range");
771
772 // handle extrapolation
773 if (extrapolate) {
774 if (x < xAt(0)) {
775 Scalar m = evalDerivative_(xAt(0), /*segmentIdx=*/0);
776 Scalar y0 = y_(0);
777 return y0 + m*(x - xAt(0));
778 }
779 else if (x > xAt(static_cast<size_t>(static_cast<long int>(numSamples()) - 1))) {
780 Scalar m = evalDerivative_(xAt(static_cast<size_t>(numSamples() - 1)),
781 /*segmentIdx=*/static_cast<size_t>(numSamples()-2));
782 Scalar y0 = y_(static_cast<size_t>(numSamples() - 1));
783 return y0 + m*(x - xAt(static_cast<size_t>(numSamples() - 1)));
784 }
785 }
786
787 return eval_(x, segmentIdx_(scalarValue(x)));
788 }
789
802 template <class Evaluation>
803 Evaluation evalDerivative(const Evaluation& x, bool extrapolate = false) const
804 {
805 if (!extrapolate && !applies(x))
806 throw NumericalProblem("Tried to evaluate the derivative of a spline outside of its range");
807
808 // handle extrapolation
809 if (extrapolate) {
810 if (x < xAt(0))
811 return evalDerivative_(xAt(0), /*segmentIdx=*/0);
812 else if (x > xAt(numSamples() - 1))
813 return evalDerivative_(xAt(numSamples() - 1), /*segmentIdx=*/numSamples() - 2);
814 }
815
816 return evalDerivative_(x, segmentIdx_(scalarValue(x)));
817 }
818
831 template <class Evaluation>
832 Evaluation evalSecondDerivative(const Evaluation& x, bool extrapolate = false) const
833 {
834 if (!extrapolate && !applies(x))
835 throw NumericalProblem("Tried to evaluate the second derivative of a spline outside of its range");
836 else if (extrapolate)
837 return 0.0;
838
839 return evalDerivative2_(x, segmentIdx_(scalarValue(x)));
840 }
841
854 template <class Evaluation>
855 Evaluation evalThirdDerivative(const Evaluation& x, bool extrapolate = false) const
856 {
857 if (!extrapolate && !applies(x))
858 throw NumericalProblem("Tried to evaluate the third derivative of a spline outside of its range");
859 else if (extrapolate)
860 return 0.0;
861
862 return evalDerivative3_(x, segmentIdx_(scalarValue(x)));
863 }
864
871 template <class Evaluation>
872 Evaluation intersect(const Evaluation& a,
873 const Evaluation& b,
874 const Evaluation& c,
875 const Evaluation& d) const
876 {
877 return intersectInterval(xAt(0), xAt(numSamples() - 1), a, b, c, d);
878 }
879
886 template <class Evaluation>
887 Evaluation intersectInterval(Scalar x0, Scalar x1,
888 const Evaluation& a,
889 const Evaluation& b,
890 const Evaluation& c,
891 const Evaluation& d) const
892 {
893 assert(applies(x0) && applies(x1));
894
895 Evaluation tmpSol[3], sol = 0;
896 size_t nSol = 0;
897 size_t iFirst = segmentIdx_(x0);
898 size_t iLast = segmentIdx_(x1);
899 for (size_t i = iFirst; i <= iLast; ++i)
900 {
901 size_t nCur = intersectSegment_(tmpSol, i, a, b, c, d, x0, x1);
902 if (nCur == 1)
903 sol = tmpSol[0];
904
905 nSol += nCur;
906 if (nSol > 1) {
907 throw std::runtime_error("Spline has more than one intersection"); //<<a<<"x^3 + "<<b<"x^2 + "<<c<"x + "<<d);
908 }
909 }
910
911 if (nSol != 1)
912 throw std::runtime_error("Spline has no intersection"); //<<a<"x^3 + " <<b<"x^2 + "<<c<"x + "<<d<<"!");
913
914 return sol;
915 }
916
925 int monotonic(Scalar x0, Scalar x1,
926 [[maybe_unused]] bool extrapolate = false) const
927 {
928 assert(std::abs(x0 - x1) > 1e-30);
929
930 // make sure that x0 is smaller than x1
931 if (x0 > x1)
932 std::swap(x0, x1);
933
934 assert(x0 < x1);
935
936 int r = 3;
937 if (x0 < xAt(0)) {
938 assert(extrapolate);
939 Scalar m = evalDerivative_(xAt(0), /*segmentIdx=*/0);
940 if (std::abs(m) < 1e-20)
941 r = (m < 0)?-1:1;
942 x0 = xAt(0);
943 };
944
945 size_t i = segmentIdx_(x0);
946 if (x_(i + 1) >= x1) {
947 // interval is fully contained within a single spline
948 // segment
949 monotonic_(i, x0, x1, r);
950 return r;
951 }
952
953 // the first segment overlaps with the specified interval
954 // partially
955 monotonic_(i, x0, x_(i+1), r);
956 ++ i;
957
958 // make sure that the segments which are completly in the
959 // interval [x0, x1] all exhibit the same monotonicity.
960 size_t iEnd = segmentIdx_(x1);
961 for (; i < iEnd - 1; ++i) {
962 monotonic_(i, x_(i), x_(i + 1), r);
963 if (!r)
964 return 0;
965 }
966
967 // if the user asked for a part of the spline which is
968 // extrapolated, we need to check the slope at the spline's
969 // endpoint
970 if (x1 > xAt(numSamples() - 1)) {
971 assert(extrapolate);
972
973 Scalar m = evalDerivative_(xAt(numSamples() - 1), /*segmentIdx=*/numSamples() - 2);
974 if (m < 0)
975 return (r < 0 || r==3) ? -1 : 0;
976 else if (m > 0)
977 return r > 0 ? 1 : 0;
978
979 return r;
980 }
981
982 // check for the last segment
983 monotonic_(iEnd, x_(iEnd), x1, r);
984
985 return r;
986 }
987
992 int monotonic() const
993 { return monotonic(xAt(0), xAt(numSamples() - 1)); }
994
995protected:
999 struct ComparatorX_
1000 {
1001 explicit ComparatorX_(const std::vector<Scalar>& x)
1002 : x_(x)
1003 {}
1004
1005 bool operator ()(unsigned idxA, unsigned idxB) const
1006 { return x_.at(idxA) < x_.at(idxB); }
1007
1008 const std::vector<Scalar>& x_;
1009 };
1010
1015 {
1016 size_t n = numSamples();
1017
1018 // create a vector containing 0...n-1
1019 std::vector<unsigned> idxVector(n);
1020 for (unsigned i = 0; i < n; ++i)
1021 idxVector[i] = i;
1022
1023 // sort the indices according to the x values of the sample
1024 // points
1025 ComparatorX_ cmp(xPos_);
1026 std::ranges::sort(idxVector, cmp);
1027
1028 // reorder the sample points
1029 std::vector<Scalar> tmpX(n), tmpY(n);
1030 for (size_t i = 0; i < idxVector.size(); ++ i) {
1031 tmpX[i] = xPos_[idxVector[i]];
1032 tmpY[i] = yPos_[idxVector[i]];
1033 }
1034 xPos_ = tmpX;
1035 yPos_ = tmpY;
1036 }
1037
1043 {
1044 // reverse the arrays
1045 size_t n = numSamples();
1046 for (unsigned i = 0; i <= (n - 1)/2; ++i) {
1047 std::swap(xPos_[i], xPos_[n - i - 1]);
1048 std::swap(yPos_[i], yPos_[n - i - 1]);
1049 }
1050 }
1051
1055 void setNumSamples_(size_t nSamples)
1056 {
1057 xPos_.resize(nSamples);
1058 yPos_.resize(nSamples);
1059 slopeVec_.resize(nSamples);
1060 }
1061
1067 void makeFullSpline_(Scalar m0, Scalar m1)
1068 {
1069 Matrix M(numSamples());
1070 std::vector<Scalar> d(numSamples());
1071 std::vector<Scalar> moments(numSamples());
1072
1073 // create linear system of equations
1074 this->makeFullSystem_(M, d, m0, m1);
1075
1076 // solve for the moments (-> second derivatives)
1077 M.solve(moments, d);
1078
1079 // convert the moments to slopes at the sample points
1080 this->setSlopesFromMoments_(slopeVec_, moments);
1081 }
1082
1089 {
1090 Matrix M(numSamples(), numSamples());
1091 Vector d(numSamples());
1092 Vector moments(numSamples());
1093
1094 // create linear system of equations
1095 this->makeNaturalSystem_(M, d);
1096
1097 // solve for the moments (-> second derivatives)
1098 M.solve(moments, d);
1099
1100 // convert the moments to slopes at the sample points
1101 this->setSlopesFromMoments_(slopeVec_, moments);
1102 }
1103
1110 {
1111 Matrix M(numSamples() - 1);
1112 Vector d(numSamples() - 1);
1113 Vector moments(numSamples() - 1);
1114
1115 // create linear system of equations. This is a bit hacky,
1116 // because it assumes that std::vector internally stores its
1117 // data as a big C-style array, but it saves us from yet
1118 // another copy operation
1119 this->makePeriodicSystem_(M, d);
1120
1121 // solve for the moments (-> second derivatives)
1122 M.solve(moments, d);
1123
1124 moments.resize(numSamples());
1125 for (int i = static_cast<int>(numSamples()) - 2; i >= 0; --i) {
1126 unsigned ui = static_cast<unsigned>(i);
1127 moments[ui+1] = moments[ui];
1128 }
1129 moments[0] = moments[numSamples() - 1];
1130
1131 // convert the moments to slopes at the sample points
1132 this->setSlopesFromMoments_(slopeVec_, moments);
1133 }
1134
1141 template <class DestVector, class SourceVector>
1142 void assignSamplingPoints_(DestVector& destX,
1143 DestVector& destY,
1144 const SourceVector& srcX,
1145 const SourceVector& srcY,
1146 unsigned nSamples)
1147 {
1148 assert(nSamples >= 2);
1149
1150 // copy sample points, make sure that the first x value is
1151 // smaller than the last one
1152 for (unsigned i = 0; i < nSamples; ++i) {
1153 unsigned idx = i;
1154 if (srcX[0] > srcX[nSamples - 1])
1155 idx = nSamples - i - 1;
1156 destX[i] = srcX[idx];
1157 destY[i] = srcY[idx];
1158 }
1159 }
1160
1161 template <class DestVector, class ListIterator>
1162 void assignFromArrayList_(DestVector& destX,
1163 DestVector& destY,
1164 const ListIterator& srcBegin,
1165 const ListIterator& srcEnd,
1166 unsigned nSamples)
1167 {
1168 assert(nSamples >= 2);
1169
1170 // find out wether the x values are in reverse order
1171 ListIterator it = srcBegin;
1172 ++it;
1173 bool reverse = false;
1174 if ((*srcBegin)[0] > (*it)[0])
1175 reverse = true;
1176 --it;
1177
1178 // loop over all sampling points
1179 for (unsigned i = 0; it != srcEnd; ++i, ++it) {
1180 unsigned idx = i;
1181 if (reverse)
1182 idx = nSamples - i - 1;
1183 destX[idx] = (*it)[0];
1184 destY[idx] = (*it)[1];
1185 }
1186 }
1187
1195 template <class DestVector, class ListIterator>
1196 void assignFromTupleList_(DestVector& destX,
1197 DestVector& destY,
1198 ListIterator srcBegin,
1199 ListIterator srcEnd,
1200 unsigned nSamples)
1201 {
1202 assert(nSamples >= 2);
1203
1204 // copy sample points, make sure that the first x value is
1205 // smaller than the last one
1206
1207 // find out wether the x values are in reverse order
1208 ListIterator it = srcBegin;
1209 ++it;
1210 bool reverse = false;
1211 if (std::get<0>(*srcBegin) > std::get<0>(*it))
1212 reverse = true;
1213 --it;
1214
1215 // loop over all sampling points
1216 for (unsigned i = 0; it != srcEnd; ++i, ++it) {
1217 unsigned idx = i;
1218 if (reverse)
1219 idx = nSamples - i - 1;
1220 destX[idx] = std::get<0>(*it);
1221 destY[idx] = std::get<1>(*it);
1222 }
1223 }
1224
1229 template <class Vector, class Matrix>
1230 void makeFullSystem_(Matrix& M, Vector& d, Scalar m0, Scalar m1)
1231 {
1232 makeNaturalSystem_(M, d);
1233
1234 size_t n = numSamples() - 1;
1235 // first row
1236 M[0][1] = 1;
1237 d[0] = 6/h_(1) * ( (y_(1) - y_(0))/h_(1) - m0);
1238
1239 // last row
1240 M[n][n - 1] = 1;
1241
1242 // right hand side
1243 d[n] =
1244 6/h_(n)
1245 *
1246 (m1 - (y_(n) - y_(n - 1))/h_(n));
1247 }
1248
1253 template <class Vector, class Matrix>
1254 void makeNaturalSystem_(Matrix& M, Vector& d)
1255 {
1256 M = 0.0;
1257
1258 // See: J. Stoer: "Numerische Mathematik 1", 9th edition,
1259 // Springer, 2005, p. 111
1260 size_t n = numSamples() - 1;
1261
1262 // second to next to last rows
1263 for (size_t i = 1; i < n; ++i) {
1264 Scalar lambda_i = h_(i + 1) / (h_(i) + h_(i + 1));
1265 Scalar mu_i = 1 - lambda_i;
1266 Scalar d_i =
1267 6 / (h_(i) + h_(i + 1))
1268 *
1269 ( (y_(i + 1) - y_(i))/h_(i + 1) - (y_(i) - y_(i - 1))/h_(i));
1270
1271 M[i][i-1] = mu_i;
1272 M[i][i] = 2;
1273 M[i][i + 1] = lambda_i;
1274 d[i] = d_i;
1275 };
1276
1277 // See Stroer, equation (2.5.2.7)
1278 Scalar lambda_0 = 0;
1279 Scalar d_0 = 0;
1280
1281 Scalar mu_n = 0;
1282 Scalar d_n = 0;
1283
1284 // first row
1285 M[0][0] = 2;
1286 M[0][1] = lambda_0;
1287 d[0] = d_0;
1288
1289 // last row
1290 M[n][n-1] = mu_n;
1291 M[n][n] = 2;
1292 d[n] = d_n;
1293 }
1294
1299 template <class Matrix, class Vector>
1300 void makePeriodicSystem_(Matrix& M, Vector& d)
1301 {
1302 M = 0.0;
1303
1304 // See: J. Stoer: "Numerische Mathematik 1", 9th edition,
1305 // Springer, 2005, p. 111
1306 size_t n = numSamples() - 1;
1307
1308 assert(M.rows() == n);
1309
1310 // second to next to last rows
1311 for (size_t i = 2; i < n; ++i) {
1312 Scalar lambda_i = h_(i + 1) / (h_(i) + h_(i + 1));
1313 Scalar mu_i = 1 - lambda_i;
1314 Scalar d_i =
1315 6 / (h_(i) + h_(i + 1))
1316 *
1317 ( (y_(i + 1) - y_(i))/h_(i + 1) - (y_(i) - y_(i - 1))/h_(i));
1318
1319 M[i-1][i-2] = mu_i;
1320 M[i-1][i-1] = 2;
1321 M[i-1][i] = lambda_i;
1322 d[i-1] = d_i;
1323 };
1324
1325 Scalar lambda_n = h_(1) / (h_(n) + h_(1));
1326 Scalar lambda_1 = h_(2) / (h_(1) + h_(2));;
1327 Scalar mu_1 = 1 - lambda_1;
1328 Scalar mu_n = 1 - lambda_n;
1329
1330 Scalar d_1 =
1331 6 / (h_(1) + h_(2))
1332 *
1333 ( (y_(2) - y_(1))/h_(2) - (y_(1) - y_(0))/h_(1));
1334 Scalar d_n =
1335 6 / (h_(n) + h_(1))
1336 *
1337 ( (y_(1) - y_(n))/h_(1) - (y_(n) - y_(n-1))/h_(n));
1338
1339
1340 // first row
1341 M[0][0] = 2;
1342 M[0][1] = lambda_1;
1343 M[0][n-1] = mu_1;
1344 d[0] = d_1;
1345
1346 // last row
1347 M[n-1][0] = lambda_n;
1348 M[n-1][n-2] = mu_n;
1349 M[n-1][n-1] = 2;
1350 d[n-1] = d_n;
1351 }
1352
1361 template <class Vector>
1362 void makeMonotonicSpline_(Vector& slopes)
1363 {
1364 auto n = numSamples();
1365
1366 // calculate the slopes of the secant lines
1367 std::vector<Scalar> delta(n);
1368 for (size_t k = 0; k < n - 1; ++k)
1369 delta[k] = (y_(k + 1) - y_(k))/(x_(k + 1) - x_(k));
1370
1371 // calculate the "raw" slopes at the sample points
1372 for (size_t k = 1; k < n - 1; ++k)
1373 slopes[k] = (delta[k - 1] + delta[k])/2;
1374 slopes[0] = delta[0];
1375 slopes[n - 1] = delta[n - 2];
1376
1377 // post-process the "raw" slopes at the sample points
1378 for (size_t k = 0; k < n - 1; ++k) {
1379 if (std::abs(delta[k]) < 1e-50) {
1380 // make the spline flat if the inputs are equal
1381 slopes[k] = 0;
1382 slopes[k + 1] = 0;
1383 ++ k;
1384 continue;
1385 }
1386 else {
1387 Scalar alpha = slopes[k] / delta[k];
1388 Scalar beta = slopes[k + 1] / delta[k];
1389
1390 if (alpha < 0 || (k > 0 && slopes[k] / delta[k - 1] < 0)) {
1391 slopes[k] = 0;
1392 }
1393 // limit (alpha, beta) to a circle of radius 3
1394 else if (alpha*alpha + beta*beta > 3*3) {
1395 Scalar tau = 3.0/std::sqrt(alpha*alpha + beta*beta);
1396 slopes[k] = tau*alpha*delta[k];
1397 slopes[k + 1] = tau*beta*delta[k];
1398 }
1399 }
1400 }
1401 }
1402
1410 template <class MomentsVector, class SlopeVector>
1411 void setSlopesFromMoments_(SlopeVector& slopes, const MomentsVector& moments)
1412 {
1413 size_t n = numSamples();
1414
1415 // evaluate slope at the rightmost point.
1416 // See: J. Stoer: "Numerische Mathematik 1", 9th edition,
1417 // Springer, 2005, p. 109
1418 Scalar mRight;
1419
1420 {
1421 Scalar h = this->h_(n - 1);
1422 Scalar x = h;
1423 //Scalar x_1 = 0;
1424
1425 Scalar A =
1426 (y_(n - 1) - y_(n - 2))/h
1427 -
1428 h/6*(moments[n-1] - moments[n - 2]);
1429
1430 mRight =
1431 //- moments[n - 2] * x_1*x_1 / (2 * h)
1432 //+
1433 moments[n - 1] * x*x / (2 * h)
1434 +
1435 A;
1436 }
1437
1438 // evaluate the slope for the first n-1 sample points
1439 for (size_t i = 0; i < n - 1; ++ i) {
1440 // See: J. Stoer: "Numerische Mathematik 1", 9th edition,
1441 // Springer, 2005, p. 109
1442 Scalar h_i = this->h_(i + 1);
1443 //Scalar x_i = 0;
1444 Scalar x_i1 = h_i;
1445
1446 Scalar A_i =
1447 (y_(i+1) - y_(i))/h_i
1448 -
1449 h_i/6*(moments[i+1] - moments[i]);
1450
1451 slopes[i] =
1452 - moments[i] * x_i1*x_i1 / (2 * h_i)
1453 +
1454 //moments[i + 1] * x_i*x_i / (2 * h_i)
1455 //+
1456 A_i;
1457
1458 }
1459 slopes[n - 1] = mRight;
1460 }
1461
1462
1463 // evaluate the spline at a given the position and given the
1464 // segment index
1465 template <class Evaluation>
1466 Evaluation eval_(const Evaluation& x, size_t i) const
1467 {
1468 // See http://en.wikipedia.org/wiki/Cubic_Hermite_spline
1469 Scalar delta = h_(i + 1);
1470 Evaluation t = (x - x_(i))/delta;
1471
1472 return
1473 h00_(t) * y_(i)
1474 + h10_(t) * slope_(i)*delta
1475 + h01_(t) * y_(i + 1)
1476 + h11_(t) * slope_(i + 1)*delta;
1477 }
1478
1479 // evaluate the derivative of a spline given the actual position
1480 // and the segment index
1481 template <class Evaluation>
1482 Evaluation evalDerivative_(const Evaluation& x, size_t i) const
1483 {
1484 // See http://en.wikipedia.org/wiki/Cubic_Hermite_spline
1485 Scalar delta = h_(i + 1);
1486 Evaluation t = (x - x_(i))/delta;
1487 Evaluation alpha = 1 / delta;
1488
1489 return
1490 alpha *
1491 (h00_prime_(t) * y_(i)
1492 + h10_prime_(t) * slope_(i)*delta
1493 + h01_prime_(t) * y_(i + 1)
1494 + h11_prime_(t) * slope_(i + 1)*delta);
1495 }
1496
1497 // evaluate the second derivative of a spline given the actual
1498 // position and the segment index
1499 template <class Evaluation>
1500 Evaluation evalDerivative2_(const Evaluation& x, size_t i) const
1501 {
1502 // See http://en.wikipedia.org/wiki/Cubic_Hermite_spline
1503 Scalar delta = h_(i + 1);
1504 Evaluation t = (x - x_(i))/delta;
1505 Evaluation alpha = 1 / delta;
1506
1507 return
1508 alpha*alpha
1509 *(h00_prime2_(t) * y_(i)
1510 + h10_prime2_(t) * slope_(i)*delta
1511 + h01_prime2_(t) * y_(i + 1)
1512 + h11_prime2_(t) * slope_(i + 1)*delta);
1513 }
1514
1515 // evaluate the third derivative of a spline given the actual
1516 // position and the segment index
1517 template <class Evaluation>
1518 Evaluation evalDerivative3_(const Evaluation& x, size_t i) const
1519 {
1520 // See http://en.wikipedia.org/wiki/Cubic_Hermite_spline
1521 Scalar delta = h_(i + 1);
1522 Evaluation t = (x - x_(i))/delta;
1523 Evaluation alpha = 1 / delta;
1524
1525 return
1526 alpha*alpha*alpha
1527 *(h00_prime3_(t)*y_(i)
1528 + h10_prime3_(t)*slope_(i)*delta
1529 + h01_prime3_(t)*y_(i + 1)
1530 + h11_prime3_(t)*slope_(i + 1)*delta);
1531 }
1532
1533 // hermite basis functions
1534 template <class Evaluation>
1535 Evaluation h00_(const Evaluation& t) const
1536 { return (2*t - 3)*t*t + 1; }
1537
1538 template <class Evaluation>
1539 Evaluation h10_(const Evaluation& t) const
1540 { return ((t - 2)*t + 1)*t; }
1541
1542 template <class Evaluation>
1543 Evaluation h01_(const Evaluation& t) const
1544 { return (-2*t + 3)*t*t; }
1545
1546 template <class Evaluation>
1547 Evaluation h11_(const Evaluation& t) const
1548 { return (t - 1)*t*t; }
1549
1550 // first derivative of the hermite basis functions
1551 template <class Evaluation>
1552 Evaluation h00_prime_(const Evaluation& t) const
1553 { return (3*2*t - 2*3)*t; }
1554
1555 template <class Evaluation>
1556 Evaluation h10_prime_(const Evaluation& t) const
1557 { return (3*t - 2*2)*t + 1; }
1558
1559 template <class Evaluation>
1560 Evaluation h01_prime_(const Evaluation& t) const
1561 { return (-3*2*t + 2*3)*t; }
1562
1563 template <class Evaluation>
1564 Evaluation h11_prime_(const Evaluation& t) const
1565 { return (3*t - 2)*t; }
1566
1567 // second derivative of the hermite basis functions
1568 template <class Evaluation>
1569 Evaluation h00_prime2_(const Evaluation& t) const
1570 { return 2*3*2*t - 2*3; }
1571
1572 template <class Evaluation>
1573 Evaluation h10_prime2_(const Evaluation& t) const
1574 { return 2*3*t - 2*2; }
1575
1576 template <class Evaluation>
1577 Evaluation h01_prime2_(const Evaluation& t) const
1578 { return -2*3*2*t + 2*3; }
1579
1580 template <class Evaluation>
1581 Evaluation h11_prime2_(const Evaluation& t) const
1582 { return 2*3*t - 2; }
1583
1584 // third derivative of the hermite basis functions
1585 template <class Evaluation>
1586 Scalar h00_prime3_(const Evaluation&) const
1587 { return 2*3*2; }
1588
1589 template <class Evaluation>
1590 Scalar h10_prime3_(const Evaluation&) const
1591 { return 2*3; }
1592
1593 template <class Evaluation>
1594 Scalar h01_prime3_(const Evaluation&) const
1595 { return -2*3*2; }
1596
1597 template <class Evaluation>
1598 Scalar h11_prime3_(const Evaluation&) const
1599 { return 2*3; }
1600
1601 // returns the monotonicality of an interval of a spline segment
1602 //
1603 // The return value have the following meaning:
1604 //
1605 // 3: spline is constant within interval [x0, x1]
1606 // 1: spline is monotonously increasing in the specified interval
1607 // 0: spline is not monotonic (or constant) in the specified interval
1608 // -1: spline is monotonously decreasing in the specified interval
1609 int monotonic_(size_t i, Scalar x0, Scalar x1, int& r) const
1610 {
1611 // coefficients of derivative in monomial basis
1612 Scalar a = 3*a_(i);
1613 Scalar b = 2*b_(i);
1614 Scalar c = c_(i);
1615
1616 if (std::abs(a) < 1e-20 && std::abs(b) < 1e-20 && std::abs(c) < 1e-20)
1617 return 3; // constant in interval, r stays unchanged!
1618
1619 Scalar disc = b*b - 4*a*c;
1620 if (disc < 0) {
1621 // discriminant of derivative is smaller than 0, i.e. the
1622 // segment's derivative does not exhibit any extrema.
1623 if (x0*(x0*a + b) + c > 0) {
1624 r = (r==3 || r == 1)?1:0;
1625 return 1;
1626 }
1627 else {
1628 r = (r==3 || r == -1)?-1:0;
1629 return -1;
1630 }
1631 }
1632 disc = std::sqrt(disc);
1633 Scalar xE1 = (-b + disc)/(2*a);
1634 Scalar xE2 = (-b - disc)/(2*a);
1635
1636 if (std::abs(disc) < 1e-30) {
1637 // saddle point -> no extrema
1638 if (std::abs(xE1 - x0) < 1e-30)
1639 // make sure that we're not picking the saddle point
1640 // to determine whether we're monotonically increasing
1641 // or decreasing
1642 x0 = x1;
1643 if (x0*(x0*a + b) + c > 0) {
1644 r = (r==3 || r == 1)?1:0;
1645 return 1;
1646 }
1647 else {
1648 r = (r==3 || r == -1)?-1:0;
1649 return -1;
1650 }
1651 };
1652 if ((x0 < xE1 && xE1 < x1) ||
1653 (x0 < xE2 && xE2 < x1))
1654 {
1655 // there is an extremum in the range (x0, x1)
1656 r = 0;
1657 return 0;
1658 }
1659 // no extremum in range (x0, x1)
1660 x0 = (x0 + x1)/2; // pick point in the middle of the interval
1661 // to avoid extrema on the boundaries
1662 if (x0*(x0*a + b) + c > 0) {
1663 r = (r==3 || r == 1)?1:0;
1664 return 1;
1665 }
1666 else {
1667 r = (r==3 || r == -1)?-1:0;
1668 return -1;
1669 }
1670 }
1671
1676 template <class Evaluation>
1677 size_t intersectSegment_(Evaluation* sol,
1678 size_t segIdx,
1679 const Evaluation& a,
1680 const Evaluation& b,
1681 const Evaluation& c,
1682 const Evaluation& d,
1683 Scalar x0 = -1e30, Scalar x1 = 1e30) const
1684 {
1685 unsigned n =
1687 a_(segIdx) - a,
1688 b_(segIdx) - b,
1689 c_(segIdx) - c,
1690 d_(segIdx) - d);
1691 x0 = std::max(x_(segIdx), x0);
1692 x1 = std::min(x_(segIdx+1), x1);
1693
1694 // filter the intersections outside of the specified interval
1695 size_t k = 0;
1696 for (unsigned j = 0; j < n; ++j) {
1697 if (x0 <= sol[j] && sol[j] <= x1) {
1698 sol[k] = sol[j];
1699 ++k;
1700 }
1701 }
1702 return k;
1703 }
1704
1705 // find the segment index for a given x coordinate
1706 size_t segmentIdx_(Scalar x) const
1707 {
1708 // bisection
1709 size_t iLow = 0;
1710 size_t iHigh = numSamples() - 1;
1711
1712 while (iLow + 1 < iHigh) {
1713 size_t i = (iLow + iHigh) / 2;
1714 if (x_(i) > x)
1715 iHigh = i;
1716 else
1717 iLow = i;
1718 };
1719 return iLow;
1720 }
1721
1725 Scalar h_(size_t i) const
1726 {
1727 assert(x_(i) > x_(i-1)); // the sampling points must be given
1728 // in ascending order
1729 return x_(i) - x_(i - 1);
1730 }
1731
1735 Scalar x_(size_t i) const
1736 { return xPos_[i]; }
1737
1741 Scalar y_(size_t i) const
1742 { return yPos_[i]; }
1743
1748 Scalar slope_(size_t i) const
1749 { return slopeVec_[i]; }
1750
1751 // returns the coefficient in front of the x^3 term. In Stoer this
1752 // is delta.
1753 Scalar a_(size_t i) const
1754 { return evalDerivative3_(/*x=*/Scalar(0.0), i)/6.0; }
1755
1756 // returns the coefficient in front of the x^2 term In Stoer this
1757 // is gamma.
1758 Scalar b_(size_t i) const
1759 { return evalDerivative2_(/*x=*/Scalar(0.0), i)/2.0; }
1760
1761 // returns the coefficient in front of the x^1 term. In Stoer this
1762 // is beta.
1763 Scalar c_(size_t i) const
1764 { return evalDerivative_(/*x=*/Scalar(0.0), i); }
1765
1766 // returns the coefficient in front of the x^0 term. In Stoer this
1767 // is alpha.
1768 Scalar d_(size_t i) const
1769 { return eval_(/*x=*/Scalar(0.0), i); }
1770
1771 Vector xPos_;
1772 Vector yPos_;
1773 Vector slopeVec_;
1774};
1775}
1776
1777#endif
Provides the OPM specific exception classes.
Provides free functions to invert polynomials of degree 1, 2 and 3.
Provides a tridiagonal matrix that also supports non-zero entries in the upper right and lower left.
Definition Exceptions.hpp:40
size_t numSamples() const
Return the number of (x, y) values.
Definition Spline.hpp:721
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 Spline.hpp:1230
void setArrayOfPoints(size_t nSamples, const PointArray &points, SplineType splineType=Natural, bool sortInputs=true)
Set the sampling points of a natural spline using a C-style array.
Definition Spline.hpp:595
void setXYArrays(size_t nSamples, const ScalarArrayX &x, const ScalarArrayY &y, Scalar m0, Scalar m1, bool sortInputs=true)
Set the sampling points and the boundary slopes of a full spline using C-style arrays.
Definition Spline.hpp:326
Scalar valueAt(size_t sampleIdx) const
Return the x value of a given sampling point.
Definition Spline.hpp:733
void makeNaturalSpline_()
Create a natural spline from the already set sampling points.
Definition Spline.hpp:1088
void sortInput_()
Sort the sample points in ascending order of their x value.
Definition Spline.hpp:1014
Spline(size_t nSamples, const ScalarArrayX &x, const ScalarArrayY &y, SplineType splineType=Natural, bool sortInputs=true)
Convenience constructor for a natural or a periodic spline.
Definition Spline.hpp:141
void printCSV(Scalar xi0, Scalar xi1, size_t k, std::ostream &os) const
Prints k tuples of the format (x, y, dx/dy, isMonotonic) to stdout.
Definition Spline.cpp:33
Spline(const PointContainer &points, Scalar m0, Scalar m1, bool sortInputs=true)
Convenience constructor for a full spline.
Definition Spline.hpp:253
void setSlopesFromMoments_(SlopeVector &slopes, const MomentsVector &moments)
Convert the moments at the sample points to slopes.
Definition Spline.hpp:1411
void setContainerOfTuples(const XYContainer &points, Scalar m0, Scalar m1, bool sortInputs=true)
Set the sampling points and the boundary slopes of a full spline using a STL-compatible container of ...
Definition Spline.hpp:473
void makePeriodicSystem_(Matrix &M, Vector &d)
Make the linear system of equations Mx = d which results in the moments of the periodic spline.
Definition Spline.hpp:1300
int monotonic(Scalar x0, Scalar x1, bool extrapolate=false) const
Returns 1 if the spline is monotonically increasing, -1 if the spline is mononously decreasing and 0 ...
Definition Spline.hpp:925
bool applies(const Evaluation &x) const
Return true iff the given x is in range [x1, xn].
Definition Spline.hpp:715
Spline(const ScalarContainer &x, const ScalarContainer &y, SplineType splineType=Natural, bool sortInputs=true)
Convenience constructor for a natural or a periodic spline.
Definition Spline.hpp:172
Spline(Scalar x0, Scalar x1, Scalar y0, Scalar y1, Scalar m0, Scalar m1)
Convenience constructor for a full spline with just two sampling points.
Definition Spline.hpp:126
Evaluation intersect(const Evaluation &a, const Evaluation &b, const Evaluation &c, const Evaluation &d) const
Find the intersections of the spline with a cubic polynomial in the whole interval,...
Definition Spline.hpp:872
Scalar y_(size_t i) const
Returns the y coordinate of the i-th sampling point.
Definition Spline.hpp:1741
Evaluation evalDerivative(const Evaluation &x, bool extrapolate=false) const
Evaluate the spline's derivative at a given position.
Definition Spline.hpp:803
SplineType
The type of the spline to be created.
Definition Spline.hpp:101
Spline(const PointContainer &points, SplineType splineType=Natural, bool sortInputs=true)
Convenience constructor for a natural or a periodic spline.
Definition Spline.hpp:186
void makeNaturalSystem_(Matrix &M, Vector &d)
Make the linear system of equations Mx = d which results in the moments of the natural spline.
Definition Spline.hpp:1254
Scalar h_(size_t i) const
Returns x[i] - x[i - 1].
Definition Spline.hpp:1725
void set(Scalar x0, Scalar x1, Scalar y0, Scalar y1, Scalar m0, Scalar m1)
Set the sampling points and the boundary slopes of the spline with two sampling points.
Definition Spline.hpp:270
Spline(const ScalarContainerX &x, const ScalarContainerY &y, Scalar m0, Scalar m1, bool sortInputs=true)
Convenience constructor for a full spline.
Definition Spline.hpp:237
void setContainerOfPoints(const XYContainer &points, SplineType splineType=Natural, bool sortInputs=true)
Set the sampling points of a natural spline using a STL-compatible container of array-like objects.
Definition Spline.hpp:637
Evaluation evalThirdDerivative(const Evaluation &x, bool extrapolate=false) const
Evaluate the spline's third derivative at a given position.
Definition Spline.hpp:855
Spline(size_t nSamples, const ScalarArray &x, const ScalarArray &y, Scalar m0, Scalar m1, bool sortInputs=true)
Convenience constructor for a full spline.
Definition Spline.hpp:202
void assignSamplingPoints_(DestVector &destX, DestVector &destY, const SourceVector &srcX, const SourceVector &srcY, unsigned nSamples)
Set the sampling point vectors.
Definition Spline.hpp:1142
Scalar xAt(size_t sampleIdx) const
Return the x value of a given sampling point.
Definition Spline.hpp:727
int monotonic() const
Same as monotonic(x0, x1), but with the entire range of the spline as interval.
Definition Spline.hpp:992
void assignFromTupleList_(DestVector &destX, DestVector &destY, ListIterator srcBegin, ListIterator srcEnd, unsigned nSamples)
Set the sampling points.
Definition Spline.hpp:1196
Spline(size_t nSamples, const PointArray &points, Scalar m0, Scalar m1, bool sortInputs=true)
Convenience constructor for a full spline.
Definition Spline.hpp:220
Scalar slope_(size_t i) const
Returns the slope (i.e.
Definition Spline.hpp:1748
Evaluation eval(const Evaluation &x, bool extrapolate=false) const
Evaluate the spline at a given position.
Definition Spline.hpp:767
size_t intersectSegment_(Evaluation *sol, size_t segIdx, const Evaluation &a, const Evaluation &b, const Evaluation &c, const Evaluation &d, Scalar x0=-1e30, Scalar x1=1e30) const
Find all the intersections of a segment of the spline with a cubic polynomial within a specified inte...
Definition Spline.hpp:1677
Scalar x_(size_t i) const
Returns the y coordinate of the i-th sampling point.
Definition Spline.hpp:1735
void setArrayOfPoints(size_t nSamples, const PointArray &points, Scalar m0, Scalar m1, bool sortInputs=true)
Set the sampling points and the boundary slopes of a full spline using a C-style array.
Definition Spline.hpp:394
void makeFullSpline_(Scalar m0, Scalar m1)
Create a natural spline from the already set sampling points.
Definition Spline.hpp:1067
Evaluation evalSecondDerivative(const Evaluation &x, bool extrapolate=false) const
Evaluate the spline's second derivative at a given position.
Definition Spline.hpp:832
void makeMonotonicSpline_(Vector &slopes)
Create a monotonic spline from the already set sampling points.
Definition Spline.hpp:1362
void setNumSamples_(size_t nSamples)
Resizes the internal vectors to store the sample points.
Definition Spline.hpp:1055
void setXYArrays(size_t nSamples, const ScalarArrayX &x, const ScalarArrayY &y, SplineType splineType=Natural, bool sortInputs=true)
Set the sampling points natural spline using C-style arrays.
Definition Spline.hpp:514
Spline()
Default constructor for a spline.
Definition Spline.hpp:113
void reverseSamplingPoints_()
Reverse order of the elements in the arrays which contain the sampling points.
Definition Spline.hpp:1042
void setContainerOfTuples(const XYContainer &points, SplineType splineType=Natural, bool sortInputs=true)
Set the sampling points of a natural spline using a STL-compatible container of tuple-like objects.
Definition Spline.hpp:683
void setXYContainers(const ScalarContainerX &x, const ScalarContainerY &y, SplineType splineType=Natural, bool sortInputs=true)
Set the sampling points of a natural spline using STL-compatible containers.
Definition Spline.hpp:555
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:360
void makePeriodicSpline_()
Create a periodic spline from the already set sampling points.
Definition Spline.hpp:1109
void setContainerOfPoints(const XYContainer &points, Scalar m0, Scalar m1, bool sortInputs=true)
Set the sampling points and the boundary slopes of a full spline using a STL-compatible container of ...
Definition Spline.hpp:431
Evaluation intersectInterval(Scalar x0, Scalar x1, const Evaluation &a, const Evaluation &b, const Evaluation &c, const Evaluation &d) const
Find the intersections of the spline with a cubic polynomial in a sub-interval of the spline,...
Definition Spline.hpp:887
Spline(size_t nSamples, const PointArray &points, SplineType splineType=Natural, bool sortInputs=true)
Convenience constructor for a natural or a periodic spline.
Definition Spline.hpp:157
Provides a tridiagonal matrix that also supports non-zero entries in the upper right and lower left.
Definition TridiagonalMatrix.hpp:50
In the namespace cmp are implemented functions for approximate comparison of double values based on a...
Definition cmp.hpp:64
This class implements a small container which holds the transmissibility mulitpliers for all the face...
Definition Exceptions.hpp:30
unsigned invertCubicPolynomial(SolContainer *sol, Scalar a, Scalar b, Scalar c, Scalar d)
Invert a cubic polynomial analytically.
Definition PolynomialUtils.hpp:149
Helper class needed to sort the input sampling points.
Definition Spline.hpp:1000