93 typedef std::vector<Scalar> Vector;
127 Scalar y0, Scalar y1,
128 Scalar m0, Scalar m1)
129 {
set(x0, x1, y0, y1, m0, m1); }
140 template <
class ScalarArrayX,
class ScalarArrayY>
142 const ScalarArrayX& x,
143 const ScalarArrayY& y,
145 bool sortInputs =
true)
146 { this->
setXYArrays(nSamples, x, y, splineType, sortInputs); }
156 template <
class Po
intArray>
158 const PointArray& points,
160 bool sortInputs =
true)
171 template <
class ScalarContainer>
173 const ScalarContainer& y,
175 bool sortInputs =
true)
185 template <
class Po
intContainer>
186 explicit Spline(
const PointContainer& points,
188 bool sortInputs =
true)
201 template <
class ScalarArray>
203 const ScalarArray& x,
204 const ScalarArray& y,
207 bool sortInputs =
true)
208 { this->
setXYArrays(nSamples, x, y, m0, m1, sortInputs); }
219 template <
class Po
intArray>
221 const PointArray& points,
224 bool sortInputs =
true)
236 template <
class ScalarContainerX,
class ScalarContainerY>
238 const ScalarContainerY& y,
241 bool sortInputs =
true)
252 template <
class Po
intContainer>
256 bool sortInputs =
true)
270 void set(Scalar x0, Scalar x1,
271 Scalar y0, Scalar y1,
272 Scalar m0, Scalar m1)
325 template <
class ScalarArrayX,
class ScalarArrayY>
327 const ScalarArrayX& x,
328 const ScalarArrayY& y,
329 Scalar m0, Scalar m1,
330 bool sortInputs =
true)
332 assert(nSamples > 1);
335 for (
size_t i = 0; i < nSamples; ++i) {
359 template <
class ScalarContainerX,
class ScalarContainerY>
361 const ScalarContainerY& y,
362 Scalar m0, Scalar m1,
363 bool sortInputs =
true)
365 assert(x.size() == y.size());
366 assert(x.size() > 1);
370 std::copy(x.begin(), x.end(), xPos_.begin());
371 std::copy(y.begin(), y.end(), yPos_.begin());
393 template <
class Po
intArray>
395 const PointArray& points,
398 bool sortInputs =
true)
402 assert(nSamples > 1);
405 for (
size_t i = 0; i < nSamples; ++i) {
406 xPos_[i] = points[i][0];
407 yPos_[i] = points[i][1];
430 template <
class XYContainer>
434 bool sortInputs =
true)
438 assert(points.size() > 1);
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) {
472 template <
class XYContainer>
476 bool sortInputs =
true)
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);
513 template <
class ScalarArrayX,
class ScalarArrayY>
515 const ScalarArrayX& x,
516 const ScalarArrayY& y,
518 bool sortInputs =
true)
520 assert(nSamples > 1);
523 for (
size_t i = 0; i < nSamples; ++i) {
533 if (splineType == Periodic)
535 else if (splineType == Natural)
537 else if (splineType == Monotonic)
540 throw std::runtime_error(
"Spline type "+std::to_string(
int(splineType))+
" not supported at this place");
554 template <
class ScalarContainerX,
class ScalarContainerY>
556 const ScalarContainerY& y,
558 bool sortInputs =
true)
560 assert(x.size() == y.size());
561 assert(x.size() > 1);
564 std::copy(x.begin(), x.end(), xPos_.begin());
565 std::copy(y.begin(), y.end(), yPos_.begin());
572 if (splineType == Periodic)
574 else if (splineType == Natural)
576 else if (splineType == Monotonic)
579 throw std::runtime_error(
"Spline type "+std::to_string(
int(splineType))+
" not supported at this place");
594 template <
class Po
intArray>
596 const PointArray& points,
598 bool sortInputs =
true)
602 assert(nSamples > 1);
605 for (
size_t i = 0; i < nSamples; ++i) {
606 xPos_[i] = points[i][0];
607 yPos_[i] = points[i][1];
615 if (splineType == Periodic)
617 else if (splineType == Natural)
619 else if (splineType == Monotonic)
622 throw std::runtime_error(
"Spline type "+std::to_string(
int(splineType))+
" not supported at this place");
636 template <
class XYContainer>
639 bool sortInputs =
true)
643 assert(points.size() > 1);
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) {
658 if (splineType == Periodic)
660 else if (splineType == Natural)
662 else if (splineType == Monotonic)
665 throw std::runtime_error(
"Spline type "+std::to_string(
int(splineType))+
" not supported at this place");
682 template <
class XYContainer>
685 bool sortInputs =
true)
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);
701 if (splineType == Periodic)
703 else if (splineType == Natural)
705 else if (splineType == Monotonic)
708 throw std::runtime_error(
"Spline type "+std::to_string(
int(splineType))+
" not supported at this place");
714 template <
class Evaluation>
722 {
return xPos_.size(); }
727 Scalar
xAt(
size_t sampleIdx)
const
728 {
return x_(sampleIdx); }
734 {
return y_(sampleIdx); }
753 void printCSV(Scalar xi0, Scalar xi1,
size_t k, std::ostream& os)
const;
766 template <
class Evaluation>
767 Evaluation
eval(
const Evaluation& x,
bool extrapolate =
false)
const
769 if (!extrapolate && !
applies(x))
775 Scalar m = evalDerivative_(
xAt(0), 0);
777 return y0 + m*(x -
xAt(0));
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)),
783 return y0 + m*(x -
xAt(
static_cast<size_t>(
numSamples() - 1)));
787 return eval_(x, segmentIdx_(scalarValue(x)));
802 template <
class Evaluation>
805 if (!extrapolate && !
applies(x))
806 throw NumericalProblem(
"Tried to evaluate the derivative of a spline outside of its range");
811 return evalDerivative_(
xAt(0), 0);
816 return evalDerivative_(x, segmentIdx_(scalarValue(x)));
831 template <
class Evaluation>
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)
839 return evalDerivative2_(x, segmentIdx_(scalarValue(x)));
854 template <
class Evaluation>
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)
862 return evalDerivative3_(x, segmentIdx_(scalarValue(x)));
871 template <
class Evaluation>
875 const Evaluation& d)
const
886 template <
class Evaluation>
891 const Evaluation& d)
const
895 Evaluation tmpSol[3], sol = 0;
897 size_t iFirst = segmentIdx_(x0);
898 size_t iLast = segmentIdx_(x1);
899 for (
size_t i = iFirst; i <= iLast; ++i)
907 throw std::runtime_error(
"Spline has more than one intersection");
912 throw std::runtime_error(
"Spline has no intersection");
926 [[maybe_unused]]
bool extrapolate =
false)
const
928 assert(std::abs(x0 - x1) > 1e-30);
939 Scalar m = evalDerivative_(
xAt(0), 0);
940 if (std::abs(m) < 1e-20)
945 size_t i = segmentIdx_(x0);
946 if (
x_(i + 1) >= x1) {
949 monotonic_(i, x0, x1, r);
955 monotonic_(i, x0,
x_(i+1), r);
960 size_t iEnd = segmentIdx_(x1);
961 for (; i < iEnd - 1; ++i) {
962 monotonic_(i,
x_(i),
x_(i + 1), r);
975 return (r < 0 || r==3) ? -1 : 0;
977 return r > 0 ? 1 : 0;
983 monotonic_(iEnd,
x_(iEnd), x1, r);
1001 explicit ComparatorX_(
const std::vector<Scalar>& x)
1005 bool operator ()(
unsigned idxA,
unsigned idxB)
const
1006 {
return x_.at(idxA) < x_.at(idxB); }
1008 const std::vector<Scalar>& x_;
1019 std::vector<unsigned> idxVector(n);
1020 for (
unsigned i = 0; i < n; ++i)
1026 std::sort(idxVector.begin(), idxVector.end(),
cmp);
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]];
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]);
1057 xPos_.resize(nSamples);
1058 yPos_.resize(nSamples);
1059 slopeVec_.resize(nSamples);
1077 M.solve(moments, d);
1098 M.solve(moments, d);
1122 M.solve(moments, d);
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];
1141 template <
class DestVector,
class SourceVector>
1144 const SourceVector& srcX,
1145 const SourceVector& srcY,
1148 assert(nSamples >= 2);
1152 for (
unsigned i = 0; i < nSamples; ++i) {
1154 if (srcX[0] > srcX[nSamples - 1])
1155 idx = nSamples - i - 1;
1156 destX[i] = srcX[idx];
1157 destY[i] = srcY[idx];
1161 template <
class DestVector,
class ListIterator>
1162 void assignFromArrayList_(DestVector& destX,
1164 const ListIterator& srcBegin,
1165 const ListIterator& srcEnd,
1168 assert(nSamples >= 2);
1171 ListIterator it = srcBegin;
1173 bool reverse =
false;
1174 if ((*srcBegin)[0] > (*it)[0])
1179 for (
unsigned i = 0; it != srcEnd; ++i, ++it) {
1182 idx = nSamples - i - 1;
1183 destX[idx] = (*it)[0];
1184 destY[idx] = (*it)[1];
1195 template <
class DestVector,
class ListIterator>
1198 ListIterator srcBegin,
1199 ListIterator srcEnd,
1202 assert(nSamples >= 2);
1208 ListIterator it = srcBegin;
1210 bool reverse =
false;
1211 if (std::get<0>(*srcBegin) > std::get<0>(*it))
1216 for (
unsigned i = 0; it != srcEnd; ++i, ++it) {
1219 idx = nSamples - i - 1;
1220 destX[idx] = std::get<0>(*it);
1221 destY[idx] = std::get<1>(*it);
1229 template <
class Vector,
class Matrix>
1237 d[0] = 6/
h_(1) * ( (
y_(1) -
y_(0))/
h_(1) - m0);
1246 (m1 - (
y_(n) -
y_(n - 1))/
h_(n));
1253 template <
class Vector,
class Matrix>
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;
1267 6 / (
h_(i) +
h_(i + 1))
1269 ( (
y_(i + 1) -
y_(i))/
h_(i + 1) - (
y_(i) -
y_(i - 1))/
h_(i));
1273 M[i][i + 1] = lambda_i;
1278 Scalar lambda_0 = 0;
1299 template <
class Matrix,
class Vector>
1308 assert(M.rows() == n);
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;
1315 6 / (
h_(i) +
h_(i + 1))
1317 ( (
y_(i + 1) -
y_(i))/
h_(i + 1) - (
y_(i) -
y_(i - 1))/
h_(i));
1321 M[i-1][i] = lambda_i;
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;
1347 M[n-1][0] = lambda_n;
1361 template <
class Vector>
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));
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];
1378 for (
size_t k = 0; k < n - 1; ++k) {
1379 if (std::abs(delta[k]) < 1e-50) {
1387 Scalar alpha = slopes[k] / delta[k];
1388 Scalar beta = slopes[k + 1] / delta[k];
1390 if (alpha < 0 || (k > 0 && slopes[k] / delta[k - 1] < 0)) {
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];
1410 template <
class MomentsVector,
class SlopeVector>
1421 Scalar h = this->
h_(n - 1);
1426 (
y_(n - 1) -
y_(n - 2))/h
1428 h/6*(moments[n-1] - moments[n - 2]);
1433 moments[n - 1] * x*x / (2 * h)
1439 for (
size_t i = 0; i < n - 1; ++ i) {
1442 Scalar h_i = this->
h_(i + 1);
1447 (
y_(i+1) -
y_(i))/h_i
1449 h_i/6*(moments[i+1] - moments[i]);
1452 - moments[i] * x_i1*x_i1 / (2 * h_i)
1459 slopes[n - 1] = mRight;
1465 template <
class Evaluation>
1466 Evaluation eval_(
const Evaluation& x,
size_t i)
const
1469 Scalar delta =
h_(i + 1);
1470 Evaluation t = (x -
x_(i))/delta;
1474 + h10_(t) *
slope_(i)*delta
1475 + h01_(t) *
y_(i + 1)
1476 + h11_(t) *
slope_(i + 1)*delta;
1481 template <
class Evaluation>
1482 Evaluation evalDerivative_(
const Evaluation& x,
size_t i)
const
1485 Scalar delta =
h_(i + 1);
1486 Evaluation t = (x -
x_(i))/delta;
1487 Evaluation alpha = 1 / delta;
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);
1499 template <
class Evaluation>
1500 Evaluation evalDerivative2_(
const Evaluation& x,
size_t i)
const
1503 Scalar delta =
h_(i + 1);
1504 Evaluation t = (x -
x_(i))/delta;
1505 Evaluation alpha = 1 / delta;
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);
1517 template <
class Evaluation>
1518 Evaluation evalDerivative3_(
const Evaluation& x,
size_t i)
const
1521 Scalar delta =
h_(i + 1);
1522 Evaluation t = (x -
x_(i))/delta;
1523 Evaluation alpha = 1 / delta;
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);
1534 template <
class Evaluation>
1535 Evaluation h00_(
const Evaluation& t)
const
1536 {
return (2*t - 3)*t*t + 1; }
1538 template <
class Evaluation>
1539 Evaluation h10_(
const Evaluation& t)
const
1540 {
return ((t - 2)*t + 1)*t; }
1542 template <
class Evaluation>
1543 Evaluation h01_(
const Evaluation& t)
const
1544 {
return (-2*t + 3)*t*t; }
1546 template <
class Evaluation>
1547 Evaluation h11_(
const Evaluation& t)
const
1548 {
return (t - 1)*t*t; }
1551 template <
class Evaluation>
1552 Evaluation h00_prime_(
const Evaluation& t)
const
1553 {
return (3*2*t - 2*3)*t; }
1555 template <
class Evaluation>
1556 Evaluation h10_prime_(
const Evaluation& t)
const
1557 {
return (3*t - 2*2)*t + 1; }
1559 template <
class Evaluation>
1560 Evaluation h01_prime_(
const Evaluation& t)
const
1561 {
return (-3*2*t + 2*3)*t; }
1563 template <
class Evaluation>
1564 Evaluation h11_prime_(
const Evaluation& t)
const
1565 {
return (3*t - 2)*t; }
1568 template <
class Evaluation>
1569 Evaluation h00_prime2_(
const Evaluation& t)
const
1570 {
return 2*3*2*t - 2*3; }
1572 template <
class Evaluation>
1573 Evaluation h10_prime2_(
const Evaluation& t)
const
1574 {
return 2*3*t - 2*2; }
1576 template <
class Evaluation>
1577 Evaluation h01_prime2_(
const Evaluation& t)
const
1578 {
return -2*3*2*t + 2*3; }
1580 template <
class Evaluation>
1581 Evaluation h11_prime2_(
const Evaluation& t)
const
1582 {
return 2*3*t - 2; }
1585 template <
class Evaluation>
1586 Scalar h00_prime3_(
const Evaluation&)
const
1589 template <
class Evaluation>
1590 Scalar h10_prime3_(
const Evaluation&)
const
1593 template <
class Evaluation>
1594 Scalar h01_prime3_(
const Evaluation&)
const
1597 template <
class Evaluation>
1598 Scalar h11_prime3_(
const Evaluation&)
const
1609 int monotonic_(
size_t i, Scalar x0, Scalar x1,
int& r)
const
1616 if (std::abs(a) < 1e-20 && std::abs(b) < 1e-20 && std::abs(c) < 1e-20)
1619 Scalar disc = b*b - 4*a*c;
1623 if (x0*(x0*a + b) + c > 0) {
1624 r = (r==3 || r == 1)?1:0;
1628 r = (r==3 || r == -1)?-1:0;
1632 disc = std::sqrt(disc);
1633 Scalar xE1 = (-b + disc)/(2*a);
1634 Scalar xE2 = (-b - disc)/(2*a);
1636 if (std::abs(disc) < 1e-30) {
1638 if (std::abs(xE1 - x0) < 1e-30)
1643 if (x0*(x0*a + b) + c > 0) {
1644 r = (r==3 || r == 1)?1:0;
1648 r = (r==3 || r == -1)?-1:0;
1652 if ((x0 < xE1 && xE1 < x1) ||
1653 (x0 < xE2 && xE2 < x1))
1662 if (x0*(x0*a + b) + c > 0) {
1663 r = (r==3 || r == 1)?1:0;
1667 r = (r==3 || r == -1)?-1:0;
1676 template <
class Evaluation>
1679 const Evaluation& a,
1680 const Evaluation& b,
1681 const Evaluation& c,
1682 const Evaluation& d,
1683 Scalar x0 = -1e30, Scalar x1 = 1e30)
const
1691 x0 = std::max(
x_(segIdx), x0);
1692 x1 = std::min(
x_(segIdx+1), x1);
1696 for (
unsigned j = 0; j < n; ++j) {
1697 if (x0 <= sol[j] && sol[j] <= x1) {
1706 size_t segmentIdx_(Scalar x)
const
1712 while (iLow + 1 < iHigh) {
1713 size_t i = (iLow + iHigh) / 2;
1725 Scalar
h_(
size_t i)
const
1727 assert(
x_(i) >
x_(i-1));
1729 return x_(i) -
x_(i - 1);
1735 Scalar
x_(
size_t i)
const
1736 {
return xPos_[i]; }
1741 Scalar
y_(
size_t i)
const
1742 {
return yPos_[i]; }
1749 {
return slopeVec_[i]; }
1753 Scalar a_(
size_t i)
const
1754 {
return evalDerivative3_(Scalar(0.0), i)/6.0; }
1758 Scalar b_(
size_t i)
const
1759 {
return evalDerivative2_(Scalar(0.0), i)/2.0; }
1763 Scalar c_(
size_t i)
const
1764 {
return evalDerivative_(Scalar(0.0), i); }
1768 Scalar d_(
size_t i)
const
1769 {
return eval_(Scalar(0.0), i); }
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:148
Helper class needed to sort the input sampling points.
Definition Spline.hpp:1000