Public Types | Public Member Functions | Static Public Member Functions | Public Attributes | Protected Attributes
vpgl_ba_fixed_k_lsqr Class Reference

a class for bundle adjustment with fixed intrinsic parameters. More...

#include <vpgl_ba_fixed_k_lsqr.h>

Inheritance diagram for vpgl_ba_fixed_k_lsqr:
Inheritance graph
[legend]

List of all members.

Public Types

enum  UseGradient
enum  UseWeights

Public Member Functions

 vpgl_ba_fixed_k_lsqr (const vcl_vector< vpgl_calibration_matrix< double > > &K, const vcl_vector< vgl_point_2d< double > > &image_points, const vcl_vector< vcl_vector< bool > > &mask)
 Constructor.
 vpgl_ba_fixed_k_lsqr (const vcl_vector< vpgl_calibration_matrix< double > > &K, const vcl_vector< vgl_point_2d< double > > &image_points, const vcl_vector< vnl_matrix< double > > &inv_covars, const vcl_vector< vcl_vector< bool > > &mask)
 Constructor.
virtual ~vpgl_ba_fixed_k_lsqr ()
virtual void jac_Aij (unsigned int i, unsigned int j, vnl_double_3x4 const &Pi, vnl_vector< double > const &ai, vnl_vector< double > const &bj, vnl_vector< double > const &c, vnl_matrix< double > &Aij)
 compute the Jacobian Aij.
virtual void jac_Bij (unsigned int i, unsigned int j, vnl_double_3x4 const &Pi, vnl_vector< double > const &ai, vnl_vector< double > const &bj, vnl_vector< double > const &c, vnl_matrix< double > &Bij)
 compute the Jacobian Bij.
virtual void jac_Cij (unsigned int i, unsigned int j, vnl_double_3x4 const &Pi, vnl_vector< double > const &ai, vnl_vector< double > const &bj, vnl_vector< double > const &c, vnl_matrix< double > &Cij)
 compute the Jacobian Cij.
virtual vgl_homg_point_3d< double > param_to_point (int j, const double *bj, const vnl_vector< double > &c) const
 construct the
virtual vnl_vector_fixed
< double, 4 > 
param_to_pt_vector (int j, const double *bj, const vnl_vector< double > &c) const
 construct the
virtual
vpgl_perspective_camera
< double > 
param_to_cam (int i, const double *ai, const vnl_vector< double > &c) const
 construct the
virtual vnl_double_3x4 param_to_cam_matrix (int i, const double *ai, const vnl_vector< double > &c) const
 compute a 3x4 camera matrix of camera
virtual void f (vnl_vector< double > const &a, vnl_vector< double > const &b, vnl_vector< double > const &c, vnl_vector< double > &e)
 Compute all the reprojection errors.
virtual void fij (int i, int j, vnl_vector< double > const &ai, vnl_vector< double > const &bj, vnl_vector< double > const &c, vnl_vector< double > &fij)
 Compute the residuals from the ith component of a, the jth component of b, and all of c.
virtual void jac_blocks (vnl_vector< double > const &a, vnl_vector< double > const &b, vnl_vector< double > const &c, vcl_vector< vnl_matrix< double > > &A, vcl_vector< vnl_matrix< double > > &B, vcl_vector< vnl_matrix< double > > &C)
 Compute the sparse Jacobian in block form.
virtual void jac_Aij (int i, int j, vnl_vector< double > const &ai, vnl_vector< double > const &bj, vnl_vector< double > const &c, vnl_matrix< double > &Aij)
virtual void jac_Bij (int i, int j, vnl_vector< double > const &ai, vnl_vector< double > const &bj, vnl_vector< double > const &c, vnl_matrix< double > &Bij)
virtual void jac_Cij (int i, int j, vnl_vector< double > const &ai, vnl_vector< double > const &bj, vnl_vector< double > const &c, vnl_matrix< double > &Cij)
void compute_weight_ij (int i, int j, vnl_vector< double > const &ai, vnl_vector< double > const &bj, vnl_vector< double > const &c, vnl_vector< double > const &fij, double &weight)
 Use an M-estimator to compute weights.
void set_residual_scale (double scale)
 set the residual scale for the robust estimation.
vgl_homg_point_3d< double > param_to_point (int j, const vnl_vector< double > &b, const vnl_vector< double > &c) const
 construct the
vnl_vector_fixed< double, 4 > param_to_pt_vector (int j, const vnl_vector< double > &b, const vnl_vector< double > &c) const
 construct the
vpgl_perspective_camera< double > param_to_cam (int i, const vnl_vector< double > &a, const vnl_vector< double > &c) const
 construct the
vnl_double_3x4 param_to_cam_matrix (int i, const vnl_vector< double > &a, const vnl_vector< double > &c) const
 construct the
void reset ()
 reset the weights.
void throw_failure ()
void clear_failure ()
virtual void fd_jac_blocks (vnl_vector< double > const &a, vnl_vector< double > const &b, vnl_vector< double > const &c, vcl_vector< vnl_matrix< double > > &A, vcl_vector< vnl_matrix< double > > &B, vcl_vector< vnl_matrix< double > > &C, double stepsize)
virtual void compute_weights (vnl_vector< double > const &a, vnl_vector< double > const &b, vnl_vector< double > const &c, vnl_vector< double > const &f, vnl_vector< double > &weights)
virtual void apply_weights (vnl_vector< double > const &weights, vnl_vector< double > &f)
virtual void apply_weights (vnl_vector< double > const &weights, vcl_vector< vnl_matrix< double > > &A, vcl_vector< vnl_matrix< double > > &B, vcl_vector< vnl_matrix< double > > &C)
void fd_jac_Aij (int i, int j, vnl_vector< double > const &ai, vnl_vector< double > const &bj, vnl_vector< double > const &c, vnl_matrix< double > &Aij, double stepsize)
void fd_jac_Bij (int i, int j, vnl_vector< double > const &ai, vnl_vector< double > const &bj, vnl_vector< double > const &c, vnl_matrix< double > &Bij, double stepsize)
void fd_jac_Cij (int i, int j, vnl_vector< double > const &ai, vnl_vector< double > const &bj, vnl_vector< double > const &c, vnl_matrix< double > &Cij, double stepsize)
virtual void apply_weight_ij (int i, int j, double const &weight, vnl_vector< double > &fij)
virtual void apply_weight_ij (int i, int j, double const &weight, vnl_matrix< double > &Aij, vnl_matrix< double > &Bij, vnl_matrix< double > &Cij)
virtual void trace (int iteration, vnl_vector< double > const &a, vnl_vector< double > const &b, vnl_vector< double > const &c, vnl_vector< double > const &e)
unsigned int number_of_params_a (int i) const
unsigned int number_of_params_b (int j) const
unsigned int number_of_params_c () const
unsigned int number_of_residuals (int k) const
unsigned int number_of_residuals (int i, int j) const
unsigned int index_a (int i) const
unsigned int index_b (int j) const
unsigned int index_e (int k) const
unsigned int number_of_a () const
unsigned int number_of_b () const
unsigned int number_of_e () const
bool has_gradient () const
bool has_weights () const
const vnl_crs_indexresidual_indices () const

Static Public Member Functions

static vnl_vector< double > create_param_vector (const vcl_vector< vpgl_perspective_camera< double > > &cameras)
 Create the parameter vector a from a vector of cameras.
static vnl_vector< double > create_param_vector (const vcl_vector< vgl_point_3d< double > > &world_points)
 Create the parameter vector b from a vector of 3D points.
static vnl_double_3x3 rod_to_matrix (vnl_vector< double > const &r)
 Fast conversion of rotation from Rodrigues vector to matrix.
static void jac_inhomg_3d_point (vnl_double_3x4 const &P, vnl_vector< double > const &pt, vnl_matrix< double > &J)
 compute the 2x3 Jacobian of camera projection with respect to point location df/dpt where $f(pt) = P*pt$.
static void jac_camera_center (vnl_double_3x3 const &M, vnl_vector< double > const &C, vnl_vector< double > const &pt, vnl_matrix< double > &J)
 compute the 2x3 Jacobian of camera projection with respect to camera center df/dC where $f(C) = [M | -M*C]*pt$.
static void jac_camera_rotation (vnl_double_3x3 const &K, vnl_vector< double > const &C, vnl_vector< double > const &r, vnl_vector< double > const &pt, vnl_matrix< double > &J)
 compute the 2x3 Jacobian of camera projection with respect to camera rotation df/dr where $f(r) = K*rod_to_matrix(r)*[I | -C]*pt$.

Public Attributes

 no_gradient
 use_gradient
 no_weights
 use_weights
bool failure

Protected Attributes

vcl_vector
< vpgl_calibration_matrix
< double > > 
K_
 The fixed internal camera calibration.
vcl_vector< vnl_double_3x3Km_
 The fixed internal camera calibration in matrix form.
vcl_vector< vgl_point_2d
< double > > 
image_points_
 The corresponding points in the image.
vcl_vector< vnl_matrix< double > > factored_inv_covars_
 The Cholesky factored inverse covariances for each image point.
bool use_covars_
 Flag to enable covariance weighted errors.
double scale2_
 The square of the scale of the robust estimator.
int iteration_count_
vnl_crs_index residual_indices_
vcl_vector< unsigned int > indices_a_
vcl_vector< unsigned int > indices_b_
unsigned int num_params_c_
vcl_vector< unsigned int > indices_e_
bool use_gradient_
bool use_weights_

Detailed Description

a class for bundle adjustment with fixed intrinsic parameters.

Definition at line 18 of file vpgl_ba_fixed_k_lsqr.h.


Constructor & Destructor Documentation

vpgl_ba_fixed_k_lsqr::vpgl_ba_fixed_k_lsqr ( const vcl_vector< vpgl_calibration_matrix< double > > &  K,
const vcl_vector< vgl_point_2d< double > > &  image_points,
const vcl_vector< vcl_vector< bool > > &  mask 
)

Constructor.

Note:
image points are not homogeneous because they require finite points to measure projection error

Definition at line 12 of file vpgl_ba_fixed_k_lsqr.cxx.

vpgl_ba_fixed_k_lsqr::vpgl_ba_fixed_k_lsqr ( const vcl_vector< vpgl_calibration_matrix< double > > &  K,
const vcl_vector< vgl_point_2d< double > > &  image_points,
const vcl_vector< vnl_matrix< double > > &  inv_covars,
const vcl_vector< vcl_vector< bool > > &  mask 
)

Constructor.

Each image point is assigned an inverse covariance (error projector) matrix

Note:
image points are not homogeneous because they require finite points to measure projection error

Definition at line 27 of file vpgl_ba_fixed_k_lsqr.cxx.

virtual vpgl_ba_fixed_k_lsqr::~vpgl_ba_fixed_k_lsqr ( ) [inline, virtual]

Definition at line 38 of file vpgl_ba_fixed_k_lsqr.h.


Member Function Documentation

void vpgl_bundle_adjust_lsqr::compute_weight_ij ( int  i,
int  j,
vnl_vector< double > const &  ai,
vnl_vector< double > const &  bj,
vnl_vector< double > const &  c,
vnl_vector< double > const &  fij,
double &  weight 
) [virtual, inherited]

Use an M-estimator to compute weights.

Reimplemented from vnl_sparse_lst_sqr_function.

Definition at line 198 of file vpgl_bundle_adjust_lsqr.cxx.

vnl_vector< double > vpgl_ba_fixed_k_lsqr::create_param_vector ( const vcl_vector< vpgl_perspective_camera< double > > &  cameras) [static]

Create the parameter vector a from a vector of cameras.

Definition at line 136 of file vpgl_ba_fixed_k_lsqr.cxx.

vnl_vector< double > vpgl_ba_fixed_k_lsqr::create_param_vector ( const vcl_vector< vgl_point_3d< double > > &  world_points) [static]

Create the parameter vector b from a vector of 3D points.

Definition at line 158 of file vpgl_ba_fixed_k_lsqr.cxx.

void vpgl_bundle_adjust_lsqr::f ( vnl_vector< double > const &  a,
vnl_vector< double > const &  b,
vnl_vector< double > const &  c,
vnl_vector< double > &  e 
) [virtual, inherited]

Compute all the reprojection errors.

Given the parameter vectors a, b, and c, compute the vector of residuals e. e has been sized appropriately before the call. The parameters in for each camera are in a. The 3D point parameters are in b. Global parameters are in c.

Given the parameter vectors a, b, and c, compute the vector of residuals e. e has been sized appropriately before the call. The parameters in a for each camera are {wx, wy, wz, tx, ty, tz} where w is the Rodrigues vector of the rotation and t is the translation. The parameters in b for each 3D point are {px, py, pz} the non-homogeneous position.

Construct the ith camera.

Reimplemented from vnl_sparse_lst_sqr_function.

Definition at line 85 of file vpgl_bundle_adjust_lsqr.cxx.

void vpgl_bundle_adjust_lsqr::fij ( int  i,
int  j,
vnl_vector< double > const &  ai,
vnl_vector< double > const &  bj,
vnl_vector< double > const &  c,
vnl_vector< double > &  fij 
) [virtual, inherited]

Compute the residuals from the ith component of a, the jth component of b, and all of c.

Compute the residuals from the ith component of a and the jth component of b.

This function is not normally used because f() has a self-contained efficient implementation. It is used for finite-differencing if the gradient is marked as unavailable

This is not normally used because f() has a self-contained efficient implementation It is used for finite-differencing if the gradient is marked as unavailable

Construct the ith camera.

Reimplemented from vnl_sparse_lst_sqr_function.

Definition at line 127 of file vpgl_bundle_adjust_lsqr.cxx.

void vpgl_ba_fixed_k_lsqr::jac_Aij ( unsigned int  i,
unsigned int  j,
vnl_double_3x4 const &  Pi,
vnl_vector< double > const &  ai,
vnl_vector< double > const &  bj,
vnl_vector< double > const &  c,
vnl_matrix< double > &  Aij 
) [virtual]

compute the Jacobian Aij.

Implements vpgl_bundle_adjust_lsqr.

Definition at line 40 of file vpgl_ba_fixed_k_lsqr.cxx.

void vpgl_ba_fixed_k_lsqr::jac_Bij ( unsigned int  i,
unsigned int  j,
vnl_double_3x4 const &  Pi,
vnl_vector< double > const &  ai,
vnl_vector< double > const &  bj,
vnl_vector< double > const &  c,
vnl_matrix< double > &  Bij 
) [virtual]

compute the Jacobian Bij.

Implements vpgl_bundle_adjust_lsqr.

Definition at line 66 of file vpgl_ba_fixed_k_lsqr.cxx.

void vpgl_bundle_adjust_lsqr::jac_blocks ( vnl_vector< double > const &  a,
vnl_vector< double > const &  b,
vnl_vector< double > const &  c,
vcl_vector< vnl_matrix< double > > &  A,
vcl_vector< vnl_matrix< double > > &  B,
vcl_vector< vnl_matrix< double > > &  C 
) [virtual, inherited]

Compute the sparse Jacobian in block form.

Construct the ith camera.

Reimplemented from vnl_sparse_lst_sqr_function.

Definition at line 157 of file vpgl_bundle_adjust_lsqr.cxx.

void vpgl_bundle_adjust_lsqr::jac_camera_center ( vnl_double_3x3 const &  M,
vnl_vector< double > const &  C,
vnl_vector< double > const &  pt,
vnl_matrix< double > &  J 
) [static, inherited]

compute the 2x3 Jacobian of camera projection with respect to camera center df/dC where $f(C) = [M | -M*C]*pt$.

Definition at line 249 of file vpgl_bundle_adjust_lsqr.cxx.

void vpgl_bundle_adjust_lsqr::jac_camera_rotation ( vnl_double_3x3 const &  K,
vnl_vector< double > const &  C,
vnl_vector< double > const &  r,
vnl_vector< double > const &  pt,
vnl_matrix< double > &  J 
) [static, inherited]

compute the 2x3 Jacobian of camera projection with respect to camera rotation df/dr where $f(r) = K*rod_to_matrix(r)*[I | -C]*pt$.

Here r is a Rodrigues vector, K is an upper triangular calibration matrix

Definition at line 265 of file vpgl_bundle_adjust_lsqr.cxx.

void vpgl_ba_fixed_k_lsqr::jac_Cij ( unsigned int  i,
unsigned int  j,
vnl_double_3x4 const &  Pi,
vnl_vector< double > const &  ai,
vnl_vector< double > const &  bj,
vnl_vector< double > const &  c,
vnl_matrix< double > &  Cij 
) [virtual]

compute the Jacobian Cij.

Implements vpgl_bundle_adjust_lsqr.

Definition at line 78 of file vpgl_ba_fixed_k_lsqr.cxx.

void vpgl_bundle_adjust_lsqr::jac_inhomg_3d_point ( vnl_double_3x4 const &  P,
vnl_vector< double > const &  pt,
vnl_matrix< double > &  J 
) [static, inherited]

compute the 2x3 Jacobian of camera projection with respect to point location df/dpt where $f(pt) = P*pt$.

Definition at line 217 of file vpgl_bundle_adjust_lsqr.cxx.

vpgl_perspective_camera< double > vpgl_ba_fixed_k_lsqr::param_to_cam ( int  i,
const double *  ai,
const vnl_vector< double > &  c 
) const [virtual]

construct the

Parameters:
i-thperspective camera from a pointer to the i-th parameters of
aand parameters
c.

Implements vpgl_bundle_adjust_lsqr.

Definition at line 109 of file vpgl_ba_fixed_k_lsqr.cxx.

vpgl_perspective_camera<double> vpgl_bundle_adjust_lsqr::param_to_cam ( int  i,
const vnl_vector< double > &  a,
const vnl_vector< double > &  c 
) const [inline, inherited]

construct the

Parameters:
i-thperspective camera from parameter vector
a.

Definition at line 153 of file vpgl_bundle_adjust_lsqr.h.

vnl_double_3x4 vpgl_ba_fixed_k_lsqr::param_to_cam_matrix ( int  i,
const double *  ai,
const vnl_vector< double > &  c 
) const [virtual]

compute a 3x4 camera matrix of camera

Parameters:
ifrom a pointer to the i-th parameters of
aand parameters
c.

Implements vpgl_bundle_adjust_lsqr.

Definition at line 120 of file vpgl_ba_fixed_k_lsqr.cxx.

vnl_double_3x4 vpgl_bundle_adjust_lsqr::param_to_cam_matrix ( int  i,
const vnl_vector< double > &  a,
const vnl_vector< double > &  c 
) const [inline, inherited]

construct the

Parameters:
i-thperspective camera matrix from parameter vectors
aand
c.

Definition at line 167 of file vpgl_bundle_adjust_lsqr.h.

vgl_homg_point_3d< double > vpgl_ba_fixed_k_lsqr::param_to_point ( int  j,
const double *  bj,
const vnl_vector< double > &  c 
) const [virtual]

construct the

Parameters:
j-thperspective camera from a pointer to the j-th parameters of
band parameters
c.

Implements vpgl_bundle_adjust_lsqr.

Definition at line 91 of file vpgl_ba_fixed_k_lsqr.cxx.

vgl_homg_point_3d<double> vpgl_bundle_adjust_lsqr::param_to_point ( int  j,
const vnl_vector< double > &  b,
const vnl_vector< double > &  c 
) const [inline, inherited]

construct the

Parameters:
j-th3D point from parameter vector
band
c.

Definition at line 123 of file vpgl_bundle_adjust_lsqr.h.

vnl_vector_fixed< double, 4 > vpgl_ba_fixed_k_lsqr::param_to_pt_vector ( int  j,
const double *  bj,
const vnl_vector< double > &  c 
) const [virtual]

construct the

Parameters:
j-thperspective camera from a pointer to the j-th parameters of
band parameters
c.

Implements vpgl_bundle_adjust_lsqr.

Definition at line 100 of file vpgl_ba_fixed_k_lsqr.cxx.

vnl_vector_fixed<double,4> vpgl_bundle_adjust_lsqr::param_to_pt_vector ( int  j,
const vnl_vector< double > &  b,
const vnl_vector< double > &  c 
) const [inline, inherited]

construct the

Parameters:
j-th3D point from parameter vector
band
c.

Definition at line 138 of file vpgl_bundle_adjust_lsqr.h.

void vpgl_bundle_adjust_lsqr::reset ( ) [inline, inherited]

reset the weights.

Definition at line 180 of file vpgl_bundle_adjust_lsqr.h.

vnl_double_3x3 vpgl_bundle_adjust_lsqr::rod_to_matrix ( vnl_vector< double > const &  r) [static, inherited]

Fast conversion of rotation from Rodrigues vector to matrix.

Definition at line 359 of file vpgl_bundle_adjust_lsqr.cxx.

void vpgl_bundle_adjust_lsqr::set_residual_scale ( double  scale) [inline, inherited]

set the residual scale for the robust estimation.

Definition at line 118 of file vpgl_bundle_adjust_lsqr.h.


Member Data Documentation

vcl_vector<vnl_matrix<double> > vpgl_bundle_adjust_lsqr::factored_inv_covars_ [protected, inherited]

The Cholesky factored inverse covariances for each image point.

Definition at line 216 of file vpgl_bundle_adjust_lsqr.h.

vcl_vector<vgl_point_2d<double> > vpgl_bundle_adjust_lsqr::image_points_ [protected, inherited]

The corresponding points in the image.

Definition at line 214 of file vpgl_bundle_adjust_lsqr.h.

int vpgl_bundle_adjust_lsqr::iteration_count_ [protected, inherited]

Definition at line 222 of file vpgl_bundle_adjust_lsqr.h.

vcl_vector<vpgl_calibration_matrix<double> > vpgl_ba_fixed_k_lsqr::K_ [protected]

The fixed internal camera calibration.

Definition at line 105 of file vpgl_ba_fixed_k_lsqr.h.

vcl_vector<vnl_double_3x3> vpgl_ba_fixed_k_lsqr::Km_ [protected]

The fixed internal camera calibration in matrix form.

Definition at line 107 of file vpgl_ba_fixed_k_lsqr.h.

double vpgl_bundle_adjust_lsqr::scale2_ [protected, inherited]

The square of the scale of the robust estimator.

Definition at line 220 of file vpgl_bundle_adjust_lsqr.h.

bool vpgl_bundle_adjust_lsqr::use_covars_ [protected, inherited]

Flag to enable covariance weighted errors.

Definition at line 218 of file vpgl_bundle_adjust_lsqr.h.


The documentation for this class was generated from the following files: