22 #include <Eigen/Dense> 55 for (
const auto& m : m_modes)
57 const double w = exp(m.log_w);
58 se_averager.
append(m.mean, w);
70 const size_t N = m_modes.size();
75 this->getMean(estMean2D);
86 for (
const auto& m : m_modes)
89 sumW += w = exp(m.log_w);
92 estMean_i -= estMeanMat;
101 if (sumW != 0)
cov *= (1.0 / sumW);
103 return {
cov, estMean2D};
109 uint32_t N = m_modes.size();
112 for (
const auto m : m_modes)
114 out << m.log_w << m.mean;
130 for (
auto& m : m_modes)
135 if (version == 0) m.log_w = log(max(1e-300, m.log_w));
161 if (
this == &o)
return;
165 m_modes =
dynamic_cast<const CPosePDFSOG*
>(&o)->m_modes;
171 m_modes[0].log_w = 0;
185 if (!f)
return false;
187 for (
const auto& m : m_modes)
189 f,
"%e %e %e %e %e %e %e %e %e %e\n", exp(m.log_w), m.mean.x(),
190 m.mean.y(), m.mean.phi(), m.cov(0, 0), m.cov(1, 1), m.cov(2, 2),
191 m.cov(0, 1), m.cov(0, 2), m.cov(1, 2));
216 for (
auto& m : m_modes)
219 m.mean.composeFrom(newReferenceBase, m.mean);
227 enforceCovSymmetry();
236 rot(0, 0) = rot(1, 1) = cos(ang);
237 rot(0, 1) = -sin(ang);
238 rot(1, 0) = sin(ang);
241 for (
auto& m : m_modes)
259 [[maybe_unused]]
size_t N,
260 [[maybe_unused]] std::vector<CVectorDouble>& outSamples)
const 274 [[maybe_unused]]
const double minMahalanobisDistToDrop)
283 const auto* p1 =
dynamic_cast<const CPosePDFSOG*
>(&p1_);
298 -0.5 * (3 * log(
M_2PI) - log(covInv.
det()) +
299 (eta.transpose() * p2->cov.asEigen() * eta.asEigen())(0, 0));
301 this->m_modes.clear();
302 for (
const auto& m : p1->m_modes)
304 auxSOG_Kernel_i.
mean = m.mean;
318 newKernel.
mean = auxGaussianProduct.
mean;
319 newKernel.
cov = auxGaussianProduct.
cov;
324 eta_i = covInv_i * eta_i;
329 new_eta_i = new_covInv_i * new_eta_i;
332 -0.5 * (3 * log(
M_2PI) - log(new_covInv_i.
det()) +
333 (eta_i.transpose() * auxSOG_Kernel_i.
cov.
asEigen() *
334 eta_i.asEigen())(0, 0));
336 -0.5 * (3 * log(
M_2PI) - log(new_covInv_i.
det()) +
337 (new_eta_i.transpose() * newKernel.
cov.
asEigen() *
338 new_eta_i.asEigen())(0, 0));
341 newKernel.
log_w = m.log_w + a + a_i - new_a_i;
344 this->m_modes.push_back(newKernel);
363 out->m_modes.resize(m_modes.size());
365 for (itSrc = m_modes.begin(), itDest =
out->m_modes.begin();
366 itSrc != m_modes.end(); ++itSrc, ++itDest)
369 (itDest)->
mean = -(itSrc)->mean;
372 (itDest)->
cov = (itSrc)->cov;
381 for (
auto& m : m_modes) m.mean = m.mean + Ap;
383 this->rotateAllCovariances(Ap.
phi());
398 for (
const auto& m : m_modes)
415 for (
const auto& m : m_modes)
417 MU(0, 0) = m.mean.x();
418 MU(1, 0) = m.mean.y();
420 COV(0, 0) = m.cov(0, 0);
421 COV(1, 1) = m.cov(1, 1);
422 COV(0, 1) = COV(1, 0) = m.cov(0, 1);
440 for (
const auto& m : m_modes)
457 for (
auto& m : m_modes)
459 m.cov(0, 1) = m.cov(1, 0);
460 m.cov(0, 2) = m.cov(2, 0);
461 m.cov(1, 2) = m.cov(2, 1);
472 if (!m_modes.size())
return;
474 double maxW = m_modes[0].log_w;
475 for (
auto& m : m_modes) maxW = max(maxW, m.log_w);
477 for (
auto& m : m_modes) m.log_w -= maxW;
486 double x_min,
double x_max,
double y_min,
double y_max,
double resolutionXY,
495 const auto Nx = (size_t)ceil((x_max - x_min) / resolutionXY);
496 const auto Ny = (size_t)ceil((y_max - y_min) / resolutionXY);
500 for (
size_t i = 0; i < Ny; i++)
502 double y = y_min + i * resolutionXY;
503 for (
size_t j = 0; j < Nx; j++)
505 double x = x_min + j * resolutionXY;
506 outMatrix(i, j) = evaluatePDF(
CPose2D(x, y, phi), sumOverAllPhis);
522 size_t N = m_modes.size();
530 for (
size_t i = 0; i < (N - 1);)
537 for (
size_t j = 0; j < N; j++) sumW += exp(m_modes[j].log_w);
540 const double Wi = exp(m_modes[i].log_w) / sumW;
542 double min_Bij = std::numeric_limits<double>::max();
550 for (
size_t j = 0; j < N; j++)
553 const double Wj = exp(m_modes[j].log_w) / sumW;
554 const double Wij_ = 1.0 / (Wi + Wj);
558 Pij.
asEigen() += m_modes[j].cov.asEigen() * Wj * Wij_;
568 AUX *= Wi * Wj * Wij_ * Wij_;
571 double Bij = (Wi + Wj) * log(Pij.det()) -
572 Wi * log(m_modes[i].
cov.
det()) -
573 Wj * log(m_modes[j].
cov.
det());
576 cout <<
"try merge[" << i <<
", " << j
577 <<
"] -> Bij: " << Bij << endl;
581 cout <<
"Pij: " << Pij << endl
582 <<
" Pi: " << m_modes[i].cov << endl
583 <<
" Pj: " << m_modes[j].cov << endl;
596 cout <<
"merge[" << i <<
", " << best_j
597 <<
"] Tempting merge: KLd = " << min_Bij;
599 if (min_Bij < max_KLd)
601 if (
verbose) cout <<
" Accepted." << endl;
612 const double Wj = exp(Mj.
log_w) / sumW;
613 const double Wij_ = 1.0 / (Wi + Wj);
614 const double Wi_ = Wi * Wij_;
615 const double Wj_ = Wj * Wij_;
623 Mij.
cov = min_Bij_COV;
627 m_modes.erase(m_modes.begin() + best_j);
631 if (
verbose) cout <<
" Nope." << endl;
648 auto it_best =
end();
649 double best_log_w = -std::numeric_limits<double>::max();
651 for (
auto i =
begin(); i !=
end(); ++i)
653 if (i->log_w > best_log_w)
655 best_log_w = i->log_w;
660 if (it_best !=
end())
662 mean_point = it_best->mean;
Computes weighted and un-weighted averages of SE(2) poses.
void normalizeWeights()
Normalize the weights in m_modes such as the maximum log-weight is 0.
void clear()
Clear the list of modes.
void inverse(CPosePDF &o) const override
Returns a new PDF such as: NEW_PDF = (0,0,0) - THIS_PDF.
void mergeModes(double max_KLd=0.5, bool verbose=false)
Merge very close modes so the overall number of modes is reduced while preserving the total distribut...
void serializeSymmetricMatrixTo(MAT &m, mrpt::serialization::CArchive &out)
Binary serialization of symmetric matrices, saving the space of duplicated values.
CPose2D mean
The mean value.
uint8_t serializeGetVersion() const override
Must return the current versioning number of the object.
The struct for each mode:
#define THROW_EXCEPTION(msg)
Declares a class that represents a Probability Density function (PDF) of a 2D pose ...
int void fclose(FILE *f)
An OS-independent version of fclose.
void rotateAllCovariances(double ang)
Rotate all the covariance matrixes by replacing them by , where .
CListGaussianModes::const_iterator const_iterator
#define IMPLEMENTS_SERIALIZABLE(class_name, base, NameSpace)
To be added to all CSerializable-classes implementation files.
void serializeTo(mrpt::serialization::CArchive &out) const override
Pure virtual method for writing (serializing) to an abstract archive.
double evaluatePDF(const mrpt::poses::CPose2D &x, bool sumOverAllPhis=false) const
Evaluates the PDF at a given point.
void bayesianFusion(const CPosePDF &p1, const CPosePDF &p2, const double minMahalanobisDistToDrop=0) override
Bayesian fusion of two points gauss.
void getMean(CPose2D &mean_pose) const override
void bayesianFusion(const CPosePDF &p1, const CPosePDF &p2, const double minMahalanobisDistToDrop=0) override
Bayesian fusion of two pose distributions, then save the result in this object (WARNING: Currently p1...
void enforceCovSymmetry()
Ensures the symmetry of the covariance matrix (eventually certain operations in the math-coprocessor ...
void getCovariance(mrpt::math::CMatrixDouble &cov) const
Returns the estimate of the covariance matrix (STATE_LEN x STATE_LEN covariance matrix) ...
virtual void getMean(type_value &mean_point) const =0
Returns the mean, or mathematical expectation of the probability density distribution (PDF)...
void append(const mrpt::poses::CPose2D &p)
Adds a new pose to the computation.
mrpt::math::CMatrixDouble33 cov
The 3x3 covariance matrix.
void matProductOf_AAt(const MAT_A &A)
this = A * AT
std::tuple< cov_mat_t, type_value > getCovarianceAndMean() const override
Returns an estimate of the pose covariance matrix (STATE_LENxSTATE_LEN cov matrix) and the mean...
#define MRPT_THROW_UNKNOWN_SERIALIZATION_VERSION(__V)
For use in CSerializable implementations.
#define ASSERT_(f)
Defines an assertion mechanism.
CListGaussianModes::iterator iterator
CMatrixFixed< double, 3, 3 > CMatrixDouble33
void wrapToPiInPlace(T &a)
Modifies the given angle to translate it into the ]-pi,pi] range.
This base provides a set of functions for maths stuff.
#define CLASS_ID(T)
Access to runtime class ID for a defined class name.
CMatrixFixed< double, 3, 1 > CMatrixDouble31
Declares a class that represents a Probability Density function (PDF) of a 2D pose ...
double phi() const
Get the phi angle of the 2D pose (in radians)
virtual const mrpt::rtti::TRuntimeClassId * GetRuntimeClass() const override
Returns information about the class of an object in runtime.
auto block(int start_row, int start_col)
non-const block(): Returns an Eigen::Block reference to the block
bool saveToTextFile(const std::string &file) const override
Save the density to a text file, with the following format: There is one row per Gaussian "mode"...
mrpt::math::CMatrixDouble33 cov
void deserializeSymmetricMatrixFrom(MAT &m, mrpt::serialization::CArchive &in)
Binary serialization of symmetric matrices, saving the space of duplicated values.
Scalar det() const
Determinant of matrix.
double x() const
Common members of all points & poses classes.
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...
Declares a class that represents a probability density function (pdf) of a 2D pose (x...
Classes for 2D/3D geometry representation, both of single values and probability density distribution...
double evaluateNormalizedPDF(const mrpt::poses::CPose2D &x) const
Evaluates the ratio PDF(x) / max_PDF(x*), that is, the normalized PDF in the range [0...
int fprintf(FILE *fil, const char *format,...) noexcept MRPT_printf_format_check(2
An OS-independent version of fprintf.
const_iterator end() const
This is the global namespace for all Mobile Robot Programming Toolkit (MRPT) libraries.
const_iterator begin() const
Virtual base class for "archives": classes abstracting I/O streams.
void serializeFrom(mrpt::serialization::CArchive &in, uint8_t serial_version) override
Pure virtual method for reading (deserializing) from an abstract archive.
void getMostLikelyCovarianceAndMean(mrpt::math::CMatrixDouble33 &cov, CPose2D &mean_point) const
For the most likely Gaussian mode in the SOG, returns the pose covariance matrix (3x3 cov matrix) and...
A class used to store a 2D pose, including the 2D coordinate point and a heading (phi) angle...
A class used to store a 3D pose (a 3D translation + a rotation in 3D).
mrpt::vision::TStereoCalibResults out
This file implements matrix/vector text and binary serialization.
void drawManySamples(size_t N, std::vector< mrpt::math::CVectorDouble > &outSamples) const override
Draws a number of samples from the distribution, and saves as a list of 1x3 vectors, where each row contains a (x,y,phi) datum.
void getHomogeneousMatrix(mrpt::math::CMatrixDouble44 &out_HM) const
Returns the corresponding 4x4 homogeneous transformation matrix for the point(translation) or pose (t...
void setSize(size_t row, size_t col, bool zeroNewElements=false)
Changes the size of matrix, maintaining the previous contents.
void changeCoordinatesReference(const CPose3D &newReferenceBase) override
this = p (+) this.
double mean(const CONTAINER &v)
Computes the mean value of a vector.
EIGEN_MAP asEigen()
Get as an Eigen-compatible Eigen::Map object.
void copyFrom(const CPosePDF &o) override
Copy operator, translating if necesary (for example, between particles and gaussian representations) ...
FILE * fopen(const char *fileName, const char *mode) noexcept
An OS-independent version of fopen.
void operator+=(const mrpt::poses::CPose2D &Ap)
Makes: thisPDF = thisPDF + Ap, where "+" is pose composition (both the mean, and the covariance matri...
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.
double log_w
The log-weight.
void evaluatePDFInArea(double x_min, double x_max, double y_min, double y_max, double resolutionXY, double phi, mrpt::math::CMatrixDouble &outMatrix, bool sumOverAllPhis=false)
Evaluates the PDF within a rectangular grid (and a fixed orientation) and saves the result in a matri...
CMatrixFixed< double, ROWS, COLS > cast_double() const
void resize(const size_t N)
Resize the number of SOG modes.
void drawSingleSample(CPose2D &outPart) const override
Draws a single sample from the distribution.
CMatrixDynamic< double > CMatrixDouble
Declares a matrix of double numbers (non serializable).
void multiply_HCHt(const MAT_H &H, const MAT_C &C, MAT_R &R, bool accumResultInOutput=false)
R = H * C * H^t.
void get_average(mrpt::poses::CPose2D &out_mean) const
Returns the average pose.