vnl_qr.h
Go to the documentation of this file.
1 // This is core/vnl/algo/vnl_qr.h
2 #ifndef vnl_qr_h_
3 #define vnl_qr_h_
4 //:
5 // \file
6 // \brief Calculate inverse of a matrix using QR
7 // \author Andrew W. Fitzgibbon, Oxford RRG
8 // \date 08 Dec 1996
9 //
10 // \verbatim
11 // Modifications
12 // 081296 AWF Temporarily abandoned as I realized my problem was symmetric...
13 // 080697 AWF Recovered, implemented solve().
14 // 200897 AWF Added determinant().
15 // 071097 AWF Added Q(), R().
16 // Christian Stoecklin, ETH Zurich, added QtB(v)
17 // 31-mar-2000 fsm: templated
18 // 28/03/2001 - dac (Manchester) - tidied up documentation
19 // 13 Jan.2003 - Peter Vanroose - added missing implementation for inverse(),
20 // tinverse(), solve(matrix), extract_q_and_r().
21 // \endverbatim
22 
23 #include <iosfwd>
24 #include <vnl/vnl_vector.h>
25 #include <vnl/vnl_matrix.h>
26 #include <vnl/algo/vnl_algo_export.h>
27 #ifdef _MSC_VER
28 # include <vcl_msvc_warnings.h>
29 #endif
30 
31 //: Extract the Q*R decomposition of matrix M.
32 // The decomposition is stored in a compact and time-efficient
33 // packed form, which is most easily used via the "solve" and
34 // "determinant" methods.
35 
36 template <class T>
37 class vnl_qr
38 {
39  public:
40  vnl_qr(vnl_matrix<T> const & M);
41  ~vnl_qr();
42 
43  //: return the inverse matrix of M
44  vnl_matrix<T> inverse () const;
45  //: return the transpose of the inverse matrix of M
46  vnl_matrix<T> tinverse () const;
47  //: return the original matrix M
48  vnl_matrix<T> recompose () const;
49 
50  //: Solve equation M x = rhs for x using the computed decomposition.
51  vnl_matrix<T> solve (const vnl_matrix<T>& rhs) const;
52  //: Solve equation M x = rhs for x using the computed decomposition.
53  vnl_vector<T> solve (const vnl_vector<T>& rhs) const;
54 
55  //: Return the determinant of M. This is computed from M = Q R as follows:
56  // |M| = |Q| |R|.
57  // |R| is the product of the diagonal elements.
58  // |Q| is (-1)^n as it is a product of Householder reflections.
59  // So det = -prod(-r_ii).
60  T determinant() const;
61  //: Unpack and return unitary part Q.
62  vnl_matrix<T> const& Q() const;
63  //: Unpack and return R.
64  vnl_matrix<T> const& R() const;
65  //: Return residual vector d of M x = b -> d = Q'b
66  vnl_vector<T> QtB(const vnl_vector<T>& b) const;
67 
68  void extract_q_and_r(vnl_matrix<T>* q, vnl_matrix<T>* r) const { *q = Q(); *r = R(); }
69 
70  private:
74  mutable vnl_matrix<T>* Q_;
75  mutable vnl_matrix<T>* R_;
76 
77  // Disallow assignment.
78  vnl_qr(const vnl_qr<T> &) { }
79  vnl_qr<T>& operator=(const vnl_qr<T> &) { return *this; }
80 };
81 
82 //: Compute determinant of matrix "M" using QR.
83 template <class T>
85 {
86  return vnl_qr<T>(m).determinant();
87 }
88 
89 template <class T>
90 std::ostream& operator<<(std::ostream&, vnl_qr<T> const & qr);
91 
92 #endif // vnl_qr_h_
vnl_matrix< T > * Q_
Definition: vnl_qr.h:74
vnl_matrix< T > inverse() const
return the inverse matrix of M.
Definition: vnl_qr.hxx:251
vnl_matrix< T > solve(const vnl_matrix< T > &rhs) const
Solve equation M x = rhs for x using the computed decomposition.
Definition: vnl_qr.hxx:289
vnl_matrix< T > tinverse() const
return the transpose of the inverse matrix of M.
Definition: vnl_qr.hxx:270
An ordinary mathematical matrix.
void extract_q_and_r(vnl_matrix< T > *q, vnl_matrix< T > *r) const
Definition: vnl_qr.h:68
#define m
Definition: vnl_vector.h:43
vnl_matrix< T > qrdc_out_
Definition: vnl_qr.h:71
vnl_matrix< T > const & Q() const
Unpack and return unitary part Q.
Definition: vnl_qr.hxx:91
vnl_matrix< T > const & R() const
Unpack and return R.
Definition: vnl_qr.hxx:156
~vnl_qr()
Definition: vnl_qr.hxx:66
std::ostream & operator<<(std::ostream &s, vnl_decnum const &r)
decimal output.
Definition: vnl_decnum.h:393
vnl_matrix< T > * R_
Definition: vnl_qr.h:75
vnl_qr(const vnl_qr< T > &)
Definition: vnl_qr.h:78
Extract the Q*R decomposition of matrix M.
Definition: vnl_algo_fwd.h:8
T determinant() const
Return the determinant of M. This is computed from M = Q R as follows:.
Definition: vnl_qr.hxx:78
An ordinary mathematical matrix.
Definition: vnl_adjugate.h:22
vnl_qr< T > & operator=(const vnl_qr< T > &)
Definition: vnl_qr.h:79
vnl_vector< T > qraux_
Definition: vnl_qr.h:72
vnl_vector< T > QtB(const vnl_vector< T > &b) const
Return residual vector d of M x = b -> d = Q'b.
Definition: vnl_qr.hxx:220
T vnl_qr_determinant(vnl_matrix< T > const &m)
Compute determinant of matrix "M" using QR.
Definition: vnl_qr.h:84
Mathematical vector class, templated by type of element.
Definition: vnl_fwd.h:16
vnl_qr(vnl_matrix< T > const &M)
Definition: vnl_qr.hxx:37
vnl_vector< long > jpvt_
Definition: vnl_qr.h:73
vnl_matrix< T > recompose() const
return the original matrix M.
Definition: vnl_qr.hxx:176