preCICE v3.2.0
Loading...
Searching...
No Matches
Bspline.cpp
Go to the documentation of this file.
1#include "math/Bspline.hpp"
3#include "utils/assertion.hpp"
4
5#include <Eigen/Core>
6#include <Eigen/Sparse>
7#include <algorithm>
8#include <cstdlib>
9#include <unsupported/Eigen/Splines>
10#include "math/differences.hpp"
11#include "utils/assertion.hpp"
12
13namespace precice::math {
14
15Bspline::Bspline(Eigen::VectorXd ts, const Eigen::MatrixXd &xs, int splineDegree)
16{
17
18 PRECICE_ASSERT(ts.size() >= 2, "Interpolation requires at least 2 samples");
19 PRECICE_ASSERT(std::is_sorted(ts.begin(), ts.end()), "Timestamps must be sorted");
20
21 // organize data in columns. Each column represents one sample in time.
22 PRECICE_ASSERT(xs.cols() == ts.size());
23 _ndofs = xs.rows(); // number of dofs. Each dof needs its own interpolant.
24 _tsMin = ts(0);
25 _tsMax = ts(ts.size() - 1);
26 auto relativeTime = [tsMin = _tsMin, tsMax = _tsMax](double t) -> double { return (t - tsMin) / (tsMax - tsMin); };
27 ts = ts.unaryExpr(relativeTime);
28
29 // The code for computing the knots and the control points is copied from Eigens bspline interpolation with some modifications
30 // https://gitlab.com/libeigen/eigen/-/blob/master/unsupported/Eigen/src/Splines/SplineFitting.h
31
32 // 1. Compute the knot vector
33 Eigen::KnotAveraging(ts, splineDegree, _knots);
34
35 // 2. Compute the control points
36 // We use a nxn sparse matrix with 2 + (n-2) * (d+1) entries and thus a fill-factor < 0.5.
37 Eigen::DenseIndex n = xs.cols();
39 matrixEntries.reserve(2 + (n - 2) * (splineDegree + 1));
40
41 matrixEntries.emplace_back(0, 0, 1.0);
42 for (Eigen::DenseIndex i = 1; i < n - 1; ++i) {
43 const Eigen::DenseIndex span = Eigen::Spline<double, 1>::Span(ts[i], splineDegree, _knots);
44 auto basisFunc = Eigen::Spline<double, 1>::BasisFunctions(ts[i], splineDegree, _knots);
45
46 for (Eigen::DenseIndex j = 0; j < splineDegree + 1; ++j) {
47 matrixEntries.emplace_back(i, span - splineDegree + j, basisFunc(j));
48 }
49 }
50 matrixEntries.emplace_back(n - 1, n - 1, 1.0);
51 PRECICE_ASSERT(matrixEntries.capacity() == matrixEntries.size(), matrixEntries.capacity(), matrixEntries.size(), n, splineDegree);
52
53 Eigen::SparseMatrix<double> A(n, n);
54 A.setFromTriplets(matrixEntries.begin(), matrixEntries.end());
55 A.makeCompressed();
56
57 Eigen::SparseQR<Eigen::SparseMatrix<double>, Eigen::COLAMDOrdering<int>> qr;
58 qr.analyzePattern(A);
59 qr.factorize(A);
60
61 _ctrls = qr.solve(xs.transpose());
62}
63
64Eigen::VectorXd Bspline::interpolateAt(double t) const
65{
66 // transform t to the relative interval [0; 1]
67 const double tRelative = std::clamp((t - _tsMin) / (_tsMax - _tsMin), 0.0, 1.0);
68
69 Eigen::VectorXd interpolated(_ndofs);
70 constexpr int splineDimension = 1;
71
72 for (int i = 0; i < _ndofs; i++) {
73 interpolated[i] = Eigen::Spline<double, splineDimension>(_knots, _ctrls.col(i))(tRelative)[0];
74 }
75
76 return interpolated;
77}
78} // namespace precice::math
Eigen::Vector2d ts
#define PRECICE_ASSERT(...)
Definition assertion.hpp:85
T begin(T... args)
T capacity(T... args)
T clamp(T... args)
Eigen::MatrixXd _ctrls
Definition Bspline.hpp:31
Eigen::VectorXd _knots
Definition Bspline.hpp:30
Eigen::VectorXd interpolateAt(double t) const
Samples the B-Spline interpolation.
Definition Bspline.cpp:64
Bspline(Eigen::VectorXd ts, const Eigen::MatrixXd &xs, int splineDegree)
Initialises the B-Spline interpolation with the given data (x0,t0), (x1,t1), ..., (xn,...
Definition Bspline.cpp:15
A C++ 11 implementation of the non-owning C++20 std::span type.
Definition span.hpp:284
T emplace_back(T... args)
T end(T... args)
T is_sorted(T... args)
provides general mathematical constants and functions.
Definition barycenter.cpp:9
T reserve(T... args)
T size(T... args)