contrib/rpl/rgrl/rgrl_est_dis_homo2d_lm.cxx
Go to the documentation of this file.
00001 #include "rgrl_est_dis_homo2d_lm.h"
00002 //:
00003 // \file
00004 #include <rgrl/rgrl_est_homography2d.h>
00005 #include <rgrl/rgrl_trans_homography2d.h>
00006 #include <rgrl/rgrl_trans_rad_dis_homo2d.h>
00007 #include <rgrl/rgrl_trans_rigid.h>
00008 #include <rgrl/rgrl_trans_translation.h>
00009 #include <rgrl/rgrl_trans_similarity.h>
00010 #include <rgrl/rgrl_trans_affine.h>
00011 #include <rgrl/rgrl_match_set.h>
00012 #include <rgrl/rgrl_internal_util.h>
00013 
00014 #include <vnl/vnl_double_2.h>
00015 #include <vnl/vnl_double_3.h>
00016 #include <vnl/vnl_double_2x2.h>
00017 #include <vnl/vnl_double_3x3.h>
00018 #include <vnl/vnl_matrix_fixed.h>
00019 #include <vnl/vnl_vector_fixed.h>
00020 #include <vnl/vnl_math.h>
00021 #include <vnl/vnl_fastops.h>
00022 #include <vnl/vnl_transpose.h>
00023 #include <vnl/vnl_least_squares_function.h>
00024 #include <vnl/algo/vnl_orthogonal_complement.h>
00025 #include <vnl/algo/vnl_levenberg_marquardt.h>
00026 #include <vnl/algo/vnl_svd.h>
00027 
00028 #include <vcl_cassert.h>
00029 
00030 static
00031 inline
00032 void
00033 H2h( vnl_matrix_fixed<double,3,3> const& H, vnl_vector<double>& h )
00034 {
00035   // h is already size 11
00036   // h.set_size( 9 );
00037   for ( unsigned i=0; i<H.rows(); ++i )
00038     for ( unsigned j=0; j<H.cols(); ++j )
00039       h( i*3+j ) = H(i,j);
00040 }
00041 
00042 static
00043 inline
00044 void
00045 h2H( vnl_vector<double> const& h, vnl_matrix_fixed<double,3,3>& H )
00046 {
00047   for ( unsigned i=0; i<3; ++i )
00048     for ( unsigned j=0; j<3; ++j )
00049       H(i,j) = h( i*3+j );
00050 }
00051 
00052 // map homography
00053 static
00054 inline
00055 void
00056 map_inhomo_point( vnl_double_2& mapped, vnl_matrix_fixed<double, 3, 3> const& H, vnl_vector<double> const& loc )
00057 {
00058   vnl_double_3 homo_from( loc[0], loc[1], 1 );
00059   vnl_double_3 homo_to = H * homo_from;
00060   mapped[0] = homo_to[0]/homo_to[2];
00061   mapped[1] = homo_to[1]/homo_to[2];
00062 }
00063 
00064 //: Return the jacobian of the transform.
00065 static
00066 void
00067 homo_wrt_loc( vnl_matrix_fixed<double, 2, 2> &      jac_loc,
00068               vnl_matrix_fixed<double, 3, 3> const& H,
00069               vnl_vector_fixed<double, 2>    const& from_loc )
00070 {
00071   // The jacobian is a 2x2 matrix with entries
00072   // [d(f_0)/dx   d(f_0)/dy;
00073   //  d(f_1)/dx   d(f_1)/dy]
00074   //
00075   const double mapped_w = H(2,0)*from_loc[0] + H(2,1)*from_loc[1] + H(2,2);
00076 
00077   // w/ respect to x
00078   jac_loc(0,0) = H(0,0)*( H(2,1)*from_loc[1]+H(2,2) ) - H(2,0)*( H(0,1)*from_loc[1] + H(0,2) );
00079   jac_loc(1,0) = H(1,0)*( H(2,1)*from_loc[1]+H(2,2) ) - H(2,0)*( H(1,1)*from_loc[1] + H(1,2) );
00080   // w/ respect to y
00081   jac_loc(0,1) = H(0,1)*( H(2,0)*from_loc[0]+H(2,2) ) - H(2,1)*( H(0,0)*from_loc[0] + H(0,2) );
00082   jac_loc(1,1) = H(1,1)*( H(2,0)*from_loc[0]+H(2,2) ) - H(2,1)*( H(1,0)*from_loc[0] + H(1,2) );
00083 
00084   jac_loc *= (1/(mapped_w*mapped_w));
00085 }
00086 
00087 
00088 static
00089 void
00090 homo_wrt_h( vnl_matrix_fixed<double, 2, 9> &       jac_h,
00091             vnl_matrix_fixed<double, 3, 3> const&  H,
00092             vnl_vector_fixed<double, 2>    const&  from_loc )
00093 {
00094   vnl_matrix_fixed<double, 3, 9 > jf(0.0); // homogeneous coordinate
00095   vnl_matrix_fixed<double, 2, 3 > jg(0.0); // inhomo, [u/w, v/w]^T
00096 
00097   // transform coordinate
00098   vnl_double_3 from_homo( from_loc[0], from_loc[1], 1 );
00099   vnl_double_3 mapped_homo = H * from_homo;
00100 
00101   // homogeneous coordinate w.r.t homography parameters
00102   jf(0,0) = jf(1,3) = jf(2,6) = from_homo[0]; // x
00103   jf(0,1) = jf(1,4) = jf(2,7) = from_homo[1]; // y
00104   jf(0,2) = jf(1,5) = jf(2,8) = 1.0;
00105 
00106   // derivatives w.r.t division
00107   jg(0,0) = 1.0/mapped_homo[2];
00108   jg(0,2) = -mapped_homo[0]/vnl_math_sqr(mapped_homo[2]);
00109   jg(1,1) = 1.0/mapped_homo[2];
00110   jg(1,2) = -mapped_homo[1]/vnl_math_sqr(mapped_homo[2]);
00111 
00112   // Apply chain rule: Jab_g(f(p)) = Jac_g * Jac_f
00113   jac_h = jg * jf;
00114 }
00115 
00116 // distort image coordinate
00117 static
00118 inline
00119 void
00120 distort( vnl_double_2& dis_loc, vnl_double_2 const& true_loc, double k1 )
00121 {
00122   const double c = 1 + k1 * true_loc.squared_magnitude();
00123   dis_loc = c * true_loc;
00124 }
00125 
00126 // jacobian w.r.t k1 parameter
00127 static
00128 inline
00129 void
00130 distort_wrt_k1( vnl_double_2& jac_k1, vnl_double_2 const& true_loc )
00131 {
00132   const double c = true_loc.squared_magnitude();
00133   jac_k1 = c * true_loc;
00134 }
00135 
00136 
00137 // jacobian w.r.t location
00138 static
00139 inline
00140 void
00141 distort_wrt_loc( vnl_double_2x2& jac_loc, vnl_double_2 const& true_loc, double k1 )
00142 {
00143   const double c = 1 + k1 * true_loc.squared_magnitude();
00144 
00145   jac_loc(0,0) = c + 2*k1*vnl_math_sqr(true_loc[0]);
00146   jac_loc(1,1) = c + 2*k1*vnl_math_sqr(true_loc[1]);
00147   jac_loc(0,1) = jac_loc(1,0) = 2 * k1 * true_loc[0] * true_loc[1];
00148 }
00149 
00150 
00151 class rgrl_rad_dis_homo2d_func
00152 : public vnl_least_squares_function
00153 {
00154  public:
00155   //: ctor
00156   rgrl_rad_dis_homo2d_func( rgrl_set_of<rgrl_match_set_sptr> const& matches,
00157                             int num_res, bool with_grad = true )
00158   : vnl_least_squares_function( 11, num_res, with_grad ? use_gradient : no_gradient ),
00159     matches_ptr_( &matches ),
00160     from_centre_(2, 0.0), to_centre_(2, 0.0)
00161   {      }
00162 
00163   void set_centres( vnl_vector<double> const& fc, vnl_vector<double> const& tc )
00164   {
00165     assert( fc.size() == 2 && tc.size() == 2 );
00166     from_centre_ = fc;
00167     to_centre_ = tc;
00168   }
00169 
00170   //: obj func value
00171   void f(vnl_vector<double> const& x, vnl_vector<double>& fx);
00172 
00173   //: Jacobian
00174   void gradf(vnl_vector<double> const& x, vnl_matrix<double>& jacobian);
00175 
00176  protected:
00177   typedef rgrl_match_set::const_from_iterator FIter;
00178   typedef FIter::to_iterator TIter;
00179 
00180   rgrl_set_of<rgrl_match_set_sptr> const* matches_ptr_;
00181   vnl_double_2                            from_centre_, to_centre_;
00182 };
00183 
00184 void
00185 rgrl_rad_dis_homo2d_func::
00186 f(vnl_vector<double> const& x, vnl_vector<double>& fx)
00187 {
00188   vnl_double_2 true_mapped, true_from, from, dis_mapped;
00189   vnl_matrix_fixed<double,2,2> error_proj_sqrt;
00190   double k1_from = x[9];
00191   double k1_to   = x[10];
00192   vnl_double_3x3  H;
00193   h2H( x, H );
00194 
00195   unsigned int ind = 0;
00196   for ( unsigned ms = 0; ms<matches_ptr_->size(); ++ms )
00197     if ( (*matches_ptr_)[ms] != 0 ) { // if pointer is valid
00198 
00199       rgrl_match_set const& one_set = *((*matches_ptr_)[ms]);
00200 
00201       for ( FIter fi=one_set.from_begin(); fi!=one_set.from_end(); ++fi )
00202       {
00203         // map from point
00204         from = fi.from_feature()->location();
00205         from -= from_centre_;
00206 
00207         // Step 1.
00208         distort( true_from, from, k1_from );
00209         // Step 2.
00210         map_inhomo_point( true_mapped, H, true_from.as_ref() );
00211         // Step 3.
00212         distort( dis_mapped, true_mapped, k1_to );
00213 
00214         for ( TIter ti=fi.begin(); ti!=fi.end(); ++ti )
00215         {
00216           vnl_double_2 to = ti.to_feature()->location();
00217           to -= to_centre_;
00218           error_proj_sqrt = ti.to_feature()->error_projector_sqrt();
00219           double const wgt = vcl_sqrt(ti.cumulative_weight());
00220           vnl_double_2 diff = error_proj_sqrt * (dis_mapped - to);
00221 
00222           // fill in
00223           fx(ind) = wgt*diff[0];
00224           fx(ind+1) = wgt*diff[1];
00225           ind+=2;
00226         }
00227       }
00228   }
00229 
00230   // check
00231   assert( ind == get_number_of_residuals() );
00232 }
00233 
00234 void
00235 rgrl_rad_dis_homo2d_func::
00236 gradf(vnl_vector<double> const& x, vnl_matrix<double>& jacobian)
00237 {
00238   assert( jacobian.rows() == get_number_of_residuals() && jacobian.cols() == 11 );
00239 
00240   double k1_from = x[9];
00241   double k1_to   = x[10];
00242   vnl_double_3x3  H;
00243   h2H( x, H );
00244 
00245   vnl_double_2x2 pu_pd;
00246   vnl_double_2 pu_k1_from;
00247   vnl_double_2x2 qu_pu;
00248   vnl_matrix_fixed<double, 2, 9> qu_h;
00249   vnl_double_2x2 qd_qu;
00250   vnl_double_2 qd_k1_to;
00251   vnl_matrix_fixed<double, 2, 9> qd_h;
00252   vnl_double_2  qd_k1_from;
00253 
00254   vnl_matrix_fixed<double,2,2> error_proj_sqrt;
00255 
00256   unsigned int ind = 0;
00257   for ( unsigned ms = 0; ms<matches_ptr_->size(); ++ms )
00258 
00259     if ( (*matches_ptr_)[ms] ) { // if pointer is valid
00260 
00261       rgrl_match_set const& one_set = *((*matches_ptr_)[ms]);
00262 
00263       for ( FIter fi=one_set.from_begin(); fi!=one_set.from_end(); ++fi )
00264       {
00265         // Step 1. undistorted from coordinate and compute apu/apd
00266         vnl_double_2 dis_from_loc = fi.from_feature()->location();
00267         dis_from_loc -= from_centre_;
00268 
00269         vnl_double_2 true_from_loc;
00270         // make the trick: *distort*
00271         distort( true_from_loc, dis_from_loc, k1_from );
00272         distort_wrt_loc( pu_pd, dis_from_loc, k1_from );
00273         distort_wrt_k1( pu_k1_from, dis_from_loc );
00274 
00275         // Step 2. homography transformation
00276         vnl_double_2 true_to_loc;
00277         map_inhomo_point( true_to_loc, H, true_from_loc.as_ref() );
00278         homo_wrt_loc( qu_pu, H, true_from_loc );
00279         homo_wrt_h( qu_h, H, true_from_loc );
00280 
00281         // Step 3. distorted To coodinates
00282         distort_wrt_loc( qd_qu, true_to_loc, k1_to );
00283         distort_wrt_k1( qd_k1_to, true_to_loc );
00284 
00285         // Steop 4. apply chain rule
00286         qd_h = qd_qu * qu_h;
00287         qd_k1_from = qd_qu * qu_pu * pu_k1_from;
00288 
00289         for ( TIter ti=fi.begin(); ti!=fi.end(); ++ti ) {
00290           //vnl_double_2 to = ti.to_feature()->location();
00291           error_proj_sqrt = ti.to_feature()->error_projector_sqrt();
00292           double const wgt = vcl_sqrt(ti.cumulative_weight());
00293 
00294           qd_k1_from = wgt * error_proj_sqrt * qd_k1_from;
00295           qd_k1_to   = wgt * error_proj_sqrt * qd_k1_to;
00296           qd_h       = wgt * error_proj_sqrt * qd_h;
00297 
00298           // fill in
00299           for ( unsigned i=0; i<9; i++ ) {
00300             jacobian(ind, i)   = qd_h(0, i);
00301             jacobian(ind+1, i) = qd_h(1, i);
00302           }
00303           // k1_from
00304           jacobian(ind, 9)   = qd_k1_from[0];
00305           jacobian(ind+1, 9) = qd_k1_from[1];
00306           // k1_to
00307           jacobian(ind, 10)   = qd_k1_to[0];
00308           jacobian(ind+1, 10) = qd_k1_to[1];
00309 
00310           ind+=2;
00311         }
00312       }
00313   }
00314 }
00315 
00316 
00317 // --------------------------------------------------------------------
00318 
00319 rgrl_est_dis_homo2d_lm::
00320 rgrl_est_dis_homo2d_lm( vnl_vector<double> const& from_centre,
00321                         vnl_vector<double> const& to_centre,
00322                         bool with_grad )
00323   : from_centre_( from_centre),
00324     to_centre_( to_centre ),
00325     with_grad_( with_grad )
00326 {
00327   assert( from_centre.size() == 2 && to_centre.size() == 2 );
00328 
00329   rgrl_estimator::set_param_dof( 10 );
00330 
00331   // default value
00332   rgrl_nonlinear_estimator::set_max_num_iter( 50 );
00333   rgrl_nonlinear_estimator::set_rel_thres( 1e-5 );
00334 }
00335 
00336 rgrl_transformation_sptr
00337 rgrl_est_dis_homo2d_lm::
00338 estimate( rgrl_set_of<rgrl_match_set_sptr> const& matches,
00339           rgrl_transformation const& cur_transform ) const
00340 {
00341   // get initialization
00342   vnl_matrix_fixed<double, 3, 3> init_H;
00343   double k1_from = 0, k1_to = 0;
00344 
00345   if ( cur_transform.is_type( rgrl_trans_rad_dis_homo2d::type_id() ) )
00346   {
00347     rgrl_trans_rad_dis_homo2d const& trans = static_cast<rgrl_trans_rad_dis_homo2d const&>( cur_transform );
00348     init_H = trans.H();
00349     k1_from = trans.k1_from();
00350     k1_to   = trans.k1_to();
00351 
00352     // check centre
00353     assert( from_centre_ == trans.from_centre() );
00354     assert( to_centre_ == trans.to_centre() );
00355   }
00356   else {
00357 
00358     if ( !rgrl_internal_util_upgrade_to_homography2D( init_H, cur_transform ) ) {
00359 
00360       // use normalized DLT to initialize
00361       DebugMacro( 0, "Use normalized DLT to initialize" );
00362       rgrl_est_homography2d est_homo;
00363       rgrl_transformation_sptr tmp_trans= est_homo.estimate( matches, cur_transform );
00364       if ( !tmp_trans )
00365         return 0;
00366       rgrl_trans_homography2d const& trans = static_cast<rgrl_trans_homography2d const&>( *tmp_trans );
00367       init_H = trans.H();
00368     }
00369 
00370     // make the init homography as a CENTERED one
00371     // centered H_ = to_matrix * H * from_matrix^-1
00372     //
00373     vnl_matrix_fixed<double, 3, 3> to_trans;
00374     to_trans.set_identity();
00375     to_trans(0,2) = -to_centre_[0];
00376     to_trans(1,2) = -to_centre_[1];
00377 
00378     vnl_matrix_fixed<double, 3, 3> from_inv;
00379     from_inv.set_identity();
00380     from_inv(0,2) = from_centre_[0];
00381     from_inv(1,2) = from_centre_[1];
00382 
00383     init_H = to_trans * init_H * from_inv;
00384   }
00385 
00386 
00387   // count the number of constraints/residuals
00388   typedef rgrl_match_set::const_from_iterator FIter;
00389   typedef FIter::to_iterator TIter;
00390   unsigned int tot_num = 0;
00391   for ( unsigned ms = 0; ms<matches.size(); ++ms )
00392     if ( matches[ms] ) { // if pointer is valid
00393       rgrl_match_set const& one_set = *(matches[ms]);
00394       for ( FIter fi=one_set.from_begin(); fi!=one_set.from_end(); ++fi ) {
00395         tot_num += fi.size() * 2;  // each point provides two constraints
00396       }
00397     }
00398 
00399   DebugMacro( 3, "From center: " << from_centre_
00400               << "  To center: " << to_centre_ << vcl_endl );
00401 
00402   // normalize H
00403   init_H /= init_H.array_two_norm();
00404 
00405   // convert to vector form
00406   vnl_vector<double> p( 11, 0.0 );
00407   H2h( init_H, p );
00408   p[9] = k1_from;
00409   p[10] = k1_to;
00410 
00411   // construct least square cost function
00412   rgrl_rad_dis_homo2d_func dis_homo_func( matches, tot_num, with_grad_ );
00413   dis_homo_func.set_centres( from_centre_.as_ref(), to_centre_.as_ref() );
00414 
00415   vnl_levenberg_marquardt lm( dis_homo_func );
00416   //lm.set_trace( true );
00417   //lm.set_check_derivatives( 3 );
00418   // we don't need it to be super accurate
00419   lm.set_f_tolerance( relative_threshold_ );
00420   lm.set_max_function_evals( max_num_iterations_ );
00421 
00422   bool ret;
00423   if ( with_grad_ )
00424     ret = lm.minimize_using_gradient(p);
00425   else
00426     ret = lm.minimize_without_gradient(p);
00427   if ( !ret ) {
00428     WarningMacro( "Levenberg-Marquardt failed" );
00429     return 0;
00430   }
00431   // lm.diagnose_outcome(vcl_cout);
00432 
00433   // convert parameters back into matrix form
00434   h2H( p, init_H );
00435   k1_from = p[9];
00436   k1_to   = p[10];
00437   // normalize H
00438   init_H /= init_H.array_two_norm();
00439 
00440   // compute covariance
00441   // JtJ is INVERSE of jacobian
00442   // vnl_svd<double> svd( lm.get_JtJ(), 1e-4 );
00443   // Cannot use get_JtJ() because it is affected by the
00444   // scale in homography parameter vector
00445   // Thus, use the normalized p vector to compute Jacobian again
00446   vnl_matrix<double> jac(tot_num, 11), jtj(11, 11);
00447   dis_homo_func.gradf( p, jac );
00448   //
00449   vnl_fastops::AtA( jtj, jac );
00450   // vcl_cout << "Jacobian:\n" << jac << "\n\nJtJ:\n" << jtj << "\n\nReal JtJ:\n" << jac.transpose() * jac << vcl_endl;
00451 
00452 
00453   // compute the inverse of scatter matrix
00454 #if 1
00455   vnl_vector<double> homo_vector( 9 );
00456   H2h( init_H, homo_vector );
00457   vnl_matrix<double> block_compliment = vnl_orthogonal_complement( homo_vector );
00458 
00459   vnl_matrix<double> compliment( 11, 10, 0.0 );
00460 
00461   for ( unsigned i=0; i<9; ++i )
00462     for ( unsigned j=0; j<8; ++j )
00463       compliment(i,j) = block_compliment(i,j);
00464 
00465   // distortion parameters
00466   compliment(9,8) = 1.0;
00467   compliment(10,9) = 1.0;
00468 
00469   vnl_svd<double> svd( vnl_transpose(compliment) * jtj *compliment, 1e-6 );
00470   if ( svd.rank() < 10 ) {
00471     WarningMacro( "The covariance of homography ranks less than 8! ");
00472     return 0;
00473   }
00474 
00475   vnl_matrix<double>covar = compliment * svd.inverse() * compliment.transpose();
00476 
00477 
00478 #else
00479   // compute inverse
00480   //
00481   vnl_svd<double> svd( jtj, 1e-6 );
00482   DebugMacro(3, "SVD of JtJ: " << svd << vcl_endl);
00483   // the second least singular value shall be greater than 0
00484   // or Rank 11-1 = 10
00485   if ( svd.rank() < 10 ) {
00486     WarningMacro( "The covariance of homography ranks less than 8! ");
00487   }
00488   // pseudo inverse only use first 10 singular values
00489   vnl_matrix<double> covar ( svd.pinverse(10) );
00490   DebugMacro(3, "Covariance: " << covar << vcl_endl );
00491 
00492   DebugMacro(2, "null vector: " << svd.nullvector() << "   estimate: " << p << vcl_endl );
00493 #endif
00494 
00495   return new rgrl_trans_rad_dis_homo2d( init_H.as_ref(), k1_from, k1_to, covar, from_centre_.as_ref(), to_centre_.as_ref() );
00496 }
00497 
00498 
00499 rgrl_transformation_sptr
00500 rgrl_est_dis_homo2d_lm::
00501 estimate( rgrl_match_set_sptr matches,
00502           rgrl_transformation const& cur_transform ) const
00503 {
00504   // use base class implementation
00505   return rgrl_estimator::estimate( matches, cur_transform );
00506 }
00507 
00508 //: change frome_centre
00509 void
00510 rgrl_est_dis_homo2d_lm::
00511 set_centres( vnl_vector<double> const& from_centre,
00512              vnl_vector<double> const& to_centre )
00513 {
00514   assert( from_centre.size() == 2 && to_centre.size() == 2 );
00515 
00516   from_centre_ = from_centre;
00517   to_centre_ = to_centre;
00518 }
00519 
00520 const vcl_type_info&
00521 rgrl_est_dis_homo2d_lm::
00522 transformation_type() const
00523 {
00524   return rgrl_trans_rad_dis_homo2d::type_id();
00525 }