10#include "bsplinebasis.h"
11#include "mykroneckerproduct.h"
12#include <unsupported/Eigen/KroneckerProduct>
19BSplineBasis::BSplineBasis()
23BSplineBasis::BSplineBasis(std::vector< std::vector<double>>& knotVectors,
24 std::vector<unsigned int> basisDegrees)
25 : numVariables(knotVectors.size())
27 if (knotVectors.size() != basisDegrees.size())
29 throw Exception(
"BSplineBasis::BSplineBasis: Incompatible sizes. Number of knot vectors is not equal to size of degree vector.");
35 for (
unsigned int i = 0;
i < numVariables;
i++)
37 bases.push_back(BSplineBasis1D(knotVectors.at(
i), basisDegrees.at(
i)));
43 bases.at(
i).setNumBasisFunctionsTarget((basisDegrees.at(
49SparseVector BSplineBasis::eval(
const DenseVector& x)
const
52 std::vector<SparseVector> basisFunctionValues;
54 for (
int var = 0; var < x.size(); var++)
56 basisFunctionValues.push_back(bases.at(var).eval(x(var)));
63DenseMatrix BSplineBasis::evalBasisJacobianOld(DenseVector& x)
const
67 J.setZero(getNumBasisFunctions(), numVariables);
70 for (
unsigned int i = 0;
i < numVariables;
i++)
76 for (
unsigned int j = 0; j < numVariables; j++)
78 DenseVector temp = bi;
84 xi = bases.at(j).evalFirstDerivative(x(j));
89 xi = bases.at(j).eval(x(j));
92 bi = kroneckerProduct(temp, xi);
96 J.block(0,
i, bi.rows(), 1) = bi.block(0, 0, bi.rows(), 1);
103SparseMatrix BSplineBasis::evalBasisJacobian(DenseVector& x)
const
106 SparseMatrix J(getNumBasisFunctions(), numVariables);
110 for (
unsigned int i = 0;
i < numVariables; ++
i)
113 std::vector<SparseVector> values(numVariables);
115 for (
unsigned int j = 0; j < numVariables; ++j)
120 values.at(j) = bases.at(j).evalDerivative(x(j), 1);
125 values.at(j) = bases.at(j).eval(x(j));
132 for (
int k = 0; k < Ji.outerSize(); ++k)
133 for (SparseVector::InnerIterator it(Ji, k); it; ++it)
137 J.insert(it.row(),
i) = it.value();
148SparseMatrix BSplineBasis::evalBasisJacobian2(DenseVector& x)
const
151 SparseMatrix J(getNumBasisFunctions(), numVariables);
153 std::vector<SparseVector> funcValues(numVariables);
154 std::vector<SparseVector> gradValues(numVariables);
156 for (
unsigned int i = 0;
i < numVariables; ++
i)
158 funcValues[
i] = bases.at(
i).eval(x(
i));
159 gradValues[
i] = bases.at(
i).evalFirstDerivative(x(
i));
163 for (
unsigned int i = 0;
i < numVariables;
i++)
165 std::vector<SparseVector> values(numVariables);
167 for (
unsigned int j = 0; j < numVariables; j++)
171 values.at(j) = gradValues.at(j);
175 values.at(j) = funcValues.at(j);
182 for (SparseVector::InnerIterator it(Ji); it; ++it)
184 J.insert(it.row(),
i) = it.value();
191SparseMatrix BSplineBasis::evalBasisHessian(DenseVector& x)
const
203 SparseMatrix H(getNumBasisFunctions()*numVariables, numVariables);
209 for (
unsigned int i = 0;
i < numVariables;
i++)
211 for (
unsigned int j = 0; j <=
i; j++)
214 SparseMatrix Hi(1, 1);
217 for (
unsigned int k = 0; k < numVariables; k++)
219 SparseMatrix temp = Hi;
222 if (
i == j && k ==
i)
225 Bk = bases.at(k).evalDerivative(x(k), 2);
227 else if (k ==
i || k == j)
229 Bk = bases.at(k).evalDerivative(x(k), 1);
233 Bk = bases.at(k).eval(x(k));
236 Hi = kroneckerProduct(temp, Bk);
240 for (
int k = 0; k < Hi.outerSize(); ++k)
241 for (SparseMatrix::InnerIterator it(Hi, k); it; ++it)
245 int row =
i * getNumBasisFunctions() + it.row();
247 H.insert(row,
col) = it.value();
257SparseMatrix BSplineBasis::insertKnots(
double tau,
unsigned int dim,
258 unsigned int multiplicity)
260 SparseMatrix
A(1, 1);
265 for (
unsigned int i = 0;
i < numVariables;
i++)
267 SparseMatrix temp =
A;
273 Ai = bases.at(
i).insertKnots(
tau, multiplicity);
278 int m = bases.at(
i).getNumBasisFunctions();
291SparseMatrix BSplineBasis::refineKnots()
293 SparseMatrix
A(1, 1);
296 for (
unsigned int i = 0;
i < numVariables;
i++)
298 SparseMatrix temp =
A;
299 SparseMatrix Ai = bases.at(
i).refineKnots();
308SparseMatrix BSplineBasis::refineKnotsLocally(DenseVector x)
310 SparseMatrix
A(1, 1);
313 for (
unsigned int i = 0;
i < numVariables;
i++)
315 SparseMatrix temp =
A;
316 SparseMatrix Ai = bases.at(
i).refineKnotsLocally(x(
i));
325SparseMatrix BSplineBasis::decomposeToBezierForm()
327 SparseMatrix
A(1, 1);
330 for (
unsigned int i = 0;
i < numVariables;
i++)
332 SparseMatrix temp =
A;
333 SparseMatrix Ai = bases.at(
i).decomposeToBezierForm();
342SparseMatrix BSplineBasis::reduceSupport(std::vector<double>& lb,
343 std::vector<double>& ub)
345 if (lb.size() != ub.size() || lb.size() != numVariables)
347 throw Exception(
"BSplineBasis::reduceSupport: Incompatible dimension of domain bounds.");
350 SparseMatrix
A(1, 1);
353 for (
unsigned int i = 0;
i < numVariables;
i++)
355 SparseMatrix temp =
A;
357 Ai = bases.at(
i).reduceSupport(lb.at(
i), ub.at(
i));
366std::vector<unsigned int> BSplineBasis::getBasisDegrees()
const
368 std::vector<unsigned int> degrees;
370 for (
const auto& basis : bases)
372 degrees.push_back(basis.getBasisDegree());
378unsigned int BSplineBasis::getBasisDegree(
unsigned int dim)
const
380 return bases.at(dim).getBasisDegree();
383unsigned int BSplineBasis::getNumBasisFunctions(
unsigned int dim)
const
385 return bases.at(dim).getNumBasisFunctions();
388unsigned int BSplineBasis::getNumBasisFunctions()
const
390 unsigned int prod = 1;
392 for (
unsigned int dim = 0; dim < numVariables; dim++)
394 prod *= bases.at(dim).getNumBasisFunctions();
400BSplineBasis1D BSplineBasis::getSingleBasis(
int dim)
402 return bases.at(dim);
405std::vector<double> BSplineBasis::getKnotVector(
int dim)
const
407 return bases.at(dim).getKnotVector();
410std::vector< std::vector<double>> BSplineBasis::getKnotVectors()
const
412 std::vector< std::vector<double>> knots;
414 for (
unsigned int i = 0;
i < numVariables;
i++)
416 knots.push_back(bases.at(
i).getKnotVector());
422unsigned int BSplineBasis::getKnotMultiplicity(
unsigned int dim,
425 return bases.at(dim).knotMultiplicity(
tau);
428double BSplineBasis::getKnotValue(
int dim,
int index)
const
430 return bases.at(dim).getKnotValue(index);
433unsigned int BSplineBasis::getLargestKnotInterval(
unsigned int dim)
const
435 return bases.at(dim).indexLongestInterval();
438std::vector<unsigned int> BSplineBasis::getNumBasisFunctionsTarget()
const
440 std::vector<unsigned int> ret;
442 for (
unsigned int dim = 0; dim < numVariables; dim++)
444 ret.push_back(bases.at(dim).getNumBasisFunctionsTarget() );
450int BSplineBasis::supportedPrInterval()
const
454 for (
unsigned int dim = 0; dim < numVariables; dim++)
456 ret *= (bases.at(dim).getBasisDegree() + 1);
462bool BSplineBasis::insideSupport(DenseVector& x)
const
464 for (
unsigned int dim = 0; dim < numVariables; dim++)
466 if (!bases.at(dim).insideSupport(x(dim)))
475std::vector<double> BSplineBasis::getSupportLowerBound()
const
477 std::vector<double> lb;
479 for (
unsigned int dim = 0; dim < numVariables; dim++)
481 std::vector<double> knots = bases.at(dim).getKnotVector();
482 lb.push_back(knots.front());
488std::vector<double> BSplineBasis::getSupportUpperBound()
const
490 std::vector<double> ub;
492 for (
unsigned int dim = 0; dim < numVariables; dim++)
494 std::vector<double> knots = bases.at(dim).getKnotVector();
495 ub.push_back(knots.back());
SparseMatrix myKroneckerProduct(const SparseMatrix &A, const SparseMatrix &B)
SparseVector kroneckerProductVectors(const std::vector< SparseVector > &vectors)