OpenMS
Loading...
Searching...
No Matches
BSplineSmoothingSpline Class Reference

Smoothing spline implementation using B-spline basis. More...

#include <OpenMS/MATH/MISC/BSplineSmoothingSpline.h>

Collaboration diagram for BSplineSmoothingSpline:
[legend]

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)
 
BSplineSmoothingSplineoperator= (const BSplineSmoothingSpline &)=delete
 Deleted copy assignment (cannot copy owned resources)
 
 BSplineSmoothingSpline (BSplineSmoothingSpline &&)=default
 Defaulted move constructor.
 
BSplineSmoothingSplineoperator= (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< BSpline2dspline_
 

Detailed Description

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:

  • Uses automatic knot selection based on smoothing parameter
  • Default s = m - sqrt(2*m) where m is number of data points
  • Allows curve to deviate from data points (controlled smoothing)

Implementation strategy: Since we cannot modify the eol-bspline library directly, we implement the smoothing spline optimization by:

  1. Testing different wavelength/node combinations in BSpline2d
  2. Finding the combination that best satisfies the smoothing constraint
  3. Using the resulting spline for evaluation

Member Enumeration Documentation

◆ FitType

enum FitType
private
Enumerator
POLYNOMIAL 
BSPLINE 

Constructor & Destructor Documentation

◆ BSplineSmoothingSpline() [1/3]

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.

Parameters
xData points (must be sorted)
yValues at data points
sSmoothing parameter (RSS target). If negative, uses automatic: s = m - sqrt(2*m)
kSpline degree (default: 3 for cubic)

The smoothing parameter s controls the trade-off between fit and smoothness:

  • s = 0: Interpolating spline (passes through all points)
  • s > 0: Allows deviation from data points
  • s = m - sqrt(2*m): Default used by scipy.UnivariateSpline
  • s -> infinity: Approaches simple polynomial fit
Note
For small datasets (n < 4), falls back to interpolation

◆ ~BSplineSmoothingSpline()

Destructor - cleans up internal BSpline2d if used.

◆ BSplineSmoothingSpline() [2/3]

Deleted copy constructor (cannot copy owned resources)

◆ BSplineSmoothingSpline() [3/3]

Defaulted move constructor.

Member Function Documentation

◆ compute_polynomial_rss()

double compute_polynomial_rss ( const std::vector< double > &  x,
const std::vector< double > &  y 
) const
private

Compute RSS for polynomial fit.

◆ compute_rss()

double compute_rss ( const BSpline2d spl,
const std::vector< double > &  x,
const std::vector< double > &  y 
) const
private

Compute RSS for given spline.

◆ eval()

double eval ( double  x) const

Evaluate the smoothing spline at x.

Parameters
xPoint at which to evaluate
Returns
Smoothed value, or NaN if fit failed

◆ eval_polynomial()

double eval_polynomial ( double  x) const
private

Evaluate polynomial at x.

◆ fit_smoothing_spline()

void fit_smoothing_spline ( const std::vector< double > &  x,
const std::vector< double > &  y,
double  s_target 
)
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:

  1. Start with minimal nodes (single cubic segment)
  2. Gradually increase nodes if RSS > s (underfitting)
  3. Find the configuration where RSS ~= s

◆ num_interior_knots()

int num_interior_knots ( ) const
inline

Get the number of interior knots selected.

◆ ok()

bool ok ( ) const
inline

Check if spline fit was successful.

◆ operator=() [1/2]

BSplineSmoothingSpline & operator= ( BSplineSmoothingSpline &&  )
default

Defaulted move assignment.

◆ operator=() [2/2]

BSplineSmoothingSpline & operator= ( const BSplineSmoothingSpline )
delete

Deleted copy assignment (cannot copy owned resources)

◆ rss()

double rss ( ) const
inline

Get the achieved RSS.

◆ smoothing_param()

double smoothing_param ( ) const
inline

Get the smoothing parameter used.

◆ try_polynomial_fit()

bool try_polynomial_fit ( const std::vector< double > &  x,
const std::vector< double > &  y,
double  s_target,
int  degree = 3 
)
private

Try polynomial fit (matches scipy when 0 interior knots)

Returns
true if polynomial fit is acceptable (RSS <= s_target)

Member Data Documentation

◆ fit_type_

FitType fit_type_ = BSPLINE
private

◆ k_

int k_ = 3
private

◆ num_interior_knots_

int num_interior_knots_ = 0
private

◆ ok_

bool ok_ = false
private

◆ poly_coeffs_

std::vector<double> poly_coeffs_
private

◆ rss_

double rss_ = 0.0
private

◆ s_

double s_ = 0.0
private

◆ spline_

std::unique_ptr<BSpline2d> spline_
private

◆ x_

std::vector<double> x_
private

◆ y_

std::vector<double> y_
private