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

Linear system solver for Mx = b using LU decomposition of a sparse matrix. More...

#include <vnl_sparse_lu.h>

List of all members.

Public Types

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

Public Member Functions

 vnl_sparse_lu (vnl_sparse_matrix< double > const &M, operation mode=quiet)
 Make sparse_lu decomposition of M optionally computing the reciprocal condition number.
 ~vnl_sparse_lu ()
void set_pivot_thresh (double pivot_thresh)
 set the relative pivot threshold should be between 0 and 1.
void set_absolute_thresh (double absolute_thresh)
 set the threshold on absolute element magnitude for pivoting.
void set_diagonal_pivoting (int diag_pivoting)
 set diagonal pivoting mode, normally 1 which gives priority to diagonal elements.
vnl_vector< double > solve (vnl_vector< double > const &b)
 Solve problem M x = b.
void solve (vnl_vector< double > const &b, vnl_vector< double > *x)
 Solve problem M x = b.
vnl_vector< double > solve_transpose (vnl_vector< double > const &b)
 Solve problem M^t x = b.
void solve_transpose (vnl_vector< double > const &b, vnl_vector< double > *x)
 Solve problem M^t x = b.
double determinant ()
 Compute determinant.
double rcond ()
 Return reciprocal condition number (smallest/largest singular values).
double max_error_bound ()
 An estimate of maximum error in solution.

Protected Member Functions

void decompose_matrix ()
bool est_condition ()
 estimate the condition number.

Protected Attributes

vnl_sparse_matrix< double > A_
bool factored_
bool condition_computed_
operation mode_
double norm_
double rcond_
double largest_
double pivot_thresh_
double absolute_thresh_
int diag_pivoting_

Private Member Functions

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

Private Attributes

void * pmatrix_
 The internal matrix representation.

Detailed Description

Linear system solver for Mx = b using LU decomposition of a sparse matrix.

Encapsulating Sparse 1.3 by Kenneth S. Kundert. The matrix is factored before solution. Any number of b vectors can be applied after the matrix is factored with out recomputation of the LU form. It is advised to construct with mode==estimate_condition and check that rcond()>sqrt(machine precision). If this is not the case the solution may be suspect. An upper bound on error is provided. The solution of M^t x = b is also available.

Definition at line 30 of file vnl_sparse_lu.h.


Member Enumeration Documentation

Modes of computation. See constructor for details.

Enumerator:
quiet 
verbose 
estimate_condition 
estimate_condition_verbose 

Definition at line 34 of file vnl_sparse_lu.h.


Constructor & Destructor Documentation

vnl_sparse_lu::vnl_sparse_lu ( vnl_sparse_matrix< double > const &  M,
operation  mode = quiet 
)

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

constructor - controls if condition information is computed.

Definition at line 21 of file vnl_sparse_lu.cxx.

vnl_sparse_lu::~vnl_sparse_lu ( )

Definition at line 15 of file vnl_sparse_lu.cxx.

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

Copy constructor - privatised to avoid it being used.


Member Function Documentation

void vnl_sparse_lu::decompose_matrix ( ) [protected]
double vnl_sparse_lu::determinant ( )

Compute determinant.

Definition at line 245 of file vnl_sparse_lu.cxx.

bool vnl_sparse_lu::est_condition ( ) [protected]

estimate the condition number.

Definition at line 58 of file vnl_sparse_lu.cxx.

double vnl_sparse_lu::max_error_bound ( )

An estimate of maximum error in solution.

Estimated upper bound of error in solution.

Not calculated unless operation mode at construction includes estimate_condition.

Definition at line 292 of file vnl_sparse_lu.cxx.

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

Assignment operator - privatised to avoid it being used.

double vnl_sparse_lu::rcond ( )

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

the reciprocal of the condition number.

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

Definition at line 278 of file vnl_sparse_lu.cxx.

void vnl_sparse_lu::set_absolute_thresh ( double  absolute_thresh) [inline]

set the threshold on absolute element magnitude for pivoting.

Should be either zero or significantly smaller than the absolute value of the smallest diagonal element.

Definition at line 54 of file vnl_sparse_lu.h.

void vnl_sparse_lu::set_diagonal_pivoting ( int  diag_pivoting) [inline]

set diagonal pivoting mode, normally 1 which gives priority to diagonal elements.

Definition at line 56 of file vnl_sparse_lu.h.

void vnl_sparse_lu::set_pivot_thresh ( double  pivot_thresh) [inline]

set the relative pivot threshold should be between 0 and 1.

If set to one then pivoting is complete and slow If near zero then roundoff error may be prohibitive but computation is fast Typical values are between 0.01 and 0.1.

Definition at line 49 of file vnl_sparse_lu.h.

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

Solve problem M x = b.

Solve least squares problem M x = b.

Definition at line 123 of file vnl_sparse_lu.cxx.

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

Solve problem M x = b.

Solve least squares problem M x = b.

Definition at line 72 of file vnl_sparse_lu.cxx.

vnl_vector< double > vnl_sparse_lu::solve_transpose ( vnl_vector< double > const &  b)

Solve problem M^t x = b.

Definition at line 182 of file vnl_sparse_lu.cxx.

void vnl_sparse_lu::solve_transpose ( vnl_vector< double > const &  b,
vnl_vector< double > *  x 
)

Solve problem M^t x = b.

Definition at line 131 of file vnl_sparse_lu.cxx.


Member Data Documentation

vnl_sparse_matrix<double> vnl_sparse_lu::A_ [protected]

Definition at line 93 of file vnl_sparse_lu.h.

double vnl_sparse_lu::absolute_thresh_ [protected]

Definition at line 101 of file vnl_sparse_lu.h.

Definition at line 95 of file vnl_sparse_lu.h.

Definition at line 102 of file vnl_sparse_lu.h.

bool vnl_sparse_lu::factored_ [protected]

Definition at line 94 of file vnl_sparse_lu.h.

double vnl_sparse_lu::largest_ [protected]

Definition at line 99 of file vnl_sparse_lu.h.

Definition at line 96 of file vnl_sparse_lu.h.

double vnl_sparse_lu::norm_ [protected]

Definition at line 97 of file vnl_sparse_lu.h.

double vnl_sparse_lu::pivot_thresh_ [protected]

Definition at line 100 of file vnl_sparse_lu.h.

void* vnl_sparse_lu::pmatrix_ [private]

The internal matrix representation.

We don't use the typedef spMatrix directly here to avoid exposing the implementation detail (sparse/spMatrix.h) to the user.

Definition at line 112 of file vnl_sparse_lu.h.

double vnl_sparse_lu::rcond_ [protected]

Definition at line 98 of file vnl_sparse_lu.h.


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