vnl_sparse_lu.h
Go to the documentation of this file.
1 // This is core/vnl/algo/vnl_sparse_lu.h
2 #ifndef vnl_sparse_lu_h_
3 #define vnl_sparse_lu_h_
4 //:
5 // \file
6 // \brief Solution of the linear system Mx = b for sparse matrices
7 // \author J. L. Mundy
8 // \date 27 Dec 2006
9 //
10 // \verbatim
11 // Modifications
12 // None
13 // \endverbatim
14 
15 #include <vnl/vnl_vector.h>
16 #include <vnl/vnl_sparse_matrix.h>
17 #include <vnl/algo/vnl_algo_export.h>
18 
19 //: Linear system solver for Mx = b using LU decomposition of a sparse matrix
20 // Encapsulating Sparse 1.3 by Kenneth S. Kundert.
21 // The matrix is factored before solution. Any number of b vectors can
22 // be applied after the matrix is factored with out recomputation of the
23 // LU form. It is advised to construct with mode==estimate_condition and
24 // check that rcond()>sqrt(machine precision). If this is not the case
25 // the solution may be suspect. An upper bound on error is provided.
26 // The solution of M^t x = b is also available.
27 
28 class VNL_ALGO_EXPORT vnl_sparse_lu
29 {
30  public:
31  //: Modes of computation. See constructor for details.
32  enum operation {
36  estimate_condition_verbose
37  };
38 
39  //: Make sparse_lu decomposition of M optionally computing the reciprocal condition number.
40  vnl_sparse_lu(vnl_sparse_matrix<double> const& M, operation mode = quiet);
41  ~vnl_sparse_lu();
42 
43  //: set the relative pivot threshold should be between 0 and 1
44  // If set to one then pivoting is complete and slow
45  // If near zero then roundoff error may be prohibitive but computation is fast
46  // Typical values are between 0.01 and 0.1.
47  void set_pivot_thresh(double pivot_thresh){pivot_thresh_=pivot_thresh;}
48 
49  //: set the threshold on absolute element magnitude for pivoting
50  // Should be either zero or significantly smaller than the absolute
51  // value of the smallest diagonal element.
52  void set_absolute_thresh(double absolute_thresh){absolute_thresh_=absolute_thresh;}
53  //: set diagonal pivoting mode, normally 1 which gives priority to diagonal elements.
54  void set_diagonal_pivoting(int diag_pivoting){diag_pivoting_=diag_pivoting;}
55 
56  //: Solve problem M x = b
57  vnl_vector<double> solve(vnl_vector<double> const& b);
58 
59  //: Solve problem M x = b
60  void solve(vnl_vector<double> const& b, vnl_vector<double>* x);
61 
62  //: Solve problem M^t x = b
63  vnl_vector<double> solve_transpose(vnl_vector<double> const& b);
64 
65  //: Solve problem M^t x = b
66  void solve_transpose(vnl_vector<double> const& b, vnl_vector<double>* x);
67 
68  //: Compute determinant
69  double determinant();
70 
71  //: Return reciprocal condition number (smallest/largest singular values).
72  // As long as rcond()>sqrt(precision) the decomposition can be used for
73  // solving equations safely.
74  // Not calculated unless operation mode at construction includes estimate_condition.
75  double rcond();
76 
77  //: An estimate of maximum error in solution.
78  // Not calculated unless operation mode at construction includes estimate_condition.
79  double max_error_bound();
80 
81  protected:
82  // Internals
83  void decompose_matrix();
84  bool est_condition();
85  // Data Members--------------------------------------------------------------
87  bool factored_;
90  double norm_;
91  double rcond_;
92  double largest_;
93  double pivot_thresh_;
96  private:
97  //: Copy constructor - privatised to avoid it being used
98  vnl_sparse_lu(vnl_sparse_lu const & that);
99  //: Assignment operator - privatised to avoid it being used
100  vnl_sparse_lu& operator=(vnl_sparse_lu const & that);
101  //: The internal matrix representation
102  //
103  // We don't use the typedef spMatrix directly here to avoid exposing
104  // the implementation detail (sparse/spMatrix.h) to the user.
105  void* pmatrix_;
106 };
107 
108 #endif // vnl_sparse_lu_h_
void set_absolute_thresh(double absolute_thresh)
set the threshold on absolute element magnitude for pivoting.
Definition: vnl_sparse_lu.h:52
void * pmatrix_
The internal matrix representation.
double absolute_thresh_
Definition: vnl_sparse_lu.h:94
vnl_sparse_matrix< double > A_
Definition: vnl_sparse_lu.h:86
void set_diagonal_pivoting(int diag_pivoting)
set diagonal pivoting mode, normally 1 which gives priority to diagonal elements.
Definition: vnl_sparse_lu.h:54
void set_pivot_thresh(double pivot_thresh)
set the relative pivot threshold should be between 0 and 1.
Definition: vnl_sparse_lu.h:47
operation
Modes of computation. See constructor for details.
Definition: vnl_sparse_lu.h:32
double pivot_thresh_
Definition: vnl_sparse_lu.h:93
bool condition_computed_
Definition: vnl_sparse_lu.h:88
operation mode_
Definition: vnl_sparse_lu.h:89
Linear system solver for Mx = b using LU decomposition of a sparse matrix.
Definition: vnl_sparse_lu.h:28
Simple sparse matrix.