OpenMS  2.7.0
Classes | Functions
OpenMS::Math Namespace Reference

Math namespace. More...

Classes

class  BilinearInterpolation
 Provides access to bilinearly interpolated values (and derivatives) from discrete data points. Values beyond the given range of data points are implicitly taken as zero. More...
 
class  LinearInterpolation
 Provides access to linearly interpolated values (and derivatives) from discrete data points. Values beyond the given range of data points are implicitly taken as zero. More...
 
class  RandomShuffler
 
struct  RANSACParam
 A simple struct to carry all the parameters required for a RANSAC run. More...
 
class  RANSAC
 This class provides a generic implementation of the RANSAC outlier detection algorithm. Is implemented and tested after the SciPy reference: http://wiki.scipy.org/Cookbook/RANSAC. More...
 
class  RansacModel
 Generic plug-in template base class using 'Curiously recurring template pattern' (CRTP) to allow for arbitrary RANSAC models (e.g. linear or quadratic fits). More...
 
class  RansacModelLinear
 Implementation of a linear RANSAC model fit. More...
 
class  RansacModelQuadratic
 Implementation of a quadratic RANSAC model fit. More...
 
class  AsymmetricStatistics
 Internal class for asymmetric distributions. More...
 
class  AveragePosition
 Maintain an average position by summing up positions with weights. More...
 
class  BasicStatistics
 Calculates some basic statistical parameters of a distribution: sum, mean, variance, and provides the normal approximation. More...
 
class  GammaDistributionFitter
 Implements a fitter for the Gamma distribution. More...
 
class  GaussFitter
 Implements a fitter for Gaussian functions. More...
 
class  GumbelDistributionFitter
 Implements a fitter for the Gumbel distribution. More...
 
class  GumbelMaxLikelihoodFitter
 Implements a fitter for the Gumbel distribution. More...
 
class  Histogram
 Representation of a histogram. More...
 
class  LinearRegression
 This class offers functions to perform least-squares fits to a straight line model, $ Y(c,x) = c_0 + c_1 x $. More...
 
class  LinearRegressionWithoutIntercept
 This class offers functions to perform least-squares fits to a straight line model, $ Y(c,x) = c_0 + c_1 x $. More...
 
class  PosteriorErrorProbabilityModel
 Implements a mixture model of the inverse gumbel and the gauss distribution or a gaussian mixture. More...
 
class  QuadraticRegression
 
class  ROCCurve
 ROCCurves show the trade-off in sensitivity and specificity for binary classifiers using different cutoff values. More...
 
struct  SummaryStatistics
 Helper class to gather (and dump) some statistics from a e.g. vector<double>. More...
 

Functions

static double ceilDecimal (double x, int decPow)
 rounds x up to the next decimal power 10 ^ decPow More...
 
static double roundDecimal (double x, int decPow)
 rounds x to the next decimal power 10 ^ decPow More...
 
static double intervalTransformation (double x, double left1, double right1, double left2, double right2)
 transforms point x of interval [left1,right1] into interval [left2,right2] More...
 
double linear2log (double x)
 Transforms a number from linear to log10 scale. Avoids negative logarithms by adding 1. More...
 
double log2linear (double x)
 Transforms a number from log10 to to linear scale. Subtracts the 1 added by linear2log(double) More...
 
bool isOdd (UInt x)
 Returns true if the given integer is odd. More...
 
template<typename T >
round (T x)
 Rounds the value. More...
 
static bool approximatelyEqual (double a, double b, double tol)
 Returns if a is approximately equal b , allowing a tolerance of tol. More...
 
template<typename T >
gcd (T a, T b)
 Returns the greatest common divisor (gcd) of two numbers by applying the Euclidean algorithm. More...
 
template<typename T >
gcd (T a, T b, T &u1, T &u2)
 Returns the greatest common divisor by applying the extended Euclidean algorithm (Knuth TAoCP vol. 2, p342). Calculates u1, u2 and u3 (which is returned) so that a * u1 + b * u2 = u3 = gcd(a, b, u1, u2) More...
 
template<typename T >
getPPM (T mz_obs, T mz_ref)
 Compute parts-per-million of two m/z values. More...
 
template<typename T >
getPPMAbs (T mz_obs, T mz_ref)
 Compute absolute parts-per-million of two m/z values. More...
 
template<typename T >
ppmToMass (T ppm, T mz_ref)
 Compute the mass diff in [Th], given a ppm value and a reference point. More...
 
template<typename T >
ppmToMassAbs (T ppm, T mz_ref)
 
static std::pair< double, double > getTolWindow (double val, double tol, bool ppm)
 Return tolerance window around val given tolerance tol. More...
 
double factLn (UInt x)
 Return the ln(x!) of a value. More...
 
template<typename T1 >
T1::value_type quantile (const T1 &x, double q)
 Returns the value of the q th quantile (0-1) in a sorted non-empty vector @x. More...
 
template<class T >
void spline_bisection (const T &peak_spline, double const left_neighbor_mz, double const right_neighbor_mz, double &max_peak_mz, double &max_peak_int, double const threshold=1e-6)
 
template<typename ValueType , typename BinSizeType >
std::ostream & operator<< (std::ostream &os, const Histogram< ValueType, BinSizeType > &hist)
 Print the contents to a stream. More...
 
template<typename Iterator >
std::vector< Wm5::Vector2d > iteratorRange2Wm5Vectors (Iterator x_begin, Iterator x_end, Iterator y_begin)
 Copies the distance(x_begin,x_end) elements starting at x_begin and y_begin into the Wm5::Vector. More...
 
template<typename IteratorType >
static void checkIteratorsNotNULL (IteratorType begin, IteratorType end)
 Helper function checking if two iterators are not equal. More...
 
template<typename IteratorType >
static void checkIteratorsEqual (IteratorType begin, IteratorType end)
 Helper function checking if two iterators are equal. More...
 
template<typename IteratorType1 , typename IteratorType2 >
static void checkIteratorsAreValid (IteratorType1 begin_b, IteratorType1 end_b, IteratorType2 begin_a, IteratorType2 end_a)
 Helper function checking if an iterator and a co-iterator both have a next element. More...
 
template<typename IteratorType >
static double sum (IteratorType begin, IteratorType end)
 Calculates the sum of a range of values. More...
 
template<typename IteratorType >
static double mean (IteratorType begin, IteratorType end)
 Calculates the mean of a range of values. More...
 
template<typename IteratorType >
static double median (IteratorType begin, IteratorType end, bool sorted=false)
 Calculates the median of a range of values. More...
 
template<typename IteratorType >
double MAD (IteratorType begin, IteratorType end, double median_of_numbers)
 median absolute deviation (MAD) More...
 
template<typename IteratorType >
double MeanAbsoluteDeviation (IteratorType begin, IteratorType end, double mean_of_numbers)
 mean absolute deviation (MeanAbsoluteDeviation) More...
 
template<typename IteratorType >
static double quantile1st (IteratorType begin, IteratorType end, bool sorted=false)
 Calculates the first quantile of a range of values. More...
 
template<typename IteratorType >
static double quantile3rd (IteratorType begin, IteratorType end, bool sorted=false)
 Calculates the third quantile of a range of values. More...
 
template<typename IteratorType >
static double variance (IteratorType begin, IteratorType end, double mean=std::numeric_limits< double >::max())
 Calculates the variance of a range of values. More...
 
template<typename IteratorType >
static double sd (IteratorType begin, IteratorType end, double mean=std::numeric_limits< double >::max())
 Calculates the standard deviation of a range of values. More...
 
template<typename IteratorType >
static double absdev (IteratorType begin, IteratorType end, double mean=std::numeric_limits< double >::max())
 Calculates the absolute deviation of a range of values. More...
 
template<typename IteratorType1 , typename IteratorType2 >
static double covariance (IteratorType1 begin_a, IteratorType1 end_a, IteratorType2 begin_b, IteratorType2 end_b)
 Calculates the covariance of two ranges of values. More...
 
template<typename IteratorType1 , typename IteratorType2 >
static double meanSquareError (IteratorType1 begin_a, IteratorType1 end_a, IteratorType2 begin_b, IteratorType2 end_b)
 Calculates the mean square error for the values in [begin_a, end_a) and [begin_b, end_b) More...
 
template<typename IteratorType1 , typename IteratorType2 >
static double classificationRate (IteratorType1 begin_a, IteratorType1 end_a, IteratorType2 begin_b, IteratorType2 end_b)
 Calculates the classification rate for the values in [begin_a, end_a) and [begin_b, end_b) More...
 
template<typename IteratorType1 , typename IteratorType2 >
static double matthewsCorrelationCoefficient (IteratorType1 begin_a, IteratorType1 end_a, IteratorType2 begin_b, IteratorType2 end_b)
 Calculates the Matthews correlation coefficient for the values in [begin_a, end_a) and [begin_b, end_b) More...
 
template<typename IteratorType1 , typename IteratorType2 >
static double pearsonCorrelationCoefficient (IteratorType1 begin_a, IteratorType1 end_a, IteratorType2 begin_b, IteratorType2 end_b)
 Calculates the Pearson correlation coefficient for the values in [begin_a, end_a) and [begin_b, end_b) More...
 
template<typename Value >
static void computeRank (std::vector< Value > &w)
 Replaces the elements in vector w by their ranks. More...
 
template<typename IteratorType1 , typename IteratorType2 >
static double rankCorrelationCoefficient (IteratorType1 begin_a, IteratorType1 end_a, IteratorType2 begin_b, IteratorType2 end_b)
 calculates the rank correlation coefficient for the values in [begin_a, end_a) and [begin_b, end_b) More...
 

Detailed Description

Math namespace.

Uses bisection to find the maximum point of a spline.

Contains mathematical auxiliary functions.

Should work with BSpline2d and CubicSpline2d

Function Documentation

◆ computeRank()

static void OpenMS::Math::computeRank ( std::vector< Value > &  w)
static

Replaces the elements in vector w by their ranks.

Referenced by rankCorrelationCoefficient().

◆ factLn()

double OpenMS::Math::factLn ( UInt  x)
inline

Return the ln(x!) of a value.

This functions comes handy when there are large factorials in a ratio formula.

Parameters
xan integer value
Returns
natural logarithm of factorial x

◆ getPPM()

T OpenMS::Math::getPPM ( mz_obs,
mz_ref 
)

Compute parts-per-million of two m/z values.

The returned ppm value can be either positive (mz_obs > mz_ref) or negative (mz_obs < mz_ref)!

Parameters
mz_obsObserved (experimental) m/z
mz_refReference (theoretical) m/z
Returns
The ppm value

Referenced by getPPMAbs().

◆ getPPMAbs()

T OpenMS::Math::getPPMAbs ( mz_obs,
mz_ref 
)

Compute absolute parts-per-million of two m/z values.

The returned ppm value is always >= 0.

Parameters
mz_obsObserved (experimental) m/z
mz_refReference (theoretical) m/z
Returns
The absolute ppm value

References getPPM().

◆ getTolWindow()

static std::pair<double, double> OpenMS::Math::getTolWindow ( double  val,
double  tol,
bool  ppm 
)
inlinestatic

Return tolerance window around val given tolerance tol.

Note that when ppm is used, the window is not symmetric. In this case, (right - val) > (val - left), i.e., the tolerance window also includes the largest value x which still has val in *its* tolerance window for the given ppms, so the compatibility relation is symmetric.

Parameters
valValue
tolTolerance
ppmWhether tol is in ppm or absolute
Returns
Tolerance window boundaries

◆ iteratorRange2Wm5Vectors()

std::vector<Wm5::Vector2d> OpenMS::Math::iteratorRange2Wm5Vectors ( Iterator  x_begin,
Iterator  x_end,
Iterator  y_begin 
)

Copies the distance(x_begin,x_end) elements starting at x_begin and y_begin into the Wm5::Vector.

Referenced by LinearRegression::computeRegression(), LinearRegression::computeRegressionWeighted(), and QuadraticRegression::computeRegressionWeighted().

◆ operator<<()

std::ostream& OpenMS::Math::operator<< ( std::ostream &  os,
const Histogram< ValueType, BinSizeType > &  hist 
)

◆ ppmToMass()

T OpenMS::Math::ppmToMass ( ppm,
mz_ref 
)

Compute the mass diff in [Th], given a ppm value and a reference point.

The returned mass diff can be either positive (ppm > 0) or negative (ppm < 0)!

Parameters
ppmParts-per-million error
mz_refReference m/z
Returns
The mass diff in [Th]

Referenced by PpmTrait::allowedTol(), and ppmToMassAbs().

◆ ppmToMassAbs()

T OpenMS::Math::ppmToMassAbs ( ppm,
mz_ref 
)

References ppmToMass().

◆ quantile()

T1::value_type OpenMS::Math::quantile ( const T1 &  x,
double  q 
)

Returns the value of the q th quantile (0-1) in a sorted non-empty vector @x.

References OpenMS::Constants::h.

◆ spline_bisection()

void OpenMS::Math::spline_bisection ( const T &  peak_spline,
double const  left_neighbor_mz,
double const  right_neighbor_mz,
double &  max_peak_mz,
double &  max_peak_int,
double const  threshold = 1e-6 
)