OpenMS
Loading...
Searching...
No Matches
BSplineSmoothingSpline.h
Go to the documentation of this file.
1// --------------------------------------------------------------------------
2// OpenMS -- Open-Source Mass Spectrometry
3// --------------------------------------------------------------------------
4// SPDX-License-Identifier: BSD-3-Clause
5//
6// --------------------------------------------------------------------------
7// $Maintainer: Justin Sing $
8// $Authors: Justin Sing $
9// --------------------------------------------------------------------------
10
11#pragma once
12
13#include <OpenMS/config.h>
15#include <vector>
16#include <memory>
17
18namespace OpenMS
19{
20 namespace Math
21 {
55 class OPENMS_DLLAPI BSplineSmoothingSpline
56 {
57 public:
74 BSplineSmoothingSpline(const std::vector<double>& x,
75 const std::vector<double>& y,
76 double s = -1.0,
77 int k = 3);
78
83
88
93
98
103
110 double eval(double x) const;
111
115 bool ok() const { return ok_; }
116
120 int num_interior_knots() const { return num_interior_knots_; }
121
125 double rss() const { return rss_; }
126
130 double smoothing_param() const { return s_; }
131
132 private:
133 bool ok_ = false;
134 int num_interior_knots_ = 0;
135 double rss_ = 0.0;
136 double s_ = 0.0;
137 int k_ = 3;
138 std::vector<double> x_;
139 std::vector<double> y_;
140
141 // Fit type
142 enum FitType { POLYNOMIAL, BSPLINE };
143 FitType fit_type_ = BSPLINE;
144
145 // Polynomial coefficients [a0, a1, a2, a3] for a0 + a1*x + a2*x^2 + a3*x^3
146 std::vector<double> poly_coeffs_;
147
148 // The fitted BSpline (owned, only used if fit_type_ == BSPLINE)
149 std::unique_ptr<BSpline2d> spline_;
150
162 void fit_smoothing_spline(const std::vector<double>& x,
163 const std::vector<double>& y,
164 double s_target);
165
171 bool try_polynomial_fit(const std::vector<double>& x,
172 const std::vector<double>& y,
173 double s_target,
174 int degree = 3);
175
179 double compute_rss(const BSpline2d* spl,
180 const std::vector<double>& x,
181 const std::vector<double>& y) const;
182
186 double compute_polynomial_rss(const std::vector<double>& x,
187 const std::vector<double>& y) const;
188
192 double eval_polynomial(double x) const;
193 };
194
195 } // namespace Math
196} // namespace OpenMS
b spline interpolation
Definition BSpline2d.h:31
Smoothing spline implementation using B-spline basis.
Definition BSplineSmoothingSpline.h:56
bool ok() const
Check if spline fit was successful.
Definition BSplineSmoothingSpline.h:115
void fit_smoothing_spline(const std::vector< double > &x, const std::vector< double > &y, double s_target)
Find optimal wavelength/nodes combination to match smoothing target.
std::vector< double > poly_coeffs_
Definition BSplineSmoothingSpline.h:146
FitType
Definition BSplineSmoothingSpline.h:142
BSplineSmoothingSpline(BSplineSmoothingSpline &&)=default
Defaulted move constructor.
double smoothing_param() const
Get the smoothing parameter used.
Definition BSplineSmoothingSpline.h:130
double eval(double x) const
Evaluate the smoothing spline at x.
int num_interior_knots() const
Get the number of interior knots selected.
Definition BSplineSmoothingSpline.h:120
BSplineSmoothingSpline(const BSplineSmoothingSpline &)=delete
Deleted copy constructor (cannot copy owned resources)
double compute_polynomial_rss(const std::vector< double > &x, const std::vector< double > &y) const
Compute RSS for polynomial fit.
std::vector< double > y_
Definition BSplineSmoothingSpline.h:139
std::vector< double > x_
Definition BSplineSmoothingSpline.h:138
BSplineSmoothingSpline & operator=(BSplineSmoothingSpline &&)=default
Defaulted move assignment.
BSplineSmoothingSpline(const std::vector< double > &x, const std::vector< double > &y, double s=-1.0, int k=3)
Construct a smoothing spline with automatic knot selection.
~BSplineSmoothingSpline()
Destructor - cleans up internal BSpline2d if used.
bool try_polynomial_fit(const std::vector< double > &x, const std::vector< double > &y, double s_target, int degree=3)
Try polynomial fit (matches scipy when 0 interior knots)
double compute_rss(const BSpline2d *spl, const std::vector< double > &x, const std::vector< double > &y) const
Compute RSS for given spline.
std::unique_ptr< BSpline2d > spline_
Definition BSplineSmoothingSpline.h:149
double eval_polynomial(double x) const
Evaluate polynomial at x.
double rss() const
Get the achieved RSS.
Definition BSplineSmoothingSpline.h:125
BSplineSmoothingSpline & operator=(const BSplineSmoothingSpline &)=delete
Deleted copy assignment (cannot copy owned resources)
Main OpenMS namespace.
Definition openswathalgo/include/OpenMS/OPENSWATHALGO/DATAACCESS/ISpectrumAccess.h:19