00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028 #ifndef MRPT_MATH_H
00029 #define MRPT_MATH_H
00030
00031 #include <mrpt/utils/utils_defs.h>
00032 #include <mrpt/math/CMatrixTemplateNumeric.h>
00033 #include <mrpt/math/CMatrixFixedNumeric.h>
00034 #include <mrpt/math/CVectorTemplate.h>
00035 #include <mrpt/math/CHistogram.h>
00036
00037 #include <mrpt/math/ops_vectors.h>
00038 #include <mrpt/math/ops_matrices.h>
00039
00040 #include <numeric>
00041 #include <cmath>
00042
00043
00044
00045
00046 namespace mrpt
00047 {
00048
00049
00050 namespace math
00051 {
00052 using namespace mrpt::utils;
00053
00054
00055
00056
00057
00058 bool BASE_IMPEXP loadVector( utils::CFileStream &f, std::vector<int> &d);
00059
00060
00061
00062
00063
00064 bool BASE_IMPEXP loadVector( utils::CFileStream &f, std::vector<double> &d);
00065
00066
00067
00068 bool BASE_IMPEXP isNaN(float f) MRPT_NO_THROWS;
00069
00070
00071 bool BASE_IMPEXP isNaN(double f) MRPT_NO_THROWS;
00072
00073
00074 bool BASE_IMPEXP isFinite(float f) MRPT_NO_THROWS;
00075
00076
00077 bool BASE_IMPEXP isFinite(double f) MRPT_NO_THROWS;
00078
00079 #ifdef HAVE_LONG_DOUBLE
00080
00081 bool BASE_IMPEXP isNaN(long double f) MRPT_NO_THROWS;
00082
00083
00084 bool BASE_IMPEXP isFinite(long double f) MRPT_NO_THROWS;
00085 #endif
00086
00087
00088
00089 template<typename T,typename K>
00090 void linspace(T first,T last, size_t count, std::vector<K> &out_vector)
00091 {
00092 if (count<2)
00093 {
00094 out_vector.assign(1,last);
00095 return;
00096 }
00097 else
00098 {
00099 out_vector.resize(count);
00100 const T incr = (last-first)/T(count-1);
00101 T c = first;
00102 for (size_t i=0;i<count;i++,c+=incr)
00103 out_vector[i] = K(c);
00104 }
00105 }
00106
00107
00108
00109 template<class T>
00110 inline std::vector<T> linspace(T first,T last, size_t count)
00111 {
00112 std::vector<T> ret;
00113 mrpt::math::linspace(first,last,count,ret);
00114 return ret;
00115 }
00116
00117
00118 template<class T,T STEP>
00119 inline std::vector<T> sequence(T first,size_t length)
00120 {
00121 std::vector<T> ret(length);
00122 if (!length) return ret;
00123 size_t i=0;
00124 while (length--) { ret[i++]=first; first+=STEP; }
00125 return ret;
00126 }
00127
00128
00129 template<class T>
00130 std::vector<T> ones(size_t count)
00131 {
00132 return std::vector<T>(count,1);
00133 }
00134
00135
00136 template<class T>
00137 std::vector<T> zeros(size_t count)
00138 {
00139 return std::vector<T>(count,0);
00140 }
00141
00142
00143
00144
00145
00146
00147 template <class T>
00148 inline void wrapTo2PiInPlace(T &a)
00149 {
00150 bool was_neg = a<0;
00151 a = fmod(a, static_cast<T>(M_2PI) );
00152 if (was_neg) a+=static_cast<T>(M_2PI);
00153 }
00154
00155
00156
00157
00158
00159 template <class T>
00160 inline T wrapTo2Pi(T a)
00161 {
00162 wrapTo2PiInPlace(a);
00163 return a;
00164 }
00165
00166
00167
00168
00169
00170 template <class T>
00171 inline T wrapToPi(T a)
00172 {
00173 return wrapTo2Pi( a + static_cast<T>(M_PI) )-static_cast<T>(M_PI);
00174 }
00175
00176
00177
00178
00179
00180 template <class T>
00181 inline void wrapToPiInPlace(T &a)
00182 {
00183 a = wrapToPi(a);
00184 }
00185
00186
00187
00188
00189
00190 template<class T>
00191 void normalize(const std::vector<T> &v, std::vector<T> &out_v)
00192 {
00193 T total=0;
00194 typename std::vector<T>::const_iterator it;
00195 for (it=v.begin(); it!=v.end();it++)
00196 total += square(*it);
00197 total = ::sqrt(total);
00198 if (total)
00199 {
00200 out_v.resize(v.size());
00201 typename std::vector<T>::iterator q;
00202 for (it=v.begin(),q=out_v.begin(); q!=out_v.end();it++,q++)
00203 *q = *it / total;
00204 }
00205 else out_v.assign(v.size(),0);
00206 }
00207
00208
00209
00210
00211
00212
00213
00214 template<class VECTOR_OF_VECTOR, class VECTORLIKE, class MATRIXLIKE>
00215 void meanAndCovVector(
00216 const VECTOR_OF_VECTOR &v,
00217 VECTORLIKE &out_mean,
00218 MATRIXLIKE &out_cov
00219 )
00220 {
00221 const size_t N = v.size();
00222 ASSERTMSG_(N>0,"The input vector contains no elements");
00223 const double N_inv = 1.0/N;
00224
00225 const size_t M = v[0].size();
00226 ASSERTMSG_(M>0,"The input vector contains rows of length 0");
00227
00228
00229 out_mean.assign(M,0);
00230 for (size_t i=0;i<N;i++)
00231 for (size_t j=0;j<M;j++)
00232 out_mean[j]+=v[i][j];
00233 out_mean*=N_inv;
00234
00235
00236
00237
00238 out_cov.zeros(M,M);
00239 for (size_t i=0;i<N;i++)
00240 {
00241 for (size_t j=0;j<M;j++)
00242 out_cov.get_unsafe(j,j)+=square(v[i][j]-out_mean[j]);
00243
00244 for (size_t j=0;j<M;j++)
00245 for (size_t k=j+1;k<M;k++)
00246 out_cov.get_unsafe(j,k)+=(v[i][j]-out_mean[j])*(v[i][k]-out_mean[k]);
00247 }
00248 for (size_t j=0;j<M;j++)
00249 for (size_t k=j+1;k<M;k++)
00250 out_cov.get_unsafe(k,j) = out_cov.get_unsafe(j,k);
00251 out_cov*=N_inv;
00252 }
00253
00254
00255
00256
00257
00258
00259 template<class VECTOR_OF_VECTOR>
00260 inline CMatrixDouble covVector( const VECTOR_OF_VECTOR &v )
00261 {
00262 vector_double m;
00263 CMatrixDouble C;
00264 meanAndCovVector(v,m,C);
00265 return C;
00266 }
00267
00268
00269
00270
00271
00272
00273
00274
00275
00276
00277 template<class VECTOR_OF_VECTORS, class MATRIXLIKE,class VECTORLIKE,class VECTORLIKE2,class VECTORLIKE3>
00278 inline void covariancesAndMeanWeighted(
00279 const VECTOR_OF_VECTORS &elements,
00280 MATRIXLIKE &covariances,
00281 VECTORLIKE &means,
00282 const VECTORLIKE2 *weights_mean,
00283 const VECTORLIKE3 *weights_cov,
00284 const bool *elem_do_wrap2pi = NULL
00285 )
00286 {
00287 ASSERTMSG_(elements.size()!=0,"No samples provided, so there is no way to deduce the output size.")
00288 typedef typename VECTORLIKE::value_type T;
00289 const size_t DIM = elements[0].size();
00290 means.resize(DIM);
00291 covariances.setSize(DIM,DIM);
00292 const size_t nElms=elements.size();
00293 const T NORM=1.0/nElms;
00294 if (weights_mean) { ASSERTDEB_(weights_mean->size()==nElms) }
00295
00296 for (size_t i=0;i<DIM;i++)
00297 {
00298 T accum = 0;
00299 if (!elem_do_wrap2pi || !elem_do_wrap2pi[i])
00300 {
00301 if (weights_mean)
00302 {
00303 for (size_t j=0;j<nElms;j++)
00304 accum+= (*weights_mean)[j] * elements[j][i];
00305 }
00306 else
00307 {
00308 for (size_t j=0;j<nElms;j++) accum+=elements[j][i];
00309 accum*=NORM;
00310 }
00311 }
00312 else
00313 {
00314 double accum_L=0,accum_R=0;
00315 double Waccum_L=0,Waccum_R=0;
00316 for (size_t j=0;j<nElms;j++)
00317 {
00318 double ang = elements[j][i];
00319 const double w = weights_mean!=NULL ? (*weights_mean)[j] : NORM;
00320 if (fabs( ang )>0.5*M_PI)
00321 {
00322 if (ang<0) ang = (M_2PI + ang);
00323 accum_L += ang * w;
00324 Waccum_L += w;
00325 }
00326 else
00327 {
00328 accum_R += ang * w;
00329 Waccum_R += w;
00330 }
00331 }
00332 if (Waccum_L>0) accum_L /= Waccum_L;
00333 if (Waccum_R>0) accum_R /= Waccum_R;
00334 if (accum_L>M_PI) accum_L -= M_2PI;
00335 accum = (accum_L* Waccum_L + accum_R * Waccum_R );
00336 }
00337 means[i]=accum;
00338 }
00339
00340 for (size_t i=0;i<DIM;i++)
00341 for (size_t j=0;j<=i;j++)
00342 {
00343 typename MATRIXLIKE::value_type elem=0;
00344 if (weights_cov)
00345 {
00346 ASSERTDEB_(weights_cov->size()==nElms)
00347 for (size_t k=0;k<nElms;k++)
00348 {
00349 const T Ai = (elements[k][i]-means[i]);
00350 const T Aj = (elements[k][j]-means[j]);
00351 if (!elem_do_wrap2pi || !elem_do_wrap2pi[i])
00352 elem+= (*weights_cov)[k] * Ai * Aj;
00353 else elem+= (*weights_cov)[k] * mrpt::math::wrapToPi(Ai) * mrpt::math::wrapToPi(Aj);
00354 }
00355 }
00356 else
00357 {
00358 for (size_t k=0;k<nElms;k++)
00359 {
00360 const T Ai = (elements[k][i]-means[i]);
00361 const T Aj = (elements[k][j]-means[j]);
00362 if (!elem_do_wrap2pi || !elem_do_wrap2pi[i])
00363 elem+= Ai * Aj;
00364 else elem+= mrpt::math::wrapToPi(Ai) * mrpt::math::wrapToPi(Aj);
00365 }
00366 elem*=NORM;
00367 }
00368 covariances.get_unsafe(i,j) = elem;
00369 if (i!=j) covariances.get_unsafe(j,i)=elem;
00370 }
00371 }
00372
00373
00374
00375
00376
00377
00378
00379 template<class VECTOR_OF_VECTORS, class MATRIXLIKE,class VECTORLIKE>
00380 void covariancesAndMean(const VECTOR_OF_VECTORS &elements,MATRIXLIKE &covariances,VECTORLIKE &means, const bool *elem_do_wrap2pi = NULL)
00381 {
00382 covariancesAndMeanWeighted<VECTOR_OF_VECTORS,MATRIXLIKE,VECTORLIKE,vector_double,vector_double>(elements,covariances,means,NULL,NULL,elem_do_wrap2pi);
00383 }
00384
00385
00386
00387
00388
00389
00390
00391
00392
00393
00394 template<class VECTORLIKE1,class VECTORLIKE2>
00395 void weightedHistogram(
00396 const VECTORLIKE1 &values,
00397 const VECTORLIKE1 &weights,
00398 float binWidth,
00399 VECTORLIKE2 &out_binCenters,
00400 VECTORLIKE2 &out_binValues )
00401 {
00402 MRPT_START;
00403 typedef typename VECTORLIKE1::value_type TNum;
00404
00405 ASSERT_( values.size() == weights.size() );
00406 ASSERT_( binWidth > 0 );
00407 TNum minBin = minimum( values );
00408 unsigned int nBins = static_cast<unsigned>(ceil((maximum( values )-minBin) / binWidth));
00409
00410
00411 out_binCenters.resize(nBins);
00412 out_binValues.clear(); out_binValues.resize(nBins,0);
00413 TNum halfBin = TNum(0.5)*binWidth;;
00414 VECTORLIKE2 binBorders(nBins+1,minBin-halfBin);
00415 for (unsigned int i=0;i<nBins;i++)
00416 {
00417 binBorders[i+1] = binBorders[i]+binWidth;
00418 out_binCenters[i] = binBorders[i]+halfBin;
00419 }
00420
00421
00422 TNum totalSum = 0;
00423 typename VECTORLIKE1::const_iterator itVal, itW;
00424 for (itVal = values.begin(), itW = weights.begin(); itVal!=values.end(); ++itVal, ++itW )
00425 {
00426 int idx = round(((*itVal)-minBin)/binWidth);
00427 if (idx>=int(nBins)) idx=nBins-1;
00428 ASSERTDEB_(idx>=0);
00429 out_binValues[idx] += *itW;
00430 totalSum+= *itW;
00431 }
00432
00433 if (totalSum)
00434 out_binValues /= totalSum;
00435
00436
00437 MRPT_END;
00438 }
00439
00440
00441
00442
00443
00444
00445
00446
00447
00448 template<class VECTORLIKE1,class VECTORLIKE2>
00449 void weightedHistogramLog(
00450 const VECTORLIKE1 &values,
00451 const VECTORLIKE1 &log_weights,
00452 float binWidth,
00453 VECTORLIKE2 &out_binCenters,
00454 VECTORLIKE2 &out_binValues )
00455 {
00456 MRPT_START;
00457 typedef typename VECTORLIKE1::value_type TNum;
00458
00459 ASSERT_( values.size() == log_weights.size() );
00460 ASSERT_( binWidth > 0 );
00461 TNum minBin = minimum( values );
00462 unsigned int nBins = static_cast<unsigned>(ceil((maximum( values )-minBin) / binWidth));
00463
00464
00465 out_binCenters.resize(nBins);
00466 out_binValues.clear(); out_binValues.resize(nBins,0);
00467 TNum halfBin = TNum(0.5)*binWidth;;
00468 VECTORLIKE2 binBorders(nBins+1,minBin-halfBin);
00469 for (unsigned int i=0;i<nBins;i++)
00470 {
00471 binBorders[i+1] = binBorders[i]+binWidth;
00472 out_binCenters[i] = binBorders[i]+halfBin;
00473 }
00474
00475
00476 const TNum max_log_weight = maximum(log_weights);
00477 TNum totalSum = 0;
00478 typename VECTORLIKE1::const_iterator itVal, itW;
00479 for (itVal = values.begin(), itW = log_weights.begin(); itVal!=values.end(); ++itVal, ++itW )
00480 {
00481 int idx = round(((*itVal)-minBin)/binWidth);
00482 if (idx>=int(nBins)) idx=nBins-1;
00483 ASSERTDEB_(idx>=0);
00484 const TNum w = exp(*itW-max_log_weight);
00485 out_binValues[idx] += w;
00486 totalSum+= w;
00487 }
00488
00489 if (totalSum)
00490 out_binValues /= totalSum;
00491
00492 MRPT_END;
00493 }
00494
00495
00496
00497
00498
00499
00500
00501
00502 template <class VECTOR_OF_VECTORS, class VECTORLIKE>
00503 inline void extractColumnFromVectorOfVectors(const size_t colIndex, const VECTOR_OF_VECTORS &data, VECTORLIKE &out_column)
00504 {
00505 const size_t N = data.size();
00506 out_column.resize(N);
00507 for (size_t i=0;i<N;i++)
00508 out_column[i]=data[i][colIndex];
00509 }
00510
00511
00512
00513 uint64_t BASE_IMPEXP factorial64(unsigned int n);
00514
00515
00516
00517 double BASE_IMPEXP factorial(unsigned int n);
00518
00519
00520
00521 template <class T>
00522 T round2up(T val)
00523 {
00524 T n = 1;
00525 while (n < val)
00526 {
00527 n <<= 1;
00528 if (n<=1)
00529 THROW_EXCEPTION("Overflow!");
00530 }
00531 return n;
00532 }
00533
00534
00535
00536
00537 template <class T>
00538 T round_10power(T val, int power10)
00539 {
00540 long double F = ::pow((long double)10.0,-(long double)power10);
00541 long int t = round_long( val * F );
00542 return T(t/F);
00543 }
00544
00545
00546
00547
00548 template<class T>
00549 double correlate_matrix(const CMatrixTemplateNumeric<T> &a1, const CMatrixTemplateNumeric<T> &a2)
00550 {
00551 if ((a1.getColCount()!=a2.getColCount())|(a1.getRowCount()!=a2.getRowCount()))
00552 THROW_EXCEPTION("Correlation Error!, images with no same size");
00553
00554 int i,j;
00555 T x1,x2;
00556 T syy=0, sxy=0, sxx=0, m1=0, m2=0 ,n=a1.getRowCount()*a2.getColCount();
00557
00558
00559 for (i=0;i<a1.getRowCount();i++)
00560 {
00561 for (j=0;j<a1.getColCount();j++)
00562 {
00563 m1 += a1(i,j);
00564 m2 += a2(i,j);
00565 }
00566 }
00567 m1 /= n;
00568 m2 /= n;
00569
00570 for (i=0;i<a1.getRowCount();i++)
00571 {
00572 for (j=0;j<a1.getColCount();j++)
00573 {
00574 x1 = a1(i,j) - m1;
00575 x2 = a2(i,j) - m2;
00576 sxx += x1*x1;
00577 syy += x2*x2;
00578 sxy += x1*x2;
00579 }
00580 }
00581
00582 return sxy / sqrt(sxx * syy);
00583 }
00584
00585
00586
00587
00588
00589
00590
00591
00592
00593
00594
00595
00596
00597
00598 template<class T>
00599 void BASE_IMPEXP qr_decomposition(
00600 CMatrixTemplateNumeric<T> &A,
00601 CMatrixTemplateNumeric<T> &R,
00602 CMatrixTemplateNumeric<T> &Q,
00603 CVectorTemplate<T> &c,
00604 int &sing);
00605
00606
00607
00608
00609
00610 template<class T>
00611 void BASE_IMPEXP UpdateCholesky(
00612 CMatrixTemplateNumeric<T> &chol,
00613 CVectorTemplate<T> &r1Modification);
00614
00615
00616
00617
00618
00619
00620
00621 void BASE_IMPEXP computeEigenValues2x2(
00622 const CMatrixFloat &in_matrix,
00623 float &min_eigenvalue,
00624 float &max_eigenvalue );
00625
00626
00627
00628
00629
00630
00631
00632
00633 double BASE_IMPEXP averageLogLikelihood( const vector_double &logLikelihoods );
00634
00635
00636
00637 double BASE_IMPEXP averageWrap2Pi(const vector_double &angles );
00638
00639
00640
00641
00642
00643
00644
00645
00646 double BASE_IMPEXP averageLogLikelihood(
00647 const vector_double &logWeights,
00648 const vector_double &logLikelihoods );
00649
00650
00651
00652
00653
00654
00655
00656
00657 std::string BASE_IMPEXP MATLAB_plotCovariance2D(
00658 const CMatrixFloat &cov22,
00659 const CVectorFloat &mean,
00660 const float &stdCount,
00661 const std::string &style = std::string("b"),
00662 const size_t &nEllipsePoints = 30 );
00663
00664
00665
00666
00667
00668
00669
00670
00671 std::string BASE_IMPEXP MATLAB_plotCovariance2D(
00672 const CMatrixDouble &cov22,
00673 const CVectorDouble &mean,
00674 const float &stdCount,
00675 const std::string &style = std::string("b"),
00676 const size_t &nEllipsePoints = 30 );
00677
00678
00679
00680
00681
00682
00683
00684 template <class MATRIXLIKE1,class MATRIXLIKE2> RET_VOID_ASSERT_MRPTMATRICES(MATRIXLIKE1,MATRIXLIKE2)
00685 homogeneousMatrixInverse(const MATRIXLIKE1 &M, MATRIXLIKE2 &out_inverse_M)
00686 {
00687 MRPT_START;
00688 ASSERT_( M.IsSquare() && size(M,1)==4);
00689
00690
00691
00692
00693
00694
00695
00696
00697
00698
00699
00700
00701 out_inverse_M.setSize(4,4);
00702
00703
00704 out_inverse_M.set_unsafe(0,0, M.get_unsafe(0,0));
00705 out_inverse_M.set_unsafe(0,1, M.get_unsafe(1,0));
00706 out_inverse_M.set_unsafe(0,2, M.get_unsafe(2,0));
00707
00708 out_inverse_M.set_unsafe(1,0, M.get_unsafe(0,1));
00709 out_inverse_M.set_unsafe(1,1, M.get_unsafe(1,1));
00710 out_inverse_M.set_unsafe(1,2, M.get_unsafe(2,1));
00711
00712 out_inverse_M.set_unsafe(2,0, M.get_unsafe(0,2));
00713 out_inverse_M.set_unsafe(2,1, M.get_unsafe(1,2));
00714 out_inverse_M.set_unsafe(2,2, M.get_unsafe(2,2));
00715
00716 const double tx = -M.get_unsafe(0,3);
00717 const double ty = -M.get_unsafe(1,3);
00718 const double tz = -M.get_unsafe(2,3);
00719
00720 const double tx_ = tx*M.get_unsafe(0,0)+ty*M.get_unsafe(1,0)+tz*M.get_unsafe(2,0);
00721 const double ty_ = tx*M.get_unsafe(0,1)+ty*M.get_unsafe(1,1)+tz*M.get_unsafe(2,1);
00722 const double tz_ = tx*M.get_unsafe(0,2)+ty*M.get_unsafe(1,2)+tz*M.get_unsafe(2,2);
00723
00724 out_inverse_M.set_unsafe(0,3, tx_ );
00725 out_inverse_M.set_unsafe(1,3, ty_ );
00726 out_inverse_M.set_unsafe(2,3, tz_ );
00727
00728 out_inverse_M.set_unsafe(3,0, 0);
00729 out_inverse_M.set_unsafe(3,1, 0);
00730 out_inverse_M.set_unsafe(3,2, 0);
00731 out_inverse_M.set_unsafe(3,3, 1);
00732
00733 MRPT_END;
00734 }
00735
00736
00737
00738
00739
00740
00741
00742
00743
00744 template <class VECTORLIKE,class VECTORLIKE2, class VECTORLIKE3, class MATRIXLIKE, class USERPARAM >
00745 void estimateJacobian(
00746 const VECTORLIKE &x,
00747 void (*functor) (const VECTORLIKE &x,const USERPARAM &y, VECTORLIKE3 &out),
00748 const VECTORLIKE2 &increments,
00749 const USERPARAM &userParam,
00750 MATRIXLIKE &out_Jacobian )
00751 {
00752 MRPT_START
00753 ASSERT_(x.size()>0 && increments.size() == x.size());
00754
00755 size_t m = 0;
00756 const size_t n = x.size();
00757
00758 for (size_t j=0;j<n;j++) { ASSERT_( increments[j]>0 ) }
00759
00760 VECTORLIKE3 f_minus, f_plus;
00761 VECTORLIKE x_mod(x);
00762
00763
00764 for (size_t j=0;j<n;j++)
00765 {
00766
00767 x_mod[j]=x[j]+increments[j];
00768 functor(x_mod,userParam, f_plus);
00769
00770 x_mod[j]=x[j]-increments[j];
00771 functor(x_mod,userParam, f_minus);
00772
00773 x_mod[j]=x[j];
00774 const double Ax_2_inv = 0.5/increments[j];
00775
00776
00777 if (j==0)
00778 {
00779 m = f_plus.size();
00780 out_Jacobian.setSize(m,n);
00781 }
00782
00783 for (size_t i=0;i<m;i++)
00784 out_Jacobian.get_unsafe(i,j) = Ax_2_inv* (f_plus[i]-f_minus[i]);
00785
00786 }
00787
00788 MRPT_END
00789 }
00790
00791
00792
00793
00794
00795
00796
00797
00798 template <class T>
00799 T interpolate(
00800 const T &x,
00801 const std::vector<T> &ys,
00802 const T &x0,
00803 const T &x1 )
00804 {
00805 MRPT_START
00806 ASSERT_(x1>x0); ASSERT_(!ys.empty());
00807 const size_t N = ys.size();
00808 if (x<=x0) return ys[0];
00809 if (x>=x1) return ys[N-1];
00810 const T Ax = (x1-x0)/T(N);
00811 const size_t i = int( (x-x0)/Ax );
00812 if (i>=N-1) return ys[N-1];
00813 const T Ay = ys[i+1]-ys[i];
00814 return ys[i] + (x-(x0+i*Ax))*Ay/Ax;
00815 MRPT_END
00816 }
00817
00818
00819
00820
00821
00822 double BASE_IMPEXP interpolate2points(const double x, const double x0, const double y0, const double x1, const double y1, bool wrap2pi = false);
00823
00824
00825
00826
00827
00828 double BASE_IMPEXP spline(const double t, const std::vector<double> &x, const std::vector<double> &y, bool wrap2pi = false);
00829
00830
00831
00832
00833
00834
00835
00836 template <typename NUMTYPE,class VECTORLIKE>
00837 NUMTYPE leastSquareLinearFit(const NUMTYPE t, const VECTORLIKE &x, const VECTORLIKE &y, bool wrap2pi = false)
00838 {
00839 MRPT_START
00840
00841
00842 ASSERT_(x.size()==y.size());
00843 ASSERT_(x.size()>1);
00844
00845 const size_t N = x.size();
00846
00847 typedef typename VECTORLIKE::value_type NUM;
00848
00849
00850 const NUM x_min = mrpt::math::minimum(x);
00851 CMatrixTemplateNumeric<NUM> Xt(2,N);
00852 for (size_t i=0;i<N;i++)
00853 {
00854 Xt.set_unsafe(0,i, 1);
00855 Xt.set_unsafe(1,i, x[i]-x_min);
00856 }
00857
00858 CMatrixTemplateNumeric<NUM> XtX;
00859 XtX.multiply_AAt(Xt);
00860
00861 CMatrixTemplateNumeric<NUM> XtXinv;
00862 XtX.inv_fast(XtXinv);
00863
00864 CMatrixTemplateNumeric<NUM> XtXinvXt;
00865 XtXinvXt.multiply(XtXinv,Xt);
00866
00867 VECTORLIKE B;
00868 XtXinvXt.multiply_Ab(y,B);
00869
00870 ASSERT_(B.size()==2)
00871
00872 NUM ret = B[0] + B[1]*(t-x_min);
00873
00874
00875 if (!wrap2pi)
00876 return ret;
00877 else return mrpt::math::wrapToPi(ret);
00878
00879 MRPT_END
00880 }
00881
00882
00883
00884
00885
00886 template <class VECTORLIKE1,class VECTORLIKE2,class VECTORLIKE3>
00887 void leastSquareLinearFit(
00888 const VECTORLIKE1 &ts,
00889 VECTORLIKE2 &outs,
00890 const VECTORLIKE3 &x,
00891 const VECTORLIKE3 &y,
00892 bool wrap2pi = false)
00893 {
00894 MRPT_START
00895
00896
00897 ASSERT_(x.size()==y.size());
00898 ASSERT_(x.size()>1);
00899
00900 const size_t N = x.size();
00901
00902
00903 typedef typename VECTORLIKE3::value_type NUM;
00904 const NUM x_min = mrpt::math::minimum(x);
00905 CMatrixTemplateNumeric<NUM> Xt(2,N);
00906 for (size_t i=0;i<N;i++)
00907 {
00908 Xt.set_unsafe(0,i, 1);
00909 Xt.set_unsafe(1,i, x[i]-x_min);
00910 }
00911
00912 CMatrixTemplateNumeric<NUM> XtX;
00913 XtX.multiply_AAt(Xt);
00914
00915 CMatrixTemplateNumeric<NUM> XtXinv;
00916 XtX.inv_fast(XtXinv);
00917
00918 CMatrixTemplateNumeric<NUM> XtXinvXt;
00919 XtXinvXt.multiply(XtXinv,Xt);
00920
00921 VECTORLIKE3 B;
00922 XtXinvXt.multiply_Ab(y,B);
00923
00924 ASSERT_(B.size()==2)
00925
00926 outs.resize(ts.size());
00927 if (!wrap2pi)
00928 for (size_t k=0;k<ts.size();k++)
00929 outs[k] = B[0] + B[1]*(ts[k]-x_min);
00930 else
00931 for (size_t k=0;k<ts.size();k++)
00932 outs[k] = mrpt::math::wrapToPi( B[0] + B[1]*(ts[k]-x_min) );
00933
00934 MRPT_END
00935 }
00936
00937
00938
00939
00940
00941
00942
00943
00944
00945
00946
00947
00948 template <typename T, typename At, size_t N>
00949 std::vector<T>& loadVector( std::vector<T> &v, At (&theArray)[N] )
00950 {
00951 MRPT_COMPILE_TIME_ASSERT(N!=0)
00952 v.resize(N);
00953 for (size_t i=0; i < N; i++)
00954 v[i] = static_cast<T>(theArray[i]);
00955 return v;
00956 }
00957
00958
00959
00960
00961
00962 void unwrap2PiSequence(vector_double &x);
00963
00964
00965
00966
00967
00968
00969
00970 template<class VECTORLIKE1,class VECTORLIKE2,class MAT>
00971 typename VECTORLIKE1::value_type mahalanobisDistance2(
00972 const VECTORLIKE1 &X,
00973 const VECTORLIKE2 &MU,
00974 const MAT &COV )
00975 {
00976 MRPT_START
00977 #if defined(_DEBUG) || (MRPT_ALWAYS_CHECKS_DEBUG_MATRICES)
00978 ASSERT_( !X.empty() );
00979 ASSERT_( X.size()==MU.size() );
00980 ASSERT_( X.size()==size(COV,1) && COV.IsSquare() );
00981 #endif
00982 const size_t N = X.size();
00983 std::vector<typename VECTORLIKE1::value_type> X_MU(N);
00984 for (size_t i=0;i<N;i++) X_MU[i]=X[i]-MU[i];
00985 return detail::multiply_HCHt_scalar(X_MU, COV.inv() );
00986 MRPT_END
00987 }
00988
00989
00990
00991
00992
00993 template<class VECTORLIKE1,class VECTORLIKE2,class MAT>
00994 inline typename VECTORLIKE1::value_type mahalanobisDistance(
00995 const VECTORLIKE1 &X,
00996 const VECTORLIKE2 &MU,
00997 const MAT &COV )
00998 {
00999 return std::sqrt( mahalanobisDistance2(X,MU,COV) );
01000 }
01001
01002
01003
01004
01005
01006 template<class VECTORLIKE,class MAT1,class MAT2,class MAT3>
01007 typename VECTORLIKE::value_type
01008 mahalanobisDistance2(
01009 const VECTORLIKE &mean_diffs,
01010 const MAT1 &COV1,
01011 const MAT2 &COV2,
01012 const MAT3 &CROSS_COV12 )
01013 {
01014 MRPT_START
01015 #if defined(_DEBUG) || (MRPT_ALWAYS_CHECKS_DEBUG_MATRICES)
01016 ASSERT_( !mean_diffs.empty() );
01017 ASSERT_( mean_diffs.size()==size(COV1,1));
01018 ASSERT_( COV1.IsSquare() && COV2.IsSquare() );
01019 ASSERT_( size(COV1,1)==size(COV2,1));
01020 #endif
01021 const size_t N = size(COV1,1);
01022 MAT1 COV = COV1;
01023 COV+=COV2;
01024 COV.substract_An(CROSS_COV12,2);
01025 MAT1 COV_inv;
01026 COV.inv_fast(COV_inv);
01027 return detail::multiply_HCHt_scalar(mean_diffs,COV_inv);
01028 MRPT_END
01029 }
01030
01031
01032
01033
01034 template<class VECTORLIKE,class MAT1,class MAT2,class MAT3> inline typename VECTORLIKE::value_type
01035 mahalanobisDistance(
01036 const VECTORLIKE &mean_diffs,
01037 const MAT1 &COV1,
01038 const MAT2 &COV2,
01039 const MAT3 &CROSS_COV12 )
01040 {
01041 return std::sqrt( mahalanobisDistance( mean_diffs, COV1,COV2,CROSS_COV12 ));
01042 }
01043
01044
01045
01046
01047 template<class VECTORLIKE,class MATRIXLIKE>
01048 inline typename MATRIXLIKE::value_type
01049 mahalanobisDistance2(const VECTORLIKE &delta_mu,const MATRIXLIKE &cov)
01050 {
01051 ASSERTDEB_(cov.IsSquare())
01052 ASSERTDEB_(cov.getColCount()==delta_mu.size())
01053 MATRIXLIKE C_inv;
01054 cov.inv(C_inv);
01055 return detail::multiply_HCHt_scalar(delta_mu,C_inv);
01056 }
01057
01058
01059
01060
01061 template<class VECTORLIKE,class MATRIXLIKE>
01062 inline typename MATRIXLIKE::value_type
01063 mahalanobisDistance(const VECTORLIKE &delta_mu,const MATRIXLIKE &cov)
01064 {
01065 return std::sqrt(mahalanobisDistance2(delta_mu,cov));
01066 }
01067
01068
01069
01070
01071 template <typename T>
01072 T productIntegralTwoGaussians(
01073 const std::vector<T> &mean_diffs,
01074 const CMatrixTemplateNumeric<T> &COV1,
01075 const CMatrixTemplateNumeric<T> &COV2
01076 )
01077 {
01078 const size_t vector_dim = mean_diffs.size();
01079 ASSERT_(vector_dim>=1)
01080
01081 CMatrixTemplateNumeric<T> C = COV1;
01082 C+= COV2;
01083 const T cov_det = C.det();
01084 CMatrixTemplateNumeric<T> C_inv;
01085 C.inv_fast(C_inv);
01086
01087 return std::pow( M_2PI, -0.5*vector_dim ) * (1.0/std::sqrt( cov_det ))
01088 * exp( -0.5 * mean_diffs.multiply_HCHt_scalar(C_inv) );
01089 }
01090
01091
01092
01093
01094 template <typename T, size_t DIM>
01095 T productIntegralTwoGaussians(
01096 const std::vector<T> &mean_diffs,
01097 const CMatrixFixedNumeric<T,DIM,DIM> &COV1,
01098 const CMatrixFixedNumeric<T,DIM,DIM> &COV2
01099 )
01100 {
01101 ASSERT_(mean_diffs.size()==DIM);
01102
01103 CMatrixFixedNumeric<T,DIM,DIM> C = COV1;
01104 C+= COV2;
01105 const T cov_det = C.det();
01106 CMatrixFixedNumeric<T,DIM,DIM> C_inv(UNINITIALIZED_MATRIX);
01107 C.inv_fast(C_inv);
01108
01109 return std::pow( M_2PI, -0.5*DIM ) * (1.0/std::sqrt( cov_det ))
01110 * exp( -0.5 * mean_diffs.multiply_HCHt_scalar(C_inv) );
01111 }
01112
01113
01114
01115
01116 template <typename T, class VECLIKE,class MATLIKE1, class MATLIKE2>
01117 void productIntegralAndMahalanobisTwoGaussians(
01118 const VECLIKE &mean_diffs,
01119 const MATLIKE1 &COV1,
01120 const MATLIKE2 &COV2,
01121 T &maha2_out,
01122 T &intprod_out,
01123 const MATLIKE1 *CROSS_COV12=NULL
01124 )
01125 {
01126 const size_t vector_dim = mean_diffs.size();
01127 ASSERT_(vector_dim>=1)
01128
01129 MATLIKE1 C = COV1;
01130 C+= COV2;
01131 if (CROSS_COV12) { C-=*CROSS_COV12; C-=*CROSS_COV12; }
01132 const T cov_det = C.det();
01133 MATLIKE1 C_inv;
01134 C.inv_fast(C_inv);
01135
01136 maha2_out = mean_diffs.multiply_HCHt_scalar(C_inv);
01137 intprod_out = std::pow( M_2PI, -0.5*vector_dim ) * (1.0/std::sqrt( cov_det ))*exp(-0.5*maha2_out);
01138 }
01139
01140
01141
01142
01143 template <typename T, class VECLIKE,class MATRIXLIKE>
01144 void mahalanobisDistance2AndLogPDF(
01145 const VECLIKE &diff_mean,
01146 const MATRIXLIKE &cov,
01147 T &maha2_out,
01148 T &log_pdf_out)
01149 {
01150 MRPT_START
01151 ASSERTDEB_(cov.IsSquare())
01152 ASSERTDEB_(cov.getColCount()==diff_mean.size())
01153 MATRIXLIKE C_inv;
01154 cov.inv(C_inv);
01155 maha2_out = multiply_HCHt_scalar(diff_mean,C_inv);
01156 log_pdf_out = static_cast<typename MATRIXLIKE::value_type>(-0.5)* (
01157 maha2_out+
01158 static_cast<typename MATRIXLIKE::value_type>(cov.getColCount())*::log(static_cast<typename MATRIXLIKE::value_type>(M_2PI))+
01159 ::log(cov.det())
01160 );
01161 MRPT_END
01162 }
01163
01164
01165
01166
01167 template <typename T, class VECLIKE,class MATRIXLIKE>
01168 inline void mahalanobisDistance2AndPDF(
01169 const VECLIKE &diff_mean,
01170 const MATRIXLIKE &cov,
01171 T &maha2_out,
01172 T &pdf_out)
01173 {
01174 mahalanobisDistance2AndLogPDF(diff_mean,cov,maha2_out,pdf_out);
01175 pdf_out = std::exp(pdf_out);
01176 }
01177
01178
01179
01180
01181
01182
01183
01184
01185
01186
01187
01188
01189 template <size_t N, typename T>
01190 std::vector<T> make_vector(const T val1, ...)
01191 {
01192 MRPT_COMPILE_TIME_ASSERT( N>0 )
01193 std::vector<T> ret;
01194 ret.reserve(N);
01195
01196 ret.push_back(val1);
01197
01198 va_list args;
01199 va_start(args,val1);
01200 for (size_t i=0;i<N-1;i++)
01201 ret.push_back( va_arg(args,T) );
01202
01203 va_end(args);
01204 return ret;
01205 }
01206
01207 }
01208
01209 }
01210
01211 #endif