![]() |
OpenMS
|
Math namespace. More...
Classes | |
| struct | AdaptiveQuantileResult |
| Result of adaptiveQuantile computation. More... | |
| class | BasicStatistics |
| Calculates some basic statistical parameters of a distribution: sum, mean, variance, and provides the normal approximation. More... | |
| 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 | BSplineSmoothingSpline |
| Smoothing spline implementation using B-spline basis. 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... | |
| struct | KernelDensityEstimation |
| Kernel Density Estimation utilities using FFT-based methods. 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 | 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... | |
| struct | MultipleTesting |
| Statistical functions for multiple testing correction. More... | |
| struct | Pi0Result |
| Result of pi0 estimation for multiple testing correction. More... | |
| class | PosteriorErrorProbabilityModel |
| Implements a mixture model of the inverse gumbel and the gauss distribution or a gaussian mixture. More... | |
| class | QuadraticRegression |
| class | RandomShuffler |
| struct | RankData |
| Rank items (1-based) with SciPy-like tie handling. 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... | |
| struct | RANSACParam |
| A simple struct to carry all the parameters required for a RANSAC run. More... | |
| 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... | |
Typedefs | |
| using | BinContainer = std::vector< RangeBase > |
Functions | |
| template<typename T > | |
| bool | extendRange (T &min, T &max, const T &value) |
| Given an interval/range and a new value, extend the range to include the new value if needed. | |
| template<typename T > | |
| bool | contains (T value, T min, T max) |
Is a value contained in [min, max] ? | |
| std::pair< double, double > | zoomIn (const double left, const double right, const float factor, const float align) |
Zoom into an interval [left, right], decreasing its width by factor (which must be in [0,inf]). | |
| BinContainer | createBins (double min, double max, uint32_t number_of_bins, double extend_margin=0) |
Split a range [min,max] into number_of_bins (with optional overlap) and return the ranges of each bin. | |
| double | ceilDecimal (double x, int decPow) |
rounds x up to the next decimal power 10 ^ decPow | |
| double | roundDecimal (double x, int decPow) |
rounds x to the next decimal power 10 ^ decPow | |
| double | intervalTransformation (double x, double left1, double right1, double left2, double right2) |
transforms point x of interval [left1,right1] into interval [left2,right2] | |
| double | linear2log (double x) |
| Transforms a number from linear to log10 scale. Avoids negative logarithms by adding 1. | |
| double | log2linear (double x) |
| Transforms a number from log10 to to linear scale. Subtracts the 1 added by linear2log(double) | |
| bool | isOdd (UInt x) |
| Returns true if the given integer is odd. | |
| template<typename T > | |
| T | round (T x) |
| Rounds the value. | |
| template<typename T > | |
| T | roundTo (const T value, int digits) |
| template<typename T > | |
| double | percentOf (T value, T total, int digits) |
| bool | approximatelyEqual (double a, double b, double tol) |
Returns if a is approximately equal b , allowing a tolerance of tol. | |
| template<typename T > | |
| T | gcd (T a, T b) |
| Returns the greatest common divisor (gcd) of two numbers by applying the Euclidean algorithm. | |
| template<typename T > | |
| 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) | |
| template<typename T > | |
| T | getPPM (T mz_obs, T mz_ref) |
| Compute parts-per-million of two m/z values. | |
| template<typename T > | |
| T | getPPMAbs (T mz_obs, T mz_ref) |
| Compute absolute parts-per-million of two m/z values. | |
| template<typename T > | |
| T | ppmToMass (T ppm, T mz_ref) |
| Compute the mass diff in [Th], given a ppm value and a reference point. | |
| template<typename T > | |
| T | ppmToMassAbs (T ppm, T mz_ref) |
| std::pair< double, double > | getTolWindow (double val, double tol, bool ppm) |
Return tolerance window around val given tolerance tol. | |
| 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. | |
| double | log_binomial_coef (unsigned n, unsigned k) |
| Calculate logarithm of binomial coefficient C(n,k) using log-gamma function. | |
| double | log_sum_exp (double x, double y) |
| Log-sum-exp operation for numerical stability. | |
| double | binomial_cdf_complement (unsigned N, unsigned n, double p) |
| Calculate binomial cumulative distribution function P(X ≥ n) | |
| 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 IteratorType > | |
| static void | checkIteratorsNotNULL (IteratorType begin, IteratorType end) |
| Helper function checking if two iterators are not equal. | |
| template<typename IteratorType > | |
| static void | checkIteratorsEqual (IteratorType begin, IteratorType end) |
| Helper function checking if two iterators are equal. | |
| 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. | |
| template<typename IteratorType > | |
| static double | sum (IteratorType begin, IteratorType end) |
| Calculates the sum of a range of values. | |
| template<typename IteratorType > | |
| static double | mean (IteratorType begin, IteratorType end) |
| Calculates the mean of a range of values. | |
| template<typename IteratorType > | |
| static double | median (IteratorType begin, IteratorType end, bool sorted=false) |
| Calculates the median of a range of values. | |
| template<typename IteratorType > | |
| double | MAD (IteratorType begin, IteratorType end, double median_of_numbers) |
| median absolute deviation (MAD) | |
| template<typename IteratorType > | |
| double | MeanAbsoluteDeviation (IteratorType begin, IteratorType end, double mean_of_numbers) |
| mean absolute deviation (MeanAbsoluteDeviation) | |
| template<typename IteratorType > | |
| static double | quantile1st (IteratorType begin, IteratorType end, bool sorted=false) |
| Calculates the first quantile of a range of values. | |
| template<typename IteratorType > | |
| static double | quantile3rd (IteratorType begin, IteratorType end, bool sorted=false) |
| Calculates the third quantile of a range of values. | |
| template<typename IteratorType > | |
| static double | quantile (IteratorType begin, IteratorType end, double q) |
| Calculates the q-quantile (0 <= q <= 1) of a sorted range of values. | |
| template<typename IteratorType > | |
| double | tukeyUpperFence (IteratorType begin, IteratorType end, double k=1.5) |
| Tukey upper fence (UF) for outlier detection. | |
| template<typename IteratorType > | |
| double | tailFractionAbove (IteratorType begin, IteratorType end, double threshold) |
| Fraction of values above a threshold. | |
| template<typename IteratorType > | |
| double | winsorizedQuantile (IteratorType begin, IteratorType end, double q, double upper_fence) |
| Quantile after winsorizing at an upper fence. | |
| template<typename IteratorType > | |
| AdaptiveQuantileResult | adaptiveQuantile (IteratorType begin, IteratorType end, double q, double k=1.5, double r_sparse=0.01, double r_dense=0.10) |
| Adaptive quantile that blends RAW and IQR-winsorized quantiles based on tail density beyond the Tukey upper fence. | |
| template<typename IteratorType > | |
| static double | variance (IteratorType begin, IteratorType end, double mean=std::numeric_limits< double >::max()) |
| 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. | |
| 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. | |
| 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. | |
| 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) | |
| template<typename IteratorType1 , typename IteratorType2 > | |
| static double | rootMeanSquareError (IteratorType1 begin_a, IteratorType1 end_a, IteratorType2 begin_b, IteratorType2 end_b) |
| Calculates the root mean square error (RMSE) for the values in [begin_a, end_a) and [begin_b, end_b) | |
| 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) | |
| 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) | |
| 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) | |
| template<typename Value > | |
| static void | computeRank (std::vector< Value > &w) |
Replaces the elements in vector w by their ranks. | |
| 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) | |
| template<typename ValueType , typename BinSizeType > | |
| std::ostream & | operator<< (std::ostream &os, const Histogram< ValueType, BinSizeType > &hist) |
| Print the contents to a stream. | |
| std::vector< double > | qValue (const std::vector< double > &p_values, double pi0, bool pfdr=false) |
| Backward-compatible wrapper for MultipleTesting::qValue. | |
| Pi0Result | pi0Est (const std::vector< double > &p_values, const std::vector< double > &lambda_=std::vector< double >(), const std::string &pi0_method="smoother", int smooth_df=3, bool smooth_log_pi0=false) |
| Backward-compatible wrapper for MultipleTesting::pi0Est (string-based API) | |
| std::vector< double > | lfdr (const std::vector< double > &p_values, double pi0, bool trunc=true, bool monotone=true, const std::string &transf="probit", double adj=1.5, double eps=1e-8, std::size_t gridsize=512, double cut=3.0) |
| Backward-compatible wrapper for MultipleTesting::lfdr (string-based API) | |
| std::vector< double > | pNorm (const std::vector< double > &stat, const std::vector< double > &stat0) |
| Backward-compatible wrapper for MultipleTesting::pNorm. | |
| template<class T > | |
| std::vector< double > | computeModelFDR (const std::vector< T > &data_in) |
| Backward-compatible wrapper for MultipleTesting::computeModelFDR. | |
| template<class T > | |
| std::vector< double > | pEmp (const std::vector< T > &stat, const std::vector< T > &stat0) |
| Backward-compatible wrapper for MultipleTesting::pEmp. | |
| std::vector< double > | rankdata_double (std::vector< double > a, RankData::Method m=RankData::Method::Average, RankData::NaNPolicy p=RankData::NaNPolicy::Propagate) |
| Free function overloads (by-value) for double input. | |
| std::vector< double > | rankdata_float (std::vector< float > a, RankData::Method m=RankData::Method::Average, RankData::NaNPolicy p=RankData::NaNPolicy::Propagate) |
| Free function overloads (by-value) for float input. | |
| std::vector< double > | rankdata_int (std::vector< int > a, RankData::Method m=RankData::Method::Average, RankData::NaNPolicy p=RankData::NaNPolicy::Propagate) |
| Free function overloads (by-value) for int input. | |
Math namespace.
Uses bisection to find the maximum point of a spline.
Contains mathematical auxiliary functions.
Should work with BSpline2d and CubicSpline2d
| struct OpenMS::Math::AdaptiveQuantileResult |
Result of adaptiveQuantile computation.
Fields:
| Class Members | ||
|---|---|---|
| double | blended {0.0} | |
| double | half_raw {0.0} | |
| double | half_rob {0.0} | |
| double | tail_fraction {0.0} | |
| double | upper_fence {std::numeric_limits<double>::infinity()} | |
| double | weight {0.0} | |
| struct OpenMS::Math::Pi0Result |
Result of pi0 estimation for multiple testing correction.
This structure contains the estimated proportion of true null hypotheses (pi0), which is a key parameter in FDR estimation. The pi0 value represents the proportion of hypotheses that are truly null (i.e., no real effect) among all tested hypotheses.
The estimation uses the method described in: Storey JD and Tibshirani R. (2003) "Statistical significance for genome-wide experiments." PNAS 100: 9440-9445. doi: 10.1073/pnas.1530509100
Fields:
| Class Members | ||
|---|---|---|
| vector< double > | lambda_ | |
| double | pi0 = 1.0 | |
| vector< double > | pi0_lambda | |
| bool | pi0_smooth = false | |
| using BinContainer = std::vector<RangeBase> |
| AdaptiveQuantileResult adaptiveQuantile | ( | IteratorType | begin, |
| IteratorType | end, | ||
| double | q, | ||
| double | k = 1.5, |
||
| double | r_sparse = 0.01, |
||
| double | r_dense = 0.10 |
||
| ) |
Adaptive quantile that blends RAW and IQR-winsorized quantiles based on tail density beyond the Tukey upper fence.
Let UF = Q3 + k*IQR on the (finite) inputs. Compute:
Blend with weight w(r): r ≤ r_sparse -> w=0 (use robust) r ≥ r_dense -> w=1 (use raw) otherwise -> linear interpolation between 0 and 1
Returned value = (1-w)*half_rob + w*half_raw.
This keeps windows stable when outliers are sparse, while respecting genuinely broad tails (dense outliers) by leaning toward the raw quantile.
References:
| IteratorType | input iterator over arithmetic values |
| [in] | begin | start iterator |
| [in] | end | past-the-end iterator |
| [in] | q | target quantile in [0,1] (e.g., 0.99 for 99% half-width) |
| [in] | k | Tukey factor (default 1.5) |
| [in] | r_sparse | tail density below which robust wins (default 0.01 = 1%) |
| [in] | r_dense | tail density above which raw wins (default 0.10 = 10%) |
References AdaptiveQuantileResult::blended, AdaptiveQuantileResult::half_raw, AdaptiveQuantileResult::half_rob, quantile(), AdaptiveQuantileResult::tail_fraction, tailFractionAbove(), tukeyUpperFence(), AdaptiveQuantileResult::upper_fence, AdaptiveQuantileResult::weight, and winsorizedQuantile().
|
inline |
Calculate binomial cumulative distribution function P(X ≥ n)
Calculates P(X ≥ n) for a binomial distribution with parameters N and p, using numerically stable algorithms in the log domain to handle large values.
| [in] | N | Total number of trials |
| [in] | n | Minimum number of successes |
| [in] | p | Probability of success in each trial |
| std::invalid_argument | if parameters are invalid |
|
inline |
Backward-compatible wrapper for MultipleTesting::computeModelFDR.
References MultipleTesting::computeModelFDR().
Referenced by MultipleTesting::pNorm().
|
static |
Replaces the elements in vector w by their ranks.
Referenced by rankCorrelationCoefficient().
| bool contains | ( | T | value, |
| T | min, | ||
| T | max | ||
| ) |
Is a value contained in [min, max] ?
| T | Type, e.g. double |
Referenced by RangeBase::contains(), ListUtils::contains(), and QApplicationTOPP::QApplicationTOPP().
|
inline |
Split a range [min,max] into number_of_bins (with optional overlap) and return the ranges of each bin.
Optionally, bins can be made overlapping, by extending each bins' left and right margin by extend_margin. The overlap between neighboring bins will thus be 2 x extend_margin. The borders of the original interval will not be extended.
| [in] | min | The minimum of the range; must be smaller than max |
| [in] | max | The maximum of the range |
| [in] | number_of_bins | How many bins should the range be divided into? Must be 1 or larger |
| [in] | extend_margin | Overlap of neighboring bins (=0 for no overlap). Negative values will shrink the range (feature). |
number_of_bins elements, each representing the margins of one bin| OpenMS::Precondition | if min >= max or number_of_bins == 0 |
References OPENMS_PRECONDITION.
| bool extendRange | ( | T & | min, |
| T & | max, | ||
| const T & | value | ||
| ) |
Given an interval/range and a new value, extend the range to include the new value if needed.
| [in,out] | min | The current minimum of the range |
| [in,out] | max | The current maximum of the range |
| [in] | value | The new value which may extend the range |
| T getPPM | ( | T | mz_obs, |
| T | 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)!
| [in] | mz_obs | Observed (experimental) m/z |
| [in] | mz_ref | Reference (theoretical) m/z |
Referenced by getPPMAbs().
| T getPPMAbs | ( | T | mz_obs, |
| T | mz_ref | ||
| ) |
Compute absolute parts-per-million of two m/z values.
The returned ppm value is always >= 0.
| [in] | mz_obs | Observed (experimental) m/z |
| [in] | mz_ref | Reference (theoretical) m/z |
References getPPM().
|
inline |
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.
| [in] | val | Value |
| [in] | tol | Tolerance |
| [in] | ppm | Whether tol is in ppm or absolute |
| std::vector< double > lfdr | ( | const std::vector< double > & | p_values, |
| double | pi0, | ||
| bool | trunc = true, |
||
| bool | monotone = true, |
||
| const std::string & | transf = "probit", |
||
| double | adj = 1.5, |
||
| double | eps = 1e-8, |
||
| std::size_t | gridsize = 512, |
||
| double | cut = 3.0 |
||
| ) |
Backward-compatible wrapper for MultipleTesting::lfdr (string-based API)
|
inline |
|
inline |
Log-sum-exp operation for numerical stability.
| [in] | x | First logarithmic value |
| [in] | y | Second logarithmic value |
| std::ostream & operator<< | ( | std::ostream & | os, |
| const Histogram< ValueType, BinSizeType > & | hist | ||
| ) |
Print the contents to a stream.
References Histogram< ValueType, BinSizeType >::centerOfBin(), and Histogram< ValueType, BinSizeType >::size().
|
inline |
Backward-compatible wrapper for MultipleTesting::pEmp.
References MultipleTesting::pEmp().
Referenced by MultipleTesting::pNorm().
| double percentOf | ( | T | value, |
| T | total, | ||
| int | digits | ||
| ) |
Computes the percentage of value in relation to total, rounded to digits.
total is zero, the function returns 0.0 to avoid division by zero.| [in] | value | The value to compute the percentage for |
| [in] | total | The total value to compute the percentage against |
| [in] | digits | The number of digits to round the result to |
value in relation to total, rounded to digits | OpenMS::Exception::InvalidValue | if value or total is negative. |
References roundTo().
| Pi0Result pi0Est | ( | const std::vector< double > & | p_values, |
| const std::vector< double > & | lambda_ = std::vector< double >(), |
||
| const std::string & | pi0_method = "smoother", |
||
| int | smooth_df = 3, |
||
| bool | smooth_log_pi0 = false |
||
| ) |
Backward-compatible wrapper for MultipleTesting::pi0Est (string-based API)
|
inline |
Backward-compatible wrapper for MultipleTesting::pNorm.
References MultipleTesting::pNorm().
| T ppmToMass | ( | T | ppm, |
| T | 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)!
| [in] | ppm | Parts-per-million error |
| [in] | mz_ref | Reference m/z |
Referenced by PpmTrait::allowedTol(), and ppmToMassAbs().
| T ppmToMassAbs | ( | T | ppm, |
| T | mz_ref | ||
| ) |
References ppmToMass().
| 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.
Referenced by adaptiveQuantile(), tukeyUpperFence(), and winsorizedQuantile().
|
inline |
Backward-compatible wrapper for MultipleTesting::qValue.
References MultipleTesting::qValue().
|
inline |
Free function overloads (by-value) for double input.
|
inline |
Free function overloads (by-value) for float input.
|
inline |
Free function overloads (by-value) for int input.
| T roundTo | ( | const T | value, |
| int | digits | ||
| ) |
Rounds to the i'th digit after the decimal point, also works for negative digits.
e.g.
| [in] | value | The value to round |
| [in] | digits | The number of digits to round to (can be negative) |
Referenced by percentOf().
| 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 |
||
| ) |
| double tailFractionAbove | ( | IteratorType | begin, |
| IteratorType | end, | ||
| double | threshold | ||
| ) |
Fraction of values above a threshold.
| IteratorType | input iterator over arithmetic values |
| [in] | begin | start iterator |
| [in] | end | past-the-end iterator |
| [in] | threshold | threshold T |
Referenced by adaptiveQuantile().
| double tukeyUpperFence | ( | IteratorType | begin, |
| IteratorType | end, | ||
| double | k = 1.5 |
||
| ) |
Tukey upper fence (UF) for outlier detection.
Computes Q3 + k * IQR on the (finite) values in [begin,end). If there are too few values or IQR ≤ 0, returns +infinity.
References: J. W. Tukey (1977). Exploratory Data Analysis.
| IteratorType | input iterator over arithmetic values |
| [in] | begin | start iterator |
| [in] | end | past-the-end iterator |
| [in] | k | Tukey factor (default 1.5) |
References quantile().
Referenced by adaptiveQuantile().
|
static |
@brief Calculates the variance of a range of values
The mean can be provided explicitly to save computation time. If left at default, it will be computed internally.
@exception Exception::InvalidRange is thrown if the range is empty @ingroup MathFunctionsStatistics
References checkIteratorsNotNULL(), and mean().
Referenced by BasicStatistics< RealT >::normalDensity_sqrt2pi(), sd(), BasicStatistics< RealT >::setVariance(), and SummaryStatistics< T >::SummaryStatistics().
| double winsorizedQuantile | ( | IteratorType | begin, |
| IteratorType | end, | ||
| double | q, | ||
| double | upper_fence | ||
| ) |
Quantile after winsorizing at an upper fence.
Copies the (finite) values in [begin,end), caps them at upper_fence (and at 0 on the lower side, which is convenient for absolute residuals), then returns the requested quantile.
If upper_fence is not finite, this falls back to the raw quantile.
References: J. W. Tukey (1962). The Future of Data Analysis.
| IteratorType | input iterator over arithmetic values |
| [in] | begin | start iterator |
| [in] | end | past-the-end iterator |
| [in] | q | quantile in [0,1] |
| [in] | upper_fence | winsorization cap (Q3+k*IQR), or +inf to disable |
References quantile().
Referenced by adaptiveQuantile().
|
inline |
Zoom into an interval [left, right], decreasing its width by factor (which must be in [0,inf]).
To actually zoom in, the factor needs to be within [0,1]. Chosing a factor > 1 actually zooms out. align (between [0,1]) determines where the zoom happens: i.e. align = 0 leaves left the same and reduces right (or extends if factor>1) whereas align = 0.5 zooms into the center of the range etc
You can do round trips, i.e. undo a zoom in, by inverting the factor:
| left | Start of interval |
| right | End of interval |
| factor | Number between [0,1] to shrink, or >1 to extend the span (=right-left) |
| align | Where to position the smaller/shrunk interval (0 = left, 1 = right, 0.5=center etc) |
References OPENMS_PRECONDITION.