Loading...
Searching...
No Matches
bsplinebasis.C
Go to the documentation of this file.
1/*
2 * This file is part of the SPLINTER library.
3 * Copyright (C) 2012 Bjarne Grimstad (bjarne.grimstad@gmail.com).
4 *
5 * This Source Code Form is subject to the terms of the Mozilla Public
6 * License, v. 2.0. If a copy of the MPL was not distributed with this
7 * file, You can obtain one at http://mozilla.org/MPL/2.0/.
8*/
9
10#include "bsplinebasis.h"
11#include "mykroneckerproduct.h"
12#include <unsupported/Eigen/KroneckerProduct>
13
14#include <iostream>
15
16namespace SPLINTER
17{
18
19BSplineBasis::BSplineBasis()
20{
21}
22
23BSplineBasis::BSplineBasis(std::vector< std::vector<double>>& knotVectors,
24 std::vector<unsigned int> basisDegrees)
25 : numVariables(knotVectors.size())
26{
27 if (knotVectors.size() != basisDegrees.size())
28 {
29 throw Exception("BSplineBasis::BSplineBasis: Incompatible sizes. Number of knot vectors is not equal to size of degree vector.");
30 }
31
32 // Set univariate bases
33 bases.clear();
34
35 for (unsigned int i = 0; i < numVariables; i++)
36 {
37 bases.push_back(BSplineBasis1D(knotVectors.at(i), basisDegrees.at(i)));
38
39 // Adjust target number of basis functions used in e.g. refinement
40 if (numVariables > 2)
41 {
42 // One extra knot is allowed
43 bases.at(i).setNumBasisFunctionsTarget((basisDegrees.at(
44 i) + 1) + 1); // Minimum degree+1
45 }
46 }
47}
48
49SparseVector BSplineBasis::eval(const DenseVector& x) const
50{
51 // Evaluate basisfunctions for each variable i and compute the tensor product of the function values
52 std::vector<SparseVector> basisFunctionValues;
53
54 for (int var = 0; var < x.size(); var++)
55 {
56 basisFunctionValues.push_back(bases.at(var).eval(x(var)));
57 }
58
59 return kroneckerProductVectors(basisFunctionValues);
60}
61
62// Old implementation of Jacobian
63DenseMatrix BSplineBasis::evalBasisJacobianOld(DenseVector& x) const
64{
65 // Jacobian basis matrix
66 DenseMatrix J;
67 J.setZero(getNumBasisFunctions(), numVariables);
68
69 // Calculate partial derivatives
70 for (unsigned int i = 0; i < numVariables; i++)
71 {
72 // One column in basis jacobian
73 DenseVector bi;
74 bi.setOnes(1);
75
76 for (unsigned int j = 0; j < numVariables; j++)
77 {
78 DenseVector temp = bi;
79 DenseVector xi;
80
81 if (j == i)
82 {
83 // Differentiated basis
84 xi = bases.at(j).evalFirstDerivative(x(j));
85 }
86 else
87 {
88 // Normal basis
89 xi = bases.at(j).eval(x(j));
90 }
91
92 bi = kroneckerProduct(temp, xi);
93 }
94
95 // Fill out column
96 J.block(0, i, bi.rows(), 1) = bi.block(0, 0, bi.rows(), 1);
97 }
98
99 return J;
100}
101
102// NOTE: does not pass tests
103SparseMatrix BSplineBasis::evalBasisJacobian(DenseVector& x) const
104{
105 // Jacobian basis matrix
106 SparseMatrix J(getNumBasisFunctions(), numVariables);
107 //J.setZero(numBasisFunctions(), numInputs);
108
109 // Calculate partial derivatives
110 for (unsigned int i = 0; i < numVariables; ++i)
111 {
112 // One column in basis jacobian
113 std::vector<SparseVector> values(numVariables);
114
115 for (unsigned int j = 0; j < numVariables; ++j)
116 {
117 if (j == i)
118 {
119 // Differentiated basis
120 values.at(j) = bases.at(j).evalDerivative(x(j), 1);
121 }
122 else
123 {
124 // Normal basis
125 values.at(j) = bases.at(j).eval(x(j));
126 }
127 }
128
129 SparseVector Ji = kroneckerProductVectors(values);
130
131 // Fill out column
132 for (int k = 0; k < Ji.outerSize(); ++k)
133 for (SparseVector::InnerIterator it(Ji, k); it; ++it)
134 {
135 if (it.value() != 0)
136 {
137 J.insert(it.row(), i) = it.value();
138 }
139 }
140
141 //J.block(0,i,Ji.rows(),1) = bi.block(0,0,Ji.rows(),1);
142 }
143
144 J.makeCompressed();
145 return J;
146}
147
148SparseMatrix BSplineBasis::evalBasisJacobian2(DenseVector& x) const
149{
150 // Jacobian basis matrix
151 SparseMatrix J(getNumBasisFunctions(), numVariables);
152 // Evaluate B-spline basis functions before looping
153 std::vector<SparseVector> funcValues(numVariables);
154 std::vector<SparseVector> gradValues(numVariables);
155
156 for (unsigned int i = 0; i < numVariables; ++i)
157 {
158 funcValues[i] = bases.at(i).eval(x(i));
159 gradValues[i] = bases.at(i).evalFirstDerivative(x(i));
160 }
161
162 // Calculate partial derivatives
163 for (unsigned int i = 0; i < numVariables; i++)
164 {
165 std::vector<SparseVector> values(numVariables);
166
167 for (unsigned int j = 0; j < numVariables; j++)
168 {
169 if (j == i)
170 {
171 values.at(j) = gradValues.at(j); // Differentiated basis
172 }
173 else
174 {
175 values.at(j) = funcValues.at(j); // Normal basis
176 }
177 }
178
179 SparseVector Ji = kroneckerProductVectors(values);
180
181 // Fill out column
182 for (SparseVector::InnerIterator it(Ji); it; ++it)
183 {
184 J.insert(it.row(), i) = it.value();
185 }
186 }
187
188 return J;
189}
190
191SparseMatrix BSplineBasis::evalBasisHessian(DenseVector& x) const
192{
193 // Hessian basis matrix
194 /* Hij = B1 x ... x DBi x ... x DBj x ... x Bn
195 * (Hii = B1 x ... x DDBi x ... x Bn)
196 * Where B are basis functions evaluated at x,
197 * DB are the derivative of the basis functions,
198 * and x is the kronecker product.
199 * Hij is in R^(numBasisFunctions x 1)
200 * so that basis hessian H is in R^(numBasisFunctions*numInputs x numInputs)
201 * The real B-spline Hessian is calculated as (c^T x 1^(numInputs x 1))*H
202 */
203 SparseMatrix H(getNumBasisFunctions()*numVariables, numVariables);
204 //H.setZero(numBasisFunctions()*numInputs, numInputs);
205
206 // Calculate partial derivatives
207 // Utilizing that Hessian is symmetric
208 // Filling out lower left triangular
209 for (unsigned int i = 0; i < numVariables; i++) // row
210 {
211 for (unsigned int j = 0; j <= i; j++) // col
212 {
213 // One column in basis jacobian
214 SparseMatrix Hi(1, 1);
215 Hi.insert(0, 0) = 1;
216
217 for (unsigned int k = 0; k < numVariables; k++)
218 {
219 SparseMatrix temp = Hi;
220 SparseMatrix Bk;
221
222 if (i == j && k == i)
223 {
224 // Diagonal element
225 Bk = bases.at(k).evalDerivative(x(k), 2);
226 }
227 else if (k == i || k == j)
228 {
229 Bk = bases.at(k).evalDerivative(x(k), 1);
230 }
231 else
232 {
233 Bk = bases.at(k).eval(x(k));
234 }
235
236 Hi = kroneckerProduct(temp, Bk);
237 }
238
239 // Fill out column
240 for (int k = 0; k < Hi.outerSize(); ++k)
241 for (SparseMatrix::InnerIterator it(Hi, k); it; ++it)
242 {
243 if (it.value() != 0)
244 {
245 int row = i * getNumBasisFunctions() + it.row();
246 int col = j;
247 H.insert(row, col) = it.value();
248 }
249 }
250 }
251 }
252
253 H.makeCompressed();
254 return H;
255}
256
257SparseMatrix BSplineBasis::insertKnots(double tau, unsigned int dim,
258 unsigned int multiplicity)
259{
260 SparseMatrix A(1, 1);
261 // A.resize(1,1);
262 A.insert(0, 0) = 1;
263
264 // Calculate multivariate knot insertion matrix
265 for (unsigned int i = 0; i < numVariables; i++)
266 {
267 SparseMatrix temp = A;
268 SparseMatrix Ai;
269
270 if (i == dim)
271 {
272 // Build knot insertion matrix
273 Ai = bases.at(i).insertKnots(tau, multiplicity);
274 }
275 else
276 {
277 // No insertion - identity matrix
278 int m = bases.at(i).getNumBasisFunctions();
279 Ai.resize(m, m);
280 Ai.setIdentity();
281 }
282
283 //A = kroneckerProduct(temp, Ai);
284 A = myKroneckerProduct(temp, Ai);
285 }
286
287 A.makeCompressed();
288 return A;
289}
290
291SparseMatrix BSplineBasis::refineKnots()
292{
293 SparseMatrix A(1, 1);
294 A.insert(0, 0) = 1;
295
296 for (unsigned int i = 0; i < numVariables; i++)
297 {
298 SparseMatrix temp = A;
299 SparseMatrix Ai = bases.at(i).refineKnots();
300 //A = kroneckerProduct(temp, Ai);
301 A = myKroneckerProduct(temp, Ai);
302 }
303
304 A.makeCompressed();
305 return A;
306}
307
308SparseMatrix BSplineBasis::refineKnotsLocally(DenseVector x)
309{
310 SparseMatrix A(1, 1);
311 A.insert(0, 0) = 1;
312
313 for (unsigned int i = 0; i < numVariables; i++)
314 {
315 SparseMatrix temp = A;
316 SparseMatrix Ai = bases.at(i).refineKnotsLocally(x(i));
317 //A = kroneckerProduct(temp, Ai);
318 A = myKroneckerProduct(temp, Ai);
319 }
320
321 A.makeCompressed();
322 return A;
323}
324
325SparseMatrix BSplineBasis::decomposeToBezierForm()
326{
327 SparseMatrix A(1, 1);
328 A.insert(0, 0) = 1;
329
330 for (unsigned int i = 0; i < numVariables; i++)
331 {
332 SparseMatrix temp = A;
333 SparseMatrix Ai = bases.at(i).decomposeToBezierForm();
334 //A = kroneckerProduct(temp, Ai);
335 A = myKroneckerProduct(temp, Ai);
336 }
337
338 A.makeCompressed();
339 return A;
340}
341
342SparseMatrix BSplineBasis::reduceSupport(std::vector<double>& lb,
343 std::vector<double>& ub)
344{
345 if (lb.size() != ub.size() || lb.size() != numVariables)
346 {
347 throw Exception("BSplineBasis::reduceSupport: Incompatible dimension of domain bounds.");
348 }
349
350 SparseMatrix A(1, 1);
351 A.insert(0, 0) = 1;
352
353 for (unsigned int i = 0; i < numVariables; i++)
354 {
355 SparseMatrix temp = A;
356 SparseMatrix Ai;
357 Ai = bases.at(i).reduceSupport(lb.at(i), ub.at(i));
358 //A = kroneckerProduct(temp, Ai);
359 A = myKroneckerProduct(temp, Ai);
360 }
361
362 A.makeCompressed();
363 return A;
364}
365
366std::vector<unsigned int> BSplineBasis::getBasisDegrees() const
367{
368 std::vector<unsigned int> degrees;
369
370 for (const auto& basis : bases)
371 {
372 degrees.push_back(basis.getBasisDegree());
373 }
374
375 return degrees;
376}
377
378unsigned int BSplineBasis::getBasisDegree(unsigned int dim) const
379{
380 return bases.at(dim).getBasisDegree();
381}
382
383unsigned int BSplineBasis::getNumBasisFunctions(unsigned int dim) const
384{
385 return bases.at(dim).getNumBasisFunctions();
386}
387
388unsigned int BSplineBasis::getNumBasisFunctions() const
389{
390 unsigned int prod = 1;
391
392 for (unsigned int dim = 0; dim < numVariables; dim++)
393 {
394 prod *= bases.at(dim).getNumBasisFunctions();
395 }
396
397 return prod;
398}
399
400BSplineBasis1D BSplineBasis::getSingleBasis(int dim)
401{
402 return bases.at(dim);
403}
404
405std::vector<double> BSplineBasis::getKnotVector(int dim) const
406{
407 return bases.at(dim).getKnotVector();
408}
409
410std::vector< std::vector<double>> BSplineBasis::getKnotVectors() const
411{
412 std::vector< std::vector<double>> knots;
413
414 for (unsigned int i = 0; i < numVariables; i++)
415 {
416 knots.push_back(bases.at(i).getKnotVector());
417 }
418
419 return knots;
420}
421
422unsigned int BSplineBasis::getKnotMultiplicity(unsigned int dim,
423 double tau) const
424{
425 return bases.at(dim).knotMultiplicity(tau);
426}
427
428double BSplineBasis::getKnotValue(int dim, int index) const
429{
430 return bases.at(dim).getKnotValue(index);
431}
432
433unsigned int BSplineBasis::getLargestKnotInterval(unsigned int dim) const
434{
435 return bases.at(dim).indexLongestInterval();
436}
437
438std::vector<unsigned int> BSplineBasis::getNumBasisFunctionsTarget() const
439{
440 std::vector<unsigned int> ret;
441
442 for (unsigned int dim = 0; dim < numVariables; dim++)
443 {
444 ret.push_back(bases.at(dim).getNumBasisFunctionsTarget() );
445 }
446
447 return ret;
448}
449
450int BSplineBasis::supportedPrInterval() const
451{
452 int ret = 1;
453
454 for (unsigned int dim = 0; dim < numVariables; dim++)
455 {
456 ret *= (bases.at(dim).getBasisDegree() + 1);
457 }
458
459 return ret;
460}
461
462bool BSplineBasis::insideSupport(DenseVector& x) const
463{
464 for (unsigned int dim = 0; dim < numVariables; dim++)
465 {
466 if (!bases.at(dim).insideSupport(x(dim)))
467 {
468 return false;
469 }
470 }
471
472 return true;
473}
474
475std::vector<double> BSplineBasis::getSupportLowerBound() const
476{
477 std::vector<double> lb;
478
479 for (unsigned int dim = 0; dim < numVariables; dim++)
480 {
481 std::vector<double> knots = bases.at(dim).getKnotVector();
482 lb.push_back(knots.front());
483 }
484
485 return lb;
486}
487
488std::vector<double> BSplineBasis::getSupportUpperBound() const
489{
490 std::vector<double> ub;
491
492 for (unsigned int dim = 0; dim < numVariables; dim++)
493 {
494 std::vector<double> knots = bases.at(dim).getKnotVector();
495 ub.push_back(knots.back());
496 }
497
498 return ub;
499}
500
501} // namespace SPLINTER
volScalarField & A
SparseMatrix myKroneckerProduct(const SparseMatrix &A, const SparseMatrix &B)
SparseVector kroneckerProductVectors(const std::vector< SparseVector > &vectors)
singVal col(i)
label i
Definition pEqn.H:46
dimensionedScalar & tau