sablib
Loading...
Searching...
No Matches
pspline.cpp
Go to the documentation of this file.
1
6
7#include "../misc/diff.h"
8#include "pspline.h"
9
10namespace sablib {
11
12//
13// Implementation of PSpline() function
14//
15const std::vector<double> PSpline(
16 const std::vector<double> & y, const unsigned int knots_num,
17 const unsigned int degree, const unsigned int s, const double lambda
18)
19{
20 if(y.size() == 0) {
21 throw std::invalid_argument("PSpline(): the length of y is zero.");
22 }
23
24 if(knots_num == 0) {
25 throw std::invalid_argument("PSpline(): knots_num is zero.");
26 }
27
28 if(degree == 0) {
29 throw std::invalid_argument("PSpline(): degree is zero.");
30 }
31
32 if(y.size() < knots_num) {
33 throw std::invalid_argument("PSpline(): knots_num is larger than length of y.");
34 }
35
36 if(s == 0 || s > 3) {
37 throw std::invalid_argument("PSpline(): s must be 1, 2 or 3.");
38 }
39
40 if(lambda <= 0) {
41 throw std::invalid_argument("PSpline(): non-positive lambda value is given.");
42 }
43
44 Eigen::VectorXd yy = Eigen::VectorXd::Map(y.data(), y.size());
45 Eigen::VectorXd xx = Eigen::VectorXd::LinSpaced(yy.size(), 0, 1);
46
47 Eigen::VectorXd interior_knots = Eigen::VectorXd::LinSpaced(knots_num, 0, 1);
48 Eigen::VectorXd knots(knots_num + degree * 2);
49
50 for(unsigned int i = 0; i < degree; i++) {
51 knots(i) = interior_knots(0);
52 knots(knots.size() - 1 - i) = interior_knots(knots_num - 1);
53 }
54
55 knots.block(degree, 0, knots_num, 1) = interior_knots;
56
57 BSpline bs(degree, knots);
58 Eigen::SparseMatrix<double> B = bs.DesignMatrix(xx);
59 Eigen::SparseMatrix<double> I, D, lambdaDTD, BTB;
60
61 I.resize(B.cols(), B.cols());
62 I.setIdentity();
63 D = Diff(I, s);
64 lambdaDTD = lambda * (D.transpose() * D);
65 BTB = B.transpose() * B;
66
67 Eigen::SimplicialCholesky< Eigen::SparseMatrix<double> > solver;
68
69 solver.compute(BTB + lambdaDTD);
70
71 if(solver.info() != Eigen::Success) {
72 throw std::runtime_error("PSpline(): solver calculation fails.");
73 }
74
75 Eigen::VectorXd coefficients = solver.solve(B.transpose() * yy);
76 std::vector<double> result(yy.size());
77
78 for(int i = 0; i < xx.size(); i++) {
79 result[i] = bs.Interpolate(xx(i), coefficients);
80 }
81
82 return result;
83}
84
85}; // namespace sablib
A class for B-Spline operations including basis function calculation and interpolation.
Definition bspline.h:24
const Eigen::SparseMatrix< Scalar > DesignMatrix(const Eigen::VectorX< Scalar > &x) const
Constructs the design matrix (collocation matrix) for a given set of x-coordinates.
Definition bspline.h:267
const Scalar Interpolate(const Scalar x, const Eigen::VectorX< Scalar > &coefficients) const
Interpolates the value at a given x-coordinate using provided B-Spline coefficients.
Definition bspline.h:315
MATLAB-like diff() function.
const Derived::PlainObject Diff(const Eigen::MatrixBase< Derived > &m0, const int n=1, const Dir dir=Dir::RowWise)
Calculates the n-th discrete difference along the given axis.
Definition diff.h:32
const std::vector< double > PSpline(const std::vector< double > &y, const unsigned int knots_num, const unsigned int degree, const unsigned int s, const double lambda)
Smoothes the input data using P-Splines (Penalized B-Splines).
Definition pspline.cpp:15
Smoothing using penalized Spline(P-Spline).