OpenMS
|
Computes the Savitzky-Golay filter coefficients using QR decomposition. More...
#include <OpenMS/FILTERING/SMOOTHING/SavitzkyGolayFilter.h>
Public Member Functions | |
SavitzkyGolayFilter () | |
Constructor. More... | |
~SavitzkyGolayFilter () override | |
Destructor. More... | |
template<class InputIt , class OutputIt > | |
void | filter (InputIt first, InputIt last, OutputIt d_first) |
void | filter (MSSpectrum &spectrum) |
Removed the noise from an MSSpectrum containing profile data. More... | |
void | filter (MSChromatogram &chromatogram) |
Removed the noise from an MSChromatogram. More... | |
void | filterExperiment (PeakMap &map) |
Removed the noise from an MSExperiment containing profile data. More... | |
Public Member Functions inherited from ProgressLogger | |
ProgressLogger () | |
Constructor. More... | |
virtual | ~ProgressLogger () |
Destructor. More... | |
ProgressLogger (const ProgressLogger &other) | |
Copy constructor. More... | |
ProgressLogger & | operator= (const ProgressLogger &other) |
Assignment Operator. More... | |
void | setLogType (LogType type) const |
Sets the progress log that should be used. The default type is NONE! More... | |
LogType | getLogType () const |
Returns the type of progress log being used. More... | |
void | startProgress (SignedSize begin, SignedSize end, const String &label) const |
Initializes the progress display. More... | |
void | setProgress (SignedSize value) const |
Sets the current progress. More... | |
void | endProgress (UInt64 bytes_processed=0) const |
void | nextProgress () const |
increment progress by 1 (according to range begin-end) More... | |
Public Member Functions inherited from DefaultParamHandler | |
DefaultParamHandler (const String &name) | |
Constructor with name that is displayed in error messages. More... | |
DefaultParamHandler (const DefaultParamHandler &rhs) | |
Copy constructor. More... | |
virtual | ~DefaultParamHandler () |
Destructor. More... | |
DefaultParamHandler & | operator= (const DefaultParamHandler &rhs) |
Assignment operator. More... | |
virtual bool | operator== (const DefaultParamHandler &rhs) const |
Equality operator. More... | |
void | setParameters (const Param ¶m) |
Sets the parameters. More... | |
const Param & | getParameters () const |
Non-mutable access to the parameters. More... | |
const Param & | getDefaults () const |
Non-mutable access to the default parameters. More... | |
const String & | getName () const |
Non-mutable access to the name. More... | |
void | setName (const String &name) |
Mutable access to the name. More... | |
const std::vector< String > & | getSubsections () const |
Non-mutable access to the registered subsections. More... | |
Protected Member Functions | |
void | updateMembers_ () override |
This method is used to update extra member variables at the end of the setParameters() method. More... | |
Protected Member Functions inherited from DefaultParamHandler | |
void | defaultsToParam_ () |
Updates the parameters after the defaults have been set in the constructor. More... | |
Protected Attributes | |
std::vector< double > | coeffs_ |
Coefficients. More... | |
UInt | frame_size_ |
UInt of the filter kernel (number of pre-tabulated coefficients) More... | |
UInt | order_ |
The order of the smoothing polynomial. More... | |
Protected Attributes inherited from ProgressLogger | |
LogType | type_ |
time_t | last_invoke_ |
ProgressLoggerImpl * | current_logger_ |
Protected Attributes inherited from DefaultParamHandler | |
Param | param_ |
Container for current parameters. More... | |
Param | defaults_ |
Container for default parameters. This member should be filled in the constructor of derived classes! More... | |
std::vector< String > | subsections_ |
Container for registered subsections. This member should be filled in the constructor of derived classes! More... | |
String | error_name_ |
Name that is displayed in error messages during the parameter checking. More... | |
bool | check_defaults_ |
If this member is set to false no checking if parameters in done;. More... | |
bool | warn_empty_defaults_ |
If this member is set to false no warning is emitted when defaults are empty;. More... | |
Additional Inherited Members | |
Public Types inherited from ProgressLogger | |
enum | LogType { CMD , GUI , NONE } |
Possible log types. More... | |
Static Public Member Functions inherited from DefaultParamHandler | |
static void | writeParametersToMetaValues (const Param &write_this, MetaInfoInterface &write_here, const String &key_prefix="") |
Writes all parameters to meta values. More... | |
Static Protected Member Functions inherited from ProgressLogger | |
static String | logTypeToFactoryName_ (LogType type) |
Return the name of the factory product used for this log type. More... | |
Static Protected Attributes inherited from ProgressLogger | |
static int | recursion_depth_ |
Computes the Savitzky-Golay filter coefficients using QR decomposition.
This class represents a Savitzky-Golay lowpass-filter. The idea of the Savitzky-Golay filter is to find filtercoefficients that preserve higher moments, which means to approximate the underlying function within the moving window by a polynomial of higher order (typically quadratic or quartic). Therefore we least-squares fit for each data point a polynomial to all points \( f_i \) in the window and set \(g_i\) to be the value of that polynomial at position \( i \). This method is superior to adjacent averaging because it tends to preserve features of the data such as peak height and width, which are usually 'washed out' by adjacent averaging.
Because of the linearity of the problem, we can reduce the work computing by fitting in advance, for fictious data consisting of all zeros except for a singe 1 and then do the fits on the real data just by taking linear combinations. There are a particular sets of filter coefficients \( c_n \) which accomplish the process of polynomial least-squares fit inside a moving window. To get the symmetric coefficient-matrix \(C \in R^{frameSize \times frameSize}\) with
\[ C= \left[ \begin{array}{cccc} c_{0,0} & c_{0,1} & \cdots & c_{0,frameSize-1} \\ \vdots & & & \vdots \\ c_{frameSize-1,0} & c_{frameSize-1,2} & \ldots & c_{frameSize-1,frameSize-1} \end{array} \right]\]
The first (last) \( \frac{frameSize}{2} \) rows of \( C \) we need to smooth the first (last) \( frameSize \) data points of the signal. So we use for the smoothing of the first data point the data point itself and the next \( frameSize-1 \) future points. For the second point we take the first datapoint, the data point itself and \( frameSize-2 \) of rightward data points... . We compute the Matrix \( C \) by solving the underlying least-squares problems with the singular value decomposition. Here we demonstrate the computation of the first row of a coefficient-matrix \( C \) for a Savitzky-Golay Filter of order=3 and frameSize=5: The design-matrix for the least-squares fit of a linear combination of 3 basis functions to 5 data points is:
\[ A=\left[ \begin{array}{ccc} x_0^0 & x_0^1 & x_0^2 \\ \\ x_1^0 & x_1^1 & x_1^2 \\ \\ x_2^0 & x_2^1 & x_2^2 \\ \\ x_3^0 & x_3^1 & x_3^2 \\ \\ x_4^0 & x_4^1 & x_4^2 \end{array} \right]. \]
To smooth the first data point we have to create a design-matrix with \( x=[0,\ldots, frameSize-1] \). Now we have to solve the over-determined set of \( frameSize \) linear equations
\[ Ac=b \]
where \( b=[1,0,\ldots,0] \) represents the fictious data. Therefore we solve the normal equations of the least-squares problem
\[ A^TAc=A^Tb. \]
Now, it is possible to get
\[ c_n=\sum_{m=0}^8 \{(A^TA)^{-1}\}_{0,m} n^m, \]
with \( 0\le n \le 8\). Because we only need one row of the inverse matrix, it is possible to use LU decomposition with only a single backsubstitution. The vector \(c=[c_0,\ldots,c_8] \) represents the wanted coefficients. Note that the solution of a least-squares problem directly from the normal equations is faster than the singular value decomposition but rather susceptible to roundoff error!
Name | Type | Default | Restrictions | Description |
---|---|---|---|---|
frame_length | int | 11 | The number of subsequent data points used for smoothing. This number has to be uneven. If it is not, 1 will be added. |
|
polynomial_order | int | 4 | Order or the polynomial that is fitted. |
Constructor.
|
override |
Destructor.
|
inline |
|
inline |
Removed the noise from an MSChromatogram.
|
inline |
Removed the noise from an MSSpectrum containing profile data.
|
inline |
Removed the noise from an MSExperiment containing profile data.
References MSExperiment::getChromatogram(), MSExperiment::getChromatograms(), and MSExperiment::size().
|
overrideprotectedvirtual |
This method is used to update extra member variables at the end of the setParameters() method.
Also call it at the end of the derived classes' copy constructor and assignment operator.
The default implementation is empty.
Reimplemented from DefaultParamHandler.
|
protected |
Coefficients.
|
protected |
UInt of the filter kernel (number of pre-tabulated coefficients)
|
protected |
The order of the smoothing polynomial.