00001 #include "ndt_registration/ndt_matcher_d2d_2d.h"
00002 #include "ndt_map/ndt_cell.h"
00003 #include "ndt_map/lazy_grid.h"
00004 #include "ndt_map/pointcloud_utils.h"
00005 #include "Eigen/Eigen"
00006 #include <fstream>
00007 #include <omp.h>
00008 #include <sys/time.h>
00009 namespace lslgeneric
00010 {
00011
00012
00013
00014 void NDTMatcherD2D_2D::init(bool _isIrregularGrid,
00015 bool useDefaultGridResolutions, std::vector<double> _resolutions)
00016 {
00017 Jest.setZero();
00018 Jest.block<2,2>(0,0).setIdentity();
00019 Hest.setZero();
00020 Zest.setZero();
00021 ZHest.setZero();
00022
00023 isIrregularGrid = _isIrregularGrid;
00024 if(useDefaultGridResolutions)
00025 {
00026 resolutions.push_back(0.2);
00027 resolutions.push_back(0.5);
00028 resolutions.push_back(1);
00029 resolutions.push_back(2);
00030 }
00031 else
00032 {
00033 resolutions = _resolutions;
00034 }
00035
00036 current_resolution = 0.1;
00037 lfd1 = 1;
00038 lfd2 = 0.05;
00039 ITR_MAX = 50;
00040 DELTA_SCORE = 10e-3*current_resolution;
00041 step_control = true;
00042
00043 }
00044
00045 bool NDTMatcherD2D_2D::match( pcl::PointCloud<pcl::PointXYZ>& target,
00046 pcl::PointCloud<pcl::PointXYZ>& source,
00047 Eigen::Transform<double,3,Eigen::Affine,Eigen::ColMajor>& T ,
00048 bool useInitialGuess)
00049 {
00050
00051 struct timeval tv_start, tv_end;
00052 struct timeval tv_start0, tv_end0;
00053 double time_load =0, time_match=0, time_combined=0;
00054
00055 gettimeofday(&tv_start0,NULL);
00056
00057
00058 pcl::PointCloud<pcl::PointXYZ> sourceCloud = source;
00059
00060 Eigen::Transform<double,3,Eigen::Affine,Eigen::ColMajor> Temp, Tinit;
00061 Tinit.setIdentity();
00062 if(useInitialGuess)
00063 {
00064 lslgeneric::transformPointCloudInPlace(T,sourceCloud);
00065 Tinit = T;
00066 }
00067
00068 T.setIdentity();
00069 bool ret = false;
00070
00071 #if 0
00072 if(isIrregularGrid)
00073 {
00074
00075 OctTree<PointTarget> pr1;
00076 NDTMap<PointTarget> targetNDT( &pr1 );
00077 targetNDT.loadPointCloud( target );
00078 targetNDT.computeNDTCells();
00079
00080 OctTree<PointSource> pr2;
00081 NDTMap<PointSource> sourceNDT( &pr2 );
00082 sourceNDT.loadPointCloud( source );
00083 sourceNDT.computeNDTCells();
00084
00085 ret = this->match( targetNDT, sourceNDT, T );
00086
00087 }
00088 else
00089 #endif
00090 {
00091
00092
00093 for(int r_ctr = resolutions.size()-1; r_ctr >=0; r_ctr--)
00094 {
00095
00096 current_resolution = resolutions[r_ctr];
00097
00098 LazyGrid prototypeSource(current_resolution);
00099 LazyGrid prototypeTarget(current_resolution);
00100
00101 gettimeofday(&tv_start,NULL);
00102 NDTMap targetNDT( &prototypeTarget );
00103 targetNDT.loadPointCloud( target );
00104 targetNDT.computeNDTCells();
00105
00106 NDTMap sourceNDT( &prototypeSource );
00107 sourceNDT.loadPointCloud( sourceCloud );
00108 sourceNDT.computeNDTCells();
00109 gettimeofday(&tv_end,NULL);
00110
00111 time_load += (tv_end.tv_sec-tv_start.tv_sec)*1000.+(tv_end.tv_usec-tv_start.tv_usec)/1000.;
00112 Temp.setIdentity();
00113
00114 gettimeofday(&tv_start,NULL);
00115 ret = this->match( targetNDT, sourceNDT, Temp );
00116 lslgeneric::transformPointCloudInPlace(Temp,sourceCloud);
00117 gettimeofday(&tv_end,NULL);
00118
00119 time_match += (tv_end.tv_sec-tv_start.tv_sec)*1000.+(tv_end.tv_usec-tv_start.tv_usec)/1000.;
00120
00121
00122 T = Temp*T;
00123
00124 #ifdef DO_DEBUG_PROC
00125 std::cout<<"RESOLUTION: "<<current_resolution<<std::endl;
00126 std::cout<<"rotation : "<<Temp.rotation().eulerAngles(0,1,2).transpose()<<std::endl;
00127 std::cout<<"translation: "<<Temp.translation().transpose()<<std::endl;
00128 std::cout<<"--------------------------------------------------------\nOverall Transform:\n";
00129 std::cout<<"rotation : "<<T.rotation().eulerAngles(0,1,2).transpose()<<std::endl;
00130 std::cout<<"translation: "<<T.translation().transpose()<<std::endl;
00131
00132 #endif
00133 }
00134 }
00135 gettimeofday(&tv_end0,NULL);
00136 time_combined = (tv_end0.tv_sec-tv_start0.tv_sec)*1000.+(tv_end0.tv_usec-tv_start0.tv_usec)/1000.;
00137 std::cout<<"load: "<<time_load<<" match "<<time_match<<" combined "<<time_combined<<std::endl;
00138 if(useInitialGuess)
00139 {
00140 T = T*Tinit;
00141 }
00142 return ret;
00143 }
00144
00145 bool NDTMatcherD2D_2D::match( NDTMap& targetNDT,
00146 NDTMap& sourceNDT,
00147 Eigen::Transform<double,3,Eigen::Affine,Eigen::ColMajor>& T ,
00148 bool useInitialGuess)
00149 {
00150
00151
00152 bool convergence = false;
00153 int itr_ctr = 0;
00154 double step_size = 1;
00155 Eigen::Matrix<double,3,1> pose_increment_v, scg;
00156
00157 Eigen::MatrixXd H(3,3), score_gradient_2d(3,1);
00158 Eigen::Transform<double,3,Eigen::Affine,Eigen::ColMajor> TR, Tbest;
00159 Eigen::Vector3d transformed_vec, mean;
00160 bool ret = true;
00161 double score_best = INT_MAX;
00162 if(!useInitialGuess)
00163 {
00164 T.setIdentity();
00165 }
00166 Tbest = T;
00167
00168 std::vector<NDTCell*> nextNDT = sourceNDT.pseudoTransformNDT(T);
00169 while(!convergence)
00170 {
00171 TR.setIdentity();
00172 H.setZero();
00173 score_gradient_2d.setZero();
00174
00175
00176 double score_here = derivativesNDT_2d(nextNDT,targetNDT,score_gradient_2d,H,true);
00177 scg = score_gradient_2d;
00178
00179 if(score_here < score_best)
00180 {
00181 Tbest = T;
00182 score_best = score_here;
00183
00184 }
00185
00186 if (score_gradient_2d.norm()<= 10e-2*DELTA_SCORE)
00187 {
00188 std::cout<<"\%gradient vanished\n";
00189 for(unsigned int i=0; i<nextNDT.size(); i++)
00190 {
00191 if(nextNDT[i]!=NULL)
00192 delete nextNDT[i];
00193 }
00194 if(score_here > score_best)
00195 {
00196 std::cout<<"crap iterations, best was "<<score_best<<" last was "<<score_here<<std::endl;
00197 T = Tbest;
00198 }
00199 return true;
00200 }
00201
00202
00203
00204 Eigen::SelfAdjointEigenSolver<Eigen::Matrix<double,3,3> > Sol (H);
00205 Eigen::Matrix<double,3,1> evals = Sol.eigenvalues().real();
00206 double minCoeff = evals.minCoeff();
00207 double maxCoeff = evals.maxCoeff();
00208 if(minCoeff < 10e-3 )
00209 {
00210 Eigen::Matrix<double,3,3> evecs = Sol.eigenvectors().real();
00211 double regularizer = score_gradient_2d.norm();
00212 regularizer = regularizer + minCoeff > 0 ? regularizer : 0.001*maxCoeff - minCoeff;
00213
00214 Eigen::Matrix<double,3,1> reg;
00215
00216 reg<<regularizer,regularizer,regularizer;
00217 evals += reg;
00218 Eigen::Matrix<double,3,3> Lam;
00219 Lam = evals.asDiagonal();
00220 H = evecs*Lam*(evecs.transpose());
00221
00222 }
00223
00224
00225
00226 pose_increment_v = -H.ldlt().solve(scg);
00227 double dginit = pose_increment_v.dot(scg);
00228 if(dginit > 0)
00229 {
00230 std::cout<<"\%can't decrease in this direction any more, done \n";
00231 for(unsigned int i=0; i<nextNDT.size(); i++)
00232 {
00233 if(nextNDT[i]!=NULL)
00234 delete nextNDT[i];
00235 }
00236 if(score_here > score_best)
00237 {
00238
00239 T = Tbest;
00240 }
00241 return true;
00242 }
00243 if(step_control) {
00244 step_size = lineSearch2D(pose_increment_v,nextNDT,targetNDT);
00245 } else {
00246 step_size = 1;
00247 }
00248 pose_increment_v = step_size*pose_increment_v;
00249 TR.setIdentity();
00250 TR = Eigen::Translation<double,3>(pose_increment_v(0),pose_increment_v(1),0)*
00251 Eigen::AngleAxis<double>(pose_increment_v(2),Eigen::Vector3d::UnitZ()) ;
00252 T = TR*T;
00253
00254
00255
00256
00257
00258
00259
00260
00261
00262 for(unsigned int i=0; i<nextNDT.size(); i++)
00263 {
00264
00265 Eigen::Vector3d meanC = nextNDT[i]->getMean();
00266 Eigen::Matrix3d covC = nextNDT[i]->getCov();
00267 meanC = TR*meanC;
00268 covC = TR.rotation()*covC*TR.rotation().transpose();
00269 nextNDT[i]->setMean(meanC);
00270 nextNDT[i]->setCov(covC);
00271 }
00272
00273 if(itr_ctr>0)
00274 {
00275 convergence = ((pose_increment_v.norm()) < DELTA_SCORE);
00276 }
00277 if(itr_ctr>ITR_MAX)
00278 {
00279 convergence = true;
00280 ret = false;
00281 }
00282 itr_ctr++;
00283 }
00284
00285
00286 score_gradient_2d.setZero();
00287
00288 double score_here = derivativesNDT_2d(nextNDT,targetNDT,score_gradient_2d,H,false);
00289 if(score_here > score_best)
00290 {
00291
00292 T = Tbest;
00293 }
00294 for(unsigned int i=0; i<nextNDT.size(); i++)
00295 {
00296 if(nextNDT[i]!=NULL)
00297 delete nextNDT[i];
00298 }
00299
00300 return ret;
00301 }
00302
00303
00304 bool NDTMatcherD2D_2D::update_gradient_hessian_local_2d(
00305 Eigen::MatrixXd &score_gradient,
00306 Eigen::MatrixXd &Hessian,
00307 const Eigen::Vector3d & x,
00308 const Eigen::Matrix3d & B,
00309 const double &likelihood,
00310 const Eigen::Matrix<double,3,3> &_Jest,
00311 const Eigen::Matrix<double,9,3> &_Hest,
00312 const Eigen::Matrix<double,3,9> &_Zest,
00313 const Eigen::Matrix<double,9,9> &_ZHest,
00314 bool computeHessian)
00315 {
00316
00317
00318 Eigen::Matrix<double,3,1> _xtBJ, _xtBZBx, _Q;
00319
00320 Eigen::Matrix<double,3,3> _xtBZBJ, _xtBH, _xtBZBZBx, _xtBZhBx;
00321 Eigen::Matrix<double,1,3> _TMP1, _xtB;
00322
00323 _xtBJ.setZero();
00324 _xtBZBx.setZero();
00325 _Q.setZero();
00326 _xtBZBJ.setZero();
00327 _xtBH.setZero();
00328 _xtBZBZBx.setZero();
00329 _xtBZhBx.setZero();
00330 _TMP1.setZero();
00331 _xtB.setZero();
00332
00333 _xtB = x.transpose()*B;
00334 _xtBJ = _xtB*_Jest;
00335
00336
00337
00338 _TMP1 = _xtB*_Zest.block<3,3>(0,6)*B;
00339 _xtBZBx(2) = _TMP1*x;
00340 if(computeHessian)
00341 {
00342 _xtBZBJ.col(2) = (_TMP1*_Jest).transpose();
00343 for(unsigned int j=0; j<3; j++)
00344 {
00345 _xtBH(2,j) = _xtB*_Hest.block<3,1>(6,j);
00346 _xtBZBZBx(2,j) = _TMP1*_Zest.block<3,3>(0,3*j)*B*x;
00347 _xtBZhBx(2,j) = _xtB*_ZHest.block<3,3>(6,3*j)*B*x;
00348 }
00349 }
00350
00351 _Q = 2*_xtBJ-_xtBZBx;
00352 double factor = -(lfd2/2)*likelihood;
00353 score_gradient += _Q*factor;
00354
00355 if(computeHessian)
00356 {
00357 Hessian += factor*(2*_Jest.transpose()*B*_Jest+2*_xtBH -_xtBZhBx -2*_xtBZBJ.transpose()
00358 -2*_xtBZBJ +_xtBZBZBx +_xtBZBZBx.transpose() -lfd2*_Q*_Q.transpose()/2 );
00359
00360 }
00361 return true;
00362 }
00363
00364
00365 bool NDTMatcherD2D_2D::update_gradient_hessian_2d(
00366 Eigen::MatrixXd &score_gradient,
00367 Eigen::MatrixXd &Hessian,
00368 const Eigen::Vector3d & x,
00369 const Eigen::Matrix3d & B,
00370 const double &likelihood,
00371 bool computeHessian)
00372 {
00373
00374 xtBJ.setZero();
00375 xtBZBx.setZero();
00376 Q.setZero();
00377 JtBJ.setZero();
00378 xtBZBJ.setZero();
00379 xtBH.setZero();
00380 xtBZBZBx.setZero();
00381 xtBZhBx.setZero();
00382 TMP1.setZero();
00383 xtB.setZero();
00384
00385 xtB = x.transpose()*B;
00386 xtBJ = xtB*Jest;
00387
00388 TMP1 = xtB*Zest.block<3,3>(0,6)*B;
00389 xtBZBx(2) = TMP1*x;
00390 if(computeHessian)
00391 {
00392 xtBZBJ.col(2) = (TMP1*Jest).transpose();
00393 for(unsigned int j=0; j<3; j++)
00394 {
00395 xtBH(2,j) = xtB*Hest.block<3,1>(6,j);
00396 xtBZBZBx(2,j) = TMP1*Zest.block<3,3>(0,3*j)*B*x;
00397 xtBZhBx(2,j) = xtB*ZHest.block<3,3>(6,3*j)*B*x;
00398 }
00399 }
00400 Q = 2*xtBJ-xtBZBx;
00401 double factor = -(lfd2/2)*likelihood;
00402 score_gradient += Q*factor;
00403
00404 if(computeHessian)
00405 {
00406 Hessian += factor*(2*Jest.transpose()*B*Jest+2*xtBH -xtBZhBx -2*xtBZBJ.transpose()
00407 -2*xtBZBJ +xtBZBZBx +xtBZBZBx.transpose() -lfd2*Q*Q.transpose()/2 );
00408
00409 }
00410 return true;
00411 }
00412
00413
00414 void NDTMatcherD2D_2D::computeDerivativesLocal_2d(Eigen::Vector3d &x, Eigen::Matrix3d C1,
00415 Eigen::Matrix<double,3,3> &_Jest,
00416 Eigen::Matrix<double,9,3> &_Hest,
00417 Eigen::Matrix<double,3,9> &_Zest,
00418 Eigen::Matrix<double,9,9> &_ZHest,
00419 bool computeHessian)
00420 {
00421 _Jest(0,2) = -x(1);
00422 _Jest(1,2) = x(0);
00423
00424 Eigen::Matrix3d myBlock;
00425
00426 _Zest.block<3,3>(0,6)<<
00427 -2*C1(0,1), -C1(1,1) + C1(0,0), -C1(1,2),
00428 -C1(1,1) + C1(0,0), 2*C1(0,1), C1(0,2),
00429 -C1(1,2), C1(0,2), 0;
00430
00431 if(computeHessian)
00432 {
00433 _Hest.block<3,1>(6,2)<<-x(0),-x(1),0;
00434 _ZHest.block<3,3>(6,6)<<
00435 2*C1(1,1) - 2*C1(0,0), -4*C1(0,1), -C1(0,2),
00436 -4*C1(0,1), 2*C1(0,0) - 2*C1(1,1), -C1(1,2),
00437 -C1(0,2), -C1(1,2), 0;
00438 }
00439 }
00440
00441 void NDTMatcherD2D_2D::computeDerivatives_2d(Eigen::Vector3d &x, Eigen::Matrix3d C1, bool computeHessian)
00442 {
00443
00444 Jest(0,2) = -x(1);
00445 Jest(1,2) = x(0);
00446
00447 Eigen::Matrix3d myBlock;
00448
00449 Zest.block<3,3>(0,6)<<
00450 -2*C1(0,1), -C1(1,1) + C1(0,0), -C1(1,2),
00451 -C1(1,1) + C1(0,0), 2*C1(0,1), C1(0,2),
00452 -C1(1,2), C1(0,2), 0;
00453
00454 if(computeHessian)
00455 {
00456 Hest.block<3,1>(6,2)<<-x(0),-x(1),0;
00457 ZHest.block<3,3>(6,6)<<
00458 2*C1(1,1) - 2*C1(0,0), -4*C1(0,1), -C1(0,2),
00459 -4*C1(0,1), 2*C1(0,0) - 2*C1(1,1), -C1(1,2),
00460 -C1(0,2), -C1(1,2), 0;
00461 }
00462 }
00463
00467
00468 #define USE_OMP_NDT_MATCHER_D2D_2D
00469
00470 double NDTMatcherD2D_2D::derivativesNDT_2d(
00471 const std::vector<NDTCell*> &sourceNDT,
00472 const NDTMap &targetNDT,
00473 Eigen::MatrixXd &score_gradient,
00474 Eigen::MatrixXd &Hessian,
00475 bool computeHessian
00476 )
00477 {
00478
00479
00480 double score_here = 0;
00481 int n_dimensions = score_gradient.rows();
00482
00483
00484 NUMBER_OF_ACTIVE_CELLS = 0;
00485 score_gradient.setZero();
00486 Hessian.setZero();
00487
00488 #ifdef USE_OMP_NDT_MATCHER_D2D_2D
00489 Eigen::MatrixXd score_gradient_omp;
00490 Eigen::MatrixXd score_here_omp;
00491 Eigen::MatrixXd Hessian_omp;
00492
00493 #define N_THREADS_2D 2
00494
00495
00496 score_gradient_omp.resize(n_dimensions,N_THREADS_2D);
00497 score_here_omp.resize(1,N_THREADS_2D);
00498 Hessian_omp.resize(n_dimensions,n_dimensions*N_THREADS_2D);
00499
00500 score_gradient_omp.setZero();
00501 score_here_omp.setZero();
00502 Hessian_omp.setZero();
00503
00504
00505 #pragma omp parallel num_threads(N_THREADS_2D)
00506 {
00507 #pragma omp for
00508 for(unsigned int i=0; i<sourceNDT.size(); i++)
00509 {
00510 pcl::PointXYZ point;
00511 Eigen::Vector3d transformed;
00512 Eigen::Vector3d meanMoving, meanFixed;
00513 Eigen::Matrix3d CMoving, CFixed, CSum, Cinv, R;
00514 Eigen::MatrixXd score_gradient_omp_loc(n_dimensions,1);
00515 Eigen::MatrixXd Hessian_omp_loc(n_dimensions,n_dimensions);
00516 Eigen::Matrix<double,3,3> _Jest;
00517 Eigen::Matrix<double,9,3> _Hest;
00518 Eigen::Matrix<double,3,9> _Zest;
00519 Eigen::Matrix<double,9,9> _ZHest;
00520 double score_here_loc=0;
00521 int thread_id = omp_get_thread_num();
00522 NDTCell *cell;
00523 bool exists = false;
00524 double det = 0;
00525
00526
00527 score_gradient_omp_loc.setZero();
00528 Hessian_omp_loc.setZero();
00529 _Jest.setZero();
00530 _Jest.block<2,2>(0,0).setIdentity();
00531 _Hest.setZero();
00532 _Zest.setZero();
00533 _ZHest.setZero();
00534
00535 meanMoving = sourceNDT[i]->getMean();
00536 CMoving= sourceNDT[i]->getCov();
00537 computeDerivativesLocal_2d(meanMoving, CMoving, _Jest, _Hest, _Zest, _ZHest, computeHessian);
00538
00539 point.x = meanMoving(0);
00540 point.y = meanMoving(1);
00541 point.z = meanMoving(2);
00542 std::vector<NDTCell*> cells = targetNDT.getCellsForPoint(point,2);
00543 for(unsigned int j=0; j<cells.size(); j++)
00544 {
00545 cell = cells[j];
00546 if(cell == NULL)
00547 {
00548 continue;
00549 }
00550 if(cell->hasGaussian_)
00551 {
00552 transformed = meanMoving - cell->getMean();
00553 CFixed = cell->getCov();
00554 CSum = (CFixed+CMoving);
00555 CSum.computeInverseAndDetWithCheck(Cinv,det,exists);
00556 if(!exists)
00557 {
00558 continue;
00559 }
00560 double l = (transformed).dot(Cinv*(transformed));
00561 if(l*0 != 0)
00562 {
00563 continue;
00564 }
00565
00566 double sh = -lfd1*(exp(-lfd2*l/2));
00567 if(!update_gradient_hessian_local_2d(score_gradient_omp_loc,Hessian_omp_loc,transformed, Cinv, sh,
00568 _Jest, _Hest, _Zest, _ZHest, computeHessian))
00569 {
00570 continue;
00571 }
00572 score_here_loc += sh;
00573 cell = NULL;
00574 }
00575 }
00576
00577 score_gradient_omp.col(thread_id) += score_gradient_omp_loc;
00578 Hessian_omp.block(0,n_dimensions*thread_id,n_dimensions,n_dimensions) += Hessian_omp_loc;
00579 score_here_omp(0,thread_id) += score_here_loc;
00580
00581 }
00582 }
00583
00584
00585
00586 score_gradient = score_gradient_omp.rowwise().sum();
00587 score_here = score_here_omp.sum();
00588 if(computeHessian)
00589 {
00590
00591 for(int i=0; i<N_THREADS_2D; ++i)
00592 {
00593 Hessian += Hessian_omp.block(0,n_dimensions*i,n_dimensions,n_dimensions);
00594 }
00595 }
00596
00597 #else
00598 pcl::PointXYZ point;
00599 Eigen::Vector3d transformed;
00600 Eigen::Vector3d meanMoving, meanFixed;
00601 Eigen::Matrix3d CMoving, CFixed, CSum, Cinv, R;
00602 NDTCell *cell;
00603 bool exists = false;
00604 double det = 0;
00605 for(unsigned int i=0; i<sourceNDT.size(); i++)
00606 {
00607 meanMoving = sourceNDT[i]->getMean();
00608 CMoving= sourceNDT[i]->getCov();
00609 computeDerivatives_2d(meanMoving, CMoving, computeHessian);
00610
00611 point.x = meanMoving(0);
00612 point.y = meanMoving(1);
00613 point.z = meanMoving(2);
00614 std::vector<NDTCell*> cells = targetNDT.getCellsForPoint(point,2);
00615 for(int j=0; j<cells.size(); j++)
00616 {
00617 cell = cells[j];
00618 if(cell == NULL)
00619 {
00620 continue;
00621 }
00622 if(cell->hasGaussian_)
00623 {
00624 transformed = meanMoving - cell->getMean();
00625 CFixed = cell->getCov();
00626 CSum = (CFixed+CMoving);
00627 CSum.computeInverseAndDetWithCheck(Cinv,det,exists);
00628 if(!exists)
00629 {
00630
00631 continue;
00632 }
00633 double l = (transformed).dot(Cinv*(transformed));
00634 if(l*0 != 0)
00635 {
00636
00637 continue;
00638 }
00639
00640 double sh = -lfd1*(exp(-lfd2*l/2));
00641
00642
00643
00644
00645
00646 if(!update_gradient_hessian_2d(score_gradient,Hessian,transformed, Cinv, sh, computeHessian))
00647 {
00648
00649 continue;
00650 }
00651
00652 score_here += sh;
00653
00654 cell = NULL;
00655 }
00656 }
00657 }
00658 #endif
00659
00660
00661
00662
00663 return score_here;
00664 }
00665
00666 double NDTMatcherD2D_2D::lineSearch2D(
00667 Eigen::Matrix<double,3,1> &increment,
00668 std::vector<NDTCell*> &sourceNDT,
00669 NDTMap &targetNDT
00670 )
00671 {
00672
00673
00674 double stp = 1.0;
00675 double recoverystep = 0.01;
00676 double dginit = 0.0;
00677 double ftol = 0.11111;
00678 double gtol = 0.99999;
00679 double stpmax = 4.0;
00680 double stpmin = 0.0001;
00681 int maxfev = 40;
00682 double xtol = 0.01;
00683
00684
00685 std::vector<NDTCell*> sourceNDTHere;
00686 double score_init = 0.0;
00687
00688 Eigen::Transform<double,3,Eigen::Affine,Eigen::ColMajor> ps;
00689 ps.setIdentity();
00690
00691 Eigen::Matrix<double,3,1> scg_here;
00692 Eigen::MatrixXd pincr(3,1), score_gradient_here(3,1);
00693 Eigen::MatrixXd pseudoH(3,3);
00694 Eigen::Vector3d eulerAngles;
00696
00697 int info = 0;
00698 int infoc = 1;
00699
00700
00701
00702
00703
00704
00705
00706
00707
00708
00709
00710
00711
00712
00713 score_gradient_here.setZero();
00714 score_init = derivativesNDT_2d(sourceNDT,targetNDT,score_gradient_here,pseudoH,false);
00715 scg_here = score_gradient_here;
00716
00717
00718
00719
00720
00721 dginit = increment.dot(scg_here);
00722
00723
00724 if (dginit >= 0.0)
00725 {
00726 std::cout << "MoreThuente::cvsrch - wrong direction (dginit = " << dginit << ")" << std::endl;
00727
00728
00729
00730 increment = -increment;
00731 dginit = -dginit;
00732
00733 if (dginit >= 0.0)
00734 {
00735 for(unsigned int i=0; i<sourceNDTHere.size(); i++)
00736 {
00737 if(sourceNDTHere[i]!=NULL)
00738 delete sourceNDTHere[i];
00739 }
00740 return recoverystep;
00741 }
00742 }
00743 else
00744 {
00745
00746 }
00747
00748
00749
00750 bool brackt = false;
00751 bool stage1 = true;
00752 int nfev = 0;
00753 double dgtest = ftol * dginit;
00754 double width = stpmax - stpmin;
00755 double width1 = 2 * width;
00756
00757
00758
00759 double finit = 0.0;
00760 finit = score_init;
00761
00762
00763
00764
00765
00766
00767
00768
00769 double stx = 0.0;
00770 double fx = finit;
00771 double dgx = dginit;
00772 double sty = 0.0;
00773 double fy = finit;
00774 double dgy = dginit;
00775
00776
00777
00778
00779
00780
00781
00782
00783 double stmin, stmax;
00784 double fm, fxm, fym, dgm, dgxm, dgym;
00785
00786 while (1)
00787 {
00788
00789
00790 if (brackt)
00791 {
00792 stmin = MoreThuente::min(stx, sty);
00793 stmax = MoreThuente::max(stx, sty);
00794 }
00795 else
00796 {
00797 stmin = stx;
00798 stmax = stp + 4 * (stp - stx);
00799 }
00800
00801
00802 stp = MoreThuente::max(stp, stpmin);
00803 stp = MoreThuente::min(stp, stpmax);
00804
00805
00806
00807
00808 if ((brackt && ((stp <= stmin) || (stp >= stmax))) ||
00809 (nfev >= maxfev - 1) || (infoc == 0) ||
00810 (brackt && (stmax - stmin <= xtol * stmax)))
00811 {
00812 stp = stx;
00813 }
00814
00815
00816
00818
00819 pincr = stp*increment;
00820
00821 ps = Eigen::Translation<double,3>(pincr(0),pincr(1),0)*
00822 Eigen::AngleAxisd(pincr(2),Eigen::Vector3d::UnitZ());
00823
00824 for(unsigned int i=0; i<sourceNDTHere.size(); i++)
00825 {
00826 if(sourceNDTHere[i]!=NULL)
00827 delete sourceNDTHere[i];
00828 }
00829 sourceNDTHere.clear();
00830 for(unsigned int i=0; i<sourceNDT.size(); i++)
00831 {
00832 NDTCell *cell = sourceNDT[i];
00833 if(cell!=NULL)
00834 {
00835 Eigen::Vector3d mean = cell->getMean();
00836 Eigen::Matrix3d cov = cell->getCov();
00837 mean = ps*mean;
00838 cov = ps.rotation()*cov*ps.rotation().transpose();
00839 NDTCell* nd = (NDTCell*)cell->copy();
00840 nd->setMean(mean);
00841 nd->setCov(cov);
00842 sourceNDTHere.push_back(nd);
00843 }
00844 }
00845
00846 double f = 0.0;
00847 score_gradient_here.setZero();
00848
00849
00850
00851
00852
00853
00854
00855
00856 f = derivativesNDT_2d(sourceNDTHere,targetNDT,score_gradient_here,pseudoH,false);
00857
00858
00859
00860
00861
00862
00863 double dg = 0.0;
00864 scg_here = score_gradient_here;
00865
00866
00867
00868 dg = increment.dot(scg_here);
00869
00870
00871
00872
00873 nfev ++;
00874
00876
00877
00878
00879 double ftest1 = finit + stp * dgtest;
00880
00881
00882
00883
00884 if ((brackt && ((stp <= stmin) || (stp >= stmax))) || (infoc == 0))
00885 info = 6;
00886
00887 if ((stp == stpmax) && (f <= ftest1) && (dg <= dgtest))
00888 info = 5;
00889
00890 if ((stp == stpmin) && ((f > ftest1) || (dg >= dgtest)))
00891 info = 4;
00892
00893 if (nfev >= maxfev)
00894 info = 3;
00895
00896 if (brackt && (stmax-stmin <= xtol*stmax))
00897 info = 2;
00898
00899
00900 bool sufficientDecreaseTest = false;
00901 sufficientDecreaseTest = (f <= ftest1);
00902
00903
00904
00905
00906 if ((sufficientDecreaseTest) && (fabs(dg) <= gtol*(-dginit)))
00907 info = 1;
00908
00909 if (info != 0)
00910 {
00911 if (info != 1)
00912 {
00913
00914
00915
00916
00917 stp = recoverystep;
00918
00919
00920
00921
00922
00923 }
00924 else
00925 {
00926
00927 }
00928
00929
00930
00931
00932
00933 for(unsigned int i=0; i<sourceNDTHere.size(); i++)
00934 {
00935 if(sourceNDTHere[i]!=NULL)
00936 delete sourceNDTHere[i];
00937 }
00938
00939 return stp;
00940
00941 }
00942
00943
00944
00945
00946
00947
00948
00949 if (stage1 && (f <= ftest1) && (dg >= MoreThuente::min(ftol, gtol) * dginit))
00950 {
00951 stage1 = false;
00952 }
00953
00954
00955
00956
00957
00958
00959
00960 if (stage1 && (f <= fx) && (f > ftest1))
00961 {
00962
00963
00964
00965 fm = f - stp * dgtest;
00966 fxm = fx - stx * dgtest;
00967 fym = fy - sty * dgtest;
00968 dgm = dg - dgtest;
00969 dgxm = dgx - dgtest;
00970 dgym = dgy - dgtest;
00971
00972
00973
00974
00975
00976 infoc = MoreThuente::cstep(stx,fxm,dgxm,sty,fym,dgym,stp,fm,dgm,
00977 brackt,stmin,stmax);
00978
00979
00980
00981 fx = fxm + stx*dgtest;
00982 fy = fym + sty*dgtest;
00983 dgx = dgxm + dgtest;
00984 dgy = dgym + dgtest;
00985
00986 }
00987
00988 else
00989 {
00990
00991
00992
00993
00994
00995 infoc = MoreThuente::cstep(stx,fx,dgx,sty,fy,dgy,stp,f,dg,
00996 brackt,stmin,stmax);
00997
00998 }
00999
01000
01001
01002
01003 if (brackt)
01004 {
01005 if (fabs(sty - stx) >= 0.66 * width1)
01006 stp = stx + 0.5 * (sty - stx);
01007 width1 = width;
01008 width = fabs(sty-stx);
01009 }
01010
01011 }
01012
01013 }
01014
01015 int NDTMatcherD2D_2D::MoreThuente::cstep(double& stx, double& fx, double& dx,
01016 double& sty, double& fy, double& dy,
01017 double& stp, double& fp, double& dp,
01018 bool& brackt, double stmin, double stmax)
01019 {
01020 int info = 0;
01021
01022
01023
01024 if ((brackt && ((stp <= MoreThuente::min(stx, sty)) || (stp >= MoreThuente::max(stx, sty)))) ||
01025 (dx * (stp - stx) >= 0.0) || (stmax < stmin))
01026 return info;
01027
01028
01029
01030 double sgnd = dp * (dx / fabs(dx));
01031
01032
01033
01034
01035
01036
01037 bool bound;
01038 double theta;
01039 double s;
01040 double gamma;
01041 double p,q,r;
01042 double stpc, stpq, stpf;
01043
01044 if (fp > fx)
01045 {
01046 info = 1;
01047 bound = 1;
01048 theta = 3 * (fx - fp) / (stp - stx) + dx + dp;
01049
01050
01051
01052 s = MoreThuente::absmax(theta, dx, dp);
01053 gamma = s * sqrt(((theta / s) * (theta / s)) - (dx / s) * (dp / s));
01054 if (stp < stx)
01055 gamma = -gamma;
01056
01057 p = (gamma - dx) + theta;
01058 q = ((gamma - dx) + gamma) + dp;
01059 r = p / q;
01060 stpc = stx + r * (stp - stx);
01061 stpq = stx + ((dx / ((fx - fp) / (stp - stx) + dx)) / 2) * (stp - stx);
01062 if (fabs(stpc - stx) < fabs(stpq - stx))
01063 stpf = stpc;
01064 else
01065 stpf = stpc + (stpq - stpc) / 2;
01066
01067 brackt = true;
01068 }
01069
01070
01071
01072
01073
01074
01075 else if (sgnd < 0.0)
01076 {
01077 info = 2;
01078 bound = false;
01079 theta = 3 * (fx - fp) / (stp - stx) + dx + dp;
01080 s = MoreThuente::absmax(theta,dx,dp);
01081 gamma = s * sqrt(((theta/s) * (theta/s)) - (dx / s) * (dp / s));
01082 if (stp > stx)
01083 gamma = -gamma;
01084 p = (gamma - dp) + theta;
01085 q = ((gamma - dp) + gamma) + dx;
01086 r = p / q;
01087 stpc = stp + r * (stx - stp);
01088 stpq = stp + (dp / (dp - dx)) * (stx - stp);
01089 if (fabs(stpc - stp) > fabs(stpq - stp))
01090 stpf = stpc;
01091 else
01092 stpf = stpq;
01093 brackt = true;
01094 }
01095
01096
01097
01098
01099
01100
01101
01102
01103
01104
01105 else if (fabs(dp) < fabs(dx))
01106 {
01107 info = 3;
01108 bound = true;
01109 theta = 3 * (fx - fp) / (stp - stx) + dx + dp;
01110 s = MoreThuente::absmax(theta, dx, dp);
01111
01112
01113
01114
01115 gamma = s * sqrt(max(0,(theta / s) * (theta / s) - (dx / s) * (dp / s)));
01116 if (stp > stx)
01117 gamma = -gamma;
01118
01119 p = (gamma - dp) + theta;
01120 q = (gamma + (dx - dp)) + gamma;
01121 r = p / q;
01122 if ((r < 0.0) && (gamma != 0.0))
01123 stpc = stp + r * (stx - stp);
01124 else if (stp > stx)
01125 stpc = stmax;
01126 else
01127 stpc = stmin;
01128
01129 stpq = stp + (dp/ (dp - dx)) * (stx - stp);
01130 if (brackt)
01131 {
01132 if (fabs(stp - stpc) < fabs(stp - stpq))
01133 stpf = stpc;
01134 else
01135 stpf = stpq;
01136 }
01137 else
01138 {
01139 if (fabs(stp - stpc) > fabs(stp - stpq))
01140 stpf = stpc;
01141 else
01142 stpf = stpq;
01143 }
01144 }
01145
01146
01147
01148
01149
01150
01151 else
01152 {
01153 info = 4;
01154 bound = false;
01155 if (brackt)
01156 {
01157 theta = 3 * (fp - fy) / (sty - stp) + dy + dp;
01158 s = MoreThuente::absmax(theta, dy, dp);
01159 gamma = s * sqrt(((theta/s)*(theta/s)) - (dy / s) * (dp / s));
01160 if (stp > sty)
01161 gamma = -gamma;
01162 p = (gamma - dp) + theta;
01163 q = ((gamma - dp) + gamma) + dy;
01164 r = p / q;
01165 stpc = stp + r * (sty - stp);
01166 stpf = stpc;
01167 }
01168 else if (stp > stx)
01169 stpf = stmax;
01170 else
01171 stpf = stmin;
01172 }
01173
01174
01175
01176
01177 if (fp > fx)
01178 {
01179 sty = stp;
01180 fy = fp;
01181 dy = dp;
01182 }
01183 else
01184 {
01185 if (sgnd < 0.0)
01186 {
01187 sty = stx;
01188 fy = fx;
01189 dy = dx;
01190 }
01191 stx = stp;
01192 fx = fp;
01193 dx = dp;
01194 }
01195
01196
01197
01198 stpf = MoreThuente::min(stmax, stpf);
01199 stpf = MoreThuente::max(stmin, stpf);
01200 stp = stpf;
01201 if (brackt && bound)
01202 {
01203 if (sty > stx)
01204 stp = min(stx + 0.66 * (sty - stx), stp);
01205 else
01206 stp = max(stx + 0.66 * (sty - stx), stp);
01207 }
01208
01209 return info;
01210
01211 }
01212
01213 double NDTMatcherD2D_2D::MoreThuente::min(double a, double b)
01214 {
01215 return (a < b ? a : b);
01216 }
01217
01218 double NDTMatcherD2D_2D::MoreThuente::max(double a, double b)
01219 {
01220 return (a > b ? a : b);
01221 }
01222
01223 double NDTMatcherD2D_2D::MoreThuente::absmax(double a, double b, double c)
01224 {
01225 a = fabs(a);
01226 b = fabs(b);
01227 c = fabs(c);
01228
01229 if (a > b)
01230 return (a > c) ? a : c;
01231 else
01232 return (b > c) ? b : c;
01233 }
01234
01235 double NDTMatcherD2D_2D::normalizeAngle(double a)
01236 {
01237
01238 return atan2(sin(a), cos(a));
01239
01240 }
01241
01242 }