contrib/rpl/rgrl/rgrl_est_proj_rad_func.txx
Go to the documentation of this file.
00001 #ifndef rgrl_est_proj_rad_func_txx_
00002 #define rgrl_est_proj_rad_func_txx_
00003 //:
00004 // \file
00005 // \author Gehua Yang
00006 // \date   March 2007
00007 // \brief  a generic class to estimate a homogeneous projection matrix with radial distortion parameter(s)  using L-M
00008 
00009 #include "rgrl_est_proj_rad_func.h"
00010 
00011 #include <rgrl/rgrl_est_proj_func.txx>
00012 #include <rgrl/rgrl_estimator.h>
00013 #include <rgrl/rgrl_match_set.h>
00014 #include <rgrl/rgrl_set_of.h>
00015 
00016 #include <vnl/vnl_matrix_fixed.h>
00017 #include <vnl/vnl_vector_fixed.h>
00018 #include <vnl/vnl_math.h>
00019 #include <vnl/vnl_transpose.h>
00020 #include <vnl/algo/vnl_levenberg_marquardt.h>
00021 #include <vnl/algo/vnl_svd.h>
00022 
00023 #include <vcl_vector.h>
00024 #include <vcl_cassert.h>
00025 
00026 template <unsigned int Tdim, unsigned int Fdim>
00027 rgrl_est_proj_rad_func<Tdim, Fdim>::
00028 rgrl_est_proj_rad_func( rgrl_set_of<rgrl_match_set_sptr> const& matches,
00029                         unsigned int camera_dof,
00030                         bool with_grad )
00031 : rgrl_est_proj_func<Tdim, Fdim>( matches, with_grad ),
00032   camera_dof_(camera_dof),
00033   image_centre_(double(0)),
00034   centre_mag_norm_const_(0)
00035 {
00036   //modify the dof in vnl_least_squares_function
00037   vnl_least_squares_function::init(this->proj_size_-1+camera_dof_,
00038                                    this->get_number_of_residuals() );
00039 
00040   // temperary storage space
00041   temp_rad_k_.resize( camera_dof_ );
00042 }
00043 
00044 template <unsigned int Tdim, unsigned int Fdim>
00045 rgrl_est_proj_rad_func<Tdim, Fdim>::
00046 rgrl_est_proj_rad_func( unsigned int camera_dof,
00047                         bool with_grad )
00048 : rgrl_est_proj_func<Tdim, Fdim>( with_grad ),
00049   camera_dof_(camera_dof),
00050   image_centre_(double(0)),
00051   centre_mag_norm_const_(0)
00052 {
00053   //modify the dof in vnl_least_squares_function
00054   vnl_least_squares_function::init(this->proj_size_-1+camera_dof_,
00055                                    this->get_number_of_residuals() );
00056 
00057   // temperary storage space
00058   temp_rad_k_.resize( camera_dof_ );
00059 }
00060 
00061 //: convert parameters
00062 template <unsigned int Tdim, unsigned int Fdim>
00063 void
00064 rgrl_est_proj_rad_func<Tdim, Fdim>::
00065 convert_parameters( vnl_vector<double>& params,
00066                     vnl_matrix_fixed<double, Tdim+1, Fdim+1>const&  proj_matrix,
00067                     vcl_vector<double>             const& rad_dist,
00068                     vnl_vector_fixed<double, Fdim> const& fc,
00069                     vnl_vector_fixed<double, Tdim> const& tc,
00070                     vnl_vector_fixed<double, Tdim> const& camera_centre )
00071 {
00072   // get params from the projection matrix
00073   vnl_vector<double> proj_params;
00074   rgrl_est_proj_func<Tdim, Fdim>::convert_parameters( proj_params, proj_matrix, fc, tc );
00075 
00076   // check camera dof
00077   assert( rad_dist.size() == camera_dof_ );
00078 
00079   // copy to new params
00080   params.set_size( proj_params.size() + rad_dist.size() );
00081   for ( unsigned i=0; i<proj_params.size(); ++i )
00082     params[i] = proj_params[i];
00083   for ( unsigned i=0; i<rad_dist.size(); ++i )
00084     params[i+proj_params.size()] = rad_dist[i];
00085 
00086   // set camera centre
00087   set_image_centre( camera_centre );
00088 }
00089 
00090 template <unsigned int Tdim, unsigned int Fdim>
00091 void
00092 rgrl_est_proj_rad_func<Tdim, Fdim>::
00093 transfer_radial_params_into_temp_storage( vnl_vector<double> const& params ) const
00094 {
00095   assert( params.size() + 1 == this->proj_size_ + camera_dof_ );
00096   const unsigned int index_shift = this->proj_size_-1;
00097   for ( unsigned i=0; i<camera_dof_; ++i )
00098     temp_rad_k_[i] = params[index_shift+i];
00099 }
00100 
00101 //: convert parameters
00102 template <unsigned int Tdim, unsigned int Fdim>
00103 inline
00104 void
00105 rgrl_est_proj_rad_func<Tdim, Fdim>::
00106 apply_radial_distortion( vnl_vector_fixed<double, Tdim>      & mapped,
00107                          vnl_vector_fixed<double, Tdim> const& p,
00108                          vcl_vector<double> const& radk ) const
00109 {
00110   assert( radk.size() == camera_dof_ );
00111 
00112   vnl_vector_fixed<double, Tdim> centred = p-image_centre_;
00113   const double radial_dist = centred.squared_magnitude() / centre_mag_norm_const_;
00114 
00115   double base = 1;
00116   double coeff = 0;
00117   for ( unsigned i=0; i<radk.size(); ++i ) {
00118     base *= radial_dist;
00119     coeff += radk[i] * base;
00120   }
00121   mapped = p + coeff*centred;
00122 }
00123 
00124 template <unsigned int Tdim, unsigned int Fdim>
00125 void
00126 rgrl_est_proj_rad_func<Tdim, Fdim>::
00127 reduced_proj_rad_jacobian( vnl_matrix<double>                            & base_jac,
00128                            vnl_matrix_fixed<double, Tdim+1, Fdim+1> const& proj,
00129                            vcl_vector<double>                       const& rad_k,
00130                            vnl_vector_fixed<double, Fdim>           const& from ) const
00131 {
00132   assert( rad_k.size() == camera_dof_ );
00133 
00134   const unsigned proj_dof = this->proj_size_ -1;
00135   const unsigned param_size = proj_dof + camera_dof_;
00136 
00137   // 0. set size
00138   base_jac.set_size( Tdim, param_size );
00139 
00140   // 1. get projection matrix jacobian
00141   vnl_matrix_fixed<double, Tdim, (Tdim+1)*(Fdim+1)-1> dP_dp;
00142   base_type::reduced_proj_jacobian( dP_dp, proj, from );
00143 
00144   // 2. gradient w.r.t to mapped location
00145   vnl_matrix_fixed<double, Tdim, Tdim >  dD_dx;
00146   vnl_vector_fixed<double, Tdim> mapped;
00147   rgrl_est_proj_map_inhomo_point<Tdim, Fdim>( mapped, proj, from-this->from_centre_ );
00148 
00149   const vnl_vector_fixed<double, Tdim> centred = mapped-image_centre_;
00150   const double radial_dist = centred.squared_magnitude() / centre_mag_norm_const_;
00151 
00152   // compute radial distortion coefficient
00153   double base = 1;
00154   double coeff = 0;
00155   for ( unsigned i=0; i<rad_k.size(); ++i ) {
00156     base *= radial_dist;
00157     coeff += rad_k[i] * base;
00158   }
00159 
00160   // two part computation for dD_dx
00161   // first part
00162   dD_dx.set_identity();
00163   dD_dx *= (1+coeff);
00164 
00165   // second part, taking gradient on the squared radial distance
00166   base = 2 / centre_mag_norm_const_;
00167   for ( unsigned k=0; k<rad_k.size(); ++k )
00168   {
00169     const double base_coeff = double(k+1)*base*rad_k[k];
00170 
00171     //upper triangular
00172     for ( unsigned i=0; i<Tdim; ++i )
00173       for ( unsigned j=i; j<Tdim; ++j ) {
00174         dD_dx( i, j ) += base_coeff*centred[i]*centred[j];
00175       }
00176 
00177     // multiplication is at the end of loop,
00178     // because in gradient it is to the power of (k-1), not k
00179     base *= radial_dist;
00180   }
00181 
00182   // fill in lower triangular
00183   for ( unsigned i=0; i<Tdim; ++i )
00184     for ( unsigned j=i; j<Tdim; ++j )
00185       dD_dx( j, i ) = dD_dx( i, j );
00186 
00187 
00188   // 3. fill in base_jac
00189   dP_dp = dD_dx * dP_dp;
00190 
00191   for ( unsigned i=0; i<Tdim; ++i )
00192     for ( unsigned j=0; j<dP_dp.cols(); ++j )
00193       base_jac( i, j ) = dP_dp( i, j );
00194 
00195   // 3. gradient w.r.t to k
00196   base = 1;
00197   for ( unsigned k=0; k<camera_dof_; ++k ) {
00198     const unsigned index = k+proj_dof;
00199     base *= radial_dist;
00200 
00201     for ( unsigned i=0; i<Tdim; ++i )
00202       base_jac( i, index ) = centred[i] * base;
00203   }
00204 }
00205 
00206 template <unsigned int Tdim, unsigned int Fdim>
00207 void
00208 rgrl_est_proj_rad_func<Tdim, Fdim>::
00209 full_proj_rad_jacobian( vnl_matrix<double>                            & base_jac,
00210                         vnl_matrix_fixed<double, Tdim+1, Fdim+1> const& proj,
00211                         vcl_vector<double>                       const& rad_k,
00212                         vnl_vector_fixed<double, Fdim>           const& from ) const
00213 {
00214   assert( rad_k.size() == camera_dof_ );
00215 
00216   const unsigned param_size = this->proj_size_ + camera_dof_;
00217 
00218   // 0. set size
00219   base_jac.set_size( Tdim, param_size );
00220 
00221   // 1. get projection matrix jacobian
00222   vnl_matrix_fixed<double, Tdim, (Tdim+1)*(Fdim+1)> dP_dp;
00223   base_type::full_proj_jacobian( dP_dp, proj, from );
00224 
00225   // 2. gradient w.r.t to mapped location
00226   vnl_matrix_fixed<double, Tdim, Tdim >  dD_dx;
00227   vnl_vector_fixed<double, Tdim> mapped;
00228   rgrl_est_proj_map_inhomo_point<Tdim, Fdim>( mapped, proj, from-this->from_centre_ );
00229 
00230   const vnl_vector_fixed<double, Tdim> centred = mapped-image_centre_;
00231   const double radial_dist = centred.squared_magnitude() / centre_mag_norm_const_;
00232 
00233   // compute radial distortion coefficient
00234   double base = 1;
00235   double coeff = 0;
00236   for ( unsigned i=0; i<rad_k.size(); ++i ) {
00237     base *= radial_dist;
00238     coeff += rad_k[i] * base;
00239   }
00240 
00241   // two part computation for dD_dx
00242   // first part
00243   dD_dx.set_identity();
00244   dD_dx *= (1+coeff);
00245 
00246   // second part, taking gradient on the squared radial distance
00247   base = 2 / centre_mag_norm_const_;
00248   for ( unsigned k=0; k<rad_k.size(); ++k )
00249   {
00250     const double base_coeff = double(k+1)*base*rad_k[k];
00251 
00252     //upper triangular
00253     for ( unsigned i=0; i<Tdim; ++i )
00254       for ( unsigned j=i; j<Tdim; ++j ) {
00255         dD_dx( i, j ) += base_coeff*centred[i]*centred[j];
00256       }
00257 
00258     // multiplication is at the end of loop,
00259     // because in gradient it is to the power of (k-1), not k
00260     base *= radial_dist;
00261   }
00262 
00263   // fill in lower triangular
00264   for ( unsigned i=0; i<Tdim; ++i )
00265     for ( unsigned j=i; j<Tdim; ++j )
00266       dD_dx( j, i ) = dD_dx( i, j );
00267 
00268 
00269   // 3. fill in base_jac
00270   dP_dp = dD_dx * dP_dp;
00271 
00272   for ( unsigned i=0; i<Tdim; ++i )
00273     for ( unsigned j=0; j<dP_dp.cols(); ++j )
00274       base_jac( i, j ) = dP_dp( i, j );
00275 
00276   // 3. gradient w.r.t to k
00277   base = 1;
00278   for ( unsigned k=0; k<camera_dof_; ++k ) {
00279     const unsigned index = k+this->proj_size_;
00280     base *= radial_dist;
00281 
00282     for ( unsigned i=0; i<Tdim; ++i )
00283       base_jac( i, index ) = centred[i] * base;
00284   }
00285 }
00286 
00287 //: compute jacobian w.r.t. location
00288 template <unsigned int Tdim, unsigned int Fdim>
00289 void
00290 rgrl_est_proj_rad_func<Tdim, Fdim>::
00291 proj_jac_wrt_loc( vnl_matrix_fixed<double, Tdim, Fdim>          & jac_loc,
00292                   vnl_matrix_fixed<double, Tdim+1, Fdim+1> const& proj,
00293                   vcl_vector<double>                       const& rad_k,
00294                   vnl_vector_fixed<double, Fdim>           const& from ) const
00295 {
00296   // dP / dx
00297   vnl_matrix_fixed<double, Tdim, Fdim> dP_dx;
00298   rgrl_est_proj_func<Tdim, Fdim>::proj_jac_wrt_loc( dP_dx, proj, from );
00299 
00300   // 2. gradient w.r.t to mapped location
00301   vnl_matrix_fixed<double, Tdim, Tdim >  dD_dx;
00302   vnl_vector_fixed<double, Tdim> mapped;
00303   rgrl_est_proj_map_inhomo_point<Tdim, Fdim>( mapped, proj, from-this->from_centre_ );
00304 
00305   const vnl_vector_fixed<double, Tdim> centred = mapped-image_centre_;
00306   const double radial_dist = centred.squared_magnitude() / centre_mag_norm_const_;
00307 
00308   // compute radial distortion coefficient
00309   double base = 1;
00310   double coeff = 0;
00311   for ( unsigned i=0; i<rad_k.size(); ++i ) {
00312     base *= radial_dist;
00313     coeff += rad_k[i] * base;
00314   }
00315 
00316   // two part computation for dD_dx
00317   // first part
00318   dD_dx.set_identity();
00319   dD_dx *= (1+coeff);
00320 
00321   // second part, taking gradient on the squared radial distance
00322   base = 2 / centre_mag_norm_const_;
00323   for ( unsigned k=0; k<rad_k.size(); ++k )
00324   {
00325     const double base_coeff = double(k+1)*base*rad_k[k];
00326 
00327     //upper triangular
00328     for ( unsigned i=0; i<Tdim; ++i )
00329       for ( unsigned j=i; j<Tdim; ++j ) {
00330         dD_dx( i, j ) += base_coeff*centred[i]*centred[j];
00331       }
00332 
00333     // multiplication is at the end of loop,
00334     // because in gradient it is to the power of (k-1), not k
00335     base *= radial_dist;
00336   }
00337 
00338   // fill in lower triangular
00339   for ( unsigned i=0; i<Tdim; ++i )
00340     for ( unsigned j=i; j<Tdim; ++j )
00341       dD_dx( j, i ) = dD_dx( i, j );
00342 
00343   // final jac
00344   jac_loc = dD_dx * dP_dx;
00345 }
00346 
00347 template <unsigned int Tdim, unsigned int Fdim>
00348 void
00349 rgrl_est_proj_rad_func<Tdim, Fdim>::
00350 f(vnl_vector<double> const& x, vnl_vector<double>& fx)
00351 {
00352   typedef rgrl_match_set::const_from_iterator FIter;
00353   typedef FIter::to_iterator TIter;
00354 
00355   vnl_vector_fixed<double, Tdim> distorted;
00356   vnl_matrix_fixed<double, Tdim,Tdim> error_proj_sqrt;
00357   unsigned int ind = 0;
00358 
00359   // retrieve the projection matrix
00360   vnl_matrix_fixed<double, Tdim+1, Fdim+1> proj;
00361   base_type::restored_centered_proj( proj, x );
00362 
00363   // retrieve the radial distortion parameters
00364   transfer_radial_params_into_temp_storage( x );
00365 
00366   for ( unsigned ms = 0; ms<this->matches_ptr_->size(); ++ms )
00367   {
00368     if ( (*this->matches_ptr_)[ms] != 0 ) // if pointer is valid
00369     {
00370       rgrl_match_set const& one_set = *((*this->matches_ptr_)[ms]);
00371       for ( FIter fi=one_set.from_begin(); fi!=one_set.from_end(); ++fi ) {
00372         // map from point
00373         vnl_vector_fixed<double, Fdim> from = fi.from_feature()->location();
00374         map_loc( distorted, proj, temp_rad_k_, from );
00375 
00376         for ( TIter ti=fi.begin(); ti!=fi.end(); ++ti ) {
00377           vnl_vector_fixed<double, Tdim> to = ti.to_feature()->location();
00378           error_proj_sqrt = ti.to_feature()->error_projector_sqrt();
00379           double const wgt = vcl_sqrt(ti.cumulative_weight());
00380           vnl_vector_fixed<double, Tdim> diff = error_proj_sqrt * (distorted - to);
00381 
00382           // fill in
00383           for ( unsigned i=0; i<Tdim; ++i )
00384             fx(ind++) = wgt * diff[i];
00385         }
00386       }
00387     }
00388   }
00389 
00390   // check
00391   assert( ind == this->get_number_of_residuals() );
00392 }
00393 
00394 template <unsigned int Tdim, unsigned int Fdim>
00395 void
00396 rgrl_est_proj_rad_func<Tdim, Fdim>::
00397 gradf(vnl_vector<double> const& x, vnl_matrix<double>& jacobian)
00398 {
00399   typedef rgrl_match_set::const_from_iterator FIter;
00400   typedef FIter::to_iterator TIter;
00401 
00402   assert( jacobian.rows() == this->get_number_of_residuals() );
00403   assert( jacobian.cols()+1 == this->proj_size_+camera_dof_ );
00404 
00405   vnl_matrix_fixed<double, Tdim,   Tdim>        error_proj_sqrt;
00406   vnl_matrix<double> base_jac, jac;
00407 
00408   // retrieve the projection matrix
00409   vnl_matrix_fixed<double, Tdim+1, Fdim+1> proj;
00410   base_type::restored_centered_proj( proj, x );
00411 
00412   // retrieve the radial distortion parameters
00413   transfer_radial_params_into_temp_storage( x );
00414 
00415   unsigned int ind = 0;
00416   for ( unsigned ms = 0; ms<this->matches_ptr_->size(); ++ms )
00417   {
00418     if ( (*this->matches_ptr_)[ms] ) // if pointer is valid
00419     {
00420       rgrl_match_set const& one_set = *((*this->matches_ptr_)[ms]);
00421       for ( FIter fi=one_set.from_begin(); fi!=one_set.from_end(); ++fi )
00422       {
00423         // map from point
00424         vnl_vector_fixed<double, Fdim> from = fi.from_feature()->location();
00425 
00426         // jacobian computation
00427         reduced_proj_rad_jacobian( base_jac, proj, temp_rad_k_, from );
00428 
00429         for ( TIter ti=fi.begin(); ti!=fi.end(); ++ti ) {
00430           //vnl_double_2 to = ti.to_feature()->location();
00431           error_proj_sqrt = ti.to_feature()->error_projector_sqrt();
00432           double const wgt = vcl_sqrt(ti.cumulative_weight());
00433           jac = error_proj_sqrt * base_jac;
00434           jac *= wgt;
00435 
00436           // fill in
00437           for ( unsigned i=0; i<Tdim; ++i,++ind ) {
00438             for ( unsigned j=0; j<jac.cols(); ++j )
00439               jacobian(ind,j) = jac(i, j);
00440           }
00441         }
00442       }
00443     }
00444   }
00445 }
00446 
00447 template <unsigned int Tdim, unsigned int Fdim>
00448 bool
00449 rgrl_est_proj_rad_func<Tdim, Fdim>::
00450 projective_estimate(  vnl_matrix_fixed<double, Tdim+1, Fdim+1>& proj,
00451                       vcl_vector<double>& rad_dist,
00452                       vnl_matrix<double>& full_covar,
00453                       vnl_vector_fixed<double, Fdim>& from_centre,
00454                       vnl_vector_fixed<double, Tdim>& to_centre,
00455                       vnl_vector_fixed<double, Tdim> const& camera_centre)
00456 {
00457   // compute weighted centres
00458   // this function is going to resize the vnl_vector, use temporary ones instead
00459   vnl_vector<double> fc, tc;
00460   if ( !rgrl_est_compute_weighted_centres( *(this->matches_ptr_), fc, tc ) ) {
00461     return false;
00462   }
00463   assert( fc.size() == from_centre.size() );
00464   assert( tc.size() == to_centre.size() );
00465   assert( rad_dist.size() == camera_dof_ );
00466   from_centre = fc;
00467   to_centre = tc;
00468 
00469   // convert parameters
00470   vnl_vector<double> p;
00471   this->convert_parameters( p, proj, rad_dist, from_centre, to_centre, camera_centre );
00472 
00473   vnl_levenberg_marquardt lm( *this );
00474   //lm.set_trace( true );
00475   //lm.set_check_derivatives( 1 );
00476   // we don't need it to be super accurate
00477   lm.set_f_tolerance( this->relative_threshold_ );
00478   lm.set_max_function_evals( this->max_num_iterations_ );
00479 
00480   bool ret;
00481   if ( this->has_gradient() )
00482     ret = lm.minimize_using_gradient(p);
00483   else
00484     ret = lm.minimize_without_gradient(p);
00485   if ( !ret ) {
00486     vcl_cerr <<  "Levenberg-Marquardt failed\n";
00487     lm.diagnose_outcome(vcl_cerr);
00488     return false;
00489   }
00490   // lm.diagnose_outcome(vcl_cout);
00491 
00492   // convert parameters back into matrix form
00493   this->restored_centered_proj( proj, p );
00494   //vcl_cout << "Final params=" << proj << vcl_endl;
00495 
00496   // copy distortion parameters
00497   const unsigned int index_shift = this->proj_size_-1;
00498   for ( unsigned i=0; i<camera_dof_; ++i )
00499     rad_dist[i] = p[ i+index_shift ];
00500 
00501   // compute covariance
00502   // Jac^\top * Jac is INVERSE of the covariance
00503   //
00504   const unsigned int param_num = this->proj_size_+camera_dof_;
00505   vnl_matrix<double> jac( this->get_number_of_residuals(), param_num-1 );
00506   this->gradf( p, jac );
00507   //
00508 
00509   // SVD decomposition:
00510   // Jac = U W V^\top
00511   vnl_svd<double> svd( jac, this->zero_svd_thres_ );
00512   if ( svd.rank()+1 < param_num ) {
00513     vcl_cerr <<  "The covariance of projection matrix ranks less than "
00514              << param_num-1 << "!\n"
00515              << "  The singular values are " << svd.W() << vcl_endl;
00516     return false;
00517   }
00518 
00519   //invert the W matrix and square it
00520   vnl_diag_matrix<double> invW( param_num-1 );
00521   for ( unsigned i=0; i+1<param_num; ++i )
00522     invW[i] = vnl_math_sqr( 1.0/svd.W(i) );
00523 
00524   //compute inverse of Jac^\top * Jac
00525   const vnl_matrix<double>  covar( svd.V() * invW * vnl_transpose( svd.V() ) );
00526 
00527   // convert the covar to full dof+1 matrix
00528   // fill in the row and column of the fixed element
00529   // with 0
00530   full_covar.set_size( param_num, param_num );
00531   full_covar.fill( 0.0 );
00532 
00533   const unsigned int param_index
00534     = this->index_row_*(Fdim+1)
00535     + this->index_col_;
00536   for ( unsigned i=0,ic=0; i<param_num; ++i,++ic ) {
00537     if ( i==param_index ) { --ic; continue; }
00538     for ( unsigned j=0,jc=0; j<param_num; ++j,++jc ) {
00539       if ( j==param_index ) { --jc; continue; }
00540       full_covar( i, j ) = covar( ic, jc );
00541     }
00542   }
00543 
00544   return true;
00545 }
00546 
00547 
00548 #undef  RGRL_EST_PROJ_RAD_FUNC_INSTANTIATE
00549 #define RGRL_EST_PROJ_RAD_FUNC_INSTANTIATE( tdim, fdim ) \
00550   template class rgrl_est_proj_rad_func< tdim, fdim >
00551 
00552 #endif //rgrl_est_proj_rad_func_txx_