Public Types | Public Member Functions | Protected Attributes | Private Member Functions | List of all members
vnl_cholesky Class Reference

Decomposition of symmetric matrix. More...

#include <vnl_cholesky.h>

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. More...
 
 ~vnl_cholesky ()=default
 
vnl_vector< double > solve (vnl_vector< double > const &b) const
 Solve LS problem M x = b. More...
 
void solve (vnl_vector< double > const &b, vnl_vector< double > *x) const
 Solve LS problem M x = b. More...
 
double determinant () const
 Compute determinant. More...
 
vnl_matrix< double > inverse () const
 
vnl_matrix< double > lower_triangle () const
 Return lower-triangular factor. More...
 
vnl_matrix< double > upper_triangle () const
 Return upper-triangular factor. More...
 
vnl_matrix< double > const & L_badly_named_method () const
 Return the decomposition matrix. More...
 
int rank_deficiency () const
 A Success/failure flag. More...
 
double rcond () const
 Return reciprocal condition number (smallest/largest singular values). More...
 
vnl_vector< double > & nullvector ()
 Return computed nullvector. More...
 
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)=delete
 Copy constructor - privatised to avoid it being used. More...
 
vnl_choleskyoperator= (vnl_cholesky const &that)=delete
 Assignment operator - privatised to avoid it being used. More...
 

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 31 of file vnl_cholesky.h.

Member Enumeration Documentation

◆ Operation

Modes of computation. See constructor for details.

Enumerator
quiet 
verbose 
estimate_condition 

Definition at line 35 of file vnl_cholesky.h.

Constructor & Destructor Documentation

◆ vnl_cholesky() [1/2]

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 25 of file vnl_cholesky.cxx.

◆ ~vnl_cholesky()

vnl_cholesky::~vnl_cholesky ( )
default

◆ vnl_cholesky() [2/2]

vnl_cholesky::vnl_cholesky ( vnl_cholesky const &  that)
privatedelete

Copy constructor - privatised to avoid it being used.

Member Function Documentation

◆ determinant()

double vnl_cholesky::determinant ( ) const

Compute determinant.

Definition at line 73 of file vnl_cholesky.cxx.

◆ inverse()

vnl_matrix< double > vnl_cholesky::inverse ( ) const

Definition at line 84 of file vnl_cholesky.cxx.

◆ L_badly_named_method()

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

Return the decomposition matrix.

Definition at line 67 of file vnl_cholesky.h.

◆ lower_triangle()

vnl_matrix< double > vnl_cholesky::lower_triangle ( ) const

Return lower-triangular factor.

Definition at line 105 of file vnl_cholesky.cxx.

◆ nullvector() [1/2]

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

Return computed nullvector.

Not calculated unless Operation mode at construction was estimate_condition.

Definition at line 80 of file vnl_cholesky.h.

◆ nullvector() [2/2]

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

Definition at line 81 of file vnl_cholesky.h.

◆ operator=()

vnl_cholesky& vnl_cholesky::operator= ( vnl_cholesky const &  that)
privatedelete

Assignment operator - privatised to avoid it being used.

◆ rank_deficiency()

int vnl_cholesky::rank_deficiency ( ) const
inline

A Success/failure flag.

Definition at line 70 of file vnl_cholesky.h.

◆ rcond()

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 76 of file vnl_cholesky.h.

◆ solve() [1/2]

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 62 of file vnl_cholesky.cxx.

◆ solve() [2/2]

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 std::vector x may be b, which will give a fractional increase in speed.

Definition at line 52 of file vnl_cholesky.cxx.

◆ upper_triangle()

vnl_matrix< double > vnl_cholesky::upper_triangle ( ) const

Return upper-triangular factor.

Definition at line 122 of file vnl_cholesky.cxx.

Member Data Documentation

◆ A_

vnl_matrix<double> vnl_cholesky::A_
protected

Definition at line 85 of file vnl_cholesky.h.

◆ nullvector_

vnl_vector<double> vnl_cholesky::nullvector_
protected

Definition at line 88 of file vnl_cholesky.h.

◆ num_dims_rank_def_

long vnl_cholesky::num_dims_rank_def_
protected

Definition at line 87 of file vnl_cholesky.h.

◆ rcond_

double vnl_cholesky::rcond_
protected

Definition at line 86 of file vnl_cholesky.h.


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