Loading...
Searching...
No Matches
bsplinebuilder.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 "bsplinebuilder.h"
11#include "mykroneckerproduct.h"
12#include <unsupported/Eigen/KroneckerProduct>
13#include <linearsolvers.h>
14#include <serializer.h>
15#include <iostream>
16#include <utilities.h>
17
18namespace SPLINTER
19{
20// Default constructor
21BSpline::Builder::Builder(const DataTable& data)
22 :
23 _data(data),
24 _degrees(getBSplineDegrees(data.getNumVariables(), 3)),
25 _numBasisFunctions(std::vector<unsigned int>(data.getNumVariables(), 0)),
26 _knotSpacing(KnotSpacing::AS_SAMPLED),
27 _smoothing(Smoothing::NONE),
28 _alpha(0.1)
29{
30}
31
32/*
33 * Build B-spline
34 */
35BSpline BSpline::Builder::build() const
36{
37 // Check data
38 // TODO: Remove this test
39 if (!_data.isGridComplete())
40 {
41 throw Exception("BSpline::Builder::build: Cannot create B-spline from irregular (incomplete) grid.");
42 }
43
44 // Build knot vectors
45 auto knotVectors = computeKnotVectors();
46 // Build B-spline (with default coefficients)
47 auto bspline = BSpline(knotVectors, _degrees);
48 // Compute coefficients from samples and update B-spline
49 auto coefficients = computeCoefficients(bspline);
50 bspline.setCoefficients(coefficients);
51 return bspline;
52}
53
54/*
55 * Find coefficients of B-spline by solving:
56 * min ||A*x - b||^2 + alpha*||R||^2,
57 * where
58 * A = mxn matrix of n basis functions evaluated at m sample points,
59 * b = vector of m sample points y-values (or x-values when calculating knot averages),
60 * x = B-spline coefficients (or knot averages),
61 * R = Regularization matrix,
62 * alpha = regularization parameter.
63 */
64DenseVector BSpline::Builder::computeCoefficients(const BSpline& bspline) const
65{
66 SparseMatrix B = computeBasisFunctionMatrix(bspline);
67 SparseMatrix A = B;
68 DenseVector b = getSamplePointValues();
69
70 if (_smoothing == Smoothing::IDENTITY)
71 {
72 /*
73 * Computing B-spline coefficients with a regularization term
74 * ||Ax-b||^2 + alpha*x^T*x
75 *
76 * NOTE: This corresponds to a Tikhonov regularization (or ridge regression) with the Identity matrix.
77 * See: https://en.wikipedia.org/wiki/Tikhonov_regularization
78 *
79 * NOTE2: consider changing regularization factor to (alpha/numSample)
80 */
81 SparseMatrix Bt = B.transpose();
82 A = Bt * B;
83 b = Bt * b;
84 auto I = SparseMatrix(A.cols(), A.cols());
85 I.setIdentity();
86 A += _alpha * I;
87 }
88 else if (_smoothing == Smoothing::PSPLINE)
89 {
90 /*
91 * The P-Spline is a smooting B-spline which relaxes the interpolation constraints on the control points to allow
92 * smoother spline curves. It minimizes an objective which penalizes both deviation from sample points (to lower bias)
93 * and the magnitude of second derivatives (to lower variance).
94 *
95 * Setup and solve equations Ax = b,
96 * A = B'*W*B + l*D'*D
97 * b = B'*W*y
98 * x = control coefficients or knot averages.
99 * B = basis functions at sample x-values,
100 * W = weighting matrix for interpolating specific points
101 * D = second-order finite difference matrix
102 * l = penalizing parameter (increase for more smoothing)
103 * y = sample y-values when calculating control coefficients,
104 * y = sample x-values when calculating knot averages
105 */
106 // Assuming regular grid
107 unsigned int numSamples = _data.getNumSamples();
108 SparseMatrix Bt = B.transpose();
109 // Weight matrix
110 SparseMatrix W;
111 W.resize(numSamples, numSamples);
112 W.setIdentity();
113 // Second order finite difference matrix
114 SparseMatrix D = getSecondOrderFiniteDifferenceMatrix(bspline);
115 // Left-hand side matrix
116 A = Bt * W * B + _alpha * D.transpose() * D;
117 // Compute right-hand side matrices
118 b = Bt * W * b;
119 }
120
121 DenseVector x;
122 int numEquations = A.rows();
123 int maxNumEquations = 100;
124 bool solveAsDense = (numEquations < maxNumEquations);
125
126 if (!solveAsDense)
127 {
128#ifndef NDEBUG
129 std::cout <<
130 "BSpline::Builder::computeBSplineCoefficients: Computing B-spline control points using sparse solver."
131 << std::endl;
132#endif // NDEBUG
133 SparseLU<> s;
134 //bool successfulSolve = (s.solve(A,Bx,Cx) && s.solve(A,By,Cy));
135 solveAsDense = !s.solve(A, b, x);
136 }
137
138 if (solveAsDense)
139 {
140#ifndef NDEBUG
141 std::cout <<
142 "BSpline::Builder::computeBSplineCoefficients: Computing B-spline control points using dense solver."
143 << std::endl;
144#endif // NDEBUG
145 DenseMatrix Ad = A.toDense();
146 DenseQR<DenseVector> s;
147
148 // DenseSVD<DenseVector> s;
149 //bool successfulSolve = (s.solve(Ad,Bx,Cx) && s.solve(Ad,By,Cy));
150 if (!s.solve(Ad, b, x))
151 {
152 throw Exception("BSpline::Builder::computeBSplineCoefficients: Failed to solve for B-spline coefficients.");
153 }
154 }
155
156 return x;
157}
158
159SparseMatrix BSpline::Builder::computeBasisFunctionMatrix(
160 const BSpline& bspline) const
161{
162 unsigned int numVariables = _data.getNumVariables();
163 unsigned int numSamples = _data.getNumSamples();
164 // TODO: Reserve nnz per row (degree+1)
165 //int nnzPrCol = bspline.basis.supportedPrInterval();
166 SparseMatrix A(numSamples, bspline.getNumBasisFunctions());
167 //A.reserve(DenseVector::Constant(numSamples, nnzPrCol)); // TODO: should reserve nnz per row!
168 int i = 0;
169
170 for (auto it = _data.cbegin(); it != _data.cend(); ++it, ++i)
171 {
172 DenseVector xi(numVariables);
173 xi.setZero();
174 std::vector<double> xv = it->getX();
175
176 for (unsigned int j = 0; j < numVariables; ++j)
177 {
178 xi(j) = xv.at(j);
179 }
180
181 SparseVector basisValues = bspline.evalBasis(xi);
182
183 for (SparseVector::InnerIterator it2(basisValues); it2; ++it2)
184 {
185 A.insert(i, it2.index()) = it2.value();
186 }
187 }
188
189 A.makeCompressed();
190 return A;
191}
192
193DenseVector BSpline::Builder::getSamplePointValues() const
194{
195 DenseVector B = DenseVector::Zero(_data.getNumSamples());
196 int i = 0;
197
198 for (auto it = _data.cbegin(); it != _data.cend(); ++it, ++i)
199 {
200 B(i) = it->getY();
201 }
202
203 return B;
204}
205
206/*
207* Function for generating second order finite-difference matrix, which is used for penalizing the
208* (approximate) second derivative in control point calculation for P-splines.
209*/
210SparseMatrix BSpline::Builder::getSecondOrderFiniteDifferenceMatrix(
211 const BSpline& bspline) const
212{
213 unsigned int numVariables = bspline.getNumVariables();
214 // Number of (total) basis functions - defines the number of columns in D
215 unsigned int numCols = bspline.getNumBasisFunctions();
216 std::vector<unsigned int> numBasisFunctions =
217 bspline.getNumBasisFunctionsPerVariable();
218 // Number of basis functions (and coefficients) in each variable
219 std::vector<unsigned int> dims;
220
221 for (unsigned int i = 0; i < numVariables; i++)
222 {
223 dims.push_back(numBasisFunctions.at(i));
224 }
225
226 std::reverse(dims.begin(), dims.end());
227
228 for (unsigned int i = 0; i < numVariables; ++i)
229 if (numBasisFunctions.at(i) < 3)
230 {
231 throw Exception("BSpline::Builder::getSecondOrderDifferenceMatrix: Need at least three coefficients/basis function per variable.");
232 }
233
234 // Number of rows in D and in each block
235 int numRows = 0;
236 std::vector< int > numBlkRows;
237
238 for (unsigned int i = 0; i < numVariables; i++)
239 {
240 int prod = 1;
241
242 for (unsigned int j = 0; j < numVariables; j++)
243 {
244 if (i == j)
245 {
246 prod *= (dims[j] - 2);
247 }
248 else
249 {
250 prod *= dims[j];
251 }
252 }
253
254 numRows += prod;
255 numBlkRows.push_back(prod);
256 }
257
258 // Resize and initialize D
259 SparseMatrix D(numRows, numCols);
260 D.reserve(DenseVector::Constant(numCols,
261 2 * numVariables)); // D has no more than two elems per col per dim
262 int i = 0; // Row index
263
264 // Loop though each dimension (each dimension has its own block)
265 for (unsigned int d = 0; d < numVariables; d++)
266 {
267 // Calculate left and right products
268 int leftProd = 1;
269 int rightProd = 1;
270
271 for (unsigned int k = 0; k < d; k++)
272 {
273 leftProd *= dims[k];
274 }
275
276 for (unsigned int k = d + 1; k < numVariables; k++)
277 {
278 rightProd *= dims[k];
279 }
280
281 // Loop through subblocks on the block diagonal
282 for (int j = 0; j < rightProd; j++)
283 {
284 // Start column of current subblock
285 int blkBaseCol = j * leftProd * dims[d];
286
287 // Block rows [I -2I I] of subblock
288 for (unsigned int l = 0; l < (dims[d] - 2); l++)
289 {
290 // Special case for first dimension
291 if (d == 0)
292 {
293 int k = j * leftProd * dims[d] + l;
294 D.insert(i, k) = 1;
295 k += leftProd;
296 D.insert(i, k) = -2;
297 k += leftProd;
298 D.insert(i, k) = 1;
299 i++;
300 }
301 else
302 {
303 // Loop for identity matrix
304 for (int n = 0; n < leftProd; n++)
305 {
306 int k = blkBaseCol + l * leftProd + n;
307 D.insert(i, k) = 1;
308 k += leftProd;
309 D.insert(i, k) = -2;
310 k += leftProd;
311 D.insert(i, k) = 1;
312 i++;
313 }
314 }
315 }
316 }
317 }
318
319 D.makeCompressed();
320 return D;
321}
322
323// Compute all knot vectors from sample data
324std::vector<std::vector<double>> BSpline::Builder::computeKnotVectors() const
325{
326 if (_data.getNumVariables() != _degrees.size())
327 {
328 throw Exception("BSpline::Builder::computeKnotVectors: Inconsistent sizes on input vectors.");
329 }
330
331 std::vector<std::vector<double>> grid = _data.getTableX();
332 std::vector<std::vector<double>> knotVectors;
333
334 for (unsigned int i = 0; i < _data.getNumVariables(); ++i)
335 {
336 // Compute knot vector
337 auto knotVec = computeKnotVector(grid.at(i), _degrees.at(i),
338 _numBasisFunctions.at(i));
339 knotVectors.push_back(knotVec);
340 }
341
342 return knotVectors;
343}
344
345// Compute a single knot vector from sample grid and degree
346std::vector<double> BSpline::Builder::computeKnotVector(
347 const std::vector<double>& values,
348 unsigned int degree,
349 unsigned int numBasisFunctions) const
350{
351 switch (_knotSpacing)
352 {
353 case KnotSpacing::AS_SAMPLED:
354 return knotVectorMovingAverage(values, degree);
355
356 case KnotSpacing::EQUIDISTANT:
357 return knotVectorEquidistant(values, degree, numBasisFunctions);
358
359 case KnotSpacing::EXPERIMENTAL:
360 return knotVectorBuckets(values, degree);
361
362 default:
363 return knotVectorMovingAverage(values, degree);
364 }
365}
366
367/*
368* Automatic construction of (p+1)-regular knot vector
369* using moving average.
370*
371* Requirement:
372* Knot vector should be of size n+p+1.
373* End knots are should be repeated p+1 times.
374*
375* Computed sizes:
376* n+2*(p) = n + p + 1 + (p - 1)
377* k = (p - 1) values must be removed from sample vector.
378* w = k + 3 window size in moving average
379*
380* Algorithm:
381* 1) compute n - k values using moving average with window size w
382* 2) repeat first and last value p + 1 times
383*
384* The resulting knot vector has n - k + 2*p = n + p + 1 knots.
385*
386* NOTE:
387* For equidistant samples, the resulting knot vector is identically to
388* the free end conditions knot vector used in cubic interpolation.
389* That is, samples (a,b,c,d,e,f) produces the knot vector (a,a,a,a,c,d,f,f,f,f) for p = 3.
390* For p = 1, (a,b,c,d,e,f) becomes (a,a,b,c,d,e,f,f).
391*
392* TODO:
393* Does not work well when number of knots is << number of samples! For such cases
394* almost all knots will lie close to the left samples. Try a bucket approach, where the
395* samples are added to buckets and the knots computed as the average of these.
396*/
397std::vector<double> BSpline::Builder::knotVectorMovingAverage(
398 const std::vector<double>& values,
399 unsigned int degree) const
400{
401 // Sort and remove duplicates
402 std::vector<double> unique = extractUniqueSorted(values);
403 // Compute sizes
404 unsigned int n = unique.size();
405 unsigned int k = degree - 1; // knots to remove
406 unsigned int w = k + 3; // Window size
407
408 // The minimum number of samples from which a free knot vector can be created
409 if (n < degree + 1)
410 {
411 std::ostringstream e;
412 e << "knotVectorMovingAverage: Only " << n
413 << " unique interpolation points are given. A minimum of degree+1 = " << degree
414 + 1
415 << " unique points are required to build a B-spline basis of degree " << degree
416 << ".";
417 throw Exception(e.str());
418 }
419
420 std::vector<double> knots(n - k - 2, 0);
421
422 // Compute (n-k-2) interior knots using moving average
423 for (unsigned int i = 0; i < n - k - 2; ++i)
424 {
425 double ma = 0;
426
427 for (unsigned int j = 0; j < w; ++j)
428 {
429 ma += unique.at(i + j);
430 }
431
432 knots.at(i) = ma / w;
433 }
434
435 // Repeat first knot p + 1 times (for interpolation of start point)
436 for (unsigned int i = 0; i < degree + 1; ++i)
437 {
438 knots.insert(knots.begin(), unique.front());
439 }
440
441 // Repeat last knot p + 1 times (for interpolation of end point)
442 for (unsigned int i = 0; i < degree + 1; ++i)
443 {
444 knots.insert(knots.end(), unique.back());
445 }
446
447 // Number of knots in a (p+1)-regular knot vector
448 //assert(knots.size() == uniqueX.size() + degree + 1);
449 return knots;
450}
451
452std::vector<double> BSpline::Builder::knotVectorEquidistant(
453 const std::vector<double>& values,
454 unsigned int degree,
455 unsigned int numBasisFunctions = 0) const
456{
457 // Sort and remove duplicates
458 std::vector<double> unique = extractUniqueSorted(values);
459 // Compute sizes
460 unsigned int n = unique.size();
461
462 if (numBasisFunctions > 0)
463 {
464 n = numBasisFunctions;
465 }
466
467 unsigned int k = degree - 1; // knots to remove
468
469 // The minimum number of samples from which a free knot vector can be created
470 if (n < degree + 1)
471 {
472 std::ostringstream e;
473 e << "knotVectorMovingAverage: Only " << n
474 << " unique interpolation points are given. A minimum of degree+1 = " << degree
475 + 1
476 << " unique points are required to build a B-spline basis of degree " << degree
477 << ".";
478 throw Exception(e.str());
479 }
480
481 // Compute (n-k-2) equidistant interior knots
482 unsigned int numIntKnots = std::max(n - k - 2, (unsigned int)0);
483 numIntKnots = std::min((unsigned int)10, numIntKnots);
484 std::vector<double> knots = linspace(unique.front(), unique.back(),
485 numIntKnots);
486
487 // Repeat first knot p + 1 times (for interpolation of start point)
488 for (unsigned int i = 0; i < degree; ++i)
489 {
490 knots.insert(knots.begin(), unique.front());
491 }
492
493 // Repeat last knot p + 1 times (for interpolation of end point)
494 for (unsigned int i = 0; i < degree; ++i)
495 {
496 knots.insert(knots.end(), unique.back());
497 }
498
499 // Number of knots in a (p+1)-regular knot vector
500 //assert(knots.size() == uniqueX.size() + degree + 1);
501 return knots;
502}
503
504std::vector<double> BSpline::Builder::knotVectorBuckets(
505 const std::vector<double>& values, unsigned int degree,
506 unsigned int maxSegments) const
507{
508 // Sort and remove duplicates
509 std::vector<double> unique = extractUniqueSorted(values);
510
511 // The minimum number of samples from which a free knot vector can be created
512 if (unique.size() < degree + 1)
513 {
514 std::ostringstream e;
515 e << "BSpline::Builder::knotVectorBuckets: Only " << unique.size()
516 << " unique sample points are given. A minimum of degree+1 = " << degree + 1
517 << " unique points are required to build a B-spline basis of degree " << degree
518 << ".";
519 throw Exception(e.str());
520 }
521
522 // Num internal knots (0 <= ni <= unique.size() - degree - 1)
523 unsigned int ni = unique.size() - degree - 1;
524 // Num segments
525 unsigned int ns = ni + degree + 1;
526
527 // Limit number of segments
528 if (ns > maxSegments && maxSegments >= degree + 1)
529 {
530 ns = maxSegments;
531 ni = ns - degree - 1;
532 }
533
534 // Num knots
535 // unsigned int nk = ns + degree + 1;
536
537 // Check numbers
538 if (ni > unique.size() - degree - 1)
539 {
540 throw Exception("BSpline::Builder::knotVectorBuckets: Invalid number of internal knots!");
541 }
542
543 // Compute window sizes
544 unsigned int w = 0;
545
546 if (ni > 0)
547 {
548 w = std::floor(unique.size() / ni);
549 }
550
551 // Residual
552 unsigned int res = unique.size() - w * ni;
553 // Create array with window sizes
554 std::vector<unsigned int> windows(ni, w);
555
556 // Add residual
557 for (unsigned int i = 0; i < res; ++i)
558 {
559 windows.at(i) += 1;
560 }
561
562 // Compute internal knots
563 std::vector<double> knots(ni, 0);
564 // Compute (n-k-2) interior knots using moving average
565 unsigned int index = 0;
566
567 for (unsigned int i = 0; i < ni; ++i)
568 {
569 for (unsigned int j = 0; j < windows.at(i); ++j)
570 {
571 knots.at(i) += unique.at(index + j);
572 }
573
574 knots.at(i) /= windows.at(i);
575 index += windows.at(i);
576 }
577
578 // Repeat first knot p + 1 times (for interpolation of start point)
579 for (unsigned int i = 0; i < degree + 1; ++i)
580 {
581 knots.insert(knots.begin(), unique.front());
582 }
583
584 // Repeat last knot p + 1 times (for interpolation of end point)
585 for (unsigned int i = 0; i < degree + 1; ++i)
586 {
587 knots.insert(knots.end(), unique.back());
588 }
589
590 return knots;
591}
592
593std::vector<double> BSpline::Builder::extractUniqueSorted(
594 const std::vector<double>& values) const
595{
596 // Sort and remove duplicates
597 std::vector<double> unique(values);
598 std::sort(unique.begin(), unique.end());
599 std::vector<double>::iterator it = unique_copy(unique.begin(), unique.end(),
600 unique.begin());
601 unique.resize(distance(unique.begin(), it));
602 return unique;
603}
604
605} // namespace SPLINTER
volScalarField & A
volScalarField & D
std::vector< double > linspace(double start, double stop, unsigned int num)
Definition utilities.C:73
label i
Definition pEqn.H:46