10#include <bsplinebasis1d.h>
19BSplineBasis1D::BSplineBasis1D()
23BSplineBasis1D::BSplineBasis1D(
const std::vector<double>& knots,
27 targetNumBasisfunctions((degree + 1) + 2 * degree + 1)
33 throw Exception(
"BSplineBasis1D::BSplineBasis1D: Knot vector is not regular.");
37SparseVector BSplineBasis1D::eval(
double x)
const
39 SparseVector values(getNumBasisFunctions());
41 if (!insideSupport(x))
47 std::vector<int> indexSupported = indexSupportedBasisfunctions(x);
48 values.reserve(indexSupported.size());
51 for (
auto it = indexSupported.begin(); it != indexSupported.end(); ++it)
53 double val = deBoorCox(x, * it, degree);
55 if (fabs(val) > 1e-12)
57 values.insert(* it) = val;
75SparseVector BSplineBasis1D::evalDerivative(
double x,
int r)
const
87 SparseVector DB(getNumBasisFunctions());
93 int knotIndex = indexHalfopenInterval(x);
98 for (
int i = 1;
i <=
p - r;
i++)
100 SparseMatrix R = buildBasisMatrix(x, knotIndex,
i);
104 for (
int i =
p - r + 1;
i <=
p;
i++)
106 SparseMatrix DR = buildBasisMatrix(x, knotIndex,
i,
true);
110 double factorial = std::tgamma(
p + 1) / std::tgamma(
p - r + 1);
113 if (B.cols() !=
p + 1)
115 throw Exception(
"BSplineBasis1D::evalDerivative: Wrong number of columns of B matrix.");
119 SparseVector DB(getNumBasisFunctions());
121 int i = knotIndex -
p;
123 for (
int k = 0; k < B.outerSize(); ++k)
124 for (SparseMatrix::InnerIterator it(B, k); it; ++it)
126 DB.insert(
i + it.col()) = it.value();
133SparseVector BSplineBasis1D::evalFirstDerivative(
double x)
const
135 SparseVector values(getNumBasisFunctions());
137 std::vector<int> supportedBasisFunctions = indexSupportedBasisfunctions(x);
139 for (
int i : supportedBasisFunctions)
143 double b1 = deBoorCox(x,
i, degree - 1);
144 double b2 = deBoorCox(x,
i + 1, degree - 1);
145 double t11 = knots.at(
i);
146 double t12 = knots.at(
i + degree);
147 double t21 = knots.at(
i + 1);
148 double t22 = knots.at(
i + degree + 1);
149 (t12 == t11) ? b1 = 0 : b1 = b1 / (t12 - t11);
150 (t22 == t21) ? b2 = 0 : b2 = b2 / (t22 - t21);
151 values.insert(
i) = degree * (b1 - b2);
158SparseMatrix BSplineBasis1D::buildBasisMatrix(
double x,
unsigned int u,
159 unsigned int k,
bool diff)
const
166 if (!(k >= 1 && k <= getBasisDegree()))
168 throw Exception(
"BSplineBasis1D::buildBasisMatrix: Incorrect input paramaters!");
173 unsigned int rows = k;
174 unsigned int cols = k + 1;
175 SparseMatrix R(rows, cols);
176 R.reserve(Eigen::VectorXi::Constant(cols, 2));
178 for (
unsigned int i = 0;
i < rows;
i++)
180 double dk = knots.at(u + 1 +
i) - knots.at(u + 1 +
i - k);
191 R.insert(
i,
i) = -1 / dk;
193 R.insert(
i,
i + 1) = 1 / dk;
198 double a = (knots.at(u + 1 +
i) - x) / dk;
206 double b = (x - knots.at(u + 1 +
i - k)) / dk;
210 R.insert(
i,
i + 1) = b;
220double BSplineBasis1D::deBoorCox(
double x,
int i,
int k)
const
224 if (inHalfopenInterval(x, knots.at(
i), knots.at(
i + 1)))
235 double s1, s2, r1, r2;
236 s1 = deBoorCoxCoeff(x, knots.at(
i), knots.at(
i + k));
237 s2 = deBoorCoxCoeff(x, knots.at(
i + 1), knots.at(
i + k + 1));
238 r1 = deBoorCox(x,
i, k - 1);
239 r2 = deBoorCox(x,
i + 1, k - 1);
240 return s1 * r1 + (1 - s2) * r2;
244double BSplineBasis1D::deBoorCoxCoeff(
double x,
double x_min,
247 if (x_min < x_max && x_min <= x && x <= x_max)
249 return (x - x_min) / (x_max - x_min);
256SparseMatrix BSplineBasis1D::insertKnots(
double tau,
unsigned int multiplicity)
258 if (!insideSupport(
tau))
260 throw Exception(
"BSplineBasis1D::insertKnots: Cannot insert knot outside domain!");
263 if (knotMultiplicity(
tau) + multiplicity > degree + 1)
265 throw Exception(
"BSplineBasis1D::insertKnots: Knot multiplicity is too high!");
269 int index = indexHalfopenInterval(
tau);
270 std::vector<double> extKnots = knots;
272 for (
unsigned int i = 0;
i < multiplicity;
i++)
274 extKnots.insert(extKnots.begin() + index + 1,
tau);
279 throw Exception(
"BSplineBasis1D::insertKnots: New knot vector is not regular!");
283 SparseMatrix
A = buildKnotInsertionMatrix(extKnots);
289SparseMatrix BSplineBasis1D::refineKnots()
292 std::vector<double> refinedKnots = knots;
293 unsigned int targetNumKnots = targetNumBasisfunctions + degree + 1;
295 while (refinedKnots.size() < targetNumKnots)
297 int index = indexLongestInterval(refinedKnots);
298 double newKnot = (refinedKnots.at(index) + refinedKnots.at(index + 1)) / 2.0;
299 refinedKnots.insert(std::lower_bound(refinedKnots.begin(), refinedKnots.end(),
305 throw Exception(
"BSplineBasis1D::refineKnots: New knot vector is not regular!");
310 throw Exception(
"BSplineBasis1D::refineKnots: New knot vector is not a proper refinement!");
314 SparseMatrix
A = buildKnotInsertionMatrix(refinedKnots);
316 knots = refinedKnots;
320SparseMatrix BSplineBasis1D::refineKnotsLocally(
double x)
322 if (!insideSupport(x))
324 throw Exception(
"BSplineBasis1D::refineKnotsLocally: Cannot refine outside support!");
327 if (getNumBasisFunctions() >= getNumBasisFunctionsTarget()
328 || assertNear(knots.front(), knots.back()))
330 unsigned int n = getNumBasisFunctions();
331 DenseMatrix
A = DenseMatrix::Identity(n, n);
332 return A.sparseView();
336 std::vector<double> refinedKnots = knots;
337 auto upper = std::lower_bound(refinedKnots.begin(), refinedKnots.end(), x);
340 if (upper == refinedKnots.begin())
342 std::advance(upper, degree + 1);
346 auto lower = std::prev(upper);
349 if (assertNear(* upper, * lower))
351 unsigned int n = getNumBasisFunctions();
352 DenseMatrix
A = DenseMatrix::Identity(n, n);
353 return A.sparseView();
357 double insertVal = x;
360 if (knotMultiplicity(x) > 0
361 || assertNear(* upper, x, 1e-6, 1e-6)
362 || assertNear(* lower, x, 1e-6, 1e-6))
364 insertVal = (* upper + * lower) / 2.0;
368 refinedKnots.insert(upper, insertVal);
372 throw Exception(
"BSplineBasis1D::refineKnotsLocally: New knot vector is not regular!");
377 throw Exception(
"BSplineBasis1D::refineKnotsLocally: New knot vector is not a proper refinement!");
381 SparseMatrix
A = buildKnotInsertionMatrix(refinedKnots);
383 knots = refinedKnots;
387SparseMatrix BSplineBasis1D::decomposeToBezierForm()
390 std::vector<double> refinedKnots = knots;
392 std::vector<double>::iterator knoti = refinedKnots.begin();
394 while (knoti != refinedKnots.end())
397 int mult = degree + 1 - knotMultiplicity(* knoti);
401 std::vector<double> newKnots(mult, * knoti);
402 refinedKnots.insert(knoti, newKnots.begin(), newKnots.end());
406 knoti = std::upper_bound(refinedKnots.begin(), refinedKnots.end(), * knoti);
411 throw Exception(
"BSplineBasis1D::refineKnots: New knot vector is not regular!");
416 throw Exception(
"BSplineBasis1D::refineKnots: New knot vector is not a proper refinement!");
420 SparseMatrix
A = buildKnotInsertionMatrix(refinedKnots);
422 knots = refinedKnots;
426SparseMatrix BSplineBasis1D::buildKnotInsertionMatrix(
const std::vector<double>
432 throw Exception(
"BSplineBasis1D::buildKnotInsertionMatrix: New knot vector is not regular!");
437 throw Exception(
"BSplineBasis1D::buildKnotInsertionMatrix: New knot vector is not a proper refinement!");
440 std::vector<double> knotsAug = refinedKnots;
441 unsigned int n = knots.size() - degree - 1;
442 unsigned int m = knotsAug.size() - degree - 1;
443 SparseMatrix
A(m, n);
445 A.reserve(Eigen::VectorXi::Constant(n, degree + 1));
448 for (
unsigned int i = 0;
i < m;
i++)
450 int u = indexHalfopenInterval(knotsAug.at(
i));
451 SparseMatrix R(1, 1);
455 for (
unsigned int j = 1; j <= degree; j++)
457 SparseMatrix Ri = buildBasisMatrix(knotsAug.at(
i + j), u, j);
462 if (R.rows() != 1 || R.cols() != (
int)degree + 1)
464 throw Exception(
"BSplineBasis1D::buildKnotInsertionMatrix: Incorrect matrix dimensions!");
470 for (
int k = 0; k < R.outerSize(); ++k)
471 for (SparseMatrix::InnerIterator it(R, k); it; ++it)
473 A.insert(
i, j + it.col()) = it.value();
486void BSplineBasis1D::supportHack(
double& x)
const
488 if (x == knots.back())
490 x = std::nextafter(x, std::numeric_limits<double>::lowest());
498int BSplineBasis1D::indexHalfopenInterval(
double x)
const
500 if (x < knots.front() || x > knots.back())
502 throw Exception(
"BSplineBasis1D::indexHalfopenInterval: x outside knot interval!");
506 std::vector<double>::const_iterator it = std::upper_bound(knots.begin(),
509 int index = it - knots.begin();
513SparseMatrix BSplineBasis1D::reduceSupport(
double lb,
double ub)
516 if (lb < knots.front() || ub > knots.back())
518 throw Exception(
"BSplineBasis1D::reduceSupport: Cannot increase support!");
521 unsigned int k = degree + 1;
522 int index_lower = indexSupportedBasisfunctions(lb).front();
523 int index_upper = indexSupportedBasisfunctions(ub).back();
526 if (k != knotMultiplicity(knots.at(index_lower)))
528 int suggested_index = index_lower - 1;
530 if (0 <= suggested_index)
532 index_lower = suggested_index;
536 throw Exception(
"BSplineBasis1D::reduceSupport: Suggested index is negative!");
541 if (knotMultiplicity(ub) == k && knots.at(index_upper) == ub)
547 std::vector<double> si;
548 si.insert(si.begin(), knots.begin() + index_lower,
549 knots.begin() + index_upper + k + 1);
551 int numOld = knots.size() - k;
552 int numNew = si.size() - k;
556 throw Exception(
"BSplineBasis1D::reduceSupport: Number of basis functions is increased instead of reduced!");
559 DenseMatrix Ad = DenseMatrix::Zero(numOld, numNew);
560 Ad.block(index_lower, 0, numNew, numNew) = DenseMatrix::Identity(numNew,
562 SparseMatrix
A = Ad.sparseView();
568double BSplineBasis1D::getKnotValue(
unsigned int index)
const
570 return knots.at(index);
573unsigned int BSplineBasis1D::knotMultiplicity(
double tau)
const
575 return std::count(knots.begin(), knots.end(),
tau);
578bool BSplineBasis1D::inHalfopenInterval(
double x,
double x_min,
581 return (x_min <= x) && (x < x_max);
584bool BSplineBasis1D::insideSupport(
double x)
const
586 return (knots.front() <= x) && (x <= knots.back());
589unsigned int BSplineBasis1D::getNumBasisFunctions()
const
591 return knots.size() - (degree + 1);
594unsigned int BSplineBasis1D::getNumBasisFunctionsTarget()
const
596 return targetNumBasisfunctions;
600std::vector<int> BSplineBasis1D::indexSupportedBasisfunctions(
double x)
const
602 std::vector<int> ret;
604 if (insideSupport(x))
606 int last = indexHalfopenInterval(x);
611 last = knots.size() - 1 - (degree + 1);
614 int first = std::max((
int)(last - degree), 0);
616 for (
int i = first;
i <= last;
i++)
625unsigned int BSplineBasis1D::indexLongestInterval()
const
627 return indexLongestInterval(knots);
630unsigned int BSplineBasis1D::indexLongestInterval(
const std::vector<double>&
635 unsigned int index = 0;
637 for (
unsigned int i = 0;
i < vec.size() - 1;
i++)
639 interval = vec.at(
i + 1) - vec.at(
i);
641 if (longest < interval)
bool isKnotVectorRegular(const std::vector< double > &knots, unsigned int degree)
bool isKnotVectorRefinement(const std::vector< double > &knots, const std::vector< double > &refinedKnots)