Main MRPT website > C++ reference for MRPT 1.9.9
distributions.h
Go to the documentation of this file.
1 /* +------------------------------------------------------------------------+
2  | Mobile Robot Programming Toolkit (MRPT) |
3  | http://www.mrpt.org/ |
4  | |
5  | Copyright (c) 2005-2017, Individual contributors, see AUTHORS file |
6  | See: http://www.mrpt.org/Authors - All rights reserved. |
7  | Released under BSD License. See details in http://www.mrpt.org/License |
8  +------------------------------------------------------------------------+ */
9 #ifndef mrpt_math_distributions_H
10 #define mrpt_math_distributions_H
11 
12 #include <mrpt/utils/utils_defs.h>
13 #include <mrpt/math/math_frwds.h>
15 
16 #include <mrpt/math/ops_matrices.h>
17 #include <mrpt/math/ops_vectors.h>
18 
19 /*---------------------------------------------------------------
20  Namespace
21  ---------------------------------------------------------------*/
22 namespace mrpt
23 {
24 namespace math
25 {
26 /** \addtogroup stats_grp Statistics functions, probability distributions
27  * \ingroup mrpt_base_grp
28  * @{ */
29 
30 /** Evaluates the univariate normal (Gaussian) distribution at a given point
31  * "x".
32  */
33 double normalPDF(double x, double mu, double std);
34 
35 /** Evaluates the multivariate normal (Gaussian) distribution at a given point
36  * "x".
37  * \param x A vector or column or row matrix with the point at which to
38  * evaluate the pdf.
39  * \param mu A vector or column or row matrix with the Gaussian mean.
40  * \param cov_inv The inverse covariance (information) matrix of the
41  * Gaussian.
42  * \param scaled_pdf If set to true, the PDF will be scaled to be in the
43  * range [0,1], in contrast to its integral from [-inf,+inf] being 1.
44  */
45 template <class VECTORLIKE1, class VECTORLIKE2, class MATRIXLIKE>
47  const VECTORLIKE1& x, const VECTORLIKE2& mu, const MATRIXLIKE& cov_inv,
48  const bool scaled_pdf = false)
49 {
51  typedef typename MATRIXLIKE::Scalar T;
52  ASSERTDEB_(cov_inv.isSquare())
53  ASSERTDEB_(
54  size_t(cov_inv.getColCount()) == size_t(x.size()) &&
55  size_t(cov_inv.getColCount()) == size_t(mu.size()))
56  T ret = ::exp(
57  static_cast<T>(-0.5) *
58  mrpt::math::multiply_HCHt_scalar((x - mu).eval(), cov_inv));
59  return scaled_pdf ? ret : ret * ::sqrt(
60  cov_inv.det() /
61  ::pow(
62  static_cast<T>(M_2PI),
63  static_cast<T>(size(cov_inv, 1))));
64  MRPT_END
65 }
66 
67 /** Evaluates the multivariate normal (Gaussian) distribution at a given point
68  * "x".
69  * \param x A vector or column or row matrix with the point at which to
70  * evaluate the pdf.
71  * \param mu A vector or column or row matrix with the Gaussian mean.
72  * \param cov The covariance matrix of the Gaussian.
73  * \param scaled_pdf If set to true, the PDF will be scaled to be in the
74  * range [0,1], in contrast to its integral from [-inf,+inf] being 1.
75  */
76 template <class VECTORLIKE1, class VECTORLIKE2, class MATRIXLIKE>
77 inline typename MATRIXLIKE::Scalar normalPDF(
78  const VECTORLIKE1& x, const VECTORLIKE2& mu, const MATRIXLIKE& cov,
79  const bool scaled_pdf = false)
80 {
81  return normalPDFInf(x, mu, cov.inverse(), scaled_pdf);
82 }
83 
84 /** Evaluates the multivariate normal (Gaussian) distribution at a given point
85  * given its distance vector "d" from the Gaussian mean.
86  */
87 template <typename VECTORLIKE, typename MATRIXLIKE>
89  const VECTORLIKE& d, const MATRIXLIKE& cov)
90 {
92  ASSERTDEB_(cov.isSquare())
93  ASSERTDEB_(size_t(cov.getColCount()) == size_t(d.size()))
94  return std::exp(
95  static_cast<typename MATRIXLIKE::Scalar>(-0.5) *
96  mrpt::math::multiply_HCHt_scalar(d, cov.inverse())) /
97  (::pow(
98  static_cast<typename MATRIXLIKE::Scalar>(M_2PI),
99  static_cast<typename MATRIXLIKE::Scalar>(
100  0.5 * cov.getColCount())) *
101  ::sqrt(cov.det()));
102  MRPT_END
103 }
104 
105 /** Kullback-Leibler divergence (KLD) between two independent multivariate
106  * Gaussians.
107  *
108  * \f$ D_\mathrm{KL}(\mathcal{N}_0 \| \mathcal{N}_1) = { 1 \over 2 } ( \log_e (
109  * { \det \Sigma_1 \over \det \Sigma_0 } ) + \mathrm{tr} ( \Sigma_1^{-1}
110  * \Sigma_0 ) + ( \mu_1 - \mu_0 )^\top \Sigma_1^{-1} ( \mu_1 - \mu_0 ) - N ) \f$
111  */
112 template <typename VECTORLIKE1, typename MATRIXLIKE1, typename VECTORLIKE2,
113  typename MATRIXLIKE2>
115  const VECTORLIKE1& mu0, const MATRIXLIKE1& cov0, const VECTORLIKE2& mu1,
116  const MATRIXLIKE2& cov1)
117 {
118  MRPT_START
119  ASSERT_(
120  size_t(mu0.size()) == size_t(mu1.size()) &&
121  size_t(mu0.size()) == size_t(size(cov0, 1)) &&
122  size_t(mu0.size()) == size_t(size(cov1, 1)) && cov0.isSquare() &&
123  cov1.isSquare())
124  const size_t N = mu0.size();
125  MATRIXLIKE2 cov1_inv;
126  cov1.inv(cov1_inv);
127  const VECTORLIKE1 mu_difs = mu0 - mu1;
128  return 0.5 * (log(cov1.det() / cov0.det()) + (cov1_inv * cov0).trace() +
129  multiply_HCHt_scalar(mu_difs, cov1_inv) - N);
130  MRPT_END
131 }
132 
133 /** Evaluates the Gaussian distribution quantile for the probability value
134  * p=[0,1].
135  * The employed approximation is that from Peter J. Acklam
136  * (pjacklam@online.no),
137  * freely available in http://home.online.no/~pjacklam.
138  */
139 double normalQuantile(double p);
140 
141 /** Evaluates the Gaussian cumulative density function.
142  * The employed approximation is that from W. J. Cody
143  * freely available in http://www.netlib.org/specfun/erf
144  * \note Equivalent to MATLAB normcdf(x,mu,s) with p=(x-mu)/s
145  */
146 double normalCDF(double p);
147 
148 /** The "quantile" of the Chi-Square distribution, for dimension "dim" and
149  * probability 0<P<1 (the inverse of chi2CDF)
150  * An aproximation from the Wilson-Hilferty transformation is used.
151  * \note Equivalent to MATLAB chi2inv(), but note that this is just an
152  * approximation, which becomes very poor for small values of "P".
153  */
154 double chi2inv(double P, unsigned int dim = 1);
155 
156 /*! Cumulative non-central chi square distribution (approximate).
157 
158  Computes approximate values of the cumulative density of a chi square
159 distribution with \a degreesOfFreedom,
160  and noncentrality parameter \a noncentrality at the given argument
161  \a arg, i.e. the probability that a random number drawn from the
162 distribution is below \a arg
163  It uses the approximate transform into a normal distribution due to Wilson
164 and Hilferty
165  (see Abramovitz, Stegun: "Handbook of Mathematical Functions", formula
166 26.3.32).
167  The algorithm's running time is independent of the inputs. The accuracy is
168 only
169  about 0.1 for few degrees of freedom, but reaches about 0.001 above dof = 5.
170 
171  \note Function code from the Vigra project
172 (http://hci.iwr.uni-heidelberg.de/vigra/); code under "MIT X11 License", GNU
173 GPL-compatible.
174 * \sa noncentralChi2PDF_CDF
175 */
176 double noncentralChi2CDF(
177  unsigned int degreesOfFreedom, double noncentrality, double arg);
178 
179 /*! Cumulative chi square distribution.
180 
181  Computes the cumulative density of a chi square distribution with \a
182  degreesOfFreedom
183  and tolerance \a accuracy at the given argument \a arg, i.e. the probability
184  that
185  a random number drawn from the distribution is below \a arg
186  by calling <tt>noncentralChi2CDF(degreesOfFreedom, 0.0, arg, accuracy)</tt>.
187 
188  \note Function code from the Vigra project
189  (http://hci.iwr.uni-heidelberg.de/vigra/); code under "MIT X11 License", GNU
190  GPL-compatible.
191 */
192 double chi2CDF(unsigned int degreesOfFreedom, double arg);
193 
194 /*! Chi square distribution PDF.
195  * Computes the density of a chi square distribution with \a degreesOfFreedom
196  * and tolerance \a accuracy at the given argument \a arg
197  * by calling <tt>noncentralChi2(degreesOfFreedom, 0.0, arg, accuracy)</tt>.
198  * \note Function code from the Vigra project
199  *(http://hci.iwr.uni-heidelberg.de/vigra/); code under "MIT X11 License", GNU
200  *GPL-compatible.
201  *
202  * \note Equivalent to MATLAB's chi2pdf(arg,degreesOfFreedom)
203  */
204 double chi2PDF(
205  unsigned int degreesOfFreedom, double arg, double accuracy = 1e-7);
206 
207 /** Returns the 'exact' PDF (first) and CDF (second) of a Non-central
208  * chi-squared probability distribution, using an iterative method.
209  * \note Equivalent to MATLAB's ncx2cdf(arg,degreesOfFreedom,noncentrality)
210  */
211 std::pair<double, double> noncentralChi2PDF_CDF(
212  unsigned int degreesOfFreedom, double noncentrality, double arg,
213  double eps = 1e-7);
214 
215 /** Return the mean and the 10%-90% confidence points (or with any other
216  * confidence value) of a set of samples by building the cummulative CDF of all
217  * the elements of the container.
218  * The container can be any MRPT container (CArray, matrices, vectors).
219  * \param confidenceInterval A number in the range (0,1) such as the confidence
220  * interval will be [100*confidenceInterval, 100*(1-confidenceInterval)].
221  */
222 template <typename CONTAINER>
224  const CONTAINER& data,
227  out_lower_conf_interval,
229  out_upper_conf_interval,
230  const double confidenceInterval = 0.1, const size_t histogramNumBins = 1000)
231 {
232  MRPT_START
233  ASSERT_(
234  data.size() != 0) // don't use .empty() here to allow using matrices
235  ASSERT_(confidenceInterval > 0 && confidenceInterval < 1)
236 
237  out_mean = mean(data);
239  minimum_maximum(data, x_min, x_max);
240 
241  const typename mrpt::math::ContainerType<CONTAINER>::element_t binWidth =
242  (x_max - x_min) / histogramNumBins;
243 
244  const std::vector<double> H =
245  mrpt::math::histogram(data, x_min, x_max, histogramNumBins);
246  std::vector<double> Hc;
247  cumsum(H, Hc); // CDF
248  Hc *= 1.0 / mrpt::math::maximum(Hc);
249 
251  std::lower_bound(Hc.begin(), Hc.end(), confidenceInterval);
252  ASSERT_(it_low != Hc.end())
254  std::upper_bound(Hc.begin(), Hc.end(), 1 - confidenceInterval);
255  ASSERT_(it_high != Hc.end())
256  const size_t idx_low = std::distance(Hc.begin(), it_low);
257  const size_t idx_high = std::distance(Hc.begin(), it_high);
258  out_lower_conf_interval = x_min + idx_low * binWidth;
259  out_upper_conf_interval = x_min + idx_high * binWidth;
260 
261  MRPT_END
262 }
263 
264 /** @} */
265 
266 } // End of MATH namespace
267 
268 } // End of namespace
269 
270 #endif
double chi2CDF(unsigned int degreesOfFreedom, double arg)
Definition: math.cpp:2043
size_t size(const MATRIXLIKE &m, const int dim)
Definition: bits.h:41
double normalCDF(double p)
Evaluates the Gaussian cumulative density function.
Definition: math.cpp:1132
Scalar * iterator
Definition: eigen_plugins.h:26
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 (a scalar) = H * C * H^t (with a vector H and a symmetric matrix C)
Definition: ops_matrices.h:69
STL namespace.
double KLD_Gaussians(const VECTORLIKE1 &mu0, const MATRIXLIKE1 &cov0, const VECTORLIKE2 &mu1, const MATRIXLIKE2 &cov1)
Kullback-Leibler divergence (KLD) between two independent multivariate Gaussians. ...
#define M_2PI
Definition: mrpt_macros.h:437
double normalQuantile(double p)
Evaluates the Gaussian distribution quantile for the probability value p=[0,1].
Definition: math.cpp:1079
#define MRPT_END
const double eps
CONTAINER::Scalar maximum(const CONTAINER &v)
void minimum_maximum(const std::vector< T > &V, T &curMin, T &curMax)
Return the maximum and minimum values of a std::vector.
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".
Definition: distributions.h:46
void cumsum(const CONTAINER1 &in_data, CONTAINER2 &out_cumsum)
#define MRPT_START
This is the global namespace for all Mobile Robot Programming Toolkit (MRPT) libraries.
#define ASSERTDEB_(f)
Defines an assertion mechanism - only when compiled in debug.
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...
Definition: math.cpp:998
std::pair< double, double > noncentralChi2PDF_CDF(unsigned int degreesOfFreedom, double noncentrality, double arg, double eps=1e-7)
Returns the &#39;exact&#39; PDF (first) and CDF (second) of a Non-central chi-squared probability distributio...
Definition: math.cpp:2066
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.
#define ASSERT_(f)
CONTAINER::value_type element_t
Definition: math_frwds.h:103
Eigen::Matrix< typename MATRIX::Scalar, MATRIX::ColsAtCompileTime, MATRIX::ColsAtCompileTime > cov(const MATRIX &v)
Computes the covariance matrix from a list of samples in an NxM matrix, where each row is a sample...
Definition: ops_matrices.h:148
void confidenceIntervals(const CONTAINER &data, typename mrpt::math::ContainerType< CONTAINER >::element_t &out_mean, typename mrpt::math::ContainerType< CONTAINER >::element_t &out_lower_conf_interval, typename mrpt::math::ContainerType< CONTAINER >::element_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...
GLenum GLint x
Definition: glext.h:3538
double normalPDF(double x, double mu, double std)
Evaluates the univariate normal (Gaussian) distribution at a given point "x".
Definition: math.cpp:989
GLsizei GLsizei GLenum GLenum const GLvoid * data
Definition: glext.h:3546
double Scalar
Definition: KmUtils.h:44
GLfloat GLfloat p
Definition: glext.h:6305
double chi2PDF(unsigned int degreesOfFreedom, double arg, double accuracy=1e-7)
Definition: math.cpp:2133
double distance(const TPoint2D &p1, const TPoint2D &p2)
Gets the distance between two points in a 2D space.
Definition: geometry.cpp:1885
double noncentralChi2CDF(unsigned int degreesOfFreedom, double noncentrality, double arg)
Definition: math.cpp:2141



Page generated by Doxygen 1.8.14 for MRPT 1.9.9 Git: ae4571287 Thu Nov 23 00:06:53 2017 +0100 at dom oct 27 23:51:55 CET 2019