OpenMS  2.8.0
LinearRegression.h
Go to the documentation of this file.
1 // --------------------------------------------------------------------------
2 // OpenMS -- Open-Source Mass Spectrometry
3 // --------------------------------------------------------------------------
4 // Copyright The OpenMS Team -- Eberhard Karls University Tuebingen,
5 // ETH Zurich, and Freie Universitaet Berlin 2002-2021.
6 //
7 // This software is released under a three-clause BSD license:
8 // * Redistributions of source code must retain the above copyright
9 // notice, this list of conditions and the following disclaimer.
10 // * Redistributions in binary form must reproduce the above copyright
11 // notice, this list of conditions and the following disclaimer in the
12 // documentation and/or other materials provided with the distribution.
13 // * Neither the name of any author or any participating institution
14 // may be used to endorse or promote products derived from this software
15 // without specific prior written permission.
16 // For a full list of authors, refer to the file AUTHORS.
17 // --------------------------------------------------------------------------
18 // THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
19 // AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
20 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
21 // ARE DISCLAIMED. IN NO EVENT SHALL ANY OF THE AUTHORS OR THE CONTRIBUTING
22 // INSTITUTIONS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
23 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
24 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS;
25 // OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY,
26 // WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR
27 // OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF
28 // ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
29 //
30 // --------------------------------------------------------------------------
31 // $Maintainer: Timo Sachsenberg $
32 // $Authors: $
33 // --------------------------------------------------------------------------
34 
35 #pragma once
36 
37 #include <OpenMS/CONCEPT/Types.h>
40 
41 #include <cmath>
42 #include <vector>
43 
44 // forward declare WildMagic types so WMs headers don't pollute OpenMS library
45 namespace Wm5
46 {
47  template <typename Real>
48  class Vector2;
49 
51 }
52 
53 namespace OpenMS
54 {
55  namespace Math
56  {
73  class OPENMS_DLLAPI LinearRegression
74  {
75 public:
76 
79  intercept_(0),
80  slope_(0),
81  x_intercept_(0),
82  lower_(0),
83  upper_(0),
84  t_star_(0),
85  r_squared_(0),
86  stand_dev_residuals_(0),
87  mean_residuals_(0),
88  stand_error_slope_(0),
89  chi_squared_(0),
90  rsd_(0)
91  {
92  }
93 
95  virtual ~LinearRegression() = default;
96 
118  void computeRegression(double confidence_interval_P,
119  std::vector<double>::const_iterator x_begin,
120  std::vector<double>::const_iterator x_end,
121  std::vector<double>::const_iterator y_begin,
122  bool compute_goodness = true);
123 
146  void computeRegressionWeighted(double confidence_interval_P,
147  std::vector<double>::const_iterator x_begin,
148  std::vector<double>::const_iterator x_end,
149  std::vector<double>::const_iterator y_begin,
150  std::vector<double>::const_iterator w_begin,
151  bool compute_goodness = true);
152 
154  double getIntercept() const;
156  double getSlope() const;
158  double getXIntercept() const;
160  double getLower() const;
162  double getUpper() const;
164  double getTValue() const;
166  double getRSquared() const;
168  double getStandDevRes() const;
170  double getMeanRes() const;
172  double getStandErrSlope() const;
174  double getChiSquared() const;
176  double getRSD() const;
177 
179  static inline double computePointY(double x, double slope, double intercept)
180  {
181  return slope * x + intercept;
182  }
183 
184 protected:
185 
187  double intercept_;
189  double slope_;
191  double x_intercept_;
193  double lower_;
195  double upper_;
197  double t_star_;
199  double r_squared_;
207  double chi_squared_;
209  double rsd_;
210 
211 
213  void computeGoodness_(const std::vector<Wm5::Vector2d>& points, double confidence_interval_P);
214 
216  template <typename Iterator>
217  double computeChiSquare(Iterator x_begin, Iterator x_end, Iterator y_begin, double slope, double intercept);
218 
220  template <typename Iterator>
221  double computeWeightedChiSquare(Iterator x_begin, Iterator x_end, Iterator y_begin, Iterator w_begin, double slope, double intercept);
222 
223 private:
224 
229 
230  }; //class
231 
232  //x, y, w must be of same size
233  template <typename Iterator>
234  double LinearRegression::computeChiSquare(Iterator x_begin, Iterator x_end, Iterator y_begin, double slope, double intercept)
235  {
236  double chi_squared = 0.0;
237  Iterator xIter = x_begin;
238  Iterator yIter = y_begin;
239  for (; xIter != x_end; ++xIter, ++yIter)
240  {
241  chi_squared += std::pow(*yIter - computePointY(*xIter, slope, intercept), 2);
242  }
243 
244  return chi_squared;
245  }
246 
247  //x, y, w must be of same size
248  template <typename Iterator>
249  double LinearRegression::computeWeightedChiSquare(Iterator x_begin, Iterator x_end, Iterator y_begin, Iterator w_begin, double slope, double intercept)
250  {
251  double chi_squared = 0.0;
252  Iterator xIter = x_begin;
253  Iterator yIter = y_begin;
254  Iterator wIter = w_begin;
255  for (; xIter != x_end; ++xIter, ++yIter, ++wIter)
256  {
257  chi_squared += *wIter * std::pow(*yIter - computePointY(*xIter, slope, intercept), 2);
258  }
259 
260  return chi_squared;
261  }
262  } // namespace Math
263 } // namespace OpenMS
264 
265 
This class offers functions to perform least-squares fits to a straight line model,...
Definition: LinearRegression.h:74
void computeRegressionWeighted(double confidence_interval_P, std::vector< double >::const_iterator x_begin, std::vector< double >::const_iterator x_end, std::vector< double >::const_iterator y_begin, std::vector< double >::const_iterator w_begin, bool compute_goodness=true)
This function computes the best-fit linear regression coefficients of the model for the weighted da...
double r_squared_
The squared correlation coefficient (Pearson)
Definition: LinearRegression.h:199
double getRSquared() const
Non-mutable access to the squared Pearson coefficient.
void computeGoodness_(const std::vector< Wm5::Vector2d > &points, double confidence_interval_P)
Computes the goodness of the fitted regression line.
double getIntercept() const
Non-mutable access to the y-intercept of the straight line.
double x_intercept_
The intercept of the fitted line with the x-axis.
Definition: LinearRegression.h:191
double getUpper() const
Non-mutable access to the upper border of confidence interval.
LinearRegression()
Constructor.
Definition: LinearRegression.h:78
double lower_
The lower bound of the confidence interval.
Definition: LinearRegression.h:193
virtual ~LinearRegression()=default
Destructor.
void computeRegression(double confidence_interval_P, std::vector< double >::const_iterator x_begin, std::vector< double >::const_iterator x_end, std::vector< double >::const_iterator y_begin, bool compute_goodness=true)
This function computes the best-fit linear regression coefficients of the model for the dataset .
double getXIntercept() const
Non-mutable access to the x-intercept of the straight line.
static double computePointY(double x, double slope, double intercept)
given x compute y = slope * x + intercept
Definition: LinearRegression.h:179
double upper_
The upper bound of the confidence interval.
Definition: LinearRegression.h:195
LinearRegression & operator=(const LinearRegression &arg)
Not implemented.
double computeWeightedChiSquare(Iterator x_begin, Iterator x_end, Iterator y_begin, Iterator w_begin, double slope, double intercept)
Compute the chi squared of a weighted linear fit.
Definition: LinearRegression.h:249
double t_star_
The value of the t-statistic.
Definition: LinearRegression.h:197
double getRSD() const
Non-mutable access to relative standard deviation.
double computeChiSquare(Iterator x_begin, Iterator x_end, Iterator y_begin, double slope, double intercept)
Compute the chi squared of a linear fit.
Definition: LinearRegression.h:234
double getTValue() const
Non-mutable access to the value of the t-distribution.
double getStandErrSlope() const
Non-mutable access to the standard error of the slope.
double getSlope() const
Non-mutable access to the slope of the straight line.
double getChiSquared() const
Non-mutable access to the chi squared value.
double chi_squared_
The value of the Chi Squared statistic.
Definition: LinearRegression.h:207
double intercept_
The intercept of the fitted line with the y-axis.
Definition: LinearRegression.h:187
double slope_
The slope of the fitted line.
Definition: LinearRegression.h:189
double getLower() const
Non-mutable access to the lower border of confidence interval.
double mean_residuals_
Mean of residuals.
Definition: LinearRegression.h:203
double stand_dev_residuals_
The standard deviation of the residuals.
Definition: LinearRegression.h:201
LinearRegression(const LinearRegression &arg)
Not implemented.
double rsd_
the relative standard deviation
Definition: LinearRegression.h:209
double getMeanRes() const
Non-mutable access to the residual mean.
double stand_error_slope_
The standard error of the slope.
Definition: LinearRegression.h:205
double getStandDevRes() const
Non-mutable access to the standard deviation of the residuals.
Main OpenMS namespace.
Definition: FeatureDeconvolution.h:47
Definition: LinearRegression.h:46
Vector2< double > Vector2d
Definition: LinearRegression.h:48
Definition: LinearRegression.h:48