libROM  v1.0
Data-driven physical simulation library
Matrix.h
1 
11 // Description: A simple, parallel Matrix class with the utility needed to
12 // support the basis generation methods of this library. A
13 // distributed Matrix has its rows distributed across processors.
14 
15 #ifndef included_Matrix_h
16 #define included_Matrix_h
17 
18 #include "Vector.h"
19 #include <vector>
20 #include <complex>
21 #include <string>
22 
23 namespace CAROM {
24 
30 class Matrix
31 {
32 public:
34  Matrix();
35 
51  Matrix(
52  int num_rows,
53  int num_cols,
54  bool distributed,
55  bool randomized = false);
56 
75  Matrix(
76  double* mat,
77  int num_rows,
78  int num_cols,
79  bool distributed,
80  bool copy_data = true);
81 
87  Matrix(
88  const Matrix& other);
89 
93  ~Matrix();
94 
102  Matrix&
103  operator = (
104  const Matrix& rhs);
105 
113  Matrix&
114  operator = (
115  const double a);
116 
124  Matrix&
125  operator += (
126  const Matrix& rhs);
127 
135  Matrix&
136  operator -= (
137  const Matrix& rhs);
138 
146  void
148  int num_rows,
149  int num_cols)
150  {
151  int new_size = num_rows*num_cols;
152  if (new_size > d_alloc_size) {
153  if (!d_owns_data) {
154  CAROM_ERROR("Can not reallocate externally owned storage.");
155  }
156  if (d_mat) {
157  delete [] d_mat;
158  }
159 
160  // Allocate new array and initialize all values to zero.
161  d_mat = new double [new_size] {0.0};
162  d_alloc_size = new_size;
163  }
164  d_num_rows = num_rows;
165  d_num_cols = num_cols;
166  if (d_distributed) {
167  calculateNumDistributedRows();
168  }
169  }
170 
176  bool
177  distributed() const
178  {
179  return d_distributed;
180  }
181 
186  bool balanced() const;
187 
193  int
194  numRows() const
195  {
196  return d_num_rows;
197  }
198 
204  int
206  {
207  if (!d_distributed) {
208  return d_num_rows;
209  }
210  return d_num_distributed_rows;
211  }
212 
219  int
220  numColumns() const
221  {
222  return d_num_cols;
223  }
224 
234  Matrix*
235  getFirstNColumns(int n) const;
236 
247  void
249  int n,
250  Matrix*& result) const;
251 
262  void
264  int n,
265  Matrix& result) const;
266 
282  Matrix*
284  const Matrix& other) const
285  {
286  Matrix* result = 0;
287  mult(other, result);
288  return result;
289  }
290 
307  Matrix*
309  const Matrix* other) const
310  {
311  CAROM_VERIFY(other != 0);
312  return mult(*other);
313  }
314 
332  void
333  mult(
334  const Matrix& other,
335  Matrix*& result) const;
336 
353  void
354  mult(
355  const Matrix& other,
356  Matrix& result) const;
357 
373  Vector*
375  const Vector& other) const
376  {
377  Vector* result = 0;
378  mult(other, result);
379  return result;
380  }
381 
398  Vector*
400  const Vector* other) const
401  {
402  CAROM_VERIFY(other != 0);
403  return mult(*other);
404  }
405 
423  void
424  mult(
425  const Vector& other,
426  Vector*& result) const;
427 
444  void
445  mult(
446  const Vector& other,
447  Vector& result) const;
448 
466  void
468  int this_row,
469  const Vector& other,
470  Vector& result) const;
471 
489  void
491  int this_row,
492  Vector& other) const;
493 
506  Matrix*
508  const Matrix& other) const
509  {
510  Matrix* result = 0;
511  elementwise_mult(other, result);
512  return result;
513  }
514 
528  Matrix*
530  const Matrix* other) const
531  {
532  CAROM_VERIFY(other != 0);
533  return elementwise_mult(*other);
534  }
535 
548  void
550  const Matrix& other,
551  Matrix*& result) const;
552 
565  void
567  const Matrix& other,
568  Matrix& result) const;
569 
575  Matrix*
577  {
578  Matrix* result = 0;
579  elementwise_square(result);
580  return result;
581  }
582 
589  void
591  Matrix*& result) const;
592 
599  void
601  Matrix& result) const;
602 
620  void
621  multPlus(
622  Vector& a,
623  const Vector& b,
624  double c) const;
625 
641  Matrix*
643  const Matrix& other) const
644  {
645  Matrix* result = 0;
646  transposeMult(other, result);
647  return result;
648  }
649 
666  Matrix*
668  const Matrix* other) const
669  {
670  CAROM_VERIFY(other != 0);
671  return transposeMult(*other);
672  }
673 
690  void
692  const Matrix& other,
693  Matrix*& result) const;
694 
710  void
712  const Matrix& other,
713  Matrix& result) const;
714 
730  Vector*
732  const Vector& other) const
733  {
734  Vector* result = 0;
735  transposeMult(other, result);
736  return result;
737  }
738 
755  Vector*
757  const Vector* other) const
758  {
759  CAROM_VERIFY(other != 0);
760  return transposeMult(*other);
761  }
762 
779  void
781  const Vector& other,
782  Vector*& result) const;
783 
800  void
802  const Vector& other,
803  Vector& result) const;
804 
813  Matrix*
814  inverse() const
815  {
816  Matrix* result = 0;
817  inverse(result);
818  return result;
819  }
820 
835  void
836  inverse(
837  Matrix*& result) const;
838 
851  void
852  inverse(
853  Matrix& result) const;
854 
861  void
862  inverse();
863 
869  Vector*
870  getColumn(int column) const
871  {
872  Vector* result = 0;
873  getColumn(column, result);
874  return result;
875  }
876 
882  void
883  getColumn(int column,
884  Vector*& result) const;
885 
891  void
892  getColumn(int column,
893  Vector& result) const;
894 
902  void transpose();
903 
913  void transposePseudoinverse();
914 
920  Matrix*
921  qr_factorize() const;
922 
938  void
939  qrcp_pivots_transpose(int* row_pivot,
940  int* row_pivot_owner,
941  int pivots_requested) const;
942 
957  void
958  orthogonalize(bool double_pass = false, double zero_tol = 1.0e-15);
959 
979  void
980  orthogonalize_last(int ncols = -1, bool double_pass = false,
981  double zero_tol = 1.0e-15);
982 
986  void
988 
992  void
994 
1006  const double&
1008  int row,
1009  int col) const
1010  {
1011  CAROM_ASSERT((0 <= row) && (row < numRows()));
1012  CAROM_ASSERT((0 <= col) && (col < numColumns()));
1013  return d_mat[row*d_num_cols+col];
1014  }
1015 
1029  double&
1031  int row,
1032  int col)
1033  {
1034  CAROM_ASSERT((0 <= row) && (row < numRows()));
1035  CAROM_ASSERT((0 <= col) && (col < numColumns()));
1036  return d_mat[row*d_num_cols+col];
1037  }
1038 
1049  const double& operator() (int row, int col) const
1050  {
1051  return item(row, col);
1052  }
1053 
1066  double& operator() (int row, int col)
1067  {
1068  return item(row, col);
1069  }
1070 
1077  void print(const char * prefix) const;
1078 
1085  void write(const std::string& base_file_name) const;
1086 
1093  void read(const std::string& base_file_name);
1094 
1102  void local_read(const std::string& base_file_name, int rank);
1103 
1107  double *getData() const
1108  {
1109  return d_mat;
1110  }
1111 
1122  void distribute(const int &local_num_rows);
1123 
1131  void gather();
1132 
1133 private:
1139  void
1140  calculateNumDistributedRows();
1141 
1156  void
1157  qrcp_pivots_transpose_serial(int* row_pivot,
1158  int* row_pivot_owner,
1159  int pivots_requested) const;
1160 
1175  void
1176  qrcp_pivots_transpose_distributed(int* row_pivot,
1177  int* row_pivot_owner,
1178  int pivots_requested) const;
1179 
1180  void
1181  qrcp_pivots_transpose_distributed_scalapack(int* row_pivot,
1182  int* row_pivot_owner,
1183  int pivots_requested) const;
1184 
1200  void
1201  qrcp_pivots_transpose_distributed_elemental(int* row_pivot,
1202  int* row_pivot_owner,
1203  int pivots_requested) const;
1204 
1220  void
1221  qrcp_pivots_transpose_distributed_elemental_balanced
1222  (int* row_pivot, int* row_pivot_owner, int pivots_requested) const;
1223 
1239  void
1240  qrcp_pivots_transpose_distributed_elemental_unbalanced
1241  (int* row_pivot, int* row_pivot_owner, int pivots_requested) const;
1242 
1246  double* d_mat;
1247 
1251  int d_num_rows;
1252 
1256  int d_num_distributed_rows;
1257 
1264  int d_num_cols;
1265 
1271  int d_alloc_size;
1272 
1278  bool d_distributed;
1279 
1283  int d_num_procs;
1284 
1288  int d_rank;
1289 
1296  bool d_owns_data;
1297 };
1298 
1313 // NOTE(goxberry@gmail.com; oxberry1@llnl.gov): Passing by value is
1314 // supposed to be more idiomatic C++11; consider passing by value
1315 // instead of passing by reference.
1316 Matrix outerProduct(const Vector &v, const Vector &w);
1317 
1332 Matrix DiagonalMatrixFactory(const Vector &v);
1333 
1348 Matrix IdentityMatrixFactory(const Vector &v);
1349 
1354 struct EigenPair {
1355 
1360 
1364  std::vector<double> eigs;
1365 };
1366 
1372 
1377 
1382 
1386  std::vector<std::complex<double>> eigs;
1387 };
1388 
1394 {
1399 
1404 
1409 };
1410 
1422 void SerialSVD(Matrix* A,
1423  Matrix* U,
1424  Vector* S,
1425  Matrix* V);
1426 
1437 
1447 
1457 
1458 Matrix* SpaceTimeProduct(const CAROM::Matrix* As, const CAROM::Matrix* At,
1459  const CAROM::Matrix* Bs, const CAROM::Matrix* Bt,
1460  const std::vector<double> *tscale=NULL,
1461  const bool At0at0=false, const bool Bt0at0=false,
1462  const bool lagB=false, const bool skip0=false);
1463 }
1464 
1465 #endif
void local_read(const std::string &base_file_name, int rank)
read a single rank of a distributed Matrix into (a) HDF file(s).
Definition: Matrix.cpp:1009
void setSize(int num_rows, int num_cols)
Sets the number of rows and columns of the matrix and reallocates storage if needed....
Definition: Matrix.h:147
Matrix * getFirstNColumns(int n) const
Get the first N columns of a matrix.
Definition: Matrix.cpp:256
double * getData() const
Get the matrix data as a pointer.
Definition: Matrix.h:1107
bool balanced() const
Returns true if rows of matrix are load-balanced.
Definition: Matrix.cpp:203
Matrix * inverse() const
Computes and returns the inverse of this.
Definition: Matrix.h:814
const double & operator()(int row, int col) const
Const Matrix member access.
Definition: Matrix.h:1049
void transposePseudoinverse()
Computes the transposePseudoinverse of this.
Definition: Matrix.cpp:882
bool distributed() const
Returns true if the Matrix is distributed.
Definition: Matrix.h:177
void rescale_rows_max()
Rescale every matrix row by its maximum absolute value.
Definition: Matrix.cpp:1903
Matrix & operator=(const Matrix &rhs)
Assignment operator.
Definition: Matrix.cpp:172
void write(const std::string &base_file_name) const
write Matrix into (a) HDF file(s).
Definition: Matrix.cpp:937
int numColumns() const
Returns the number of columns in the Matrix. This method will return the same value from each process...
Definition: Matrix.h:220
void distribute(const int &local_num_rows)
Distribute this matrix rows among MPI processes, based on the specified local number of rows....
Definition: Matrix.cpp:1046
void gather()
Gather all the distributed rows among MPI processes. This becomes not distributed after this function...
Definition: Matrix.cpp:1073
Matrix * elementwise_mult(const Matrix &other) const
Multiplies two matrices element-wise.
Definition: Matrix.h:507
void orthogonalize_last(int ncols=-1, bool double_pass=false, double zero_tol=1.0e-15)
Orthonormalizes the matrix's last column, assuming the previous columns are already orthonormal.
Definition: Matrix.cpp:1853
Matrix & operator+=(const Matrix &rhs)
Addition operator.
Definition: Matrix.cpp:183
Matrix * mult(const Matrix *other) const
Multiplies this Matrix with other and returns the product, pointer version.
Definition: Matrix.h:308
int numDistributedRows() const
Returns the number of rows of the Matrix across all processors.
Definition: Matrix.h:205
Vector * getColumn(int column) const
Returns a column of the matrix (not owned by Matrix).
Definition: Matrix.h:870
Matrix & operator-=(const Matrix &rhs)
Subtraction operator.
Definition: Matrix.cpp:193
void orthogonalize(bool double_pass=false, double zero_tol=1.0e-15)
Orthonormalizes the matrix.
Definition: Matrix.cpp:1803
void multPlus(Vector &a, const Vector &b, double c) const
Computes a += this*b*c.
Definition: Matrix.cpp:539
Vector * transposeMult(const Vector &other) const
Multiplies the transpose of this Matrix with other and returns the product, reference version.
Definition: Matrix.h:731
double & item(int row, int col)
Non-const Matrix member access. Matrix data is stored in row-major format.
Definition: Matrix.h:1030
void print(const char *prefix) const
print Matrix into (a) ascii file(s).
Definition: Matrix.cpp:918
Matrix * transposeMult(const Matrix &other) const
Multiplies the transpose of this Matrix with other and returns the product, reference version.
Definition: Matrix.h:642
void pointwise_mult(int this_row, const Vector &other, Vector &result) const
Multiplies a specified row of this Matrix with other pointwise.
Definition: Matrix.cpp:418
Matrix * mult(const Matrix &other) const
Multiplies this Matrix with other and returns the product, reference version.
Definition: Matrix.h:283
Vector * mult(const Vector *other) const
Multiplies this Matrix with other and returns the product, pointer version.
Definition: Matrix.h:399
~Matrix()
Destructor.
Definition: Matrix.cpp:164
Matrix * elementwise_mult(const Matrix *other) const
Multiplies two matrices element-wise and returns the product, pointer version.
Definition: Matrix.h:529
const double & item(int row, int col) const
Const Matrix member access. Matrix data is stored in row-major format.
Definition: Matrix.h:1007
Vector * mult(const Vector &other) const
Multiplies this Matrix with other and returns the product, reference version.
Definition: Matrix.h:374
Matrix * transposeMult(const Matrix *other) const
Multiplies the transpose of this Matrix with other and returns the product, pointer version.
Definition: Matrix.h:667
int numRows() const
Returns the number of rows of the Matrix on this processor.
Definition: Matrix.h:194
void read(const std::string &base_file_name)
read Matrix into (a) HDF file(s).
Definition: Matrix.cpp:969
Matrix * elementwise_square() const
Square every element in the matrix.
Definition: Matrix.h:576
Vector * transposeMult(const Vector *other) const
Multiplies the transpose of this Matrix with other and returns the product, pointer version.
Definition: Matrix.h:756
void transpose()
Replaces this Matrix with its transpose (in place), in the serial square case only.
Definition: Matrix.cpp:824
Matrix * qr_factorize() const
Computes and returns the Q of the QR factorization of this.
Definition: Matrix.cpp:1125
void rescale_cols_max()
Rescale every matrix column by its maximum absolute value.
Definition: Matrix.cpp:1930
void qrcp_pivots_transpose(int *row_pivot, int *row_pivot_owner, int pivots_requested) const
Compute the leading numColumns() column pivots from a QR decomposition with column pivots (QRCP) of t...
Definition: Matrix.cpp:1198
Matrix IdentityMatrixFactory(const Vector &v)
Factory function to make an identity matrix with rows distributed in the same fashion as its Vector a...
Definition: Matrix.cpp:2149
struct ComplexEigenPair NonSymmetricRightEigenSolve(Matrix *A)
Computes the eigenvectors/eigenvalues of an NxN real nonsymmetric matrix.
Definition: Matrix.cpp:2203
struct EigenPair SymmetricRightEigenSolve(Matrix *A)
Computes the eigenvectors/eigenvalues of an NxN real symmetric matrix.
Definition: Matrix.cpp:2156
Matrix outerProduct(const Vector &v, const Vector &w)
Computes the outer product of two Vectors, v and w.
Definition: Matrix.cpp:1973
void SerialSVD(Matrix *A, Matrix *U, Vector *S, Matrix *V)
Computes the SVD of a undistributed matrix in serial.
Definition: Matrix.cpp:2279
Matrix DiagonalMatrixFactory(const Vector &v)
Factory function to make a diagonal matrix with nonzero entries as in its Vector argument....
Definition: Matrix.cpp:2073
Matrix * ev_real
The real component of the eigenvectors in matrix form.
Definition: Matrix.h:1376
std::vector< std::complex< double > > eigs
The corresponding complex eigenvalues.
Definition: Matrix.h:1386
Matrix * ev_imaginary
The imaginary component of the eigenvectors in matrix form.
Definition: Matrix.h:1381
std::vector< double > eigs
The corresponding eigenvalues.
Definition: Matrix.h:1364
Matrix * ev
The eigenvectors in matrix form.
Definition: Matrix.h:1359
Matrix * U
The left singular vectors.
Definition: Matrix.h:1398
Vector * S
The singular values.
Definition: Matrix.h:1403
Matrix * V
The right singular vectors.
Definition: Matrix.h:1408