Main Page | Modules | Class List | File List | Class Members | File Members

Levenberg-Marquardt algorithm
[Vision Package]


Classes

struct  Gan_LevMarqObs
 Observation structure for Levenberg-Marquardt minimisation. More...
struct  Gan_LevMarqStruct
 Structure for holding state of Levenberg-Marquardt algorithm. More...

Typedefs

typedef Gan_Bool(* Gan_LevMarqObsFunc_h )(Gan_Vector *x, Gan_Vector *z, void *zdata, Gan_Vector *h, Gan_Matrix *H)
 Observation function for standard h-type observations.
typedef Gan_Bool(* Gan_LevMarqObsFunc_F )(Gan_Vector *x, Gan_Vector *z, void *zdata, Gan_Vector *F, Gan_Matrix *Hx, Gan_Matrix *Hz)
 Observation function for F-type observations.
typedef Gan_LevMarqObs Gan_LevMarqObs
 Observation structure for Levenberg-Marquardt minimisation.
typedef Gan_LevMarqStruct Gan_LevMarqStruct
 Structure for holding state of Levenberg-Marquardt algorithm.
typedef Gan_Bool(* Gan_LevMarqInitFunc )(Gan_Vector *x, Gan_List *obs_list, void *data)
 Callback function for initialising Levenberg-Marquardt algorithm.

Enumerations

enum  Gan_LevMarqObsType { GAN_LEV_MARQ_OBS_H, GAN_LEV_MARQ_OBS_H_ROBUST, GAN_LEV_MARQ_OBS_F }
 Observation type for Levenberg-Marquardt minimisation. More...

Functions

Gan_LevMarqStructgan_lev_marq_form (Gan_LevMarqStruct *lm)
 Forms a Levenberg-Marquardt structure.
Gan_LevMarqObsgan_lev_marq_obs_h (Gan_LevMarqStruct *lm, Gan_Vector *z, void *zdata, Gan_SquMatrix *Ni, Gan_LevMarqObsFunc_h obs_func)
 Passes an observation to a Levenberg-Marquardt structure.
Gan_LevMarqObsgan_lev_marq_obs_h_robust (Gan_LevMarqStruct *lm, Gan_Vector *z, void *zdata, Gan_SquMatrix *Ni, Gan_LevMarqObsFunc_h obs_func, double var_scale, double chi2)
 Passes a robust observation to a Levenberg-Marquardt structure.
Gan_LevMarqObsgan_lev_marq_obs_F (Gan_LevMarqStruct *lm, Gan_Vector *z, void *zdata, Gan_SquMatrix *Ni, Gan_LevMarqObsFunc_F obs_func)
 Passes an observation to a Levenberg-Marquardt structure.
Gan_Bool gan_lev_marq_init (Gan_LevMarqStruct *lm, Gan_LevMarqInitFunc init_func, void *data, double *residualp)
 Initialise Levenberg-Marquardt algorithm.
Gan_Bool gan_lev_marq_iteration (Gan_LevMarqStruct *lm, double lambda, double *residualp)
 Applies Levenberg-Marquardt iteration.
Gan_Vectorgan_lev_marq_get_x (Gan_LevMarqStruct *lm)
 Returns state of Levenberg-Marquardt minimisation.
Gan_SquMatrixgan_lev_marq_get_P (Gan_LevMarqStruct *lm)
 Returns state covariance of Levenberg-Marquardt minimisation.
void gan_lev_marq_free (Gan_LevMarqStruct *lm)
 Frees a Levenberg-Marquardt structure.
Gan_LevMarqStructgan_lev_marq_alloc (void)
 Macro: Allocates a Levenberg-Marquardt structure.

Typedef Documentation

typedef Gan_Bool(* Gan_LevMarqObsFunc_F)(Gan_Vector *x, Gan_Vector *z, void *zdata, Gan_Vector *F, Gan_Matrix *Hx, Gan_Matrix *Hz)
 

Observation function for F-type observations.

Observation function for F-type observations

\[ F(x,z+w) = 0 \]

typedef Gan_Bool(* Gan_LevMarqObsFunc_h)(Gan_Vector *x, Gan_Vector *z, void *zdata, Gan_Vector *h, Gan_Matrix *H)
 

Observation function for standard h-type observations.

Observation function for standard h-type observations

\[ z = h(x) + w \]


Enumeration Type Documentation

enum Gan_LevMarqObsType
 

Observation type for Levenberg-Marquardt minimisation.

Enumeration values:
GAN_LEV_MARQ_OBS_H  h-type observations
GAN_LEV_MARQ_OBS_H_ROBUST  robust h-type observations
GAN_LEV_MARQ_OBS_F  F-type observations


Function Documentation

Gan_LevMarqStruct* gan_lev_marq_alloc void   ) 
 

Macro: Allocates a Levenberg-Marquardt structure.

Allocates a structure for computing Levenberg-Marquardt optimisation.

Returns:
Pointer to successfully allocated structure or NULL.
See also:
gan_lev_marq_form().

Gan_LevMarqStruct * gan_lev_marq_form Gan_LevMarqStruct lm  ) 
 

Forms a Levenberg-Marquardt structure.

Parameters:
lm Pointer to a Levenberg-Marquardt structure or NULL
Allocates/fills a structure for computing Levenberg-Marquardt optimisation.
Returns:
Pointer to structure or NULL on failure.
See also:
gan_lev_marq_alloc().

void gan_lev_marq_free Gan_LevMarqStruct lm  ) 
 

Frees a Levenberg-Marquardt structure.

Parameters:
lm Pointer to Levenberg-Marquardt structure
Frees a structure created to compute Levenberg-Marquardt minimisation.
Returns:
None
See also:
gan_lev_marq_form().

Gan_SquMatrix * gan_lev_marq_get_P Gan_LevMarqStruct lm  ) 
 

Returns state covariance of Levenberg-Marquardt minimisation.

Parameters:
lm Pointer to a Levenberg-Marquardt structure
Returns the state covariance of a Levenberg-Marquardt minimisation. The matrix is not copied, but is rather passed back directly from the state structure. There is therefore no need to free the matrix after use. It is essential also that the contents of the matrix are not modified; i.e. the matrix should be treated as read only.

Returns:
Pointer to current solution covariance or NULL on failure.
See also:
gan_lev_marq_get_x().

Gan_Vector * gan_lev_marq_get_x Gan_LevMarqStruct lm  ) 
 

Returns state of Levenberg-Marquardt minimisation.

Parameters:
lm Pointer to a Levenberg-Marquardt structure
Returns the state of a Levenberg-Marquardt minimisation. The vector is not copied, but is rather passed back directly from the state structure. There is therefore no need to free the vector after use. It is essential also that the contents of the vector are not modified; i.e. the vector should be treated as read only.

Returns:
Pointer to current solution vector or NULL on failure.
See also:
gan_lev_marq_get_P().

Gan_Bool gan_lev_marq_init Gan_LevMarqStruct lm,
Gan_LevMarqInitFunc  init_func,
void *  data,
double *  residualp
 

Initialise Levenberg-Marquardt algorithm.

Parameters:
lm Pointer to a Levenberg-Marquardt structure
init_func Function to initialise value of state vector
data User pointer to be passed into init_func
residualp Pointer to initial residual
After passing some observations to the structure using gan_lev_marq_obs_h() or gan_lev_marq_obs_F(), this function invokes the callback function init_func to set the state to initial values, before applying Levenberg-Marquardt iterations using gan_lev_marq_iteration().

Returns:
GAN_TRUE on successful initialisation of structure, GAN_FALSE on failure.
See also:
gan_lev_marq_iteration().

Gan_Bool gan_lev_marq_iteration Gan_LevMarqStruct lm,
double  lambda,
double *  residual
 

Applies Levenberg-Marquardt iteration.

Parameters:
lm Pointer to a Levenberg-Marquardt structure
lambda Adjustment factor for adjusting the information matrix
residualp Output residual after iteration
Applies a single Levenberg-Marquardt iteration. The residual after the adjustment is returned in the provided pointer. The adjustment is only actually applied to the state vector $x$ if the residual has been reduced from the existing value.

Returns:
GAN_TRUE on successful completion of iteration, GAN_FALSE on failure.

Gan_LevMarqObs * gan_lev_marq_obs_F Gan_LevMarqStruct lm,
Gan_Vector z,
void *  zdata,
Gan_SquMatrix Ni,
Gan_LevMarqObsFunc_F  obs_func
 

Passes an observation to a Levenberg-Marquardt structure.

Parameters:
lm Pointer to a Levenberg-Marquardt structure
z Observation vector
zdata Data pointer to pass along with z
Ni Observation inverse covariance
obs_func Observation callback function
Adds a new observation for Levenberg-Marquardt minimisation. It relates state to observation according to the model $F(x,z+w) = 0$ where $z$ is the observation vector, $F()$ is the observation function, $x$ is the state vector and $w$ is a zero-mean noise vector, with covariance $N$.
Returns:
Handle for the new observation or NULL on failure.
See also:
gan_lev_marq_form().

Gan_LevMarqObs * gan_lev_marq_obs_h Gan_LevMarqStruct lm,
Gan_Vector z,
void *  zdata,
Gan_SquMatrix Ni,
Gan_LevMarqObsFunc_h  obs_func
 

Passes an observation to a Levenberg-Marquardt structure.

Parameters:
lm Pointer to a structure
z Observation vector
zdata Data pointer to pass along with z
Ni Observation inverse covariance
obs_func Observation callback function
Adds a new observation for Levenberg-Marquardt minimisation. The observation is a direct observation relating state to observation $z = h(x) + w$ where $z$ is the observation vector, $h()$ is the observation function, $x$ is the state vector and $w$ is a zero-mean noise vector, with covariance $N$.

Returns:
Handle for the new observation or NULL on failure.
See also:
gan_lev_marq_form().

Gan_LevMarqObs * gan_lev_marq_obs_h_robust Gan_LevMarqStruct lm,
Gan_Vector z,
void *  zdata,
Gan_SquMatrix Ni,
Gan_LevMarqObsFunc_h  obs_func,
double  var_scale,
double  chi2
 

Passes a robust observation to a Levenberg-Marquardt structure.

Parameters:
lm Pointer to a structure
z Observation vector
zdata Data pointer to pass along with z
Ni Observation inverse covariance
obs_func Observation callback function
var_scale Scaling of outlier covariance w.r.t. inlier covariance
chi2 Cutoff point on innovation for switching distributions
Adds a new observation for Levenberg-Marquardt minimisation. The observation is a direct observation relating state to observation $z = h(x) + w$ where $z$ is the observation vector, $h()$ is the observation function, $x$ is the state vector and $w$ is a zero-mean noise vector, with covariance $N$.

Returns:
Handle for the new observation or NULL on failure.
See also:
gan_lev_marq_form().


Generated on Fri Mar 17 12:45:02 2006 by  doxygen 1.3.9.1