Public Types | Public Member Functions | Protected Attributes | Private Member Functions
vnl_cholesky Class Reference

Decomposition of symmetric matrix. More...

#include <vnl_cholesky.h>

List of all members.

Public Types

enum  Operation { quiet, verbose, estimate_condition }
 Modes of computation. See constructor for details. More...

Public Member Functions

 vnl_cholesky (vnl_matrix< double > const &M, Operation mode=verbose)
 Make cholesky decomposition of M optionally computing the reciprocal condition number.
 ~vnl_cholesky ()
vnl_vector< double > solve (vnl_vector< double > const &b) const
 Solve LS problem M x = b.
void solve (vnl_vector< double > const &b, vnl_vector< double > *x) const
 Solve LS problem M x = b.
double determinant () const
 Compute determinant.
vnl_matrix< double > inverse () const
vnl_matrix< double > lower_triangle () const
 Return lower-triangular factor.
vnl_matrix< double > upper_triangle () const
 Return upper-triangular factor.
vnl_matrix< double > const & L_badly_named_method () const
 Return the decomposition matrix.
int rank_deficiency () const
 A Success/failure flag.
double rcond () const
 Return reciprocal condition number (smallest/largest singular values).
vnl_vector< double > & nullvector ()
 Return computed nullvector.
vnl_vector< double > const & nullvector () const

Protected Attributes

vnl_matrix< double > A_
double rcond_
long num_dims_rank_def_
vnl_vector< double > nullvector_

Private Member Functions

 vnl_cholesky (vnl_cholesky const &that)
 Copy constructor - privatised to avoid it being used.
vnl_choleskyoperator= (vnl_cholesky const &that)
 Assignment operator - privatised to avoid it being used.

Detailed Description

Decomposition of symmetric matrix.

A class to hold the Cholesky decomposition of a symmetric matrix and use that to solve linear systems, compute determinants and inverses. The cholesky decomposition decomposes symmetric A = L*L.transpose() where L is lower triangular

To check that the decomposition can be used safely for solving a linear equation it is wise to construct with mode==estimate_condition and check that rcond()>sqrt(machine precision). If this is not the case it might be a good idea to use vnl_svd instead.

Definition at line 33 of file vnl_cholesky.h.


Member Enumeration Documentation

Modes of computation. See constructor for details.

Enumerator:
quiet 
verbose 
estimate_condition 

Definition at line 37 of file vnl_cholesky.h.


Constructor & Destructor Documentation

vnl_cholesky::vnl_cholesky ( vnl_matrix< double > const &  M,
Operation  mode = verbose 
)

Make cholesky decomposition of M optionally computing the reciprocal condition number.

Cholesky decomposition.

Make cholesky decomposition of M optionally computing the reciprocal condition number. If mode is estimate_condition, the condition number and an approximate nullspace are estimated, at a cost of a factor of (1 + 18/n). Here's a table of 1 + 18/n:

   n:              3      5     10     50    100    500   1000
   slowdown:     7.0    4.6    2.8    1.4   1.18   1.04   1.02
   

Definition at line 28 of file vnl_cholesky.cxx.

vnl_cholesky::~vnl_cholesky ( ) [inline]

Definition at line 45 of file vnl_cholesky.h.

vnl_cholesky::vnl_cholesky ( vnl_cholesky const &  that) [private]

Copy constructor - privatised to avoid it being used.


Member Function Documentation

double vnl_cholesky::determinant ( ) const

Compute determinant.

Definition at line 76 of file vnl_cholesky.cxx.

vnl_matrix< double > vnl_cholesky::inverse ( ) const

Definition at line 87 of file vnl_cholesky.cxx.

vnl_matrix<double> const& vnl_cholesky::L_badly_named_method ( ) const [inline]

Return the decomposition matrix.

Definition at line 69 of file vnl_cholesky.h.

vnl_matrix< double > vnl_cholesky::lower_triangle ( ) const

Return lower-triangular factor.

Definition at line 108 of file vnl_cholesky.cxx.

vnl_vector<double>& vnl_cholesky::nullvector ( ) [inline]

Return computed nullvector.

Not calculated unless Operation mode at construction was estimate_condition.

Definition at line 82 of file vnl_cholesky.h.

vnl_vector<double> const& vnl_cholesky::nullvector ( ) const [inline]

Definition at line 83 of file vnl_cholesky.h.

vnl_cholesky& vnl_cholesky::operator= ( vnl_cholesky const &  that) [private]

Assignment operator - privatised to avoid it being used.

int vnl_cholesky::rank_deficiency ( ) const [inline]

A Success/failure flag.

Definition at line 72 of file vnl_cholesky.h.

double vnl_cholesky::rcond ( ) const [inline]

Return reciprocal condition number (smallest/largest singular values).

As long as rcond()>sqrt(precision) the decomposition can be used for solving equations safely. Not calculated unless Operation mode at construction was estimate_condition.

Definition at line 78 of file vnl_cholesky.h.

vnl_vector< double > vnl_cholesky::solve ( vnl_vector< double > const &  b) const

Solve LS problem M x = b.

Solve least squares problem M x = b.

Definition at line 65 of file vnl_cholesky.cxx.

void vnl_cholesky::solve ( vnl_vector< double > const &  b,
vnl_vector< double > *  x 
) const

Solve LS problem M x = b.

Solve least squares problem M x = b.

The right-hand-side vcl_vector x may be b, which will give a fractional increase in speed.

Definition at line 55 of file vnl_cholesky.cxx.

vnl_matrix< double > vnl_cholesky::upper_triangle ( ) const

Return upper-triangular factor.

Definition at line 125 of file vnl_cholesky.cxx.


Member Data Documentation

vnl_matrix<double> vnl_cholesky::A_ [protected]

Definition at line 87 of file vnl_cholesky.h.

Definition at line 90 of file vnl_cholesky.h.

Definition at line 89 of file vnl_cholesky.h.

double vnl_cholesky::rcond_ [protected]

Definition at line 88 of file vnl_cholesky.h.


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