vnl_matrix.hxx
Go to the documentation of this file.
1 // This is core/vnl/vnl_matrix.hxx
2 #ifndef vnl_matrix_hxx_
3 #define vnl_matrix_hxx_
4 //:
5 // \file
6 //
7 // Copyright (C) 1991 Texas Instruments Incorporated.
8 // Copyright (C) 1992 General Electric Company.
9 //
10 // Permission is granted to any individual or institution to use, copy, modify,
11 // and distribute this software, provided that this complete copyright and
12 // permission notice is maintained, intact, in all copies and supporting
13 // documentation.
14 //
15 // Texas Instruments Incorporated, General Electric Company,
16 // provides this software "as is" without express or implied warranty.
17 //
18 // Created: MBN Apr 21, 1989 Initial design and implementation
19 // Updated: MBN Jun 22, 1989 Removed non-destructive methods
20 // Updated: LGO Aug 09, 1989 Inherit from Generic
21 // Updated: MBN Aug 20, 1989 Changed template usage to reflect new syntax
22 // Updated: MBN Sep 11, 1989 Added conditional exception handling and base class
23 // Updated: LGO Oct 05, 1989 Don't re-allocate data in operator= when same size
24 // Updated: LGO Oct 19, 1989 Add extra parameter to varargs constructor
25 // Updated: MBN Oct 19, 1989 Added optional argument to set_compare method
26 // Updated: LGO Dec 08, 1989 Allocate column data in one chunk
27 // Updated: LGO Dec 08, 1989 Clean-up get and put, add const everywhere.
28 // Updated: LGO Dec 19, 1989 Remove the map and reduce methods
29 // Updated: MBN Feb 22, 1990 Changed size arguments from int to unsigned int
30 // Updated: MJF Jun 30, 1990 Added base class name to constructor initializer
31 // Updated: VDN Feb 21, 1992 New lite version
32 // Updated: VDN May 05, 1992 Use envelope to avoid unnecessary copying
33 // Updated: VDN Sep 30, 1992 Matrix inversion with singular value decomposition
34 // Updated: AWF Aug 21, 1996 set_identity, normalize_rows, scale_row.
35 // Updated: AWF Sep 30, 1996 set_row/column methods. Const-correct data_block().
36 // Updated: AWF 14 Feb 1997 get_n_rows, get_n_columns.
37 // Updated: PVR 20 Mar 1997 get_row, get_column.
38 //
39 // The parameterized vnl_matrix<T> class implements two dimensional arithmetic
40 // matrices of a user specified type. The only constraint placed on the type is
41 // that it must overload the following operators: +, -, *, and /. Thus, it
42 // will be possible to have a vnl_matrix over std::complex<T>. The vnl_matrix<T>
43 // class is static in size, that is once a vnl_matrix<T> of a particular size
44 // has been created, there is no dynamic growth method available. You can
45 // resize the matrix, with the loss of any existing data using set_size().
46 //
47 // Each matrix contains a protected data section that has a T** slot that
48 // points to the physical memory allocated for the two dimensional array. In
49 // addition, two integers specify the number of rows and columns for the
50 // matrix. These values are provided in the constructors. A single protected
51 // slot contains a pointer to a compare function to be used in equality
52 // operations. The default function used is the built-in == operator.
53 //
54 // Four different constructors are provided. The first constructor takes two
55 // integer arguments specifying the row and column size. Enough memory is
56 // allocated to hold row*column elements of type Type. The second constructor
57 // takes the same two first arguments, but also accepts an additional third
58 // argument that is a reference to an object of the appropriate type whose
59 // value is used as an initial fill value. The third constructor is similar to
60 // the third, except that it accepts a variable number of initialization values
61 // for the Matrix. If there are fewer values than elements, the rest are set
62 // to zero. Finally, the last constructor takes a single argument consisting of
63 // a reference to a Matrix and duplicates its size and element values.
64 //
65 // Methods are provided for destructive scalar and Matrix addition,
66 // multiplication, check for equality and inequality, fill, reduce, and access
67 // and set individual elements. Finally, both the input and output operators
68 // are overloaded to allow for formatted input and output of matrix elements.
69 //
70 // Good matrix inversion is needed. We choose singular value decomposition,
71 // since it is general and works great for nearly singular cases. Singular
72 // value decomposition is preferred to LU decomposition, since the accuracy
73 // of the pivots is independent from the left->right top->down elimination.
74 // LU decomposition also does not give eigenvectors and eigenvalues when
75 // the matrix is symmetric.
76 //
77 // Several different constructors are provided. See .h file for brief descriptions.
78 
79 //--------------------------------------------------------------------------------
80 
81 #include <cstdio>
82 #include <cstdlib>
83 #include <cctype>
84 #include <iostream>
85 #include <vector>
86 #include <algorithm>
87 #include "vnl_matrix.h"
88 
89 #include <cassert>
90 #ifdef _MSC_VER
91 # include <vcl_msvc_warnings.h>
92 #endif
93 
94 #include <vnl/vnl_math.h>
95 #include <vnl/vnl_vector.h>
96 #include <vnl/vnl_c_vector.h>
97 #include <vnl/vnl_numeric_traits.h>
98 //--------------------------------------------------------------------------------
99 # define vnl_matrix_construct_hack()
100 
101 // This macro allocates and initializes the dynamic storage used by a vnl_matrix.
102 #define vnl_matrix_alloc_blah() \
103 do { \
104  if (this->num_rows && this->num_cols) { \
105  /* Allocate memory to hold the row pointers */ \
106  this->data = vnl_c_vector<T>::allocate_Tptr(this->num_rows); \
107  /* Allocate memory to hold the elements of the matrix */ \
108  T* elmns = vnl_c_vector<T>::allocate_T(this->num_rows * this->num_cols); \
109  /* Fill in the array of row pointers */ \
110  for (unsigned int i = 0; i < this->num_rows; ++ i) \
111  this->data[i] = elmns + i*this->num_cols; \
112  } \
113  else { \
114  /* This is to make sure .begin() and .end() work for 0xN matrices: */ \
115  (this->data = vnl_c_vector<T>::allocate_Tptr(1))[0] = 0; \
116  } \
117 } while (false)
118 
119 // This macro releases the dynamic storage used by a vnl_matrix.
120 #define vnl_matrix_free_blah \
121 do { \
122  if (this->data) { \
123  if (this->num_cols && this->num_rows) { \
124  vnl_c_vector<T>::deallocate(this->data[0], this->num_cols * this->num_rows); \
125  vnl_c_vector<T>::deallocate(this->data, this->num_rows); \
126  } \
127  else { \
128  vnl_c_vector<T>::deallocate(this->data, 1); \
129  } \
130  } \
131 } while (false)
132 
133 //: Creates a matrix with given number of rows and columns.
134 // Elements are not initialized. O(m*n).
135 
136 template <class T>
137 vnl_matrix<T>::vnl_matrix (unsigned rowz, unsigned colz)
138 : num_rows(rowz), num_cols(colz)
139 {
142 }
143 
144 //: Creates a matrix with given number of rows and columns, and initialize all elements to value. O(m*n).
145 
146 template <class T>
147 vnl_matrix<T>::vnl_matrix (unsigned rowz, unsigned colz, T const& value)
148 : num_rows(rowz), num_cols(colz)
149 {
152  std::fill_n( this->data[0], rowz * colz, value );
153 }
154 
155 //: r rows, c cols, special type. Currently implements "identity" and "null".
156 template <class T>
157 vnl_matrix<T>::vnl_matrix(unsigned r, unsigned c, vnl_matrix_type t)
158 : num_rows(r), num_cols(c)
159 {
162  switch (t) {
163  case vnl_matrix_identity:
164  assert(r == c);
165  for (unsigned int i = 0; i < r; ++ i)
166  for (unsigned int j = 0; j < c; ++ j)
167  this->data[i][j] = (i==j) ? T(1) : T(0);
168  break;
169  case vnl_matrix_null:
170  std::fill_n( this->data[0], r * c, T(0) );
171  break;
172  default:
173  assert(false);
174  break;
175  }
176 }
177 
178 #if 1 // fsm: who uses this?
179 //: Creates a matrix with given dimension (rows, cols) and initialize first n elements, row-wise, to values. O(m*n).
180 
181 template <class T>
182 vnl_matrix<T>::vnl_matrix (unsigned rowz, unsigned colz, unsigned n, T const values[])
183 : num_rows(rowz), num_cols(colz)
184 {
187  if (n > rowz*colz)
188  n = rowz*colz;
189  std::copy( values, values + n, this->data[0] );
190 }
191 #endif
192 
193 //: Creates a matrix from a block array of data, stored row-wise.
194 // O(m*n).
195 
196 template <class T>
197 vnl_matrix<T>::vnl_matrix (T const* datablck, unsigned rowz, unsigned colz)
198 : num_rows(rowz), num_cols(colz)
199 {
202  std::copy( datablck, datablck + rowz * colz, this->data[0] );
203 }
204 
205 
206 //: Creates a new matrix and copies all the elements.
207 // O(m*n).
208 
209 template <class T>
211 : num_rows(from.num_rows), num_cols(from.num_cols)
212 {
214  if (from.data && from.data[0]) {
216  T const *src = from.data[0];
217  std::copy( src, src + this->num_rows * this->num_cols, this->data[0] );
218  }
219  else {
220  num_rows = 0;
221  num_cols = 0;
222  data = nullptr;
223  }
224 }
225 
226 //------------------------------------------------------------
227 
228 template <class T>
230 : num_rows(A.num_rows), num_cols(A.num_cols)
231 {
232 #ifndef NDEBUG
233  if (A.num_rows != B.num_rows || A.num_cols != B.num_cols)
234  vnl_error_matrix_dimension ("vnl_tag_add", A.num_rows, A.num_cols, B.num_rows, B.num_cols);
235 #endif
236 
239 
240  unsigned int n = A.num_rows * A.num_cols;
241  T const *a = A.data[0];
242  T const *b = B.data[0];
243  T *dst = this->data[0];
244 
245  for (unsigned int i=0; i<n; ++i)
246  dst[i] = T(a[i] + b[i]);
247 }
248 
249 template <class T>
251 : num_rows(A.num_rows), num_cols(A.num_cols)
252 {
253 #ifndef NDEBUG
254  if (A.num_rows != B.num_rows || A.num_cols != B.num_cols)
255  vnl_error_matrix_dimension ("vnl_tag_sub", A.num_rows, A.num_cols, B.num_rows, B.num_cols);
256 #endif
257 
260 
261  unsigned int n = A.num_rows * A.num_cols;
262  T const *a = A.data[0];
263  T const *b = B.data[0];
264  T *dst = this->data[0];
265 
266  for (unsigned int i=0; i<n; ++i)
267  dst[i] = T(a[i] - b[i]);
268 }
269 
270 template <class T>
272 : num_rows(M.num_rows), num_cols(M.num_cols)
273 {
276 
277  unsigned int n = M.num_rows * M.num_cols;
278  T const *m = M.data[0];
279  T *dst = this->data[0];
280 
281  for (unsigned int i=0; i<n; ++i)
282  dst[i] = T(m[i] * s);
283 }
284 
285 template <class T>
287 : num_rows(M.num_rows), num_cols(M.num_cols)
288 {
291 
292  unsigned int n = M.num_rows * M.num_cols;
293  T const *m = M.data[0];
294  T *dst = this->data[0];
295 
296  for (unsigned int i=0; i<n; ++i)
297  dst[i] = T(m[i] / s);
298 }
299 
300 template <class T>
302 : num_rows(M.num_rows), num_cols(M.num_cols)
303 {
306 
307  unsigned int n = M.num_rows * M.num_cols;
308  T const *m = M.data[0];
309  T *dst = this->data[0];
310 
311  for (unsigned int i=0; i<n; ++i)
312  dst[i] = T(m[i] + s);
313 }
314 
315 template <class T>
317 : num_rows(M.num_rows), num_cols(M.num_cols)
318 {
321 
322  unsigned int n = M.num_rows * M.num_cols;
323  T const *m = M.data[0];
324  T *dst = this->data[0];
325 
326  for (unsigned int i=0; i<n; ++i)
327  dst[i] = T(m[i] - s);
328 }
329 
330 template <class T>
332 : num_rows(A.num_rows), num_cols(B.num_cols)
333 {
334 #ifndef NDEBUG
335  if (A.num_cols != B.num_rows)
336  vnl_error_matrix_dimension("vnl_tag_mul", A.num_rows, A.num_cols, B.num_rows, B.num_cols);
337 #endif
338 
339  unsigned int l = A.num_rows;
340  unsigned int m = A.num_cols; // == B.num_rows
341  unsigned int n = B.num_cols;
342 
345 
346  for (unsigned int i=0; i<l; ++i) {
347  for (unsigned int k=0; k<n; ++k) {
348  T sum(0);
349  for (unsigned int j=0; j<m; ++j)
350  sum += T(A.data[i][j] * B.data[j][k]);
351  this->data[i][k] = sum;
352  }
353  }
354 }
355 
356 //------------------------------------------------------------
357 
358 template <class T>
360 {
361  // save some fcalls if data is 0 (i.e. in matrix_fixed)
362  if (data) destroy();
363 }
364 
365 //: Frees up the dynamic storage used by matrix.
366 // O(m*n).
367 
368 template <class T>
370 {
372 }
373 
374 template <class T>
376 {
377  if (data) {
378  destroy();
379  num_rows = 0;
380  num_cols = 0;
381  data = nullptr;
382  }
383 }
384 
385 // Resizes the data arrays of THIS matrix to (rows x cols). O(m*n).
386 // Elements are not initialized, existing data is not preserved.
387 // Returns true if size is changed.
388 
389 template <class T>
390 bool vnl_matrix<T>::set_size (unsigned rowz, unsigned colz)
391 {
392  if (this->data) {
393  // if no change in size, do not reallocate.
394  if (this->num_rows == rowz && this->num_cols == colz)
395  return false;
396 
397  // else, simply release old storage and allocate new.
399  this->num_rows = rowz; this->num_cols = colz;
401  }
402  else {
403  // This happens if the matrix is default constructed.
404  this->num_rows = rowz; this->num_cols = colz;
406  }
407 
408  return true;
409 }
410 
411 #undef vnl_matrix_alloc_blah
412 #undef vnl_matrix_free_blah
413 
414 //------------------------------------------------------------
415 
416 //: Sets all elements of matrix to specified value. O(m*n).
417 
418 template <class T>
420 {
421  // not safe if data == NULL, due to data[0] call
422  if (data && data[0])
423  std::fill_n( this->data[0], this->num_rows * this->num_cols, value );
424  return *this;
425 }
426 
427 //: Sets all diagonal elements of matrix to specified value. O(n).
428 
429 template <class T>
431 {
432  for (unsigned int i = 0; i < this->num_rows && i < this->num_cols; ++i)
433  this->data[i][i] = value;
434  return *this;
435 }
436 
437 //: Sets the diagonal elements of this matrix to the specified list of values.
438 
439 template <class T>
441 {
442  assert(diag.size() >= this->num_rows ||
443  diag.size() >= this->num_cols);
444  // The length of the diagonal of a non-square matrix is the minimum of
445  // the matrix's width & height; that explains the "||" in the assert,
446  // and the "&&" in the upper bound for the "for".
447  for (unsigned int i = 0; i < this->num_rows && i < this->num_cols; ++i)
448  this->data[i][i] = diag[i];
449  return *this;
450 }
451 
452 //: Access an element for reading or writing
453 // There are assert style boundary checks - #define NDEBUG to turn them off.
454 template <class T>
455 T & vnl_matrix<T>::operator()(unsigned r, unsigned c)
456 {
457 #if VNL_CONFIG_CHECK_BOUNDS
458  assert(r < rows()); // Check the row index is valid
459  assert(c < cols()); // Check the column index is valid
460 #endif
461  return this->data[r][c];
462 }
463 
464 //: Access an element for reading
465 // There are assert style boundary checks - #define NDEBUG to turn them off.
466 template <class T>
467 T const & vnl_matrix<T>::operator()(unsigned r, unsigned c) const
468 {
469 #if VNL_CONFIG_CHECK_BOUNDS
470  assert(r < rows()); // Check the row index is valid
471  assert(c < cols()); // Check the column index is valid
472 #endif
473  return this->data[r][c];
474 }
475 
476 
477 //: Copies all elements of rhs matrix into lhs matrix. O(m*n).
478 // If needed, the arrays in lhs matrix are freed up, and new arrays are
479 // allocated to match the dimensions of the rhs matrix.
480 
481 template <class T>
483 {
484  if (this != &rhs) { // make sure *this != m
485  if (rhs.data) {
486  this->set_size(rhs.num_rows, rhs.num_cols);
487  if (rhs.data[0]) {
488  std::copy( rhs.data[0], rhs.data[0] + this->num_rows * this->num_cols, this->data[0] );
489  }
490  }
491  else {
492  // rhs is default-constructed.
493  clear();
494  }
495  }
496  return *this;
497 }
498 
499 template <class T>
500 void vnl_matrix<T>::print(std::ostream& os) const
501 {
502  for (unsigned int i = 0; i < this->rows(); i++) {
503  for (unsigned int j = 0; j < this->columns(); j++)
504  os << this->data[i][j] << ' ';
505  os << '\n';
506  }
507 }
508 
509 //: Prints the 2D array of elements of a matrix out to a stream.
510 // O(m*n).
511 
512 template <class T>
513 std::ostream& operator<< (std::ostream& os, vnl_matrix<T> const& m)
514 {
515  for (unsigned int i = 0; i < m.rows(); ++i) {
516  for (unsigned int j = 0; j < m.columns(); ++j)
517  os << m(i, j) << ' ';
518  os << '\n';
519  }
520  return os;
521 }
522 
523 //: Read a vnl_matrix from an ascii std::istream.
524 // Automatically determines file size if the input matrix has zero size.
525 template <class T>
526 std::istream& operator>>(std::istream& s, vnl_matrix<T>& M)
527 {
528  M.read_ascii(s);
529  return s;
530 }
531 
532 template <class T>
534 {
535  for (unsigned int i = 0; i < this->num_rows; i++)
536  for (unsigned int j = 0; j < this->num_cols; j++)
537  this->data[i][j] += value;
538  return *this;
539 }
540 
541 template <class T>
543 {
544  for (unsigned int i = 0; i < this->num_rows; i++)
545  for (unsigned int j = 0; j < this->num_cols; j++)
546  this->data[i][j] -= value;
547  return *this;
548 }
549 
550 template <class T>
552 {
553  for (unsigned int i = 0; i < this->num_rows; i++)
554  for (unsigned int j = 0; j < this->num_cols; j++)
555  this->data[i][j] *= value;
556  return *this;
557 }
558 
559 template <class T>
561 {
562  for (unsigned int i = 0; i < this->num_rows; i++)
563  for (unsigned int j = 0; j < this->num_cols; j++)
564  this->data[i][j] /= value;
565  return *this;
566 }
567 
568 
569 //: Adds lhs matrix with rhs matrix, and stores in place in lhs matrix.
570 // O(m*n). The dimensions of the two matrices must be identical.
571 
572 template <class T>
574 {
575 #ifndef NDEBUG
576  if (this->num_rows != rhs.num_rows ||
577  this->num_cols != rhs.num_cols) // Size match?
578  vnl_error_matrix_dimension ("operator+=",
579  this->num_rows, this->num_cols,
580  rhs.num_rows, rhs.num_cols);
581 #endif
582  for (unsigned int i = 0; i < this->num_rows; i++) // For each row
583  for (unsigned int j = 0; j < this->num_cols; j++) // For each element in column
584  this->data[i][j] += rhs.data[i][j]; // Add elements
585  return *this;
586 }
587 
588 
589 //: Subtract lhs matrix with rhs matrix and store in place in lhs matrix.
590 // O(m*n).
591 // The dimensions of the two matrices must be identical.
592 
593 template <class T>
595 {
596 #ifndef NDEBUG
597  if (this->num_rows != rhs.num_rows ||
598  this->num_cols != rhs.num_cols) // Size?
599  vnl_error_matrix_dimension ("operator-=",
600  this->num_rows, this->num_cols,
601  rhs.num_rows, rhs.num_cols);
602 #endif
603  for (unsigned int i = 0; i < this->num_rows; i++)
604  for (unsigned int j = 0; j < this->num_cols; j++)
605  this->data[i][j] -= rhs.data[i][j];
606  return *this;
607 }
608 
609 
610 template <class T>
611 vnl_matrix<T> operator- (T const& value, vnl_matrix<T> const& m)
612 {
613  vnl_matrix<T> result(m.rows(),m.columns());
614  for (unsigned int i = 0; i < m.rows(); i++) // For each row
615  for (unsigned int j = 0; j < m.columns(); j++) // For each element in column
616  result.put(i,j, T(value - m.get(i,j)) ); // subtract from value element.
617  return result;
618 }
619 
620 //: Returns new matrix which is the negation of THIS matrix.
621 // O(m*n).
622 
623 template <class T>
625 {
626  vnl_matrix<T> result(this->num_rows, this->num_cols);
627  for (unsigned int i = 0; i < this->num_rows; i++)
628  for (unsigned int j = 0; j < this->num_cols; j++)
629  result.data[i][j] = - this->data[i][j];
630  return result;
631 }
632 
633 //: Return the matrix made by applying "f" to each element.
634 template <class T>
635 vnl_matrix<T> vnl_matrix<T>::apply(T (*f)(T const&)) const
636 {
637  vnl_matrix<T> ret(num_rows, num_cols);
638  vnl_c_vector<T>::apply(this->data[0], num_rows * num_cols, f, ret.data_block());
639  return ret;
640 }
641 
642 //: Return the matrix made by applying "f" to each element.
643 template <class T>
645 {
646  vnl_matrix<T> ret(num_rows, num_cols);
647  vnl_c_vector<T>::apply(this->data[0], num_rows * num_cols, f, ret.data_block());
648  return ret;
649 }
650 
651 //: Make a vector by applying a function across rows.
652 template <class T>
655 ::apply_rowwise(T (*f)(vnl_vector<T> const&)) const
656 {
657  vnl_vector<T> v(this->num_rows);
658  for (unsigned int i = 0; i < this->num_rows; ++i)
659  v.put(i,f(this->get_row(i)));
660  return v;
661 }
662 
663 //: Make a vector by applying a function across columns.
664 template <class T>
667 ::apply_columnwise(T (*f)(vnl_vector<T> const&)) const
668 {
669  vnl_vector<T> v(this->num_cols);
670  for (unsigned int i = 0; i < this->num_cols; ++i)
671  v.put(i,f(this->get_column(i)));
672  return v;
673 }
674 
675 
676 ////--------------------------- Additions------------------------------------
677 
678 //: Returns new matrix with rows and columns transposed.
679 // O(m*n).
680 
681 template <class T>
683 {
684  vnl_matrix<T> result(this->num_cols, this->num_rows);
685  for (unsigned int i = 0; i < this->num_cols; i++)
686  for (unsigned int j = 0; j < this->num_rows; j++)
687  result.data[i][j] = this->data[j][i];
688  return result;
689 }
690 
691 // adjoint/hermitian transpose
692 
693 template <class T>
695 {
696  vnl_matrix<T> result(transpose());
697  vnl_c_vector<T>::conjugate(result.begin(), // src
698  result.begin(), // dst
699  result.size()); // size of block
700  return result;
701 }
702 
703 //: Replaces the submatrix of THIS matrix, starting at top left corner, by the elements of matrix m. O(m*n).
704 // This is the reverse of extract().
705 
706 template <class T>
708  unsigned top, unsigned left)
709 {
710  unsigned int bottom = top + m.num_rows;
711  unsigned int right = left + m.num_cols;
712 #ifndef NDEBUG
713  if (this->num_rows < bottom || this->num_cols < right)
714  vnl_error_matrix_dimension ("update",
715  bottom, right, m.num_rows, m.num_cols);
716 #endif
717  for (unsigned int i = top; i < bottom; i++)
718  for (unsigned int j = left; j < right; j++)
719  this->data[i][j] = m.data[i-top][j-left];
720  return *this;
721 }
722 
723 
724 //: Returns a copy of submatrix of THIS matrix, specified by the top-left corner and size in rows, cols. O(m*n).
725 // Use update() to copy new values of this submatrix back into THIS matrix.
726 
727 template <class T>
728 vnl_matrix<T> vnl_matrix<T>::extract (unsigned rowz, unsigned colz,
729  unsigned top, unsigned left) const {
730  vnl_matrix<T> result(rowz, colz);
731  this->extract( result, top, left );
732  return result;
733 }
734 
735 template <class T>
737  unsigned top, unsigned left) const {
738  unsigned const rowz = submatrix.rows();
739  unsigned const colz = submatrix.cols();
740 #ifndef NDEBUG
741  unsigned int bottom = top + rowz;
742  unsigned int right = left + colz;
743  if ((this->num_rows < bottom) || (this->num_cols < right))
744  vnl_error_matrix_dimension ("extract",
745  this->num_rows, this->num_cols, bottom, right);
746 #endif
747  for (unsigned int i = 0; i < rowz; i++) // actual copy of all elements
748  for (unsigned int j = 0; j < colz; j++) // in submatrix
749  submatrix.data[i][j] = data[top+i][left+j];
750 }
751 
752 //: Returns the dot product of the two matrices. O(m*n).
753 // This is the sum of all pairwise products of the elements m1[i,j]*m2[i,j].
754 
755 template <class T>
756 T dot_product (vnl_matrix<T> const& m1, vnl_matrix<T> const& m2)
757 {
758 #ifndef NDEBUG
759  if (m1.rows() != m2.rows() || m1.columns() != m2.columns()) // Size?
760  vnl_error_matrix_dimension ("dot_product",
761  m1.rows(), m1.columns(),
762  m2.rows(), m2.columns());
763 #endif
764  return vnl_c_vector<T>::dot_product(m1.begin(), m2.begin(), m1.rows()*m1.cols());
765 }
766 
767 //: Hermitian inner product.
768 // O(mn).
769 
770 template <class T>
771 T inner_product (vnl_matrix<T> const& m1, vnl_matrix<T> const& m2)
772 {
773 #ifndef NDEBUG
774  if (m1.rows() != m2.rows() || m1.columns() != m2.columns()) // Size?
775  vnl_error_matrix_dimension ("inner_product",
776  m1.rows(), m1.columns(),
777  m2.rows(), m2.columns());
778 #endif
779  return vnl_c_vector<T>::inner_product(m1.begin(), m2.begin(), m1.rows()*m1.cols());
780 }
781 
782 // cos_angle. O(mn).
783 
784 template <class T>
785 T cos_angle (vnl_matrix<T> const& a, vnl_matrix<T> const& b)
786 {
787  typedef typename vnl_numeric_traits<T>::abs_t Abs_t;
788  typedef typename vnl_numeric_traits<Abs_t>::real_t abs_r;
789 
790  T ab = inner_product(a,b);
791  Abs_t a_b = (Abs_t)std::sqrt( (abs_r)vnl_math::abs(inner_product(a,a) * inner_product(b,b)) );
792 
793  return T( ab / a_b);
794 }
795 
796 //: Returns new matrix whose elements are the products m1[ij]*m2[ij].
797 // O(m*n).
798 
799 template <class T>
801  vnl_matrix<T> const& m2)
802 {
803 #ifndef NDEBUG
804  if (m1.rows() != m2.rows() || m1.columns() != m2.columns()) // Size?
805  vnl_error_matrix_dimension ("element_product",
806  m1.rows(), m1.columns(), m2.rows(), m2.columns());
807 #endif
808  vnl_matrix<T> result(m1.rows(), m1.columns());
809  for (unsigned int i = 0; i < m1.rows(); i++)
810  for (unsigned int j = 0; j < m1.columns(); j++)
811  result.put(i,j, T(m1.get(i,j) * m2.get(i,j)) );
812  return result;
813 }
814 
815 //: Returns new matrix whose elements are the quotients m1[ij]/m2[ij].
816 // O(m*n).
817 
818 template <class T>
820  vnl_matrix<T> const& m2)
821 {
822 #ifndef NDEBUG
823  if (m1.rows() != m2.rows() || m1.columns() != m2.columns()) // Size?
824  vnl_error_matrix_dimension("element_quotient",
825  m1.rows(), m1.columns(), m2.rows(), m2.columns());
826 #endif
827  vnl_matrix<T> result(m1.rows(), m1.columns());
828  for (unsigned int i = 0; i < m1.rows(); i++)
829  for (unsigned int j = 0; j < m1.columns(); j++)
830  result.put(i,j, T(m1.get(i,j) / m2.get(i,j)) );
831  return result;
832 }
833 
834 //: Fill this matrix with the given data.
835 // We assume that p points to a contiguous rows*cols array, stored rowwise.
836 template <class T>
838 {
839  std::copy( p, p + this->num_rows * this->num_cols, this->data[0] );
840  return *this;
841 }
842 
843 //: Fill the given array with this matrix.
844 // We assume that p points to a contiguous rows*cols array, stored rowwise.
845 template <class T>
846 void vnl_matrix<T>::copy_out(T *p) const
847 {
848  std::copy( this->data[0], this->data[0] + this->num_rows * this->num_cols, p );
849 }
850 
851 //: Fill this matrix with a matrix having 1s on the main diagonal and 0s elsewhere.
852 template <class T>
854 {
855  for (unsigned int i = 0; i < this->num_rows; ++i) // For each row in the Matrix
856  for (unsigned int j = 0; j < this->num_cols; ++j) // For each element in column
857  this->data[i][j] = (i==j) ? T(1) : T(0);
858  return *this;
859 }
860 
861 //: Make each row of the matrix have unit norm.
862 // All-zero rows are ignored.
863 template <class T>
865 {
866  typedef typename vnl_numeric_traits<T>::abs_t Abs_t;
867  typedef typename vnl_numeric_traits<T>::real_t Real_t;
868  typedef typename vnl_numeric_traits<Real_t>::abs_t abs_real_t;
869  for (unsigned int i = 0; i < this->num_rows; ++i) { // For each row in the Matrix
870  Abs_t norm(0); // double will not do for all types.
871  for (unsigned int j = 0; j < this->num_cols; ++j) // For each element in row
872  norm += vnl_math::squared_magnitude(this->data[i][j]);
873 
874  if (norm != 0) {
875  abs_real_t scale = abs_real_t(1)/(std::sqrt((abs_real_t)norm));
876  for (unsigned int j = 0; j < this->num_cols; ++j)
877  this->data[i][j] = T(Real_t(this->data[i][j]) * scale);
878  }
879  }
880  return *this;
881 }
882 
883 //: Make each column of the matrix have unit norm.
884 // All-zero columns are ignored.
885 template <class T>
887 {
888  typedef typename vnl_numeric_traits<T>::abs_t Abs_t;
889  typedef typename vnl_numeric_traits<T>::real_t Real_t;
890  typedef typename vnl_numeric_traits<Real_t>::abs_t abs_real_t;
891  for (unsigned int j = 0; j < this->num_cols; j++) { // For each column in the Matrix
892  Abs_t norm(0); // double will not do for all types.
893  for (unsigned int i = 0; i < this->num_rows; i++)
894  norm += vnl_math::squared_magnitude(this->data[i][j]);
895 
896  if (norm != 0) {
897  abs_real_t scale = abs_real_t(1)/(std::sqrt((abs_real_t)norm));
898  for (unsigned int i = 0; i < this->num_rows; i++)
899  this->data[i][j] = T(Real_t(this->data[i][j]) * scale);
900  }
901  }
902  return *this;
903 }
904 
905 //: Multiply row[row_index] by value
906 template <class T>
907 vnl_matrix<T>& vnl_matrix<T>::scale_row(unsigned row_index, T value)
908 {
909 #ifndef NDEBUG
910  if (row_index >= this->num_rows)
911  vnl_error_matrix_row_index("scale_row", row_index);
912 #endif
913  for (unsigned int j = 0; j < this->num_cols; j++) // For each element in row
914  this->data[row_index][j] *= value;
915  return *this;
916 }
917 
918 //: Multiply column[column_index] by value
919 template <class T>
920 vnl_matrix<T>& vnl_matrix<T>::scale_column(unsigned column_index, T value)
921 {
922 #ifndef NDEBUG
923  if (column_index >= this->num_cols)
924  vnl_error_matrix_col_index("scale_column", column_index);
925 #endif
926  for (unsigned int j = 0; j < this->num_rows; j++) // For each element in column
927  this->data[j][column_index] *= value;
928  return *this;
929 }
930 
931 //: Returns a copy of n rows, starting from "row"
932 template <class T>
933 vnl_matrix<T> vnl_matrix<T>::get_n_rows (unsigned row, unsigned n) const
934 {
935 #ifndef NDEBUG
936  if (row + n > this->num_rows)
937  vnl_error_matrix_row_index ("get_n_rows", row);
938 #endif
939 
940  // Extract data rowwise.
941  return vnl_matrix<T>(data[row], n, this->num_cols);
942 }
943 
944 //: Returns a copy of n columns, starting from "column".
945 template <class T>
946 vnl_matrix<T> vnl_matrix<T>::get_n_columns (unsigned column, unsigned n) const
947 {
948 #ifndef NDEBUG
949  if (column + n > this->num_cols)
950  vnl_error_matrix_col_index ("get_n_columns", column);
951 #endif
952 
953  vnl_matrix<T> result(this->num_rows, n);
954  for (unsigned int c = 0; c < n; ++c)
955  for (unsigned int r = 0; r < this->num_rows; ++r)
956  result(r, c) = data[r][column + c];
957  return result;
958 }
959 
960 //: Create a vector out of row[row_index].
961 template <class T>
962 vnl_vector<T> vnl_matrix<T>::get_row(unsigned row_index) const
963 {
964 #ifdef ERROR_CHECKING
965  if (row_index >= this->num_rows)
966  vnl_error_matrix_row_index ("get_row", row_index);
967 #endif
968 
969  vnl_vector<T> v(this->num_cols);
970  for (unsigned int j = 0; j < this->num_cols; j++) // For each element in row
971  v[j] = this->data[row_index][j];
972  return v;
973 }
974 
975 //: Create a vector out of column[column_index].
976 template <class T>
977 vnl_vector<T> vnl_matrix<T>::get_column(unsigned column_index) const
978 {
979 #ifdef ERROR_CHECKING
980  if (column_index >= this->num_cols)
981  vnl_error_matrix_col_index ("get_column", column_index);
982 #endif
983 
984  vnl_vector<T> v(this->num_rows);
985  for (unsigned int j = 0; j < this->num_rows; j++) // For each element in row
986  v[j] = this->data[j][column_index];
987  return v;
988 }
989 
990 //: Create a vector out of row[row_index].
991 template <class T>
995 {
996  vnl_matrix<T> m(i.size(), this->num_cols);
997  for (unsigned int j = 0; j < i.size(); ++j)
998  m.set_row(j, this->get_row(i.get(j)));
999  return m;
1000 }
1001 
1002 //: Create a vector out of column[column_index].
1003 template <class T>
1007 {
1008  vnl_matrix<T> m(this->num_rows, i.size());
1009  for (unsigned int j = 0; j < i.size(); ++j)
1010  m.set_column(j, this->get_column(i.get(j)));
1011  return m;
1012 }
1013 
1014 //: Return a vector with the content of the (main) diagonal
1015 template <class T>
1017 {
1018  vnl_vector<T> v(this->num_rows < this->num_cols ? this->num_rows : this->num_cols);
1019  for (unsigned int j = 0; j < this->num_rows && j < this->num_cols; ++j)
1020  v[j] = this->data[j][j];
1021  return v;
1022 }
1023 
1024 //: Flatten row-major (C-style)
1025 template <class T>
1027 {
1028  vnl_vector<T> v(this->num_rows * this->num_cols);
1029  v.copy_in(this->data_block());
1030  return v;
1031 }
1032 
1033 //: Flatten column-major (Fortran-style)
1034 template <class T>
1036 {
1037  vnl_vector<T> v(this->num_rows * this->num_cols);
1038  for (unsigned int c = 0; c < this->num_cols; ++c)
1039  for (unsigned int r = 0; r < this->num_rows; ++r)
1040  v[c*this->num_rows+r] = this->data[r][c];
1041  return v;
1042 }
1043 
1044 //--------------------------------------------------------------------------------
1045 
1046 //: Set row[row_index] to data at given address. No bounds check.
1047 template <class T>
1048 vnl_matrix<T>& vnl_matrix<T>::set_row(unsigned row_index, T const *v)
1049 {
1050  for (unsigned int j = 0; j < this->num_cols; j++) // For each element in row
1051  this->data[row_index][j] = v[j];
1052  return *this;
1053 }
1054 
1055 //: Set row[row_index] to given vector.
1056 template <class T>
1058 {
1059 #ifndef NDEBUG
1060  if (v.size() != this->num_cols)
1061  vnl_error_vector_dimension ("vnl_matrix::set_row", v.size(), this->num_cols);
1062 #endif
1063  set_row(row_index,v.data_block());
1064  return *this;
1065 }
1066 
1067 //: Set row[row_index] to given value.
1068 template <class T>
1069 vnl_matrix<T>& vnl_matrix<T>::set_row(unsigned row_index, T v)
1070 {
1071  for (unsigned int j = 0; j < this->num_cols; j++) // For each element in row
1072  this->data[row_index][j] = v;
1073  return *this;
1074 }
1075 
1076 //--------------------------------------------------------------------------------
1077 
1078 //: Set column[column_index] to data at given address.
1079 template <class T>
1080 vnl_matrix<T>& vnl_matrix<T>::set_column(unsigned column_index, T const *v)
1081 {
1082  for (unsigned int i = 0; i < this->num_rows; i++) // For each element in row
1083  this->data[i][column_index] = v[i];
1084  return *this;
1085 }
1086 
1087 //: Set column[column_index] to given vector.
1088 template <class T>
1090 {
1091 #ifndef NDEBUG
1092  if (v.size() != this->num_rows)
1093  vnl_error_vector_dimension ("vnl_matrix::set_column", v.size(), this->num_rows);
1094 #endif
1095  set_column(column_index,v.data_block());
1096  return *this;
1097 }
1098 
1099 //: Set column[column_index] to given value.
1100 template <class T>
1101 vnl_matrix<T>& vnl_matrix<T>::set_column(unsigned column_index, T v)
1102 {
1103  for (unsigned int j = 0; j < this->num_rows; j++) // For each element in row
1104  this->data[j][column_index] = v;
1105  return *this;
1106 }
1107 
1108 
1109 //: Set columns starting at starting_column to given matrix
1110 template <class T>
1111 vnl_matrix<T>& vnl_matrix<T>::set_columns(unsigned starting_column, vnl_matrix<T> const& m)
1112 {
1113 #ifndef NDEBUG
1114  if (this->num_rows != m.num_rows ||
1115  this->num_cols < m.num_cols + starting_column) // Size match?
1116  vnl_error_matrix_dimension ("set_columns",
1117  this->num_rows, this->num_cols,
1118  m.num_rows, m.num_cols);
1119 #endif
1120 
1121  for (unsigned int j = 0; j < m.num_cols; ++j)
1122  for (unsigned int i = 0; i < this->num_rows; i++) // For each element in row
1123  this->data[i][starting_column + j] = m.data[i][j];
1124  return *this;
1125 }
1126 
1127 //--------------------------------------------------------------------------------
1128 
1129 //: Two matrices are equal if and only if they have the same dimensions and the same values.
1130 // O(m*n).
1131 // Elements are compared with operator== as default.
1132 // Change this default with set_compare() at run time or by specializing
1133 // vnl_matrix_compare at compile time.
1134 
1135 template <class T>
1137 {
1138  if (this == &rhs) // same object => equal.
1139  return true;
1140 
1141  if (this->num_rows != rhs.num_rows || this->num_cols != rhs.num_cols)
1142  return false; // different sizes => not equal.
1143 
1144  for (unsigned int i = 0; i < this->num_rows; i++) // For each row
1145  for (unsigned int j = 0; j < this->num_cols; j++) // For each column
1146  if (!(this->data[i][j] == rhs.data[i][j])) // different element ?
1147  return false; // Then not equal.
1148 
1149  return true; // Else same; return true
1150 }
1151 
1152 template <class T>
1153 bool vnl_matrix<T>::is_equal(vnl_matrix<T> const& rhs, double tol) const
1154 {
1155  if (this == &rhs) // same object => equal.
1156  return true;
1157 
1158  if (this->num_rows != rhs.num_rows || this->num_cols != rhs.num_cols)
1159  return false; // different sizes => not equal.
1160 
1161  for (unsigned int i = 0; i < this->rows(); ++i)
1162  for (unsigned int j = 0; j < this->columns(); ++j)
1163  if (vnl_math::abs(this->data[i][j] - rhs.data[i][j]) > tol)
1164  return false; // difference greater than tol
1165 
1166  return true;
1167 }
1168 
1169 
1170 template <class T>
1172 {
1173  T const zero(0);
1174  T const one(1);
1175  for (unsigned int i = 0; i < this->rows(); ++i)
1176  for (unsigned int j = 0; j < this->columns(); ++j) {
1177  T xm = (*this)(i,j);
1178  if ( !((i == j) ? (xm == one) : (xm == zero)) )
1179  return false;
1180  }
1181  return true;
1182 }
1183 
1184 //: Return true if maximum absolute deviation of M from identity is <= tol.
1185 template <class T>
1186 bool vnl_matrix<T>::is_identity(double tol) const
1187 {
1188  T one(1);
1189  for (unsigned int i = 0; i < this->rows(); ++i)
1190  for (unsigned int j = 0; j < this->columns(); ++j) {
1191  T xm = (*this)(i,j);
1192  abs_t absdev = (i == j) ? vnl_math::abs(xm - one) : vnl_math::abs(xm);
1193  if (absdev > tol)
1194  return false;
1195  }
1196  return true;
1197 }
1198 
1199 template <class T>
1201 {
1202  T const zero(0);
1203  for (unsigned int i = 0; i < this->rows(); ++i)
1204  for (unsigned int j = 0; j < this->columns(); ++j)
1205  if ( !( (*this)(i, j) == zero) )
1206  return false;
1207 
1208  return true;
1209 }
1210 
1211 //: Return true if max(abs((*this))) <= tol.
1212 template <class T>
1213 bool vnl_matrix<T>::is_zero(double tol) const
1214 {
1215  for (unsigned int i = 0; i < this->rows(); ++i)
1216  for (unsigned int j = 0; j < this->columns(); ++j)
1217  if (vnl_math::abs((*this)(i,j)) > tol)
1218  return false;
1219 
1220  return true;
1221 }
1222 
1223 //: Return true if any element of (*this) is nan
1224 template <class T>
1226 {
1227  for (unsigned int i = 0; i < this->rows(); ++i)
1228  for (unsigned int j = 0; j < this->columns(); ++j)
1229  if (vnl_math::isnan((*this)(i,j)))
1230  return true;
1231 
1232  return false;
1233 }
1234 
1235 //: Return false if any element of (*this) is inf or nan
1236 template <class T>
1238 {
1239  for (unsigned int i = 0; i < this->rows(); ++i)
1240  for (unsigned int j = 0; j < this->columns(); ++j)
1241  if (!vnl_math::isfinite((*this)(i,j)))
1242  return false;
1243 
1244  return true;
1245 }
1246 
1247 //: Abort if any element of M is inf or nan
1248 template <class T>
1250 {
1251  if (is_finite())
1252  return;
1253 
1254  std::cerr << "\n\n" __FILE__ ": " << __LINE__ << ": matrix has non-finite elements\n";
1255 
1256  if (rows() <= 20 && cols() <= 20) {
1257  std::cerr << __FILE__ ": here it is:\n" << *this;
1258  }
1259  else {
1260  std::cerr << __FILE__ ": it is quite big (" << rows() << 'x' << cols() << ")\n"
1261  << __FILE__ ": in the following picture '-' means finite and '*' means non-finite:\n";
1262 
1263  for (unsigned int i=0; i<rows(); ++i) {
1264  for (unsigned int j=0; j<cols(); ++j)
1265  std::cerr << char(vnl_math::isfinite((*this)(i, j)) ? '-' : '*');
1266  std::cerr << '\n';
1267  }
1268  }
1269  std::cerr << __FILE__ ": calling abort()\n";
1270  std::abort();
1271 }
1272 
1273 //: Abort unless M has the given size.
1274 template <class T>
1275 void vnl_matrix<T>::assert_size_internal(unsigned rs,unsigned cs) const
1276 {
1277  if (this->rows()!=rs || this->cols()!=cs) {
1278  std::cerr << __FILE__ ": size is " << this->rows() << 'x' << this->cols()
1279  << ". should be " << rs << 'x' << cs << std::endl;
1280  std::abort();
1281  }
1282 }
1283 
1284 //: Read a vnl_matrix from an ascii std::istream.
1285 // Automatically determines file size if the input matrix has zero size.
1286 template <class T>
1287 bool vnl_matrix<T>::read_ascii(std::istream& s)
1288 {
1289  if (!s.good()) {
1290  std::cerr << __FILE__ ": vnl_matrix<T>::read_ascii: Called with bad stream\n";
1291  return false;
1292  }
1293 
1294  bool size_known = (this->rows() != 0);
1295 
1296  if (size_known) {
1297  for (unsigned int i = 0; i < this->rows(); ++i)
1298  for (unsigned int j = 0; j < this->columns(); ++j)
1299  s >> this->data[i][j];
1300 
1301  return s.good() || s.eof();
1302  }
1303 
1304  bool debug = false;
1305 
1306  std::vector<T> first_row_vals;
1307  if (debug)
1308  std::cerr << __FILE__ ": vnl_matrix<T>::read_ascii: Determining file dimensions: ";
1309 
1310  for (;;) {
1311  // Clear whitespace, looking for a newline
1312  while (true)
1313  {
1314  int c = s.get();
1315  if (c == EOF)
1316  goto loademup;
1317  if (!std::isspace(c)) {
1318  if (!s.putback(char(c)).good())
1319  std::cerr << "vnl_matrix<T>::read_ascii: Could not push back '" << c << "'\n";
1320 
1321  goto readfloat;
1322  }
1323  // First newline after first number tells us the column dimension
1324  if (c == '\n' && first_row_vals.size() > 0) {
1325  goto loademup;
1326  }
1327  }
1328  readfloat:
1329  T val;
1330  s >> val;
1331  if (!s.fail())
1332  first_row_vals.push_back(val);
1333  if (s.eof())
1334  goto loademup;
1335  }
1336  loademup:
1337  std::size_t colz = first_row_vals.size();
1338 
1339  if (debug) std::cerr << colz << " cols, ";
1340 
1341  if (colz == 0)
1342  return false;
1343 
1344  // need to be careful with resizing here as will often be reading humungous files
1345  // So let's just build an array of row pointers
1346  std::vector<T*> row_vals;
1347  row_vals.reserve(1000);
1348  {
1349  // Copy first row. Can't use first_row_vals, as may be a vector of bool...
1350  T* row = vnl_c_vector<T>::allocate_T(colz);
1351  for (unsigned int k = 0; k < colz; ++k)
1352  row[k] = first_row_vals[k];
1353  row_vals.push_back(row);
1354  }
1355 
1356  while (true)
1357  {
1358  T* row = vnl_c_vector<T>::allocate_T(colz);
1359  if (row == nullptr) {
1360  std::cerr << "vnl_matrix<T>::read_ascii: Error, Out of memory on row "
1361  << row_vals.size() << std::endl;
1362  return false;
1363  }
1364  s >> row[0];
1365  if (!s.good())
1366  {
1367  vnl_c_vector<T>::deallocate(row, colz);
1368  break;
1369  }
1370  for (unsigned int k = 1; k < colz; ++k) {
1371  if (s.eof()) {
1372  std::cerr << "vnl_matrix<T>::read_ascii: Error, EOF on row "
1373  << row_vals.size() << ", column " << k << std::endl;
1374 
1375  return false;
1376  }
1377  s >> row[k];
1378  if (s.fail()) {
1379  std::cerr << "vnl_matrix<T>::read_ascii: Error, row "
1380  << row_vals.size() << " failed on column " << k << std::endl;
1381  return false;
1382  }
1383  }
1384  row_vals.push_back(row);
1385  }
1386 
1387  std::size_t rowz = row_vals.size();
1388 
1389  if (debug)
1390  std::cerr << rowz << " rows.\n";
1391 
1392  set_size((unsigned int)rowz, (unsigned int)colz);
1393 
1394  T* p = this->data[0];
1395  for (unsigned int i = 0; i < rowz; ++i) {
1396  for (unsigned int j = 0; j < colz; ++j)
1397  *p++ = row_vals[i][j];
1398  /*if (i>0)*/ vnl_c_vector<T>::deallocate(row_vals[i], colz);
1399  }
1400 
1401  return true;
1402 }
1403 
1404 //: Read a vnl_matrix from an ascii std::istream.
1405 // Automatically determines file size if the input matrix has zero size.
1406 // This is a static method so you can type
1407 // <verb>
1408 // vnl_matrix<float> M = vnl_matrix<float>::read(cin);
1409 // </verb>
1410 // which many people prefer to the ">>" alternative.
1411 template <class T>
1413 {
1414  vnl_matrix<T> M;
1415  s >> M;
1416  return M;
1417 }
1418 
1419 template <class T>
1421 {
1422  std::swap(this->num_rows, that.num_rows);
1423  std::swap(this->num_cols, that.num_cols);
1424  std::swap(this->data, that.data);
1425 }
1426 
1427 //: Reverse order of rows. Name is from Matlab, meaning "flip upside down".
1428 template <class T>
1430 {
1431  const unsigned int n = this->rows();
1432  const unsigned int colz = this->columns();
1433 
1434  const unsigned int m = n / 2;
1435  for (unsigned int r = 0; r < m; ++r) {
1436  const unsigned int r1 = r;
1437  const unsigned int r2 = n - 1 - r;
1438  for (unsigned int c = 0; c < colz; ++c) {
1439  const T tmp = (*this)(r1, c);
1440  (*this)(r1, c) = (*this)(r2, c);
1441  (*this)(r2, c) = tmp;
1442  }
1443  }
1444  return *this;
1445 }
1446 
1447 //: Reverse order of columns.
1448 template <class T>
1450 {
1451  const unsigned int n = this->cols();
1452  const unsigned int rowz = this->rows();
1453 
1454  const unsigned int m = n / 2;
1455  for (unsigned int c = 0; c < m; ++c) {
1456  const unsigned int c1 = c;
1457  const unsigned int c2 = n - 1 - c;
1458  for (unsigned int r = 0; r < rowz; ++r) {
1459  const T tmp = (*this)(r, c1);
1460  (*this)(r, c1) = (*this)(r, c2);
1461  (*this)(r, c2) = tmp;
1462  }
1463  }
1464  return *this;
1465 }
1466 
1467 // || M || = \max \sum | M |
1468 // 1 j i ij
1469 template <class T>
1471 {
1472  abs_t max = 0;
1473  for (unsigned int j=0; j<this->num_cols; ++j) {
1474  abs_t tmp = 0;
1475  for (unsigned int i=0; i<this->num_rows; ++i)
1476  tmp += vnl_math::abs(this->data[i][j]);
1477  if (tmp > max)
1478  max = tmp;
1479  }
1480  return max;
1481 }
1482 
1483 // || M || = \max \sum | M |
1484 // oo i j ij
1485 template <class T>
1487 {
1488  abs_t max = 0;
1489  for (unsigned int i=0; i<this->num_rows; ++i) {
1490  abs_t tmp = 0;
1491  for (unsigned int j=0; j<this->num_cols; ++j)
1492  tmp += vnl_math::abs(this->data[i][j]);
1493  if (tmp > max)
1494  max = tmp;
1495  }
1496  return max;
1497 }
1498 
1499 template <class doublereal> // ideally, char* should be bool* - PVr
1500 int vnl_inplace_transpose(doublereal *a, unsigned m, unsigned n, char* move, unsigned iwrk)
1501 {
1502  doublereal b, c;
1503  int k = m * n - 1;
1504  int iter, i1, i2, im, i1c, i2c, ncount, max_;
1505 
1506 // *****
1507 // ALGORITHM 380 - REVISED
1508 // *****
1509 // A IS A ONE-DIMENSIONAL ARRAY OF LENGTH MN=M*N, WHICH
1510 // CONTAINS THE MXN MATRIX TO BE TRANSPOSED (STORED
1511 // COLUMNWISE). MOVE IS A ONE-DIMENSIONAL ARRAY OF LENGTH IWRK
1512 // USED TO STORE INFORMATION TO SPEED UP THE PROCESS. THE
1513 // VALUE IWRK=(M+N)/2 IS RECOMMENDED. IOK INDICATES THE
1514 // SUCCESS OR FAILURE OF THE ROUTINE.
1515 // NORMAL RETURN IOK=0
1516 // ERRORS IOK=-2 ,IWRK NEGATIVE OR ZERO
1517 // IOK.GT.0, (SHOULD NEVER OCCUR),IN THIS CASE
1518 // WE SET IOK EQUAL TO THE FINAL VALUE OF ITER WHEN THE SEARCH
1519 // IS COMPLETED BUT SOME LOOPS HAVE NOT BEEN MOVED
1520 // NOTE * MOVE(I) WILL STAY ZERO FOR FIXED POINTS
1521 
1522  if (m < 2 || n < 2)
1523  return 0; // JUST RETURN IF MATRIX IS SINGLE ROW OR COLUMN
1524  if (iwrk < 1)
1525  return -2; // ERROR RETURN
1526  if (m == n) {
1527  // IF MATRIX IS SQUARE, EXCHANGE ELEMENTS A(I,J) AND A(J,I).
1528  for (unsigned i = 0; i < n; ++i)
1529  for (unsigned j = i+1; j < n; ++j) {
1530  i1 = i + j * n;
1531  i2 = j + i * m;
1532  b = a[i1];
1533  a[i1] = a[i2];
1534  a[i2] = b;
1535  }
1536  return 0; // NORMAL RETURN
1537  }
1538  ncount = 2;
1539  for (unsigned i = 0; i < iwrk; ++i)
1540  move[i] = char(0); // false;
1541  if (m > 2 && n > 2) {
1542  // CALCULATE THE NUMBER OF FIXED POINTS, EUCLIDS ALGORITHM FOR GCD(M-1,N-1).
1543  int ir2 = m - 1;
1544  int ir1 = n - 1;
1545  int ir0 = ir2 % ir1;
1546  while (ir0 != 0) {
1547  ir2 = ir1;
1548  ir1 = ir0;
1549  ir0 = ir2 % ir1;
1550  }
1551  ncount += ir1 - 1;
1552  }
1553 // SET INITIAL VALUES FOR SEARCH
1554  iter = 1;
1555  im = m;
1556 // AT LEAST ONE LOOP MUST BE RE-ARRANGED
1557  goto L80;
1558 // SEARCH FOR LOOPS TO REARRANGE
1559 L40:
1560  max_ = k - iter;
1561  ++iter;
1562  if (iter > max_)
1563  return iter; // error return
1564  im += m;
1565  if (im > k)
1566  im -= k;
1567  i2 = im;
1568  if (iter == i2)
1569  goto L40;
1570  if (iter <= (int)iwrk) {
1571  if (move[iter-1])
1572  goto L40;
1573  else
1574  goto L80;
1575  }
1576  while (i2 > iter && i2 < max_) {
1577  i1 = i2;
1578  i2 = m * i1 - k * (i1 / n);
1579  }
1580  if (i2 != iter)
1581  goto L40;
1582 // REARRANGE THE ELEMENTS OF A LOOP AND ITS COMPANION LOOP
1583 L80:
1584  i1 = iter;
1585  b = a[i1];
1586  i1c = k - iter;
1587  c = a[i1c];
1588  while (true) {
1589  i2 = m * i1 - k * (i1 / n);
1590  i2c = k - i2;
1591  if (i1 <= (int)iwrk)
1592  move[i1-1] = '1'; // true;
1593  if (i1c <= (int)iwrk)
1594  move[i1c-1] = '1'; // true;
1595  ncount += 2;
1596  if (i2 == iter)
1597  break;
1598  if (i2+iter == k) {
1599  doublereal d = b; b = c; c = d; // interchange b and c
1600  break;
1601  }
1602  a[i1] = a[i2];
1603  a[i1c] = a[i2c];
1604  i1 = i2;
1605  i1c = i2c;
1606  }
1607 // FINAL STORE AND TEST FOR FINISHED
1608  a[i1] = b;
1609  a[i1c] = c;
1610  if (ncount > k)
1611  return 0; // NORMAL RETURN
1612  goto L40;
1613 } /* dtrans_ */
1614 
1615 
1616 //: Transpose matrix M in place.
1617 // Works for rectangular matrices using an enormously clever algorithm from ACM TOMS.
1618 template <class T>
1620 {
1621  unsigned m = rows();
1622  unsigned n = columns();
1623  unsigned iwrk = (m+n)/2;
1624  std::vector<char> move(iwrk);
1625 
1626  int iok = ::vnl_inplace_transpose(data_block(), n, m, &move[0], iwrk);
1627  if (iok != 0)
1628  std::cerr << __FILE__ " : inplace_transpose() -- iok = " << iok << '\n';
1629 
1630  this->num_rows = n;
1631  this->num_cols = m;
1632 
1633  // row pointers. we have to reallocate even when n<=m because
1634  // vnl_c_vector<T>::deallocate needs to know n_when_allocatod.
1635  {
1636  T *tmp = data[0];
1639  for (unsigned i=0; i<n; ++i)
1640  data[i] = tmp + i * m;
1641  }
1642  return *this;
1643 }
1644 
1645 //------------------------------------------------------------------------------
1646 
1647 #define VNL_MATRIX_INSTANTIATE(T) \
1648 template class VNL_EXPORT vnl_matrix<T >; \
1649 template VNL_EXPORT vnl_matrix<T > operator-(T const &, vnl_matrix<T > const &); \
1650 /*template VNL_EXPORT vnl_matrix<T > operator+(T const &, vnl_matrix<T > const &) ; */ \
1651 /*template VNL_EXPORT vnl_matrix<T > operator*(T const &, vnl_matrix<T > const &) ; */ \
1652 template VNL_EXPORT T dot_product(vnl_matrix<T > const &, vnl_matrix<T > const &); \
1653 template VNL_EXPORT T inner_product(vnl_matrix<T > const &, vnl_matrix<T > const &); \
1654 template VNL_EXPORT T cos_angle(vnl_matrix<T > const &, vnl_matrix<T > const &); \
1655 template VNL_EXPORT vnl_matrix<T > element_product(vnl_matrix<T > const &, vnl_matrix<T > const &); \
1656 template VNL_EXPORT vnl_matrix<T > element_quotient(vnl_matrix<T > const &, vnl_matrix<T > const &); \
1657 template VNL_EXPORT int vnl_inplace_transpose(T*, unsigned, unsigned, char*, unsigned); \
1658 template VNL_EXPORT std::ostream & operator<<(std::ostream &, vnl_matrix<T > const &); \
1659 template VNL_EXPORT std::istream & operator>>(std::istream &, vnl_matrix<T > &)
1660 
1661 #endif // vnl_matrix_hxx_
vnl_matrix< T > & operator=(T const &v)
Set all elements to value v.
Definition: vnl_matrix.h:274
static void deallocate(T **, const std::size_t n_when_allocated)
vnl_c_vector< std::complex< double > >::abs_t abs_t
Type def for norms.
Definition: vnl_matrix.h:499
vnl_matrix & set_diagonal(vnl_vector< T > const &)
Sets the diagonal elements of this matrix to the specified list of values.
Definition: vnl_matrix.hxx:440
vnl_vector_fixed< T, n > element_quotient(const vnl_vector_fixed< T, n > &a, const vnl_vector_fixed< T, n > &b)
unsigned int cols() const
Return the number of columns.
Definition: vnl_matrix.h:183
vnl_matrix< T > get_columns(vnl_vector< unsigned int > i) const
Get a matrix composed of columns from the indices specified in the supplied vector.
T & operator()(unsigned r, unsigned c)
Access an element for reading or writing.
Definition: vnl_matrix.hxx:455
vnl_matrix()
Default constructor creates an empty matrix of size 0,0.
Definition: vnl_matrix.h:117
bool is_zero() const
Return true if all elements equal to zero.
static void conjugate(T const *, T *, unsigned)
vnl_matrix & copy_in(T const *)
Fills (laminates) this matrix with the given data, then returns it.
Definition: vnl_matrix.hxx:837
bool has_nans() const
Return true if matrix contains NaNs.
vnl_matrix & fill_diagonal(T const &)
Sets all diagonal elements of matrix to specified value; returns "*this".
Definition: vnl_matrix.hxx:430
void swap(double &a, double &b)
An ordinary mathematical matrix.
vnl_matrix< T > get_n_columns(unsigned colstart, unsigned n) const
Get n columns beginning at colstart.
Definition: vnl_matrix.hxx:946
vnl_matrix< T > conjugate_transpose() const
Return conjugate transpose.
Definition: vnl_matrix.hxx:694
vnl_matrix< T > & operator+=(T value)
Add rhs to each element of lhs matrix in situ.
Definition: vnl_matrix.hxx:533
vnl_vector< T > apply_rowwise(T(*f)(vnl_vector< T > const &)) const
Make a vector by applying a function across rows.
Definition: vnl_matrix.hxx:655
vnl_matrix< T > transpose() const
Return transpose.
Definition: vnl_matrix.hxx:682
vnl_vector< T > flatten_row_major() const
Flatten row-major (C-style).
vnl_matrix & scale_row(unsigned row, T value)
Scales elements in given row by a factor T, and returns "*this".
Definition: vnl_matrix.hxx:907
vnl_matrix< T > operator-() const
Negate all elements of matrix.
Definition: vnl_matrix.hxx:624
Templated zero/one/precision.
vnl_matrix< T > extract(unsigned r, unsigned c, unsigned top=0, unsigned left=0) const
Extract a sub-matrix of size r x c, starting at (top,left).
Definition: vnl_matrix.hxx:728
#define m
Definition: vnl_vector.h:43
vnl_vector_fixed< T, n > element_product(const vnl_vector_fixed< T, n > &a, const vnl_vector_fixed< T, n > &b)
size_t size() const
Return the length, number of elements, dimension of this vector.
Definition: vnl_vector.h:126
#define vnl_matrix_alloc_blah()
Definition: vnl_matrix.hxx:102
vnl_matrix< T > & operator *=(T value)
Scalar multiplication in situ of lhs matrix by rhs.
Definition: vnl_matrix.hxx:551
Namespace with standard math functions.
iterator begin()
Iterator pointing to start of data.
Definition: vnl_matrix.h:620
void vnl_error_matrix_row_index(char const *fcn, int r)
Raise exception for invalid row index.
Definition: vnl_error.cxx:51
abs_t operator_one_norm() const
void vnl_error_matrix_col_index(char const *fcn, int c)
Raise exception for invalid col index.
Definition: vnl_error.cxx:61
enum VNL_EXPORT vnl_matrix_type
Definition: vnl_matrix.h:80
vnl_bignum squared_magnitude(vnl_bignum const &x)
Definition: vnl_bignum.h:433
std::ostream & operator<<(std::ostream &s, vnl_decnum const &r)
decimal output.
Definition: vnl_decnum.h:393
static void apply(T const *, unsigned, T(*f)(T), T *v_out)
vnl_matrix & normalize_rows()
Normalizes each row so it is a unit vector, and returns "*this".
Definition: vnl_matrix.hxx:864
T get(unsigned r, unsigned c) const
get element with boundary checks if error checking is on.
Definition: vnl_matrix.h:684
vnl_matrix< T > & operator-=(T value)
Subtract rhs from each element of lhs matrix in situ.
Definition: vnl_matrix.hxx:542
#define v
Definition: vnl_vector.h:42
vnl_matrix & normalize_columns()
Normalizes each column so it is a unit vector, and returns "*this".
Definition: vnl_matrix.hxx:886
vnl_matrix & inplace_transpose()
Transposes this matrix efficiently, and returns it.
void vnl_error_vector_dimension(char const *fcn, int l1, int l2)
Raise exception for invalid dimension.
Definition: vnl_error.cxx:30
bool isnan(vnl_bignum const &)
Definition: vnl_bignum.h:435
void assert_size_internal(unsigned r, unsigned c) const
Abort unless M has the given size.
vnl_matrix & set_column(unsigned i, T const *v)
Set the elements of the i'th column to v[i] (No bounds checking).
vnl_matrix< T > get_rows(vnl_vector< unsigned int > i) const
Get a matrix composed of rows from the indices specified in the supplied vector.
Definition: vnl_matrix.hxx:994
int vnl_inplace_transpose(doublereal *a, unsigned m, unsigned n, char *move, unsigned iwrk)
void assert_finite_internal() const
Abort if any element of M is inf or nan.
vnl_matrix< T > get_n_rows(unsigned rowstart, unsigned n) const
Get n rows beginning at rowstart.
Definition: vnl_matrix.hxx:933
static T * allocate_T(const std::size_t n)
vnl_matrix< T > & operator/=(T value)
Scalar division of lhs matrix in situ by rhs.
Definition: vnl_matrix.hxx:560
unsigned num_rows
Definition: vnl_matrix.h:664
An ordinary mathematical matrix.
Definition: vnl_adjugate.h:22
bool is_identity() const
Return true if all elements equal to identity.
void clear()
Make the matrix as if it had been default-constructed.
Definition: vnl_matrix.hxx:375
vnl_matrix & fill(T const &)
Sets all elements of matrix to specified value, and returns "*this".
Definition: vnl_matrix.hxx:419
static T inner_product(T const *, T const *, unsigned)
conjugate second.
vnl_matrix & flipud()
Reverses the order of rows, and returns "*this".
unsigned num_cols
Definition: vnl_matrix.h:665
bool isfinite(vnl_bignum const &x)
Definition: vnl_bignum.h:436
T dot_product(const vnl_vector_fixed< T, n > &a, const vnl_vector_fixed< T, n > &b)
static T ** allocate_Tptr(const std::size_t n)
Memory allocation.
void swap(vnl_matrix< T > &that)
Swap this matrix with that matrix.
void print(std::ostream &os) const
Print matrix to os in some hopefully sensible format.
Definition: vnl_matrix.hxx:500
vnl_decnum max(vnl_decnum const &x, vnl_decnum const &y)
Definition: vnl_decnum.h:412
Mathematical vector class, templated by type of element.
Definition: vnl_fwd.h:16
abs_t operator_inf_norm() const
T const * data_block() const
Access the contiguous block storing the elements in the matrix row-wise. O(1).
Definition: vnl_matrix.h:601
#define vnl_matrix_construct_hack()
Definition: vnl_matrix.hxx:99
vnl_vector< T > get_column(unsigned c) const
Get a vector equal to the given column.
Definition: vnl_matrix.hxx:977
bool operator_eq(vnl_matrix< T > const &rhs) const
Return true if *this == rhs.
vnl_matrix & set_columns(unsigned starting_column, vnl_matrix< T > const &M)
Set columns to those in M, starting at starting_column, then return *this.
vnl_vector< T > flatten_column_major() const
Flatten column-major (Fortran-style).
T cos_angle(vnl_matrix< T > const &a, vnl_matrix< T > const &b)
Definition: vnl_matrix.hxx:785
vnl_vector< T > get_row(unsigned r) const
Get a vector equal to the given row.
Definition: vnl_matrix.hxx:962
static vnl_matrix< T > read(std::istream &s)
Read a vnl_matrix from an ascii std::istream, automatically determining file size if the input matrix...
vnl_matrix & set_row(unsigned i, T const *v)
Set the elements of the i'th row to v[i] (No bounds checking).
vnl_matrix_null
Definition: vnl_matrix.h:82
vnl_bignum abs(vnl_bignum const &x)
Definition: vnl_bignum.h:432
unsigned int size() const
Return the total number of elements stored by the matrix.
Definition: vnl_matrix.h:176
void vnl_error_matrix_dimension(char const *fcn, int r1, int c1, int r2, int c2)
Raise exception for invalid dimensions.
Definition: vnl_error.cxx:70
void destroy()
Delete data.
Definition: vnl_matrix.hxx:369
vnl_vector< T > apply_columnwise(T(*f)(vnl_vector< T > const &)) const
Make a vector by applying a function across columns.
Definition: vnl_matrix.hxx:667
unsigned int rows() const
Return the number of rows.
Definition: vnl_matrix.h:179
bool is_equal(vnl_matrix< T > const &rhs, double tol) const
Return true if all elements of both matrices are equal, within given tolerance.
VNL_EXPORT std::istream & operator>>(std::istream &s, vnl_decnum &r)
decimal input.
Definition: vnl_decnum.cxx:411
vnl_bignum operator-(vnl_bignum const &r1, vnl_bignum const &r2)
Returns the difference of two bignum numbers.
Definition: vnl_bignum.h:290
vnl_matrix< T > apply(T(*f)(T)) const
Make a new matrix by applying function to each element.
Definition: vnl_matrix.hxx:644
bool is_finite() const
Return true if finite.
T inner_product(vnl_matrix< T > const &m1, vnl_matrix< T > const &m2)
Hermitian inner product.
Definition: vnl_matrix.hxx:771
vnl_matrix & scale_column(unsigned col, T value)
Scales elements in given column by a factor T, and returns "*this".
Definition: vnl_matrix.hxx:920
vnl_vector< T > get_diagonal() const
Return a vector with the content of the (main) diagonal.
bool set_size(unsigned r, unsigned c)
Resize to r rows by c columns. Old data lost.
Definition: vnl_matrix.hxx:390
bool read_ascii(std::istream &s)
Read a vnl_matrix from an ascii std::istream.
static T dot_product(T const *, T const *, unsigned)
void copy_out(T *) const
Fills the given array with this matrix.
Definition: vnl_matrix.hxx:846
vnl_matrix< T > & update(vnl_matrix< T > const &, unsigned top=0, unsigned left=0)
Set values of this matrix to those of M, starting at [top,left].
Definition: vnl_matrix.hxx:707
~vnl_matrix()
Matrix destructor.
Definition: vnl_matrix.hxx:359
T get(size_t i) const
Get value at element i.
Definition: vnl_vector.h:415
unsigned int columns() const
Return the number of columns.
Definition: vnl_matrix.h:187
Math on blocks of memory.
vnl_matrix & set_identity()
Sets this matrix to an identity matrix, then returns "*this".
Definition: vnl_matrix.hxx:853
#define vnl_matrix_free_blah
Definition: vnl_matrix.hxx:120
void put(unsigned r, unsigned c, T const &)
set element with boundary checks if error checking is on.
Definition: vnl_matrix.h:700
vnl_matrix & fliplr()
Reverses the order of columns, and returns "*this".