Statistics functions, probability distributions
// global functions template <class VECTORLIKE1, class VECTORLIKE2, class MAT> MAT::Scalar mrpt::math::mahalanobisDistance2( const VECTORLIKE1& X, const VECTORLIKE2& MU, const MAT& COV ); template <class VECTORLIKE1, class VECTORLIKE2, class MAT> VECTORLIKE1::Scalar mrpt::math::mahalanobisDistance( const VECTORLIKE1& X, const VECTORLIKE2& MU, const MAT& COV ); template <class VECTORLIKE, class MAT1, class MAT2, class MAT3> MAT1::Scalar mrpt::math::mahalanobisDistance2( const VECTORLIKE& mean_diffs, const MAT1& COV1, const MAT2& COV2, const MAT3& CROSS_COV12 ); template <class VECTORLIKE, class MAT1, class MAT2, class MAT3> VECTORLIKE::Scalar mrpt::math::mahalanobisDistance( const VECTORLIKE& mean_diffs, const MAT1& COV1, const MAT2& COV2, const MAT3& CROSS_COV12 ); template <class VECTORLIKE, class MATRIXLIKE> MATRIXLIKE::Scalar mrpt::math::mahalanobisDistance2( const VECTORLIKE& delta_mu, const MATRIXLIKE& cov ); template <class VECTORLIKE, class MATRIXLIKE> MATRIXLIKE::Scalar mrpt::math::mahalanobisDistance( const VECTORLIKE& delta_mu, const MATRIXLIKE& cov ); template <typename T> T mrpt::math::productIntegralTwoGaussians( const std::vector<T>& mean_diffs, const CMatrixDynamic<T>& COV1, const CMatrixDynamic<T>& COV2 ); template <typename T, size_t DIM> T mrpt::math::productIntegralTwoGaussians( const std::vector<T>& mean_diffs, const CMatrixFixed<T, DIM, DIM>& COV1, const CMatrixFixed<T, DIM, DIM>& COV2 ); template <typename T, class VECLIKE, class MATLIKE1, class MATLIKE2> void mrpt::math::productIntegralAndMahalanobisTwoGaussians( const VECLIKE& mean_diffs, const MATLIKE1& COV1, const MATLIKE2& COV2, T& maha2_out, T& intprod_out, const MATLIKE1* CROSS_COV12 = nullptr ); template <typename T, class VECLIKE, class MATRIXLIKE> void mrpt::math::mahalanobisDistance2AndLogPDF( const VECLIKE& diff_mean, const MATRIXLIKE& cov, T& maha2_out, T& log_pdf_out ); template <typename T, class VECLIKE, class MATRIXLIKE> void mrpt::math::mahalanobisDistance2AndPDF( const VECLIKE& diff_mean, const MATRIXLIKE& cov, T& maha2_out, T& pdf_out ); template < class VECTOR_OF_VECTORS, class MATRIXLIKE, class VECTORLIKE, class VECTORLIKE2, class VECTORLIKE3 > void mrpt::math::covariancesAndMeanWeighted( const VECTOR_OF_VECTORS& elements, MATRIXLIKE& covariances, VECTORLIKE& means, const VECTORLIKE2* weights_mean, const VECTORLIKE3* weights_cov, const bool* elem_do_wrap2pi = nullptr ); template <class VECTOR_OF_VECTORS, class MATRIXLIKE, class VECTORLIKE> void mrpt::math::covariancesAndMean( const VECTOR_OF_VECTORS& elements, MATRIXLIKE& covariances, VECTORLIKE& means, const bool* elem_do_wrap2pi = nullptr ); template <class VECTORLIKE1, class VECTORLIKE2> void mrpt::math::weightedHistogram( const VECTORLIKE1& values, const VECTORLIKE1& weights, float binWidth, VECTORLIKE2& out_binCenters, VECTORLIKE2& out_binValues ); template <class VECTORLIKE1, class VECTORLIKE2> void mrpt::math::weightedHistogramLog( const VECTORLIKE1& values, const VECTORLIKE1& log_weights, float binWidth, VECTORLIKE2& out_binCenters, VECTORLIKE2& out_binValues ); double mrpt::math::averageLogLikelihood(const CVectorDouble& logLikelihoods); double mrpt::math::averageWrap2Pi(const CVectorDouble& angles); double mrpt::math::averageLogLikelihood(const CVectorDouble& logWeights, const CVectorDouble& logLikelihoods); double mrpt::math::normalPDF(double x, double mu, double std); template <class VECTORLIKE1, class VECTORLIKE2, class MATRIXLIKE> MATRIXLIKE::Scalar mrpt::math::normalPDFInf( const VECTORLIKE1& x, const VECTORLIKE2& mu, const MATRIXLIKE& cov_inv, const bool scaled_pdf = false ); template <class VECTORLIKE1, class VECTORLIKE2, class MATRIXLIKE> MATRIXLIKE::Scalar mrpt::math::normalPDF( const VECTORLIKE1& x, const VECTORLIKE2& mu, const MATRIXLIKE& cov, const bool scaled_pdf = false ); template <typename VECTORLIKE, typename MATRIXLIKE> MATRIXLIKE::Scalar mrpt::math::normalPDF( const VECTORLIKE& d, const MATRIXLIKE& cov ); template < typename VECTORLIKE1, typename MATRIXLIKE1, typename VECTORLIKE2, typename MATRIXLIKE2 > double mrpt::math::KLD_Gaussians( const VECTORLIKE1& mu0, const MATRIXLIKE1& cov0, const VECTORLIKE2& mu1, const MATRIXLIKE2& cov1 ); double mrpt::math::normalQuantile(double p); double mrpt::math::normalCDF(double p); double mrpt::math::chi2inv(double P, unsigned int dim = 1); double mrpt::math::noncentralChi2CDF(unsigned int degreesOfFreedom, double noncentrality, double arg); double mrpt::math::chi2CDF(unsigned int degreesOfFreedom, double arg); double mrpt::math::chi2PDF(unsigned int degreesOfFreedom, double arg, double accuracy = 1e-7); std::pair<double, double> mrpt::math::noncentralChi2PDF_CDF( unsigned int degreesOfFreedom, double noncentrality, double arg, double eps = 1e-7 ); template <typename CONTAINER, typename T> void mrpt::math::confidenceIntervals( const CONTAINER& data, T& out_mean, T& out_lower_conf_interval, T& out_upper_conf_interval, const double confidenceInterval = 0.1, const size_t histogramNumBins = 1000 ); template <typename CONTAINER, typename T> void mrpt::math::confidenceIntervalsFromHistogram( const CONTAINER& histogramCoords, const CONTAINER& histogramNormalizedHits, T& out_lower_conf_interval, T& out_upper_conf_interval, const double confidenceInterval = 0.1 ); std::string mrpt::math::MATLAB_plotCovariance2D( const CMatrixFloat& cov22, const CVectorFloat& mean, float stdCount, const std::string& style = std::string("b"), size_t nEllipsePoints = 30 ); std::string mrpt::math::MATLAB_plotCovariance2D( const CMatrixDouble& cov22, const CVectorDouble& mean, float stdCount, const std::string& style = std::string("b"), size_t nEllipsePoints = 30 );
Global Functions
template <class VECTORLIKE1, class VECTORLIKE2, class MAT> MAT::Scalar mrpt::math::mahalanobisDistance2( const VECTORLIKE1& X, const VECTORLIKE2& MU, const MAT& COV )
Computes the squared mahalanobis distance of a vector X given the mean MU and the covariance inverse COV_inv.
template <class VECTORLIKE1, class VECTORLIKE2, class MAT> VECTORLIKE1::Scalar mrpt::math::mahalanobisDistance( const VECTORLIKE1& X, const VECTORLIKE2& MU, const MAT& COV )
Computes the mahalanobis distance of a vector X given the mean MU and the covariance inverse COV_inv.
template <class VECTORLIKE, class MAT1, class MAT2, class MAT3> MAT1::Scalar mrpt::math::mahalanobisDistance2( const VECTORLIKE& mean_diffs, const MAT1& COV1, const MAT2& COV2, const MAT3& CROSS_COV12 )
Computes the squared mahalanobis distance between two non-independent Gaussians, given the two covariance matrices and the vector with the difference of their means.
template <class VECTORLIKE, class MAT1, class MAT2, class MAT3> VECTORLIKE::Scalar mrpt::math::mahalanobisDistance( const VECTORLIKE& mean_diffs, const MAT1& COV1, const MAT2& COV2, const MAT3& CROSS_COV12 )
Computes the mahalanobis distance between two non-independent Gaussians (or independent if CROSS_COV12=nullptr), given the two covariance matrices and the vector with the difference of their means.
template <class VECTORLIKE, class MATRIXLIKE> MATRIXLIKE::Scalar mrpt::math::mahalanobisDistance2( const VECTORLIKE& delta_mu, const MATRIXLIKE& cov )
Computes the squared mahalanobis distance between a point and a Gaussian, given the covariance matrix and the vector with the difference between the mean and the point.
template <class VECTORLIKE, class MATRIXLIKE> MATRIXLIKE::Scalar mrpt::math::mahalanobisDistance( const VECTORLIKE& delta_mu, const MATRIXLIKE& cov )
Computes the mahalanobis distance between a point and a Gaussian, given the covariance matrix and the vector with the difference between the mean and the point.
template <typename T> T mrpt::math::productIntegralTwoGaussians( const std::vector<T>& mean_diffs, const CMatrixDynamic<T>& COV1, const CMatrixDynamic<T>& COV2 )
Computes the integral of the product of two Gaussians, with means separated by “mean_diffs” and covariances “COV1” and “COV2”.
template <typename T, size_t DIM> T mrpt::math::productIntegralTwoGaussians( const std::vector<T>& mean_diffs, const CMatrixFixed<T, DIM, DIM>& COV1, const CMatrixFixed<T, DIM, DIM>& COV2 )
Computes the integral of the product of two Gaussians, with means separated by “mean_diffs” and covariances “COV1” and “COV2”.
template <typename T, class VECLIKE, class MATLIKE1, class MATLIKE2> void mrpt::math::productIntegralAndMahalanobisTwoGaussians( const VECLIKE& mean_diffs, const MATLIKE1& COV1, const MATLIKE2& COV2, T& maha2_out, T& intprod_out, const MATLIKE1* CROSS_COV12 = nullptr )
Computes both, the integral of the product of two Gaussians and their square Mahalanobis distance.
See also:
productIntegralTwoGaussians, mahalanobisDistance2
template <typename T, class VECLIKE, class MATRIXLIKE> void mrpt::math::mahalanobisDistance2AndLogPDF( const VECLIKE& diff_mean, const MATRIXLIKE& cov, T& maha2_out, T& log_pdf_out )
Computes both, the logarithm of the PDF and the square Mahalanobis distance between a point (given by its difference wrt the mean) and a Gaussian.
See also:
productIntegralTwoGaussians, mahalanobisDistance2, normalPDF, mahalanobisDistance2AndPDF
template <typename T, class VECLIKE, class MATRIXLIKE> void mrpt::math::mahalanobisDistance2AndPDF( const VECLIKE& diff_mean, const MATRIXLIKE& cov, T& maha2_out, T& pdf_out )
Computes both, the PDF and the square Mahalanobis distance between a point (given by its difference wrt the mean) and a Gaussian.
See also:
productIntegralTwoGaussians, mahalanobisDistance2, normalPDF
template < class VECTOR_OF_VECTORS, class MATRIXLIKE, class VECTORLIKE, class VECTORLIKE2, class VECTORLIKE3 > void mrpt::math::covariancesAndMeanWeighted( const VECTOR_OF_VECTORS& elements, MATRIXLIKE& covariances, VECTORLIKE& means, const VECTORLIKE2* weights_mean, const VECTORLIKE3* weights_cov, const bool* elem_do_wrap2pi = nullptr )
Computes covariances and mean of any vector of containers, given optional weights for the different samples.
Parameters:
elements |
Any kind of vector of vectors/arrays, eg. std::vector<mrpt::math::CVectorDouble>, with all the input samples, each sample in a “row”. |
covariances |
Output estimated covariance; it can be a fixed/dynamic matrix or a matrixview. |
means |
Output estimated mean; it can be CVectorDouble/CVectorFixedDouble, etc… |
weights_mean |
If !=nullptr, it must point to a vector of size()==number of elements, with normalized weights to take into account for the mean. |
weights_cov |
If !=nullptr, it must point to a vector of size()==number of elements, with normalized weights to take into account for the covariance. |
elem_do_wrap2pi |
If !=nullptr; it must point to an array of “bool” of size()==dimension of each element, stating if it’s needed to do a wrap to [-pi,pi] to each dimension. |
See also:
This method is used in mrpt::math::unscented_transform_gaussian
template <class VECTOR_OF_VECTORS, class MATRIXLIKE, class VECTORLIKE> void mrpt::math::covariancesAndMean( const VECTOR_OF_VECTORS& elements, MATRIXLIKE& covariances, VECTORLIKE& means, const bool* elem_do_wrap2pi = nullptr )
Computes covariances and mean of any vector of containers.
Parameters:
elements |
Any kind of vector of vectors/arrays, eg. std::vector<mrpt::math::CVectorDouble>, with all the input samples, each sample in a “row”. |
covariances |
Output estimated covariance; it can be a fixed/dynamic matrix or a matrixview. |
means |
Output estimated mean; it can be CVectorDouble/CVectorFixedDouble, etc… |
elem_do_wrap2pi |
If !=nullptr; it must point to an array of “bool” of size()==dimension of each element, stating if it’s needed to do a wrap to [-pi,pi] to each dimension. |
template <class VECTORLIKE1, class VECTORLIKE2> void mrpt::math::weightedHistogram( const VECTORLIKE1& values, const VECTORLIKE1& weights, float binWidth, VECTORLIKE2& out_binCenters, VECTORLIKE2& out_binValues )
Computes the weighted histogram for a vector of values and their corresponding weights.
Parameters:
values |
[IN] The N values |
weights |
[IN] The weights for the corresponding N values (don’t need to be normalized) |
binWidth |
[IN] The desired width of the bins |
out_binCenters |
[OUT] The centers of the M bins generated to cover from the minimum to the maximum value of “values” with the given “binWidth” |
out_binValues |
[OUT] The ratio of values at each given bin, such as the whole vector sums up the unity. |
See also:
template <class VECTORLIKE1, class VECTORLIKE2> void mrpt::math::weightedHistogramLog( const VECTORLIKE1& values, const VECTORLIKE1& log_weights, float binWidth, VECTORLIKE2& out_binCenters, VECTORLIKE2& out_binValues )
Computes the weighted histogram for a vector of values and their corresponding log-weights.
Parameters:
values |
[IN] The N values |
weights |
[IN] The log-weights for the corresponding N values (don’t need to be normalized) |
binWidth |
[IN] The desired width of the bins |
out_binCenters |
[OUT] The centers of the M bins generated to cover from the minimum to the maximum value of “values” with the given “binWidth” |
out_binValues |
[OUT] The ratio of values at each given bin, such as the whole vector sums up the unity. |
See also:
double mrpt::math::averageLogLikelihood(const CVectorDouble& logLikelihoods)
A numerically-stable method to compute average likelihood values with strongly different ranges (unweighted likelihoods: compute the arithmetic mean).
This method implements this equation:
See also the tutorial page.
double mrpt::math::averageWrap2Pi(const CVectorDouble& angles)
Computes the average of a sequence of angles in radians taking into account the correct wrapping in the range \(]-\pi,\pi [\), for example, the mean of (2,-2) is \(\pi\), not 0.
double mrpt::math::averageLogLikelihood( const CVectorDouble& logWeights, const CVectorDouble& logLikelihoods )
A numerically-stable method to average likelihood values with strongly different ranges (weighted likelihoods).
This method implements this equation:
See also the tutorial page.
double mrpt::math::normalPDF(double x, double mu, double std)
Evaluates the univariate normal (Gaussian) distribution at a given point “x”.
template <class VECTORLIKE1, class VECTORLIKE2, class MATRIXLIKE> MATRIXLIKE::Scalar mrpt::math::normalPDFInf( const VECTORLIKE1& x, const VECTORLIKE2& mu, const MATRIXLIKE& cov_inv, const bool scaled_pdf = false )
Evaluates the multivariate normal (Gaussian) distribution at a given point “x”.
Parameters:
x |
A vector or column or row matrix with the point at which to evaluate the pdf. |
mu |
A vector or column or row matrix with the Gaussian mean. |
cov_inv |
The inverse covariance (information) matrix of the Gaussian. |
scaled_pdf |
If set to true, the PDF will be scaled to be in the range [0,1], in contrast to its integral from [-inf,+inf] being 1. |
template <class VECTORLIKE1, class VECTORLIKE2, class MATRIXLIKE> MATRIXLIKE::Scalar mrpt::math::normalPDF( const VECTORLIKE1& x, const VECTORLIKE2& mu, const MATRIXLIKE& cov, const bool scaled_pdf = false )
Evaluates the multivariate normal (Gaussian) distribution at a given point “x”.
Parameters:
x |
A vector or column or row matrix with the point at which to evaluate the pdf. |
mu |
A vector or column or row matrix with the Gaussian mean. |
cov |
The covariance matrix of the Gaussian. |
scaled_pdf |
If set to true, the PDF will be scaled to be in the range [0,1], in contrast to its integral from [-inf,+inf] being 1. |
template <typename VECTORLIKE, typename MATRIXLIKE> MATRIXLIKE::Scalar mrpt::math::normalPDF( const VECTORLIKE& d, const MATRIXLIKE& cov )
Evaluates the multivariate normal (Gaussian) distribution at a given point given its distance vector “d” from the Gaussian mean.
template < typename VECTORLIKE1, typename MATRIXLIKE1, typename VECTORLIKE2, typename MATRIXLIKE2 > double mrpt::math::KLD_Gaussians( const VECTORLIKE1& mu0, const MATRIXLIKE1& cov0, const VECTORLIKE2& mu1, const MATRIXLIKE2& cov1 )
Kullback-Leibler divergence (KLD) between two independent multivariate Gaussians.
\(D_\mathrm{KL}(\mathcal{N}_0 \| \mathcal{N}_1) = { 1 \over 2 } ( \log_e ( { \det \Sigma_1 \over \det \Sigma_0 } ) + \mathrm{tr} ( \Sigma_1^{-1} \Sigma_0 ) + ( \mu_1 - \mu_0 )^\top \Sigma_1^{-1} ( \mu_1 - \mu_0 ) - N )\)
double mrpt::math::normalQuantile(double p)
Evaluates the Gaussian distribution quantile for the probability value p=[0,1].
The employed approximation is that from Peter J. Acklam (pjacklam@online.no), freely available in http://home.online.no/~pjacklam.
double mrpt::math::normalCDF(double p)
Evaluates the Gaussian cumulative density function.
The employed approximation is that from W. J. Cody freely available in http://www.netlib.org/specfun/erf Equivalent to MATLAB normcdf(x,mu,s) with p=(x-mu)/s
double mrpt::math::chi2inv(double P, unsigned int dim = 1)
The “quantile” of the Chi-Square distribution, for dimension “dim” and probability 0<P<1 (the inverse of chi2CDF) An aproximation from the Wilson-Hilferty transformation is used.
Equivalent to MATLAB chi2inv(), but note that this is just an approximation, which becomes very poor for small values of “P”.
double mrpt::math::noncentralChi2CDF( unsigned int degreesOfFreedom, double noncentrality, double arg )
Cumulative non-central chi square distribution (approximate).
Computes approximate values of the cumulative density of a chi square distribution with degreesOfFreedom, and noncentrality parameter noncentrality at the given argument arg, i.e. the probability that a random number drawn from the distribution is below arg It uses the approximate transform into a normal distribution due to Wilson and Hilferty (see Abramovitz, Stegun: “Handbook of Mathematical Functions”, formula 26.3.32). The algorithm’s running time is independent of the inputs. The accuracy is only about 0.1 for few degrees of freedom, but reaches about 0.001 above dof = 5.
Function code from the Vigra project (http://hci.iwr.uni-heidelberg.de/vigra/); code under “MIT X11 License”, GNU GPL-compatible.
See also:
double mrpt::math::chi2CDF(unsigned int degreesOfFreedom, double arg)
Cumulative chi square distribution.
Computes the cumulative density of a chi square distribution with degreesOfFreedom and tolerance accuracy at the given argument arg, i.e. the probability that a random number drawn from the distribution is below arg by calling noncentralChi2CDF(degreesOfFreedom, 0.0, arg, accuracy)
.
Function code from the Vigra project (http://hci.iwr.uni-heidelberg.de/vigra/); code under “MIT X11 License”, GNU GPL-compatible.
double mrpt::math::chi2PDF( unsigned int degreesOfFreedom, double arg, double accuracy = 1e-7 )
Chi square distribution PDF.
Computes the density of a chi square distribution with degreesOfFreedom and tolerance accuracy at the given argument arg by calling noncentralChi2(degreesOfFreedom, 0.0, arg, accuracy)
. Function code from the Vigra project (http://hci.iwr.uni-heidelberg.de/vigra/); code under “MIT X11 License”, GNU GPL-compatible.
Equivalent to MATLAB’s chi2pdf(arg,degreesOfFreedom)
std::pair<double, double> mrpt::math::noncentralChi2PDF_CDF( unsigned int degreesOfFreedom, double noncentrality, double arg, double eps = 1e-7 )
Returns the ‘exact’ PDF (first) and CDF (second) of a Non-central chi-squared probability distribution, using an iterative method.
Equivalent to MATLAB’s ncx2cdf(arg,degreesOfFreedom,noncentrality)
template <typename CONTAINER, typename T> void mrpt::math::confidenceIntervals( const CONTAINER& data, T& out_mean, T& out_lower_conf_interval, T& out_upper_conf_interval, const double confidenceInterval = 0.1, const size_t histogramNumBins = 1000 )
Return the mean and the 10%-90% confidence points (or with any other confidence value) of a set of samples by building the cummulative CDF of all the elements of the container.
The container can be any MRPT container (CArray, matrices, vectors).
Parameters:
confidenceInterval |
A number in the range (0,1) such as the confidence interval will be [100*confidenceInterval, 100*(1-confidenceInterval)]. |
template <typename CONTAINER, typename T> void mrpt::math::confidenceIntervalsFromHistogram( const CONTAINER& histogramCoords, const CONTAINER& histogramNormalizedHits, T& out_lower_conf_interval, T& out_upper_conf_interval, const double confidenceInterval = 0.1 )
Return the mean and the 10%-90% confidence points (or with any other confidence value) of a set of samples from their histogram.
The container can be any MRPT container (CArray, matrices, vectors).
Parameters:
confidenceInterval |
A number in the range (0,1) such as the confidence interval will be [100*confidenceInterval, 100*(1-confidenceInterval)]. |
std::string mrpt::math::MATLAB_plotCovariance2D( const CMatrixFloat& cov22, const CVectorFloat& mean, float stdCount, const std::string& style = std::string("b"), size_t nEllipsePoints = 30 )
Generates a string with the MATLAB commands required to plot an confidence interval (ellipse) for a 2D Gaussian (‘float’ version).
Parameters:
cov22 |
The 2x2 covariance matrix |
mean |
The 2-length vector with the mean |
stdCount |
How many “quantiles” to get into the area of the ellipse: 2: 95%, 3:99.97%,… |
style |
A matlab style string, for colors, line styles,… |
nEllipsePoints |
The number of points in the ellipse to generate |
std::string mrpt::math::MATLAB_plotCovariance2D( const CMatrixDouble& cov22, const CVectorDouble& mean, float stdCount, const std::string& style = std::string("b"), size_t nEllipsePoints = 30 )
Generates a string with the MATLAB commands required to plot an confidence interval (ellipse) for a 2D Gaussian (‘double’ version).
Parameters:
cov22 |
The 2x2 covariance matrix |
mean |
The 2-length vector with the mean |
stdCount |
How many “quantiles” to get into the area of the ellipse: 2: 95%, 3:99.97%,… |
style |
A matlab style string, for colors, line styles,… |
nEllipsePoints |
The number of points in the ellipse to generate |