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