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  Copyright (C) 2010-2013 by Andreas Lauser
5 
6  This file is part of the Open Porous Media project (OPM).
7 
8  OPM is free software: you can redistribute it and/or modify
9  it under the terms of the GNU General Public License as published by
10  the Free Software Foundation, either version 2 of the License, or
11  (at your option) any later version.
12 
13  OPM is distributed in the hope that it will be useful,
14  but WITHOUT ANY WARRANTY; without even the implied warranty of
15  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16  GNU General Public License for more details.
17 
18  You should have received a copy of the GNU General Public License
19  along with OPM. If not, see <http://www.gnu.org/licenses/>.
20 */
25 #ifndef OPM_SPLINE_HPP
26 #define OPM_SPLINE_HPP
27 
30 #include <opm/common/ErrorMacros.hpp>
32 
33 #include <ostream>
34 #include <vector>
35 #include <tuple>
36 
37 namespace Opm
38 {
88 template<class Scalar>
89 class Spline
90 {
92  typedef std::vector<Scalar> Vector;
93 
94 public:
100  enum SplineType {
105  };
106 
113  { }
114 
125  Spline(Scalar x0, Scalar x1,
126  Scalar y0, Scalar y1,
127  Scalar m0, Scalar m1)
128  { set(x0, x1, y0, y1, m0, m1); }
129 
138  template <class ScalarArrayX, class ScalarArrayY>
139  Spline(size_t nSamples,
140  const ScalarArrayX &x,
141  const ScalarArrayY &y,
142  SplineType splineType = Natural,
143  bool sortInputs = false)
144  { this->setXYArrays(nSamples, x, y, splineType, sortInputs); }
145 
153  template <class PointArray>
154  Spline(size_t nSamples,
155  const PointArray &points,
156  SplineType splineType = Natural,
157  bool sortInputs = false)
158  { this->setArrayOfPoints(nSamples, points, splineType, sortInputs); }
159 
167  template <class ScalarContainer>
168  Spline(const ScalarContainer &x,
169  const ScalarContainer &y,
170  SplineType splineType = Natural,
171  bool sortInputs = false)
172  { this->setXYContainers(x, y, splineType, sortInputs); }
173 
180  template <class PointContainer>
181  Spline(const PointContainer &points,
182  SplineType splineType = Natural,
183  bool sortInputs = false)
184  { this->setContainerOfPoints(points, splineType, sortInputs); }
185 
196  template <class ScalarArray>
197  Spline(size_t nSamples,
198  const ScalarArray &x,
199  const ScalarArray &y,
200  Scalar m0,
201  Scalar m1,
202  bool sortInputs = false)
203  { this->setXYArrays(nSamples, x, y, m0, m1, sortInputs); }
204 
214  template <class PointArray>
215  Spline(size_t nSamples,
216  const PointArray &points,
217  Scalar m0,
218  Scalar m1,
219  bool sortInputs = false)
220  { this->setArrayOfPoints(nSamples, points, m0, m1, sortInputs); }
221 
231  template <class ScalarContainerX, class ScalarContainerY>
232  Spline(const ScalarContainerX &x,
233  const ScalarContainerY &y,
234  Scalar m0,
235  Scalar m1,
236  bool sortInputs = false)
237  { this->setXYContainers(x, y, m0, m1, sortInputs); }
238 
247  template <class PointContainer>
248  Spline(const PointContainer &points,
249  Scalar m0,
250  Scalar m1,
251  bool sortInputs = false)
252  { this->setContainerOfPoints(points, m0, m1, sortInputs); }
253 
257  size_t numSamples() const
258  { return xPos_.size(); }
259 
271  void set(Scalar x0, Scalar x1,
272  Scalar y0, Scalar y1,
273  Scalar m0, Scalar m1)
274  {
275  slopeVec_.resize(2);
276  xPos_.resize(2);
277  yPos_.resize(2);
278 
279  if (x0 > x1) {
280  xPos_[0] = x1;
281  xPos_[1] = x0;
282  yPos_[0] = y1;
283  yPos_[1] = y0;
284  }
285  else {
286  xPos_[0] = x0;
287  xPos_[1] = x1;
288  yPos_[0] = y0;
289  yPos_[1] = y1;
290  }
291 
292  slopeVec_[0] = m0;
293  slopeVec_[1] = m1;
294 
295  Matrix M(numSamples());
296  Vector d(numSamples());
297  Vector moments(numSamples());
298  this->makeFullSystem_(M, d, m0, m1);
299 
300  // solve for the moments
301  M.solve(moments, d);
302 
303  this->setSlopesFromMoments_(slopeVec_, moments);
304  }
305 
306 
310  // Full splines //
314 
326  template <class ScalarArrayX, class ScalarArrayY>
327  void setXYArrays(size_t nSamples,
328  const ScalarArrayX &x,
329  const ScalarArrayY &y,
330  Scalar m0, Scalar m1,
331  bool sortInputs = false)
332  {
333  assert(nSamples > 1);
334 
335  setNumSamples_(nSamples);
336  for (size_t i = 0; i < nSamples; ++i) {
337  xPos_[i] = x[i];
338  yPos_[i] = y[i];
339  }
340 
341  if (sortInputs)
342  sortInput_();
343  else if (xPos_[0] > xPos_[numSamples() - 1])
345 
346  makeFullSpline_(m0, m1);
347  }
348 
360  template <class ScalarContainerX, class ScalarContainerY>
361  void setXYContainers(const ScalarContainerX &x,
362  const ScalarContainerY &y,
363  Scalar m0, Scalar m1,
364  bool sortInputs = false)
365  {
366  assert(x.size() == y.size());
367  assert(x.size() > 1);
368 
369  setNumSamples_(x.size());
370 
371  std::copy(x.begin(), x.end(), xPos_.begin());
372  std::copy(y.begin(), y.end(), yPos_.begin());
373 
374  if (sortInputs)
375  sortInput_();
376  else if (xPos_[0] > xPos_[numSamples() - 1])
378 
379  makeFullSpline_(m0, m1);
380  }
381 
394  template <class PointArray>
395  void setArrayOfPoints(size_t nSamples,
396  const PointArray &points,
397  Scalar m0,
398  Scalar m1,
399  bool sortInputs = false)
400  {
401  // a spline with no or just one sampling points? what an
402  // incredible bad idea!
403  assert(nSamples > 1);
404 
405  setNumSamples_(nSamples);
406  for (size_t i = 0; i < nSamples; ++i) {
407  xPos_[i] = points[i][0];
408  yPos_[i] = points[i][1];
409  }
410 
411  if (sortInputs)
412  sortInput_();
413  else if (xPos_[0] > xPos_[numSamples() - 1])
415 
416  makeFullSpline_(m0, m1);
417  }
418 
431  template <class XYContainer>
432  void setContainerOfPoints(const XYContainer &points,
433  Scalar m0,
434  Scalar m1,
435  bool sortInputs = false)
436  {
437  // a spline with no or just one sampling points? what an
438  // incredible bad idea!
439  assert(points.size() > 1);
440 
441  setNumSamples_(points.size());
442  typename XYContainer::const_iterator it = points.begin();
443  typename XYContainer::const_iterator endIt = points.end();
444  for (size_t i = 0; it != endIt; ++i, ++it) {
445  xPos_[i] = (*it)[0];
446  yPos_[i] = (*it)[1];
447  }
448 
449  if (sortInputs)
450  sortInput_();
451  else if (xPos_[0] > xPos_[numSamples() - 1])
453 
454  // make a full spline
455  makeFullSpline_(m0, m1);
456  }
457 
473  template <class XYContainer>
474  void setContainerOfTuples(const XYContainer &points,
475  Scalar m0,
476  Scalar m1,
477  bool sortInputs = false)
478  {
479  // resize internal arrays
480  setNumSamples_(points.size());
481  typename XYContainer::const_iterator it = points.begin();
482  typename XYContainer::const_iterator endIt = points.end();
483  for (int i = 0; it != endIt; ++i, ++it) {
484  xPos_[i] = std::get<0>(*it);
485  yPos_[i] = std::get<1>(*it);
486  }
487 
488  if (sortInputs)
489  sortInput_();
490  else if (xPos_[0] > xPos_[numSamples() - 1])
492 
493  // make a full spline
494  makeFullSpline_(m0, m1);
495  }
496 
500  // Natural/Periodic splines //
504 
514  template <class ScalarArrayX, class ScalarArrayY>
515  void setXYArrays(size_t nSamples,
516  const ScalarArrayX &x,
517  const ScalarArrayY &y,
518  SplineType splineType = Natural,
519  bool sortInputs = false)
520  {
521  assert(nSamples > 1);
522 
523  setNumSamples_(nSamples);
524  for (size_t i = 0; i < nSamples; ++i) {
525  xPos_[i] = x[i];
526  yPos_[i] = y[i];
527  }
528 
529  if (sortInputs)
530  sortInput_();
531  else if (xPos_[0] > xPos_[numSamples() - 1])
533 
534  if (splineType == Periodic)
536  else if (splineType == Natural)
538  else if (splineType == Monotonic)
540  else
541  OPM_THROW(std::runtime_error, "Spline type " << splineType << " not supported at this place");
542  }
543 
555  template <class ScalarContainerX, class ScalarContainerY>
556  void setXYContainers(const ScalarContainerX &x,
557  const ScalarContainerY &y,
558  SplineType splineType = Natural,
559  bool sortInputs = false)
560  {
561  assert(x.size() == y.size());
562  assert(x.size() > 1);
563 
564  setNumSamples_(x.size());
565  std::copy(x.begin(), x.end(), xPos_.begin());
566  std::copy(y.begin(), y.end(), yPos_.begin());
567 
568  if (sortInputs)
569  sortInput_();
570  else if (xPos_[0] > xPos_[numSamples() - 1])
572 
573  if (splineType == Periodic)
575  else if (splineType == Natural)
577  else if (splineType == Monotonic)
579  else
580  OPM_THROW(std::runtime_error, "Spline type " << splineType << " not supported at this place");
581  }
582 
595  template <class PointArray>
596  void setArrayOfPoints(size_t nSamples,
597  const PointArray &points,
598  SplineType splineType = Natural,
599  bool sortInputs = false)
600  {
601  // a spline with no or just one sampling points? what an
602  // incredible bad idea!
603  assert(nSamples > 1);
604 
605  setNumSamples_(nSamples);
606  for (size_t i = 0; i < nSamples; ++i) {
607  xPos_[i] = points[i][0];
608  yPos_[i] = points[i][1];
609  }
610 
611  if (sortInputs)
612  sortInput_();
613  else if (xPos_[0] > xPos_[numSamples() - 1])
615 
616  if (splineType == Periodic)
618  else if (splineType == Natural)
620  else if (splineType == Monotonic)
622  else
623  OPM_THROW(std::runtime_error, "Spline type " << splineType << " not supported at this place");
624  }
625 
637  template <class XYContainer>
638  void setContainerOfPoints(const XYContainer &points,
639  SplineType splineType = Natural,
640  bool sortInputs = false)
641  {
642  // a spline with no or just one sampling points? what an
643  // incredible bad idea!
644  assert(points.size() > 1);
645 
646  setNumSamples_(points.size());
647  typename XYContainer::const_iterator it = points.begin();
648  typename XYContainer::const_iterator endIt = points.end();
649  for (size_t i = 0; it != endIt; ++ i, ++it) {
650  xPos_[i] = (*it)[0];
651  yPos_[i] = (*it)[1];
652  }
653 
654  if (sortInputs)
655  sortInput_();
656  else if (xPos_[0] > xPos_[numSamples() - 1])
658 
659  if (splineType == Periodic)
661  else if (splineType == Natural)
663  else if (splineType == Monotonic)
665  else
666  OPM_THROW(std::runtime_error, "Spline type " << splineType << " not supported at this place");
667  }
668 
683  template <class XYContainer>
684  void setContainerOfTuples(const XYContainer &points,
685  SplineType splineType = Natural,
686  bool sortInputs = false)
687  {
688  // resize internal arrays
689  setNumSamples_(points.size());
690  typename XYContainer::const_iterator it = points.begin();
691  typename XYContainer::const_iterator endIt = points.end();
692  for (unsigned i = 0; it != endIt; ++i, ++it) {
693  xPos_[i] = std::get<0>(*it);
694  yPos_[i] = std::get<1>(*it);
695  }
696 
697  if (sortInputs)
698  sortInput_();
699  else if (xPos_[0] > xPos_[numSamples() - 1])
701 
702  if (splineType == Periodic)
704  else if (splineType == Natural)
706  else if (splineType == Monotonic)
708  else
709  OPM_THROW(std::runtime_error, "Spline type " << splineType << " not supported at this place");
710  }
711 
715  bool applies(Scalar x) const
716  {
717  return x_(0) <= x && x <= x_(numSamples() - 1);
718  }
719 
723  Scalar xMin() const
724  { return x_(0); }
725 
729  Scalar xMax() const
730  { return x_(numSamples() - 1); }
731 
735  Scalar yFirst() const
736  { return y_(0); }
737 
741  Scalar yLast() const
742  { return y_(numSamples() - 1); }
743 
761  void printCSV(Scalar xi0, Scalar xi1, size_t k, std::ostream &os = std::cout) const
762  {
763  Scalar x0 = std::min(xi0, xi1);
764  Scalar x1 = std::max(xi0, xi1);
765  const size_t n = numSamples() - 1;
766  for (size_t i = 0; i <= k; ++i) {
767  double x = i*(x1 - x0)/k + x0;
768  double x_p1 = x + (x1 - x0)/k;
769  double y;
770  double dy_dx;
771  double mono = 1;
772  if (!applies(x)) {
773  if (x < x_(0)) {
774  dy_dx = evalDerivative(x_(0));
775  y = (x - x_(0))*dy_dx + y_(0);
776  mono = (dy_dx>0)?1:-1;
777  }
778  else if (x > x_(n)) {
779  dy_dx = evalDerivative(x_(n));
780  y = (x - x_(n))*dy_dx + y_(n);
781  mono = (dy_dx>0)?1:-1;
782  }
783  else {
784  OPM_THROW(std::runtime_error,
785  "The sampling points given to a spline must be sorted by their x value!");
786  }
787  }
788  else {
789  y = eval(x);
790  dy_dx = evalDerivative(x);
791  mono = monotonic(x, x_p1, /*extrapolate=*/true);
792  }
793 
794  os << x << " " << y << " " << dy_dx << " " << mono << "\n";
795  }
796  }
797 
809  Scalar eval(Scalar x, bool extrapolate=false) const
810  {
811  assert(extrapolate || applies(x));
812 
813  // handle extrapolation
814  if (extrapolate) {
815  if (x < xMin()) {
816  Scalar m = evalDerivative_(xMin(), /*segmentIdx=*/0);
817  Scalar y0 = y_(0);
818  return y0 + m*(x - xMin());
819  }
820  else if (x > xMax()) {
821  Scalar m = evalDerivative_(xMax(), /*segmentIdx=*/numSamples()-2);
822  Scalar y0 = y_(numSamples() - 1);
823  return y0 + m*(x - xMax());
824  }
825  }
826 
827  return eval_(x, segmentIdx_(x));
828  }
829 
841  template <class Evaluation>
842  Evaluation eval(const Evaluation& x, bool extrapolate=false) const
843  {
844  assert(extrapolate || applies(x.value));
845 
846  // handle extrapolation
847  if (extrapolate) {
848  if (x.value < xMin()) {
849  Scalar m = evalDerivative_(xMin(), /*segmentIdx=*/0);
850  Scalar y0 = y_(0);
851  return Evaluation::createConstant(y0 + m*(x.value - xMin()));
852  }
853  else if (x > xMax()) {
854  Scalar m = evalDerivative_(xMax(), /*segmentIdx=*/numSamples()-2);
855  Scalar y0 = y_(numSamples() - 1);
856  return Evaluation::createConstant(y0 + m*(x.value - xMax()));
857  }
858  }
859 
860  return eval_(x, segmentIdx_(x.value));
861  }
862 
875  Scalar evalDerivative(Scalar x, bool extrapolate=false) const
876  {
877  assert(extrapolate || applies(x));
878  if (extrapolate) {
879  if (x < xMin())
880  return evalDerivative_(xMin(), /*segmentIdx=*/0);
881  else if (x > xMax())
882  return evalDerivative_(xMax(), /*segmentIdx=*/numSamples() - 2);
883  }
884 
885  return evalDerivative_(x, segmentIdx_(x));
886  }
887 
900  Scalar evalSecondDerivative(Scalar x, bool extrapolate=false) const
901  {
902  assert(extrapolate || applies(x));
903  if (extrapolate)
904  return 0.0;
905 
906  return evalDerivative2_(x, segmentIdx_(x));
907  }
908 
921  Scalar evalThirdDerivative(Scalar x, bool extrapolate=false) const
922  {
923  assert(extrapolate || applies(x));
924  if (extrapolate)
925  return 0.0;
926 
927  return evalDerivative3_(x, segmentIdx_(x));
928  }
929 
936  template <class Evaluation>
937  Evaluation intersect(const Evaluation& a,
938  const Evaluation& b,
939  const Evaluation& c,
940  const Evaluation& d) const
941  {
942  return intersectInterval(xMin(), xMax(), a, b, c, d);
943  }
944 
951  template <class Evaluation>
952  Evaluation intersectInterval(Scalar x0, Scalar x1,
953  const Evaluation& a,
954  const Evaluation& b,
955  const Evaluation& c,
956  const Evaluation& d) const
957  {
958  assert(applies(x0) && applies(x1));
959 
960  Evaluation tmpSol[3], sol = 0;
961  size_t nSol = 0;
962  size_t iFirst = segmentIdx_(x0);
963  size_t iLast = segmentIdx_(x1);
964  for (size_t i = iFirst; i <= iLast; ++i)
965  {
966  size_t nCur = intersectSegment_(tmpSol, i, a, b, c, d, x0, x1);
967  if (nCur == 1)
968  sol = tmpSol[0];
969 
970  nSol += nCur;
971  if (nSol > 1) {
972  OPM_THROW(std::runtime_error,
973  "Spline has more than one intersection"); //<<a<<"x^3 + "<<b<"x^2 + "<<c<"x + "<<d);
974  }
975  }
976 
977  if (nSol != 1)
978  OPM_THROW(std::runtime_error,
979  "Spline has no intersection"); //<<a<"x^3 + " <<b<"x^2 + "<<c<"x + "<<d<<"!");
980 
981  return sol;
982  }
983 
992  int monotonic(Scalar x0, Scalar x1, OPM_OPTIM_UNUSED bool extrapolate = false) const
993  {
994  assert(std::abs(x0 - x1) > 1e-30);
995 
996  // make sure that x0 is smaller than x1
997  if (x0 > x1)
998  std::swap(x0, x1);
999 
1000  assert(x0 < x1);
1001 
1002  int r = 3;
1003  if (x0 < xMin()) {
1004  assert(extrapolate);
1005  Scalar m = evalDerivative_(xMin(), /*segmentIdx=*/0);
1006  if (std::abs(m) < 1e-20)
1007  r = (m < 0)?-1:1;
1008  x0 = xMin();
1009  };
1010 
1011  size_t i = segmentIdx_(x0);
1012  if (x_(i + 1) >= x1) {
1013  // interval is fully contained within a single spline
1014  // segment
1015  monotonic_(i, x0, x1, r);
1016  return r;
1017  }
1018 
1019  // the first segment overlaps with the specified interval
1020  // partially
1021  monotonic_(i, x0, x_(i+1), r);
1022  ++ i;
1023 
1024  // make sure that the segments which are completly in the
1025  // interval [x0, x1] all exhibit the same monotonicity.
1026  size_t iEnd = segmentIdx_(x1);
1027  for (; i < iEnd - 1; ++i) {
1028  monotonic_(i, x_(i), x_(i + 1), r);
1029  if (!r)
1030  return 0;
1031  }
1032 
1033  // if the user asked for a part of the spline which is
1034  // extrapolated, we need to check the slope at the spline's
1035  // endpoint
1036  if (x1 > xMax()) {
1037  assert(extrapolate);
1038 
1039  Scalar m = evalDerivative_(xMax(), /*segmentIdx=*/numSamples() - 2);
1040  if (m < 0)
1041  return (r < 0 || r==3)?-1:0;
1042  else if (m > 0)
1043  return (r > 0 || r==3)?1:0;
1044 
1045  return r;
1046  }
1047 
1048  // check for the last segment
1049  monotonic_(iEnd, x_(iEnd), x1, r);
1050 
1051  return r;
1052  }
1053 
1058  int monotonic() const
1059  { return monotonic(xMin(), xMax()); }
1060 
1061 protected:
1066  {
1067  ComparatorX_(const std::vector<Scalar> &x)
1068  : x_(x)
1069  {}
1070 
1071  bool operator ()(unsigned idxA, unsigned idxB) const
1072  { return x_.at(idxA) < x_.at(idxB); }
1073 
1074  const std::vector<Scalar> &x_;
1075  };
1076 
1080  void sortInput_()
1081  {
1082  size_t n = numSamples();
1083 
1084  // create a vector containing 0...n-1
1085  std::vector<unsigned> idxVector(static_cast<unsigned>(n));
1086  for (unsigned i = 0; i < n; ++i)
1087  idxVector[static_cast<unsigned>(i)] = i;
1088 
1089  // sort the indices according to the x values of the sample
1090  // points
1091  ComparatorX_ cmp(xPos_);
1092  std::sort(idxVector.begin(), idxVector.end(), cmp);
1093 
1094  // reorder the sample points
1095  std::vector<Scalar> tmpX(n), tmpY(n);
1096  for (size_t i = 0; i < idxVector.size(); ++ i) {
1097  tmpX[i] = xPos_[idxVector[i]];
1098  tmpY[i] = yPos_[idxVector[i]];
1099  }
1100  xPos_ = tmpX;
1101  yPos_ = tmpY;
1102  }
1103 
1109  {
1110  // reverse the arrays
1111  size_t n = numSamples();
1112  for (unsigned i = 0; i <= (n - 1)/2; ++i) {
1113  std::swap(xPos_[i], xPos_[n - i - 1]);
1114  std::swap(yPos_[i], yPos_[n - i - 1]);
1115  }
1116  }
1117 
1121  void setNumSamples_(size_t nSamples)
1122  {
1123  xPos_.resize(nSamples);
1124  yPos_.resize(nSamples);
1125  slopeVec_.resize(nSamples);
1126  }
1127 
1133  void makeFullSpline_(Scalar m0, Scalar m1)
1134  {
1135  Matrix M(numSamples());
1136  std::vector<Scalar> d(numSamples());
1137  std::vector<Scalar> moments(numSamples());
1138 
1139  // create linear system of equations
1140  this->makeFullSystem_(M, d, m0, m1);
1141 
1142  // solve for the moments (-> second derivatives)
1143  M.solve(moments, d);
1144 
1145  // convert the moments to slopes at the sample points
1146  this->setSlopesFromMoments_(slopeVec_, moments);
1147  }
1148 
1155  {
1156  Matrix M(numSamples(), numSamples());
1157  Vector d(numSamples());
1158  Vector moments(numSamples());
1159 
1160  // create linear system of equations
1161  this->makeNaturalSystem_(M, d);
1162 
1163  // solve for the moments (-> second derivatives)
1164  M.solve(moments, d);
1165 
1166  // convert the moments to slopes at the sample points
1167  this->setSlopesFromMoments_(slopeVec_, moments);
1168  }
1169 
1176  {
1177  Matrix M(numSamples() - 1);
1178  Vector d(numSamples() - 1);
1179  Vector moments(numSamples() - 1);
1180 
1181  // create linear system of equations. This is a bit hacky,
1182  // because it assumes that std::vector internally stores its
1183  // data as a big C-style array, but it saves us from yet
1184  // another copy operation
1185  this->makePeriodicSystem_(M, d);
1186 
1187  // solve for the moments (-> second derivatives)
1188  M.solve(moments, d);
1189 
1190  moments.resize(numSamples());
1191  for (int i = static_cast<int>(numSamples()) - 2; i >= 0; --i)
1192  moments[static_cast<unsigned>(i+1)] = moments[static_cast<unsigned>(i)];
1193  moments[0] = moments[numSamples() - 1];
1194 
1195  // convert the moments to slopes at the sample points
1196  this->setSlopesFromMoments_(slopeVec_, moments);
1197  }
1198 
1205  template <class DestVector, class SourceVector>
1206  void assignSamplingPoints_(DestVector &destX,
1207  DestVector &destY,
1208  const SourceVector &srcX,
1209  const SourceVector &srcY,
1210  int numSamples)
1211  {
1212  assert(numSamples >= 2);
1213 
1214  // copy sample points, make sure that the first x value is
1215  // smaller than the last one
1216  for (int i = 0; i < numSamples; ++i) {
1217  int idx = i;
1218  if (srcX[0] > srcX[numSamples - 1])
1219  idx = numSamples - i - 1;
1220  destX[i] = srcX[idx];
1221  destY[i] = srcY[idx];
1222  }
1223  }
1224 
1225  template <class DestVector, class ListIterator>
1226  void assignFromArrayList_(DestVector &destX,
1227  DestVector &destY,
1228  const ListIterator &srcBegin,
1229  const ListIterator &srcEnd,
1230  int numSamples)
1231  {
1232  assert(numSamples >= 2);
1233 
1234  // find out wether the x values are in reverse order
1235  ListIterator it = srcBegin;
1236  ++it;
1237  bool reverse = false;
1238  if ((*srcBegin)[0] > (*it)[0])
1239  reverse = true;
1240  --it;
1241 
1242  // loop over all sampling points
1243  for (int i = 0; it != srcEnd; ++i, ++it) {
1244  int idx = i;
1245  if (reverse)
1246  idx = numSamples - i - 1;
1247  destX[i] = (*it)[0];
1248  destY[i] = (*it)[1];
1249  }
1250  }
1251 
1259  template <class DestVector, class ListIterator>
1260  void assignFromTupleList_(DestVector &destX,
1261  DestVector &destY,
1262  ListIterator srcBegin,
1263  ListIterator srcEnd,
1264  int numSamples)
1265  {
1266  assert(numSamples >= 2);
1267 
1268  // copy sample points, make sure that the first x value is
1269  // smaller than the last one
1270 
1271  // find out wether the x values are in reverse order
1272  ListIterator it = srcBegin;
1273  ++it;
1274  bool reverse = false;
1275  if (std::get<0>(*srcBegin) > std::get<0>(*it))
1276  reverse = true;
1277  --it;
1278 
1279  // loop over all sampling points
1280  for (int i = 0; it != srcEnd; ++i, ++it) {
1281  int idx = i;
1282  if (reverse)
1283  idx = numSamples - i - 1;
1284  destX[i] = std::get<0>(*it);
1285  destY[i] = std::get<1>(*it);
1286  }
1287  }
1288 
1293  template <class Vector, class Matrix>
1294  void makeFullSystem_(Matrix &M, Vector &d, Scalar m0, Scalar m1)
1295  {
1296  makeNaturalSystem_(M, d);
1297 
1298  size_t n = numSamples() - 1;
1299  // first row
1300  M[0][1] = 1;
1301  d[0] = 6/h_(1) * ( (y_(1) - y_(0))/h_(1) - m0);
1302 
1303  // last row
1304  M[n][n - 1] = 1;
1305 
1306  // right hand side
1307  d[n] =
1308  6/h_(n)
1309  *
1310  (m1 - (y_(n) - y_(n - 1))/h_(n));
1311  }
1312 
1317  template <class Vector, class Matrix>
1318  void makeNaturalSystem_(Matrix &M, Vector &d)
1319  {
1320  M = 0.0;
1321 
1322  // See: J. Stoer: "Numerische Mathematik 1", 9th edition,
1323  // Springer, 2005, p. 111
1324  size_t n = numSamples() - 1;
1325 
1326  // second to next to last rows
1327  for (size_t i = 1; i < n; ++i) {
1328  Scalar lambda_i = h_(i + 1) / (h_(i) + h_(i + 1));
1329  Scalar mu_i = 1 - lambda_i;
1330  Scalar d_i =
1331  6 / (h_(i) + h_(i + 1))
1332  *
1333  ( (y_(i + 1) - y_(i))/h_(i + 1) - (y_(i) - y_(i - 1))/h_(i));
1334 
1335  M[i][i-1] = mu_i;
1336  M[i][i] = 2;
1337  M[i][i + 1] = lambda_i;
1338  d[i] = d_i;
1339  };
1340 
1341  // See Stroer, equation (2.5.2.7)
1342  Scalar lambda_0 = 0;
1343  Scalar d_0 = 0;
1344 
1345  Scalar mu_n = 0;
1346  Scalar d_n = 0;
1347 
1348  // first row
1349  M[0][0] = 2;
1350  M[0][1] = lambda_0;
1351  d[0] = d_0;
1352 
1353  // last row
1354  M[n][n-1] = mu_n;
1355  M[n][n] = 2;
1356  d[n] = d_n;
1357  }
1358 
1363  template <class Matrix, class Vector>
1364  void makePeriodicSystem_(Matrix &M, Vector &d)
1365  {
1366  M = 0.0;
1367 
1368  // See: J. Stoer: "Numerische Mathematik 1", 9th edition,
1369  // Springer, 2005, p. 111
1370  size_t n = numSamples() - 1;
1371 
1372  assert(M.rows() == n);
1373 
1374  // second to next to last rows
1375  for (size_t i = 2; i < n; ++i) {
1376  Scalar lambda_i = h_(i + 1) / (h_(i) + h_(i + 1));
1377  Scalar mu_i = 1 - lambda_i;
1378  Scalar d_i =
1379  6 / (h_(i) + h_(i + 1))
1380  *
1381  ( (y_(i + 1) - y_(i))/h_(i + 1) - (y_(i) - y_(i - 1))/h_(i));
1382 
1383  M[i-1][i-2] = mu_i;
1384  M[i-1][i-1] = 2;
1385  M[i-1][i] = lambda_i;
1386  d[i-1] = d_i;
1387  };
1388 
1389  Scalar lambda_n = h_(1) / (h_(n) + h_(1));
1390  Scalar lambda_1 = h_(2) / (h_(1) + h_(2));;
1391  Scalar mu_1 = 1 - lambda_1;
1392  Scalar mu_n = 1 - lambda_n;
1393 
1394  Scalar d_1 =
1395  6 / (h_(1) + h_(2))
1396  *
1397  ( (y_(2) - y_(1))/h_(2) - (y_(1) - y_(0))/h_(1));
1398  Scalar d_n =
1399  6 / (h_(n) + h_(1))
1400  *
1401  ( (y_(1) - y_(n))/h_(1) - (y_(n) - y_(n-1))/h_(n));
1402 
1403 
1404  // first row
1405  M[0][0] = 2;
1406  M[0][1] = lambda_1;
1407  M[0][n-1] = mu_1;
1408  d[0] = d_1;
1409 
1410  // last row
1411  M[n-1][0] = lambda_n;
1412  M[n-1][n-2] = mu_n;
1413  M[n-1][n-1] = 2;
1414  d[n-1] = d_n;
1415  }
1416 
1425  template <class Vector>
1426  void makeMonotonicSpline_(Vector &slopes)
1427  {
1428  auto n = numSamples();
1429 
1430  // calculate the slopes of the secant lines
1431  std::vector<Scalar> delta(n);
1432  for (size_t k = 0; k < n - 1; ++k)
1433  delta[k] = (y_(k + 1) - y_(k))/(x_(k + 1) - x_(k));
1434 
1435  // calculate the "raw" slopes at the sample points
1436  for (size_t k = 1; k < n - 1; ++k)
1437  slopes[k] = (delta[k - 1] + delta[k])/2;
1438  slopes[0] = delta[0];
1439  slopes[n - 1] = delta[n - 2];
1440 
1441  // post-process the "raw" slopes at the sample points
1442  for (size_t k = 0; k < n - 1; ++k) {
1443  if (std::abs(delta[k]) < 1e-50) {
1444  // make the spline flat if the inputs are equal
1445  slopes[k] = 0;
1446  slopes[k + 1] = 0;
1447  ++ k;
1448  continue;
1449  }
1450  else {
1451  Scalar alpha = slopes[k] / delta[k];
1452  Scalar beta = slopes[k + 1] / delta[k];
1453 
1454  if (alpha < 0 || (k > 0 && slopes[k] / delta[k - 1] < 0)) {
1455  slopes[k] = 0;
1456  }
1457  // limit (alpha, beta) to a circle of radius 3
1458  else if (alpha*alpha + beta*beta > 3*3) {
1459  Scalar tau = 3.0/std::sqrt(alpha*alpha + beta*beta);
1460  slopes[k] = tau*alpha*delta[k];
1461  slopes[k + 1] = tau*beta*delta[k];
1462  }
1463  }
1464  }
1465  }
1466 
1474  template <class MomentsVector, class SlopeVector>
1475  void setSlopesFromMoments_(SlopeVector &slopes, const MomentsVector &moments)
1476  {
1477  size_t n = numSamples();
1478 
1479  // evaluate slope at the rightmost point.
1480  // See: J. Stoer: "Numerische Mathematik 1", 9th edition,
1481  // Springer, 2005, p. 109
1482  Scalar mRight;
1483 
1484  {
1485  Scalar h = this->h_(n - 1);
1486  Scalar x = h;
1487  //Scalar x_1 = 0;
1488 
1489  Scalar A =
1490  (y_(n - 1) - y_(n - 2))/h
1491  -
1492  h/6*(moments[n-1] - moments[n - 2]);
1493 
1494  mRight =
1495  //- moments[n - 2] * x_1*x_1 / (2 * h)
1496  //+
1497  moments[n - 1] * x*x / (2 * h)
1498  +
1499  A;
1500  }
1501 
1502  // evaluate the slope for the first n-1 sample points
1503  for (size_t i = 0; i < n - 1; ++ i) {
1504  // See: J. Stoer: "Numerische Mathematik 1", 9th edition,
1505  // Springer, 2005, p. 109
1506  Scalar h_i = this->h_(i + 1);
1507  //Scalar x_i = 0;
1508  Scalar x_i1 = h_i;
1509 
1510  Scalar A_i =
1511  (y_(i+1) - y_(i))/h_i
1512  -
1513  h_i/6*(moments[i+1] - moments[i]);
1514 
1515  slopes[i] =
1516  - moments[i] * x_i1*x_i1 / (2 * h_i)
1517  +
1518  //moments[i + 1] * x_i*x_i / (2 * h_i)
1519  //+
1520  A_i;
1521 
1522  }
1523  slopes[n - 1] = mRight;
1524  }
1525 
1526 
1527  // evaluate the spline at a given the position and given the
1528  // segment index
1529  Scalar eval_(Scalar x, size_t i) const
1530  {
1531  // See http://en.wikipedia.org/wiki/Cubic_Hermite_spline
1532  Scalar delta = h_(i + 1);
1533  Scalar t = (x - x_(i))/delta;
1534 
1535  return
1536  h00_(t) * y_(i)
1537  + h10_(t) * slope_(i)*delta
1538  + h01_(t) * y_(i + 1)
1539  + h11_(t) * slope_(i + 1)*delta;
1540  }
1541 
1542  // evaluate the spline at a given the position and given the
1543  // segment index
1544  template <class Evaluation>
1545  Evaluation eval_(const Evaluation& x, size_t i) const
1546  {
1547  // See http://en.wikipedia.org/wiki/Cubic_Hermite_spline
1548  Scalar delta = h_(i + 1);
1549  Scalar t = (x.value - x_(i))/delta;
1550 
1551  Evaluation result;
1552  result.value =
1553  h00_(t) * y_(i)
1554  + h10_(t) * slope_(i)*delta
1555  + h01_(t) * y_(i + 1)
1556  + h11_(t) * slope_(i + 1)*delta;
1557 
1558  Scalar df_dg = evalDerivative_(x.value, i);
1559  for (unsigned varIdx = 0; varIdx < result.derivatives.size(); ++ varIdx)
1560  result.derivatives[varIdx] = df_dg*x.derivatives[varIdx];
1561 
1562  return result;
1563  }
1564 
1565  // evaluate the derivative of a spline given the actual position
1566  // and the segment index
1567  Scalar evalDerivative_(Scalar x, size_t i) const
1568  {
1569  // See http://en.wikipedia.org/wiki/Cubic_Hermite_spline
1570  Scalar delta = h_(i + 1);
1571  Scalar t = (x - x_(i))/delta;
1572  Scalar alpha = 1 / delta;
1573 
1574  return
1575  alpha *
1576  (h00_prime_(t) * y_(i)
1577  + h10_prime_(t) * slope_(i)*delta
1578  + h01_prime_(t) * y_(i + 1)
1579  + h11_prime_(t) * slope_(i + 1)*delta);
1580  }
1581 
1582  // evaluate the second derivative of a spline given the actual
1583  // position and the segment index
1584  Scalar evalDerivative2_(Scalar x, size_t i) const
1585  {
1586  // See http://en.wikipedia.org/wiki/Cubic_Hermite_spline
1587  Scalar delta = h_(i + 1);
1588  Scalar t = (x - x_(i))/delta;
1589  Scalar alpha = 1 / delta;
1590 
1591  return
1592  alpha*alpha
1593  *(h00_prime2_(t) * y_(i)
1594  + h10_prime2_(t) * slope_(i)*delta
1595  + h01_prime2_(t) * y_(i + 1)
1596  + h11_prime2_(t) * slope_(i + 1)*delta);
1597  }
1598 
1599  // evaluate the third derivative of a spline given the actual
1600  // position and the segment index
1601  Scalar evalDerivative3_(Scalar x, size_t i) const
1602  {
1603  // See http://en.wikipedia.org/wiki/Cubic_Hermite_spline
1604  Scalar delta = h_(i + 1);
1605  Scalar t = (x - x_(i))/delta;
1606  Scalar alpha = 1 / delta;
1607 
1608  return
1609  alpha*alpha*alpha
1610  *(h00_prime3_(t)*y_(i)
1611  + h10_prime3_(t)*slope_(i)*delta
1612  + h01_prime3_(t)*y_(i + 1)
1613  + h11_prime3_(t)*slope_(i + 1)*delta);
1614  }
1615 
1616  // hermite basis functions
1617  Scalar h00_(Scalar t) const
1618  { return (2*t - 3)*t*t + 1; }
1619 
1620  Scalar h10_(Scalar t) const
1621  { return ((t - 2)*t + 1)*t; }
1622 
1623  Scalar h01_(Scalar t) const
1624  { return (-2*t + 3)*t*t; }
1625 
1626  Scalar h11_(Scalar t) const
1627  { return (t - 1)*t*t; }
1628 
1629  // first derivative of the hermite basis functions
1630  Scalar h00_prime_(Scalar t) const
1631  { return (3*2*t - 2*3)*t; }
1632 
1633  Scalar h10_prime_(Scalar t) const
1634  { return (3*t - 2*2)*t + 1; }
1635 
1636  Scalar h01_prime_(Scalar t) const
1637  { return (-3*2*t + 2*3)*t; }
1638 
1639  Scalar h11_prime_(Scalar t) const
1640  { return (3*t - 2)*t; }
1641 
1642  // second derivative of the hermite basis functions
1643  Scalar h00_prime2_(Scalar t) const
1644  { return 2*3*2*t - 2*3; }
1645 
1646  Scalar h10_prime2_(Scalar t) const
1647  { return 2*3*t - 2*2; }
1648 
1649  Scalar h01_prime2_(Scalar t) const
1650  { return -2*3*2*t + 2*3; }
1651 
1652  Scalar h11_prime2_(Scalar t) const
1653  { return 2*3*t - 2; }
1654 
1655  // third derivative of the hermite basis functions
1656  Scalar h00_prime3_(Scalar /*t*/) const
1657  { return 2*3*2; }
1658 
1659  Scalar h10_prime3_(Scalar /*t*/) const
1660  { return 2*3; }
1661 
1662  Scalar h01_prime3_(Scalar /*t*/) const
1663  { return -2*3*2; }
1664 
1665  Scalar h11_prime3_(Scalar /*t*/) const
1666  { return 2*3; }
1667 
1668  // returns the monotonicality of an interval of a spline segment
1669  //
1670  // The return value have the following meaning:
1671  //
1672  // 3: spline is constant within interval [x0, x1]
1673  // 1: spline is monotonously increasing in the specified interval
1674  // 0: spline is not monotonic (or constant) in the specified interval
1675  // -1: spline is monotonously decreasing in the specified interval
1676  int monotonic_(size_t i, Scalar x0, Scalar x1, int &r) const
1677  {
1678  // coefficients of derivative in monomial basis
1679  Scalar a = 3*a_(i);
1680  Scalar b = 2*b_(i);
1681  Scalar c = c_(i);
1682 
1683  if (std::abs(a) < 1e-20 && std::abs(b) < 1e-20 && std::abs(c) < 1e-20)
1684  return 3; // constant in interval, r stays unchanged!
1685 
1686  Scalar disc = b*b - 4*a*c;
1687  if (disc < 0) {
1688  // discriminant of derivative is smaller than 0, i.e. the
1689  // segment's derivative does not exhibit any extrema.
1690  if (x0*(x0*a + b) + c > 0) {
1691  r = (r==3 || r == 1)?1:0;
1692  return 1;
1693  }
1694  else {
1695  r = (r==3 || r == -1)?-1:0;
1696  return -1;
1697  }
1698  }
1699  disc = std::sqrt(disc);
1700  Scalar xE1 = (-b + disc)/(2*a);
1701  Scalar xE2 = (-b - disc)/(2*a);
1702 
1703  if (std::abs(disc) < 1e-30) {
1704  // saddle point -> no extrema
1705  if (std::abs(xE1 - x0) < 1e-30)
1706  // make sure that we're not picking the saddle point
1707  // to determine whether we're monotonically increasing
1708  // or decreasing
1709  x0 = x1;
1710  if (x0*(x0*a + b) + c > 0) {
1711  r = (r==3 || r == 1)?1:0;
1712  return 1;
1713  }
1714  else {
1715  r = (r==3 || r == -1)?-1:0;
1716  return -1;
1717  }
1718  };
1719  if ((x0 < xE1 && xE1 < x1) ||
1720  (x0 < xE2 && xE2 < x1))
1721  {
1722  // there is an extremum in the range (x0, x1)
1723  r = 0;
1724  return 0;
1725  }
1726  // no extremum in range (x0, x1)
1727  x0 = (x0 + x1)/2; // pick point in the middle of the interval
1728  // to avoid extrema on the boundaries
1729  if (x0*(x0*a + b) + c > 0) {
1730  r = (r==3 || r == 1)?1:0;
1731  return 1;
1732  }
1733  else {
1734  r = (r==3 || r == -1)?-1:0;
1735  return -1;
1736  }
1737  }
1738 
1743  template <class Evaluation>
1744  size_t intersectSegment_(Evaluation *sol,
1745  size_t segIdx,
1746  const Evaluation& a,
1747  const Evaluation& b,
1748  const Evaluation& c,
1749  const Evaluation& d,
1750  Scalar x0 = -1e100, Scalar x1 = 1e100) const
1751  {
1752  int n = Opm::invertCubicPolynomial(sol,
1753  a_(segIdx) - a,
1754  b_(segIdx) - b,
1755  c_(segIdx) - c,
1756  d_(segIdx) - d);
1757  x0 = std::max(x_(segIdx), x0);
1758  x1 = std::min(x_(segIdx+1), x1);
1759 
1760  // filter the intersections outside of the specified interval
1761  size_t k = 0;
1762  for (int j = 0; j < n; ++j) {
1763  if (x0 <= sol[j] && sol[j] <= x1) {
1764  sol[k] = sol[j];
1765  ++k;
1766  }
1767  }
1768  return k;
1769  }
1770 
1771  // find the segment index for a given x coordinate
1772  size_t segmentIdx_(Scalar x) const
1773  {
1774  // bisection
1775  size_t iLow = 0;
1776  size_t iHigh = numSamples() - 1;
1777 
1778  while (iLow + 1 < iHigh) {
1779  size_t i = (iLow + iHigh) / 2;
1780  if (x_(i) > x)
1781  iHigh = i;
1782  else
1783  iLow = i;
1784  };
1785  return iLow;
1786  }
1787 
1791  Scalar h_(size_t i) const
1792  {
1793  assert(x_(i) > x_(i-1)); // the sampling points must be given
1794  // in ascending order
1795  return x_(i) - x_(i - 1);
1796  }
1797 
1801  Scalar x_(size_t i) const
1802  { return xPos_[i]; }
1803 
1807  Scalar y_(size_t i) const
1808  { return yPos_[i]; }
1809 
1814  Scalar slope_(size_t i) const
1815  { return slopeVec_[i]; }
1816 
1817  // returns the coefficient in front of the x^3 term. In Stoer this
1818  // is delta.
1819  Scalar a_(size_t i) const
1820  { return evalDerivative3_(/*x=*/Scalar(0.0), i)/6.0; }
1821 
1822  // returns the coefficient in front of the x^2 term In Stoer this
1823  // is gamma.
1824  Scalar b_(size_t i) const
1825  { return evalDerivative2_(/*x=*/Scalar(0.0), i)/2.0; }
1826 
1827  // returns the coefficient in front of the x^1 term. In Stoer this
1828  // is beta.
1829  Scalar c_(size_t i) const
1830  { return evalDerivative_(/*x=*/Scalar(0.0), i); }
1831 
1832  // returns the coefficient in front of the x^0 term. In Stoer this
1833  // is alpha.
1834  Scalar d_(size_t i) const
1835  { return eval_(/*x=*/Scalar(0.0), i); }
1836 
1837  Vector xPos_;
1838  Vector yPos_;
1839  Vector slopeVec_;
1840 };
1841 }
1842 
1843 #endif
int monotonic() const
Same as monotonic(x0, x1), but with the entire range of the spline as interval.
Definition: Spline.hpp:1058
void solve(XVector &x, const BVector &b) const
Calculate the solution for a linear system of equations.
Definition: TridiagonalMatrix.hpp:712
Scalar h00_(Scalar t) const
Definition: Spline.hpp:1617
Scalar h10_(Scalar t) const
Definition: Spline.hpp:1620
Provides a tridiagonal matrix that also supports non-zero entries in the upper right and lower left...
Definition: TridiagonalMatrix.hpp:48
void assignSamplingPoints_(DestVector &destX, DestVector &destY, const SourceVector &srcX, const SourceVector &srcY, int numSamples)
Set the sampling point vectors.
Definition: Spline.hpp:1206
Definition: Spline.hpp:104
Spline()
Default constructor for a spline.
Definition: Spline.hpp:112
Scalar evalDerivative2_(Scalar x, size_t i) const
Definition: Spline.hpp:1584
int monotonic(Scalar x0, Scalar x1, OPM_OPTIM_UNUSED bool extrapolate=false) const
Returns 1 if the spline is monotonically increasing, -1 if the spline is mononously decreasing and 0 ...
Definition: Spline.hpp:992
Scalar yLast() const
Return the y value of the rightmost sampling point.
Definition: Spline.hpp:741
void setXYArrays(size_t nSamples, const ScalarArrayX &x, const ScalarArrayY &y, Scalar m0, Scalar m1, bool sortInputs=false)
Set the sampling points and the boundary slopes of a full spline using C-style arrays.
Definition: Spline.hpp:327
Scalar d_(size_t i) const
Definition: Spline.hpp:1834
Vector slopeVec_
Definition: Spline.hpp:1839
Definition: Air_Mesitylene.hpp:31
void assignFromArrayList_(DestVector &destX, DestVector &destY, const ListIterator &srcBegin, const ListIterator &srcEnd, int numSamples)
Definition: Spline.hpp:1226
ComparatorX_(const std::vector< Scalar > &x)
Definition: Spline.hpp:1067
void setContainerOfPoints(const XYContainer &points, SplineType splineType=Natural, bool sortInputs=false)
Set the sampling points of a natural spline using a STL-compatible container of array-like objects...
Definition: Spline.hpp:638
Scalar y_(size_t i) const
Returns the y coordinate of the i-th sampling point.
Definition: Spline.hpp:1807
bool applies(Scalar x) const
Return true iff the given x is in range [x1, xn].
Definition: Spline.hpp:715
Evaluation< Scalar, VarSetTag, numVars > max(const Evaluation< Scalar, VarSetTag, numVars > &x1, const Evaluation< Scalar, VarSetTag, numVars > &x2)
Definition: Math.hpp:114
Scalar h11_prime_(Scalar t) const
Definition: Spline.hpp:1639
Evaluation< Scalar, VarSetTag, numVars > sqrt(const Evaluation< Scalar, VarSetTag, numVars > &x)
Definition: Math.hpp:278
void setContainerOfTuples(const XYContainer &points, SplineType splineType=Natural, bool sortInputs=false)
Set the sampling points of a natural spline using a STL-compatible container of tuple-like objects...
Definition: Spline.hpp:684
void setContainerOfTuples(const XYContainer &points, Scalar m0, Scalar m1, bool sortInputs=false)
Set the sampling points and the boundary slopes of a full spline using a STL-compatible container of ...
Definition: Spline.hpp:474
void setContainerOfPoints(const XYContainer &points, Scalar m0, Scalar m1, bool sortInputs=false)
Set the sampling points and the boundary slopes of a full spline using a STL-compatible container of ...
Definition: Spline.hpp:432
Spline(size_t nSamples, const ScalarArray &x, const ScalarArray &y, Scalar m0, Scalar m1, bool sortInputs=false)
Convenience constructor for a full spline.
Definition: Spline.hpp:197
SplineType
The type of the spline to be created.
Definition: Spline.hpp:100
Spline(size_t nSamples, const ScalarArrayX &x, const ScalarArrayY &y, SplineType splineType=Natural, bool sortInputs=false)
Convenience constructor for a natural or a periodic spline.
Definition: Spline.hpp:139
Spline(const PointContainer &points, Scalar m0, Scalar m1, bool sortInputs=false)
Convenience constructor for a full spline.
Definition: Spline.hpp:248
Scalar evalSecondDerivative(Scalar x, bool extrapolate=false) const
Evaluate the spline's second derivative at a given position.
Definition: Spline.hpp:900
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:1318
Scalar evalDerivative3_(Scalar x, size_t i) const
Definition: Spline.hpp:1601
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:271
Evaluation eval(const Evaluation &x, bool extrapolate=false) const
Evaluate the spline for a given function evaluation.
Definition: Spline.hpp:842
Vector xPos_
Definition: Spline.hpp:1837
int invertCubicPolynomial(SolContainer *sol, Scalar a, Scalar b, Scalar c, Scalar d)
Invert a cubic polynomial analytically.
Definition: PolynomialUtils.hpp:151
void makePeriodicSpline_()
Create a periodic spline from the already set sampling points.
Definition: Spline.hpp:1175
Spline(size_t nSamples, const PointArray &points, Scalar m0, Scalar m1, bool sortInputs=false)
Convenience constructor for a full spline.
Definition: Spline.hpp:215
void makeMonotonicSpline_(Vector &slopes)
Create a monotonic spline from the already set sampling points.
Definition: Spline.hpp:1426
Scalar h11_(Scalar t) const
Definition: Spline.hpp:1626
Scalar xMax() const
Return the x value of the rightmost sampling point.
Definition: Spline.hpp:729
Scalar evalDerivative_(Scalar x, size_t i) const
Definition: Spline.hpp:1567
Scalar eval_(Scalar x, size_t i) const
Definition: Spline.hpp:1529
void assignFromTupleList_(DestVector &destX, DestVector &destY, ListIterator srcBegin, ListIterator srcEnd, int numSamples)
Set the sampling points.
Definition: Spline.hpp:1260
Scalar xMin() const
Return the x value of the leftmost sampling point.
Definition: Spline.hpp:723
void setSlopesFromMoments_(SlopeVector &slopes, const MomentsVector &moments)
Convert the moments at the sample points to slopes.
Definition: Spline.hpp:1475
Scalar h01_prime3_(Scalar) const
Definition: Spline.hpp:1662
void setArrayOfPoints(size_t nSamples, const PointArray &points, SplineType splineType=Natural, bool sortInputs=false)
Set the sampling points of a natural spline using a C-style array.
Definition: Spline.hpp:596
void setXYContainers(const ScalarContainerX &x, const ScalarContainerY &y, Scalar m0, Scalar m1, bool sortInputs=false)
Set the sampling points and the boundary slopes of a full spline using STL-compatible containers...
Definition: Spline.hpp:361
Scalar yFirst() const
Return the y value of the leftmost sampling point.
Definition: Spline.hpp:735
Definition: Spline.hpp:102
Scalar h00_prime_(Scalar t) const
Definition: Spline.hpp:1630
void sortInput_()
Sort the sample points in ascending order of their x value.
Definition: Spline.hpp:1080
Scalar h00_prime2_(Scalar t) const
Definition: Spline.hpp:1643
Vector yPos_
Definition: Spline.hpp:1838
Provides a tridiagonal matrix that also supports non-zero entries in the upper right and lower left...
Scalar h00_prime3_(Scalar) const
Definition: Spline.hpp:1656
size_t segmentIdx_(Scalar x) const
Definition: Spline.hpp:1772
Evaluation< Scalar, VarSetTag, numVars > abs(const Evaluation< Scalar, VarSetTag, numVars > &)
Definition: Math.hpp:41
Scalar h01_(Scalar t) const
Definition: Spline.hpp:1623
Definition: Spline.hpp:103
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:952
Class implementing cubic splines.
Definition: Spline.hpp:89
Provides free functions to invert polynomials of degree 1, 2 and 3.
Scalar h01_prime2_(Scalar t) const
Definition: Spline.hpp:1649
Evaluation< Scalar, VarSetTag, numVars > min(const Evaluation< Scalar, VarSetTag, numVars > &x1, const Evaluation< Scalar, VarSetTag, numVars > &x2)
Definition: Math.hpp:61
Scalar evalDerivative(Scalar x, bool extrapolate=false) const
Evaluate the spline's derivative at a given position.
Definition: Spline.hpp:875
Definition: Spline.hpp:101
Scalar h_(size_t i) const
Returns x[i] - x[i - 1].
Definition: Spline.hpp:1791
void setArrayOfPoints(size_t nSamples, const PointArray &points, Scalar m0, Scalar m1, bool sortInputs=false)
Set the sampling points and the boundary slopes of a full spline using a C-style array.
Definition: Spline.hpp:395
#define OPM_OPTIM_UNUSED
Definition: Unused.hpp:42
Spline(const PointContainer &points, SplineType splineType=Natural, bool sortInputs=false)
Convenience constructor for a natural or a periodic spline.
Definition: Spline.hpp:181
const std::vector< Scalar > & x_
Definition: Spline.hpp:1074
void reverseSamplingPoints_()
Reverse order of the elements in the arrays which contain the sampling points.
Definition: Spline.hpp:1108
Scalar eval(Scalar x, bool extrapolate=false) const
Evaluate the spline at a given position.
Definition: Spline.hpp:809
bool operator()(unsigned idxA, unsigned idxB) const
Definition: Spline.hpp:1071
Scalar a_(size_t i) const
Definition: Spline.hpp:1819
void makeNaturalSpline_()
Create a natural spline from the already set sampling points.
Definition: Spline.hpp:1154
Scalar h11_prime3_(Scalar) const
Definition: Spline.hpp:1665
void makeFullSpline_(Scalar m0, Scalar m1)
Create a natural spline from the already set sampling points.
Definition: Spline.hpp:1133
Spline(const ScalarContainer &x, const ScalarContainer &y, SplineType splineType=Natural, bool sortInputs=false)
Convenience constructor for a natural or a periodic spline.
Definition: Spline.hpp:168
Scalar x_(size_t i) const
Returns the y coordinate of the i-th sampling point.
Definition: Spline.hpp:1801
Scalar h10_prime3_(Scalar) const
Definition: Spline.hpp:1659
int monotonic_(size_t i, Scalar x0, Scalar x1, int &r) const
Definition: Spline.hpp:1676
Scalar evalThirdDerivative(Scalar x, bool extrapolate=false) const
Evaluate the spline's third derivative at a given position.
Definition: Spline.hpp:921
void printCSV(Scalar xi0, Scalar xi1, size_t k, std::ostream &os=std::cout) const
Prints k tuples of the format (x, y, dx/dy, isMonotonic) to stdout.
Definition: Spline.hpp:761
Scalar b_(size_t i) const
Definition: Spline.hpp:1824
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, throws Opm::MathError exception if there is more or less than one solution.
Definition: Spline.hpp:937
Scalar h10_prime_(Scalar t) const
Definition: Spline.hpp:1633
Helper class needed to sort the input sampling points.
Definition: Spline.hpp:1065
void setNumSamples_(size_t nSamples)
Resizes the internal vectors to store the sample points.
Definition: Spline.hpp:1121
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:1364
Scalar h01_prime_(Scalar t) const
Definition: Spline.hpp:1636
size_t rows() const
Return the number of rows of the matrix.
Definition: TridiagonalMatrix.hpp:138
Spline(size_t nSamples, const PointArray &points, SplineType splineType=Natural, bool sortInputs=false)
Convenience constructor for a natural or a periodic spline.
Definition: Spline.hpp:154
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:1294
Scalar h10_prime2_(Scalar t) const
Definition: Spline.hpp:1646
Scalar c_(size_t i) const
Definition: Spline.hpp:1829
Provides the OPM_UNUSED macro.
Evaluation eval_(const Evaluation &x, size_t i) const
Definition: Spline.hpp:1545
Spline(const ScalarContainerX &x, const ScalarContainerY &y, Scalar m0, Scalar m1, bool sortInputs=false)
Convenience constructor for a full spline.
Definition: Spline.hpp:232
Scalar slope_(size_t i) const
Returns the slope (i.e. first derivative) of the spline at the i-th sampling point.
Definition: Spline.hpp:1814
Scalar h11_prime2_(Scalar t) const
Definition: Spline.hpp:1652
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:125
size_t intersectSegment_(Evaluation *sol, size_t segIdx, const Evaluation &a, const Evaluation &b, const Evaluation &c, const Evaluation &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: Spline.hpp:1744
void setXYArrays(size_t nSamples, const ScalarArrayX &x, const ScalarArrayY &y, SplineType splineType=Natural, bool sortInputs=false)
Set the sampling points natural spline using C-style arrays.
Definition: Spline.hpp:515
void setXYContainers(const ScalarContainerX &x, const ScalarContainerY &y, SplineType splineType=Natural, bool sortInputs=false)
Set the sampling points of a natural spline using STL-compatible containers.
Definition: Spline.hpp:556
size_t numSamples() const
Returns the number of sampling points.
Definition: Spline.hpp:257