37 template <
class VECTORLIKE1,
class VECTORLIKE2,
class MATRIXLIKE>
39 const VECTORLIKE1& x,
const VECTORLIKE2& mu,
const MATRIXLIKE& cov_inv,
40 const bool scaled_pdf =
false)
46 size_t(cov_inv.cols()) ==
size_t(x.size()) &&
47 size_t(cov_inv.cols()) ==
size_t(mu.size()));
49 static_cast<T>(-0.5) *
54 cov_inv.det() / ::pow(
55 static_cast<T>(
M_2PI),
56 static_cast<T>(cov_inv.rows())));
69 template <
class VECTORLIKE1,
class VECTORLIKE2,
class MATRIXLIKE>
71 const VECTORLIKE1& x,
const VECTORLIKE2& mu,
const MATRIXLIKE&
cov,
72 const bool scaled_pdf =
false)
80 template <
typename VECTORLIKE,
typename MATRIXLIKE>
82 const VECTORLIKE& d,
const MATRIXLIKE&
cov)
88 static_cast<typename MATRIXLIKE::Scalar>(-0.5) *
91 static_cast<typename MATRIXLIKE::Scalar>(
M_2PI),
92 static_cast<typename MATRIXLIKE::Scalar>(0.5 *
cov.
cols())) *
105 typename VECTORLIKE1,
typename MATRIXLIKE1,
typename VECTORLIKE2,
106 typename MATRIXLIKE2>
108 const VECTORLIKE1& mu0,
const MATRIXLIKE1& cov0,
const VECTORLIKE2& mu1,
109 const MATRIXLIKE2& cov1)
113 size_t(mu0.size()) ==
size_t(mu1.size()) &&
114 size_t(mu0.size()) ==
size_t(cov0.rows()) &&
115 size_t(mu0.size()) ==
size_t(cov1.cols()) && cov0.isSquare() &&
117 const size_t N = mu0.size();
118 MATRIXLIKE2 cov1_inv;
119 cov1.inverse_LLt(cov1_inv);
120 const VECTORLIKE1 mu_difs = mu0 - mu1;
121 return 0.5 * (log(cov1.det() / cov0.det()) + (cov1_inv * cov0).trace() +
147 double chi2inv(
double P,
unsigned int dim = 1);
170 unsigned int degreesOfFreedom,
double noncentrality,
double arg);
185 double chi2CDF(
unsigned int degreesOfFreedom,
double arg);
198 unsigned int degreesOfFreedom,
double arg,
double accuracy = 1e-7);
205 unsigned int degreesOfFreedom,
double noncentrality,
double arg,
215 template <
typename CONTAINER,
typename T>
217 const CONTAINER&
data, T& out_mean, T& out_lower_conf_interval,
218 T& out_upper_conf_interval,
const double confidenceInterval = 0.1,
219 const size_t histogramNumBins = 1000)
224 ASSERT_(confidenceInterval > 0 && confidenceInterval < 1);
227 const auto x_min =
data.minCoeff();
228 const auto x_max =
data.maxCoeff();
229 const auto binWidth = (x_max - x_min) / histogramNumBins;
231 const std::vector<double> H =
233 std::vector<double> Hc;
237 auto it_low = std::lower_bound(Hc.begin(), Hc.end(), confidenceInterval);
240 std::upper_bound(Hc.begin(), Hc.end(), 1 - confidenceInterval);
244 out_lower_conf_interval = x_min + idx_low * binWidth;
245 out_upper_conf_interval = x_min + idx_high * binWidth;
MAT_C::Scalar multiply_HtCH_scalar(const VECTOR_H &H, const MAT_C &C)
r (scalar) = H^t*C*H (H: column vector, C: symmetric matrix)
double normalCDF(double p)
Evaluates the Gaussian cumulative density function.
This file implements miscelaneous matrix and matrix/vector operations, and internal functions in mrpt...
MAT_C::Scalar multiply_HCHt_scalar(const VECTOR_H &H, const MAT_C &C)
r (scalar) = H*C*H^t (H: row vector, C: symmetric matrix)
bool isSquare() const
returns true if matrix is NxN
double KLD_Gaussians(const VECTORLIKE1 &mu0, const MATRIXLIKE1 &cov0, const VECTORLIKE2 &mu1, const MATRIXLIKE2 &cov1)
Kullback-Leibler divergence (KLD) between two independent multivariate Gaussians. ...
double normalQuantile(double p)
Evaluates the Gaussian distribution quantile for the probability value p=[0,1].
#define ASSERT_(f)
Defines an assertion mechanism.
This base provides a set of functions for maths stuff.
std::pair< double, double > 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 distributio...
Derived inverse() const
Returns the inverse of a general matrix using LU.
void 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 sa...
double chi2PDF(unsigned int degreesOfFreedom, double arg, double accuracy=1e-7)
CONTAINER::Scalar maximum(const CONTAINER &v)
Scalar det() const
Determinant of matrix.
Derived inverse_LLt() const
Returns the inverse of a symmetric matrix using LLt.
CMatrixDouble cov(const MATRIX &v)
Computes the covariance matrix from a list of samples in an NxM matrix, where each row is a sample...
size_type cols() const
Number of columns in the matrix.
MATRIXLIKE::Scalar 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".
void cumsum(const CONTAINER1 &in_data, CONTAINER2 &out_cumsum)
double chi2CDF(unsigned int degreesOfFreedom, double arg)
double chi2inv(double P, unsigned int dim=1)
The "quantile" of the Chi-Square distribution, for dimension "dim" and probability 0<P<1 (the inverse...
#define ASSERTDEB_(f)
Defines an assertion mechanism - only when compiled in debug.
std::vector< double > histogram(const CONTAINER &v, double limit_min, double limit_max, size_t number_bins, bool do_normalization=false, std::vector< double > *out_bin_centers=nullptr)
Computes the normalized or normal histogram of a sequence of numbers given the number of bins and the...
double mean(const CONTAINER &v)
Computes the mean value of a vector.
double noncentralChi2CDF(unsigned int degreesOfFreedom, double noncentrality, double arg)
double normalPDF(double x, double mu, double std)
Evaluates the univariate normal (Gaussian) distribution at a given point "x".
double distance(const TPoint2D &p1, const TPoint2D &p2)
Gets the distance between two points in a 2D space.
static struct FontData data