vnl_qr.hxx
Go to the documentation of this file.
1 // This is core/vnl/algo/vnl_qr.hxx
2 #ifndef vnl_qr_hxx_
3 #define vnl_qr_hxx_
4 //:
5 // \file
6 // \author Andrew W. Fitzgibbon, Oxford RRG
7 // \date 08 Dec 1996
8 
9 #include <iostream>
10 #include <complex>
11 #include "vnl_qr.h"
12 #include <cassert>
13 #ifdef _MSC_VER
14 # include <vcl_msvc_warnings.h>
15 #endif
16 #include <vnl/vnl_math.h>
17 #include <vnl/vnl_complex.h> // vnl_math::squared_magnitude()
18 #include <vnl/vnl_matlab_print.h>
19 #include <vnl/vnl_complex_traits.h>
20 #include <vnl/algo/vnl_netlib.h> // dqrdc_(), dqrsl_()
21 
22 // use C++ overloading to call the right linpack routine from the template code:
23 #ifndef DOXYGEN_SHOULD_SKIP_THIS
24 #define macro(p, T) \
25 inline void vnl_linpack_qrdc(vnl_netlib_qrdc_proto(T)) \
26 { v3p_netlib_##p##qrdc_(vnl_netlib_qrdc_params); } \
27 inline void vnl_linpack_qrsl(vnl_netlib_qrsl_proto(T)) \
28 { v3p_netlib_##p##qrsl_(vnl_netlib_qrsl_params); }
29 macro(s, float);
30 macro(d, double);
31 macro(c, std::complex<float>);
32 macro(z, std::complex<double>);
33 #undef macro
34 #endif
35 
36 template <class T>
38  qrdc_out_(M.columns(), M.rows()),
39  qraux_(M.columns()),
40  jpvt_(M.rows()),
41  Q_(nullptr),
42  R_(nullptr)
43 {
44  assert(! M.empty());
45 
46  // Fill transposed O/P matrix
47  long c = M.columns();
48  long r = M.rows();
49  for (int i = 0; i < r; ++i)
50  for (int j = 0; j < c; ++j)
51  qrdc_out_(j,i) = M(i,j);
52 
53  long do_pivot = 0; // Enable[!=0]/disable[==0] pivoting.
54  jpvt_.fill(0); // Allow all columns to be pivoted if pivoting is enabled.
55 
56  vnl_vector<T> work(M.rows());
57  vnl_linpack_qrdc(qrdc_out_.data_block(), // On output, UT is R, below diag is mangled Q
58  &r, &r, &c,
59  qraux_.data_block(), // Further information required to demangle Q
60  jpvt_.data_block(),
61  work.data_block(),
62  &do_pivot);
63 }
64 
65 template <class T>
67 {
68  delete Q_;
69  delete R_;
70 }
71 
72 //: Return the determinant of M. This is computed from M = Q R as follows:
73 // |M| = |Q| |R|
74 // |R| is the product of the diagonal elements.
75 // |Q| is (-1)^n as it is a product of Householder reflections.
76 // So det = -prod(-r_ii).
77 template <class T>
79 {
80  int m = std::min((int)qrdc_out_.columns(), (int)qrdc_out_.rows());
81  T det = qrdc_out_(0,0);
82 
83  for (int i = 1; i < m; ++i)
84  det *= -qrdc_out_(i,i);
85 
86  return det;
87 }
88 
89 //: Unpack and return unitary part Q.
90 template <class T>
92 {
93  int m = qrdc_out_.columns(); // column-major storage
94  int n = qrdc_out_.rows();
95 
96  bool verbose = false;
97 
98  if (!Q_) {
99  Q_ = new vnl_matrix<T>(m,m);
100  // extract Q.
101  if (verbose) {
102  std::cerr << __FILE__ ": vnl_qr<T>::Q()\n"
103  << " m,n = " << m << ", " << n << '\n'
104  << " qr0 = [" << qrdc_out_ << "];\n"
105  << " aux = [" << qraux_ << "];\n";
106  }
107 
108  Q_->set_identity();
109  vnl_matrix<T>& matrQ = *Q_;
110 
111  vnl_vector<T> v(m, T(0));
112  vnl_vector<T> w(m, T(0));
113 
114  // Golub and vanLoan, p199. backward accumulation of householder matrices
115  // Householder vector k is [zeros(1,k-1) qraux_[k] qrdc_out_[k,:]]
116  typedef typename vnl_numeric_traits<T>::abs_t abs_t;
117  for (int k = n-1; k >= 0; --k) {
118  if (k >= m) continue;
119  // Make housevec v, and accumulate norm at the same time.
120  v[k] = qraux_[k];
121  abs_t sq = vnl_math::squared_magnitude(v[k]);
122  for (int j = k+1; j < m; ++j) {
123  v[j] = qrdc_out_(k,j);
124  sq += vnl_math::squared_magnitude(v[j]);
125  }
126  if (verbose) vnl_matlab_print(std::cerr, v, "v");
127 #ifndef DOXYGEN_SHOULD_SKIP_THIS
128 # define c vnl_complex_traits<T>::conjugate
129 #endif
130  // Premultiply emerging Q by house(v), noting that v[0..k-1] == 0.
131  // Q_new = (1 - (2/v'*v) v v')Q
132  // or Q -= (2/v'*v) v (v'Q)
133  if (sq > abs_t(0)) {
134  abs_t scale = abs_t(2)/sq;
135  // w = (2/v'*v) v' Q
136  for (int i = k; i < m; ++i) {
137  w[i] = T(0);
138  for (int j = k; j < m; ++j)
139  w[i] += scale * c(v[j]) * matrQ(j, i);
140  }
141  if (verbose) vnl_matlab_print(std::cerr, w, "w");
142 
143  // Q -= v w
144  for (int i = k; i < m; ++i)
145  for (int j = k; j < m; ++j)
146  matrQ(i,j) -= (v[i]) * (w[j]);
147  }
148 #undef c
149  }
150  }
151  return *Q_;
152 }
153 
154 //: Unpack and return R.
155 template <class T>
157 {
158  if (!R_) {
159  int m = qrdc_out_.columns(); // column-major storage
160  int n = qrdc_out_.rows();
161  R_ = new vnl_matrix<T>(m,n);
162  vnl_matrix<T> & Rmatr = *R_;
163 
164  for (int i = 0; i < m; ++i)
165  for (int j = 0; j < n; ++j)
166  if (i > j)
167  Rmatr(i, j) = T(0);
168  else
169  Rmatr(i, j) = qrdc_out_(j,i);
170  }
171 
172  return *R_;
173 }
174 
175 template <class T>
177 {
178  return Q() * R();
179 }
180 
181 // JOB: ABCDE decimal
182 // A B C D E
183 // --- --- --- --- ---
184 // Qb Q'b x norm(A*x - b) A*x
185 
186 
187 //: Solve equation M x = b for x using the computed decomposition.
188 template <class T>
190 {
191  long n = qrdc_out_.columns();
192  long p = qrdc_out_.rows();
193  const T* b_data = b.data_block();
194  vnl_vector<T> Qt_B(n);
195  vnl_vector<T> x(p);
196 
197  // see comment above
198  long JOB = 100;
199 
200  long info = 0;
201  vnl_linpack_qrsl(qrdc_out_.data_block(),
202  &n, &n, &p,
203  qraux_.data_block(),
204  b_data, (T*)nullptr, Qt_B.data_block(),
205  x.data_block(),
206  (T*)nullptr/*residual*/,
207  (T*)nullptr/*Ax*/,
208  &JOB,
209  &info);
210 
211  if (info > 0)
212  std::cerr << __FILE__ ": vnl_qr<T>::solve() : matrix is rank-deficient by "
213  << info << '\n';
214 
215  return x;
216 }
217 
218 //: Return residual vector d of M x = b -> d = Q'b
219 template <class T>
221 {
222  long n = qrdc_out_.columns();
223  long p = qrdc_out_.rows();
224  const T* b_data = b.data_block();
225  vnl_vector<T> Qt_B(n);
226 
227  // see comment above
228  long JOB = 1000;
229 
230  long info = 0;
231  vnl_linpack_qrsl(qrdc_out_.data_block(),
232  &n, &n, &p,
233  qraux_.data_block(),
234  b_data,
235  (T*)nullptr, // A: Qb
236  Qt_B.data_block(), // B: Q'b
237  (T*)nullptr, // C: x
238  (T*)nullptr, // D: residual
239  (T*)nullptr, // E: Ax
240  &JOB,
241  &info);
242 
243  if (info > 0)
244  std::cerr << __FILE__ ": vnl_qr<T>::QtB() -- matrix is rank-deficient by "
245  << info << '\n';
246 
247  return Qt_B;
248 }
249 
250 template <class T>
252 {
253  unsigned int r = qrdc_out_.columns();
254  assert(r > 0 && r == qrdc_out_.rows());
255  vnl_matrix<T> inv(r,r);
256 
257  // Use solve() to compute the inverse matrix, using (00..010..00) as rhs
258  vnl_vector<T> rhs(r,T(0));
259  for (unsigned int i=0; i<r; ++i)
260  {
261  rhs(i) = T(1);
262  vnl_vector<T> col = this->solve(rhs); // returns i-th column of inverse
263  inv.set_column(i,col);
264  rhs(i) = T(0);
265  }
266  return inv;
267 }
268 
269 template <class T>
271 {
272  unsigned int r = qrdc_out_.columns();
273  assert(r > 0 && r == qrdc_out_.rows());
274  vnl_matrix<T> tinv(r,r);
275 
276  // Use solve() to compute the inverse matrix, using (00..010..00) as rhs
277  vnl_vector<T> rhs(r,T(0));
278  for (unsigned int i=0; i<r; ++i)
279  {
280  rhs(i) = T(1);
281  vnl_vector<T> col = this->solve(rhs); // returns i-th column of inverse
282  tinv.set_row(i,col);
283  rhs(i) = T(0);
284  }
285  return tinv;
286 }
287 
288 template <class T>
290 {
291  assert(rhs.rows() == qrdc_out_.columns()); // column-major storage
292  int c = qrdc_out_.rows();
293  int n = rhs.columns();
294  vnl_matrix<T> result(c,n);
295 
296  for (int i=0; i<n; ++i)
297  {
298  vnl_vector<T> b = rhs.get_column(i);
299  vnl_vector<T> col = this->solve(b); // returns i-th column of result
300  result.set_column(i,col);
301  }
302  return result;
303 }
304 
305 //------------------------------------------------------------------------------
306 
307 #define VNL_QR_INSTANTIATE(T) \
308  template class VNL_ALGO_EXPORT vnl_qr<T >; \
309  /*template VNL_EXPORT T vnl_qr_determinant(vnl_matrix<T > const&) */
310 
311 #endif
vnl_matrix< T > inverse() const
return the inverse matrix of M.
Definition: vnl_qr.hxx:251
Print matrices and vectors in nice MATLAB format.
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
Complex additions to vnl_math.
vnl_matrix< T > tinverse() const
return the transpose of the inverse matrix of M.
Definition: vnl_qr.hxx:270
VNL_EXPORT std::ostream & vnl_matlab_print(std::ostream &, vnl_diag_matrix< T > const &, char const *variable_name=nullptr, vnl_matlab_print_format=vnl_matlab_print_format_default)
print a vnl_diagonal_matrix<T>.
To allow templated algorithms to determine appropriate actions of conjugation, complexification etc.
#define m
Definition: vnl_vector.h:43
T const * data_block() const
Access the contiguous block storing the elements in the vector. O(1).
Definition: vnl_vector.h:230
vnl_matrix< T > qrdc_out_
Definition: vnl_qr.h:71
Namespace with standard math functions.
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
Calculate inverse of a matrix using QR.
~vnl_qr()
Definition: vnl_qr.hxx:66
vnl_bignum squared_magnitude(vnl_bignum const &x)
Definition: vnl_bignum.h:433
#define v
Definition: vnl_vector.h:42
vnl_matrix & set_column(unsigned i, T const *v)
Set the elements of the i'th column to v[i] (No bounds checking).
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_vector< T > qraux_
Definition: vnl_qr.h:72
Declare in a central place the list of symbols from netlib.
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
vnl_vector & fill(T const &v)
Set all values to v.
Definition: vnl_vector.hxx:314
bool empty() const
Return true iff the size is zero.
Definition: vnl_matrix.h:549
Mathematical vector class, templated by type of element.
Definition: vnl_fwd.h:16
vnl_vector< T > get_column(unsigned c) const
Get a vector equal to the given column.
Definition: vnl_matrix.hxx:977
unsigned int rows() const
Return the number of rows.
Definition: vnl_matrix.h:179
vnl_qr(vnl_matrix< T > const &M)
Definition: vnl_qr.hxx:37
vnl_decnum min(vnl_decnum const &x, vnl_decnum const &y)
Definition: vnl_decnum.h:413
vnl_vector< long > jpvt_
Definition: vnl_qr.h:73
unsigned int columns() const
Return the number of columns.
Definition: vnl_matrix.h:187
#define macro(p, T)
Definition: vnl_svd.hxx:23
vnl_matrix< T > recompose() const
return the original matrix M.
Definition: vnl_qr.hxx:176