19 #include <Eigen/Dense> 35 return ::exp(-0.5 *
square((x - mu) /
std)) /
36 (
std * 2.506628274631000502415765284811);
48 return dim * pow(1.0 - 2.0 / (9 * dim) +
60 for (
unsigned int i = 2; i <= n; i++) ret *= i;
72 for (
unsigned int i = 2; i <= n; i++) retLog += ::log((
double)n);
84 static const double a[6] = {-3.969683028665376e+01, 2.209460984245205e+02,
85 -2.759285104469687e+02, 1.383577518672690e+02,
86 -3.066479806614716e+01, 2.506628277459239e+00};
87 static const double b[5] = {-5.447609879822406e+01, 1.615858368580409e+02,
88 -1.556989798598866e+02, 6.680131188771972e+01,
89 -1.328068155288572e+01};
90 static const double c[6] = {-7.784894002430293e-03, -3.223964580411365e-01,
91 -2.400758277161838e+00, -2.549732539343734e+00,
92 4.374664141464968e+00, 2.938163982698783e+00};
93 static const double d[4] = {7.784695709041462e-03, 3.224671290700398e-01,
94 2.445134137142996e+00, 3.754408661907416e+00};
107 (((((a[0] * t + a[1]) * t + a[2]) * t + a[3]) * t + a[4]) * t +
109 (((((b[0] * t + b[1]) * t + b[2]) * t + b[3]) * t + b[4]) * t + 1);
114 t = sqrt(-2 * ::log(q));
115 u = (((((c[0] * t + c[1]) * t + c[2]) * t + c[3]) * t + c[4]) * t +
117 ((((d[0] * t + d[1]) * t + d[2]) * t + d[3]) * t + 1);
124 t = t * 2.506628274631000502415765284811 *
126 u = u - t / (1 + u * t / 2);
128 return (p > 0.5 ? -u : u);
136 static const double a[5] = {1.161110663653770e-002, 3.951404679838207e-001,
137 2.846603853776254e+001, 1.887426188426510e+002,
138 3.209377589138469e+003};
139 static const double b[5] = {1.767766952966369e-001, 8.344316438579620e+000,
140 1.725514762600375e+002, 1.813893686502485e+003,
141 8.044716608901563e+003};
142 static const double c[9] = {
143 2.15311535474403846e-8, 5.64188496988670089e-1, 8.88314979438837594e00,
144 6.61191906371416295e01, 2.98635138197400131e02, 8.81952221241769090e02,
145 1.71204761263407058e03, 2.05107837782607147e03, 1.23033935479799725E03};
146 static const double d[9] = {
147 1.00000000000000000e00, 1.57449261107098347e01, 1.17693950891312499e02,
148 5.37181101862009858e02, 1.62138957456669019e03, 3.29079923573345963e03,
149 4.36261909014324716e03, 3.43936767414372164e03, 1.23033935480374942e03};
150 static const double p[6] = {1.63153871373020978e-2, 3.05326634961232344e-1,
151 3.60344899949804439e-1, 1.25781726111229246e-1,
152 1.60837851487422766e-2, 6.58749161529837803e-4};
153 static const double q[6] = {1.00000000000000000e00, 2.56852019228982242e00,
154 1.87295284992346047e00, 5.27905102951428412e-1,
155 6.05183413124413191e-2, 2.33520497626869185e-3};
162 if (y <= 0.46875 * 1.4142135623730950488016887242097)
166 y = u * ((((a[0] * z + a[1]) * z + a[2]) * z + a[3]) * z + a[4]) /
167 ((((b[0] * z + b[1]) * z + b[2]) * z + b[3]) * z + b[4]);
171 z = ::exp(-y * y / 2) / 2;
175 y = y / 1.4142135623730950488016887242097;
176 y = ((((((((c[0] * y + c[1]) * y + c[2]) * y + c[3]) * y + c[4]) * y +
185 / ((((((((d[0] * y + d[1]) * y + d[2]) * y + d[3]) * y + d[4]) * y +
199 z = z * 1.4142135623730950488016887242097 / y;
202 (((((p[0] * y + p[1]) * y + p[2]) * y + p[3]) * y + p[4]) * y +
204 (((((q[0] * y + q[1]) * y + q[2]) * y + q[3]) * y + q[4]) * y +
206 y = z * (0.564189583547756286948 - y);
208 return (u < 0.0 ? y : 1 - y);
217 if (!std::getline(f, str))
return false;
219 const char* s = str.c_str();
221 char *nextTok, *context;
222 const char* delim =
" \t";
226 while (nextTok !=
nullptr)
228 d.push_back(atoi(nextTok));
241 if (!std::getline(f, str))
return false;
243 const char* s = str.c_str();
245 char *nextTok, *context;
246 const char* delim =
" \t";
250 while (nextTok !=
nullptr)
252 d.push_back(atof(nextTok));
270 if (!logWeights.
size())
276 double SUM1 = 0, SUM2 = 0, tmpVal;
278 for (itLW = logWeights.
begin(), itLL = logLikelihoods.
begin();
279 itLW != logWeights.
end(); itLW++, itLL++)
281 tmpVal = *itLW - lw_max;
282 SUM1 += std::exp(tmpVal);
283 SUM2 += std::exp(tmpVal + *itLL - ll_max);
286 double res = -std::log(SUM1) + std::log(SUM2) + ll_max;
299 size_t N = logLikelihoods.
size();
305 for (
size_t i = 0; i < N; i++) SUM1 += exp(logLikelihoods[i] - ll_max);
307 double res = log(SUM1) - log(static_cast<double>(N)) + ll_max;
317 if (angles.
empty())
return 0;
319 int W_phi_R = 0, W_phi_L = 0;
320 double phi_R = 0, phi_L = 0;
326 double phi = angles[i];
327 if (std::abs(phi) > 1.5707963267948966192313216916398)
330 if (phi < 0) phi = (
M_2PI + phi);
346 if (W_phi_L) phi_L /=
static_cast<double>(W_phi_L);
347 if (W_phi_R) phi_R /=
static_cast<double>(W_phi_R);
353 return ((phi_L * W_phi_L + phi_R * W_phi_R) / (W_phi_L + W_phi_R));
358 const string& style,
size_t nEllipsePoints)
367 cov2, mean2, stdCount, style, nEllipsePoints);
374 const string& style,
size_t nEllipsePoints)
384 const auto N = nEllipsePoints;
386 const double dAng = (
M_2PI / (nEllipsePoints - 1));
387 std::vector<double> X(N), Y(N), Cos(N), Sin(N);
390 for (std::size_t idx = 0; idx < N; idx++, ang += dAng)
392 Cos[idx] = std::cos(ang);
393 Sin[idx] = std::sin(ang);
396 std::vector<double> eVals;
398 cov2x2.
eig(eigVec, eVals);
408 for (std::size_t idx = 0; idx < N; idx++)
410 X[idx] =
mean[0] + stdCount * (Cos[idx] * M(0, 0) + Sin[idx] * M(1, 0));
411 X[idx] =
mean[1] + stdCount * (Cos[idx] * M(0, 1) + Sin[idx] * M(1, 1));
417 str += string(
"plot([ ");
418 for (
size_t i = 0; i < X.size(); i++)
420 str +=
format(
"%.4f", X[i]);
421 if (i != (X.size() - 1)) str +=
",";
423 str += string(
"],[ ");
424 for (
size_t i = 0; i < X.size(); i++)
426 str +=
format(
"%.4f", Y[i]);
427 if (i != (Y.size() - 1)) str +=
",";
430 str +=
format(
"],'%s');\n", style.c_str());
439 const double x,
const double x0,
const double y0,
const double x1,
440 const double y1,
bool wrap2pi)
446 const double Ax = x1 - x0;
447 const double Ay = y1 - y0;
449 double r = y0 + Ay * (x - x0) / Ax;
463 const std::vector<double>& inV, std::vector<double>& outV,
int _winSize,
467 ASSERT_((
int)inV.size() >= _winSize);
469 size_t winSize = _winSize;
474 size_t sz = inV.size();
477 std::vector<double> aux(winSize);
478 size_t mpoint = winSize / 2;
479 for (
size_t k = 0; k < sz; ++k)
483 size_t idx_to_start =
484 std::max(
size_t(0), k - mpoint);
486 std::min(std::min(winSize, sz + mpoint - k), k + mpoint + 1);
488 aux.resize(n_elements);
489 for (
size_t m = idx_to_start, n = 0; m < idx_to_start + n_elements;
493 std::sort(aux.begin(), aux.end());
495 size_t auxSz = aux.size();
496 size_t auxMPoint = auxSz / 2;
497 outV[k] = (auxSz % 2) ? (aux[auxMPoint])
498 : (0.5 * (aux[auxMPoint - 1] +
513 T arg, T& lans, T& dans, T& pans,
unsigned int& j)
518 lans = lans + std::log(arg / j);
519 dans = std::exp(lans);
523 dans = dans * arg / j;
530 unsigned int degreesOfFreedom,
double noncentrality,
double arg,
double eps)
533 noncentrality >= 0.0 && arg >= 0.0 &&
eps > 0.0,
534 "noncentralChi2PDF_CDF(): parameters must be positive.");
536 if (arg == 0.0 && degreesOfFreedom > 0)
return std::make_pair(0.0, 0.0);
539 double b1 = 0.5 * noncentrality, ao = std::exp(-
b1), eps2 =
eps / ao,
540 lnrtpi2 = 0.22579135264473, probability, density, lans, dans, pans,
542 unsigned int maxit = 500, i, m;
543 if (degreesOfFreedom % 2)
546 lans = -0.5 * (arg + std::log(arg)) - lnrtpi2;
547 dans = std::exp(lans);
548 pans = std::erf(std::sqrt(arg / 2.0));
554 dans = std::exp(lans);
559 if (degreesOfFreedom == 0)
562 degreesOfFreedom = 2;
564 sum = 1.0 / ao - 1.0 - am;
566 probability = 1.0 + am * pans;
571 degreesOfFreedom = degreesOfFreedom - 1;
573 sum = 1.0 / ao - 1.0;
574 while (i < degreesOfFreedom)
576 degreesOfFreedom = degreesOfFreedom + 1;
581 for (++m; m < maxit; ++m)
586 density = density + am * dans;
588 probability = probability + hold;
589 if ((pans *
sum < eps2) && (hold < eps2))
break;
591 if (m == maxit)
THROW_EXCEPTION(
"noncentralChi2PDF_CDF(): no convergence.");
592 return std::make_pair(
593 0.5 * ao * density, std::min(1.0, std::max(0.0, ao * probability)));
597 unsigned int degreesOfFreedom,
double arg,
double accuracy)
600 degreesOfFreedom, 0.0, arg, accuracy)
605 unsigned int degreesOfFreedom,
double noncentrality,
double arg)
607 const double a = degreesOfFreedom + noncentrality;
608 const double b = (a + noncentrality) /
square(a);
610 (std::pow((
double)arg / a, 1.0 / 3.0) - (1.0 - 2.0 / 9.0 * b)) /
611 std::sqrt(2.0 / 9.0 * b);
612 return 0.5 * (1.0 + std::erf(t / std::sqrt(2.0)));
bool eig(Derived &eVecs, std::vector< Scalar > &eVals, bool sorted=true) const
Computes the eigenvectors and eigenvalues for a square, general matrix.
double averageWrap2Pi(const CVectorDouble &angles)
Computes the average of a sequence of angles in radians taking into account the correct wrapping in t...
A compile-time fixed-size numeric matrix container.
#define THROW_EXCEPTION(msg)
std::string std::string format(std::string_view fmt, ARGS &&... args)
CArchive & WriteAs(const TYPE_FROM_ACTUAL &value)
double normalCDF(double p)
Evaluates the Gaussian cumulative density function.
void WriteBufferFixEndianness(const T *ptr, size_t ElementCount)
Writes a sequence of elemental datatypes, taking care of reordering their bytes from the running arch...
size_type size() const
Get a 2-vector with [NROWS NCOLS] (as in MATLAB command size(x))
This file implements miscelaneous matrix and matrix/vector operations, and internal functions in mrpt...
#define MRPT_END_WITH_CLEAN_UP(stuff)
mrpt::serialization::CArchive & operator>>(mrpt::serialization::CArchive &in, CMatrixD::Ptr &pObj)
double factorial(unsigned int n)
Computes the factorial of an integer number and returns it as a double value (internally it uses loga...
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...
CVectorDynamic< double > CVectorDouble
#define MRPT_CHECK_NORMAL_NUMBER(v)
Throws an exception if the number is NaN, IND, or +/-INF, or return the same number otherwise...
double averageLogLikelihood(const CVectorDouble &logLikelihoods)
A numerically-stable method to compute average likelihood values with strongly different ranges (unwe...
STORED_TYPE ReadAs()
De-serialize a variable and returns it by value.
bool loadVector(std::istream &f, std::vector< int > &d)
Loads one row of a text file as a numerical std::vector.
double chi2PDF(unsigned int degreesOfFreedom, double arg, double accuracy=1e-7)
char * strtok(char *str, const char *strDelimit, char **context) noexcept
An OS-independent method for tokenizing a string.
CONTAINER::Scalar sum(const CONTAINER &v)
Computes the sum of all the elements.
uint64_t factorial64(unsigned int n)
Computes the factorial of an integer number and returns it as a 64-bit integer number.
CONTAINER::Scalar maximum(const CONTAINER &v)
#define ASSERTMSG_(f, __ERROR_MSG)
Defines an assertion mechanism.
typename vec_t::const_iterator const_iterator
CMatrixDouble cov(const MATRIX &v)
Computes the covariance matrix from a list of samples in an NxM matrix, where each row is a sample...
T wrapToPi(T a)
Modifies the given angle to translate it into the ]-pi,pi] range.
size_type rows() const
Number of rows in the matrix.
size_type cols() const
Number of columns in the matrix.
double chi2CDF(unsigned int degreesOfFreedom, double arg)
void noncentralChi2OneIteration(T arg, T &lans, T &dans, T &pans, unsigned int &j)
return_t square(const num_t x)
Inline function for the square of a number.
This is the global namespace for all Mobile Robot Programming Toolkit (MRPT) libraries.
Virtual base class for "archives": classes abstracting I/O streams.
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...
mrpt::serialization::CArchive & operator<<(mrpt::serialization::CArchive &s, const CVectorFloat &a)
void medianFilter(const std::vector< double > &inV, std::vector< double > &outV, int winSize, int numberOfSigmas=2)
double mean(const CONTAINER &v)
Computes the mean value of a vector.
EIGEN_MAP asEigen()
Get as an Eigen-compatible Eigen::Map object.
double noncentralChi2CDF(unsigned int degreesOfFreedom, double noncentrality, double arg)
std::string 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 2...
double interpolate2points(const double x, const double x0, const double y0, const double x1, const double y1, bool wrap2pi=false)
Linear interpolation/extrapolation: evaluates at "x" the line (x0,y0)-(x1,y1).
void resize(std::size_t N, bool zeroNewElements=false)
#define THROW_EXCEPTION_FMT(_FORMAT_STRING,...)
size_t ReadBufferFixEndianness(T *ptr, size_t ElementCount)
Reads a sequence of elemental datatypes, taking care of reordering their bytes from the MRPT stream s...
double normalPDF(double x, double mu, double std)
Evaluates the univariate normal (Gaussian) distribution at a given point "x".
EIGEN_MAP asEigen()
Get as an Eigen-compatible Eigen::Map object.
This template class provides the basic functionality for a general 2D any-size, resizable container o...
void setDiagonal(const std::size_t N, const T value)
Resize to NxN, set all entries to zero, except the main diagonal which is set to value ...
CMatrixDynamic< double > CMatrixDouble
Declares a matrix of double numbers (non serializable).
#define MRPT_UNUSED_PARAM(a)
Determines whether this is an X86 or AMD64 platform.