vnl_linear_system.h
Go to the documentation of this file.
1 // This is core/vnl/vnl_linear_system.h
2 #ifndef vnl_linear_system_h_
3 #define vnl_linear_system_h_
4 //:
5 // \file
6 // \brief Abstraction for a linear system of equations.
7 // \author David Capel, capes@robots
8 // \date July 2000
9 //
10 // \verbatim
11 // Modifications
12 // LSB (Manchester) 23/1/01 Documentation tidied
13 // \endverbatim
14 
15 #include <vnl/vnl_vector.h>
16 #include "vnl/vnl_export.h"
17 
18 //: Abstraction for a linear system of equations.
19 // vnl_linear_system provides an abstraction for a linear system
20 // of equations, Ax = b, to be solved by one of the iterative linear
21 // solvers. Access to the systems is via the pure virtual methods
22 // multiply() and transpose_multiply(). This procedural access scheme
23 // makes it possible to solve very large, sparse systems which it would
24 // be inefficient to store in matrix form.
25 //
26 // To solve the system, use an algorithm like vnl_lsqr.
27 class VNL_EXPORT vnl_linear_system
28 {
29  public:
30 
31  vnl_linear_system(unsigned int number_of_unknowns, unsigned int number_of_residuals) :
32  p_(number_of_unknowns), n_(number_of_residuals) {}
33 
34  virtual ~vnl_linear_system();
35 
36  //: Compute A*x, putting result in y
37  virtual void multiply(vnl_vector<double> const& x, vnl_vector<double>& y) const = 0;
38 
39  //: Compute A_transpose * y, putting result in x
40  virtual void transpose_multiply(vnl_vector<double> const& y, vnl_vector<double>& x) const = 0;
41 
42  //; Put the right-hand side of the system Ax = b into b
43  virtual void get_rhs(vnl_vector<double>& b) const = 0;
44 
45  //; (Optional) Apply a suitable preconditioner to x.
46  // A preconditioner is an approximation of the inverse of A.
47  // Common choices are Jacobi (1/diag(A'A)), Gauss-Seidel,
48  // and incomplete LU or Cholesky decompositions.
49  // The default implementation applies the identity.
50  virtual void apply_preconditioner(vnl_vector<double> const& x, vnl_vector<double>& px) const;
51 
52  //: Return the number of unknowns
53  unsigned int get_number_of_unknowns() const { return p_; }
54 
55  //: Return the number of residuals.
56  unsigned int get_number_of_residuals() const { return n_; }
57 
58  //: Compute rms error for parameter vector x
59  double get_rms_error(vnl_vector<double> const& x) const;
60 
61  //: Compute relative residual (|Ax - b| / |b| )for parameter vector x
62  double get_relative_residual(vnl_vector<double> const& x) const;
63 
64  protected:
65  unsigned int p_;
66  unsigned int n_;
67 };
68 
69 #endif // vnl_linear_system_h_
vnl_linear_system(unsigned int number_of_unknowns, unsigned int number_of_residuals)
Abstraction for a linear system of equations.
unsigned int get_number_of_residuals() const
Return the number of residuals.
unsigned int get_number_of_unknowns() const
Return the number of unknowns.