![]() |
OpenMS
|
Smoothing spline implementation using B-spline basis. More...
#include <OpenMS/MATH/MISC/BSplineSmoothingSpline.h>
Public Member Functions | |
| 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. | |
| BSplineSmoothingSpline (const BSplineSmoothingSpline &)=delete | |
| Deleted copy constructor (cannot copy owned resources) | |
| BSplineSmoothingSpline & | operator= (const BSplineSmoothingSpline &)=delete |
| Deleted copy assignment (cannot copy owned resources) | |
| BSplineSmoothingSpline (BSplineSmoothingSpline &&)=default | |
| Defaulted move constructor. | |
| BSplineSmoothingSpline & | operator= (BSplineSmoothingSpline &&)=default |
| Defaulted move assignment. | |
| double | eval (double x) const |
| Evaluate the smoothing spline at x. | |
| bool | ok () const |
| Check if spline fit was successful. | |
| int | num_interior_knots () const |
| Get the number of interior knots selected. | |
| double | rss () const |
| Get the achieved RSS. | |
| double | smoothing_param () const |
| Get the smoothing parameter used. | |
Private Types | |
| enum | FitType { POLYNOMIAL , BSPLINE } |
Private Member Functions | |
| 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. | |
| 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. | |
| double | compute_polynomial_rss (const std::vector< double > &x, const std::vector< double > &y) const |
| Compute RSS for polynomial fit. | |
| double | eval_polynomial (double x) const |
| Evaluate polynomial at x. | |
Private Attributes | |
| bool | ok_ = false |
| int | num_interior_knots_ = 0 |
| double | rss_ = 0.0 |
| double | s_ = 0.0 |
| int | k_ = 3 |
| std::vector< double > | x_ |
| std::vector< double > | y_ |
| FitType | fit_type_ = BSPLINE |
| std::vector< double > | poly_coeffs_ |
| std::unique_ptr< BSpline2d > | spline_ |
Smoothing spline implementation using B-spline basis.
This class implements a smoothing spline that uses B-spline basis functions with automatic knot selection based on a smoothing parameter, similar to scipy's UnivariateSpline.
Unlike the BSpline2d class which requires explicit node specification, this class implements the scipy-style optimization: minimize: RSS(f) + lambda * integral(f''(x))^2 dx
Or equivalently with smoothing parameter s: minimize: integral(f''(x))^2 dx subject to: RSS(f) <= s
The key difference from BSpline2d:
For compatibility with PyProphet's UnivariateSpline:
Implementation strategy: Since we cannot modify the eol-bspline library directly, we implement the smoothing spline optimization by:
|
private |
| 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.
| x | Data points (must be sorted) |
| y | Values at data points |
| s | Smoothing parameter (RSS target). If negative, uses automatic: s = m - sqrt(2*m) |
| k | Spline degree (default: 3 for cubic) |
The smoothing parameter s controls the trade-off between fit and smoothness:
Destructor - cleans up internal BSpline2d if used.
|
delete |
Deleted copy constructor (cannot copy owned resources)
|
default |
Defaulted move constructor.
|
private |
Compute RSS for polynomial fit.
|
private |
Compute RSS for given spline.
| double eval | ( | double | x | ) | const |
Evaluate the smoothing spline at x.
| x | Point at which to evaluate |
|
private |
Evaluate polynomial at x.
|
private |
Find optimal wavelength/nodes combination to match smoothing target.
Since BSpline2d doesn't support direct smoothing parameter s, we search for a wavelength that produces RSS close to the target s.
Strategy:
|
inline |
Get the number of interior knots selected.
|
inline |
Check if spline fit was successful.
|
default |
Defaulted move assignment.
|
delete |
Deleted copy assignment (cannot copy owned resources)
|
inline |
Get the achieved RSS.
|
inline |
Get the smoothing parameter used.
|
private |
Try polynomial fit (matches scipy when 0 interior knots)
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |