15 #ifndef included_Matrix_h
16 #define included_Matrix_h
57 bool randomized =
false);
82 bool copy_data =
true);
154 if (num_rows > INT_MAX / num_cols)
155 CAROM_ERROR(
"Matrix::setSize- new size exceeds maximum integer value!\n");
157 int new_size = num_rows*num_cols;
158 if (new_size > d_alloc_size) {
160 CAROM_ERROR(
"Can not reallocate externally owned storage.");
167 d_mat =
new double [new_size] {0.0};
168 d_alloc_size = new_size;
170 d_num_rows = num_rows;
171 d_num_cols = num_cols;
173 calculateNumDistributedRows();
185 return d_distributed;
213 if (!d_distributed) {
216 return d_num_distributed_rows;
290 const Matrix& other)
const
315 const Matrix* other)
const
317 CAROM_VERIFY(other != 0);
381 const Vector& other)
const
406 const Vector* other)
const
408 CAROM_VERIFY(other != 0);
514 const Matrix& other)
const
536 const Matrix* other)
const
538 CAROM_VERIFY(other != 0);
649 const Matrix& other)
const
674 const Matrix* other)
const
676 CAROM_VERIFY(other != 0);
738 const Vector& other)
const
763 const Vector* other)
const
765 CAROM_VERIFY(other != 0);
926 std::unique_ptr<Matrix>
934 void qr_factorize(std::vector<std::unique_ptr<Matrix>> & QR)
const;
953 int* row_pivot_owner,
954 int pivots_requested)
const;
971 orthogonalize(
bool double_pass =
false,
double zero_tol = 1.0e-15);
994 double zero_tol = 1.0e-15);
1024 CAROM_ASSERT((0 <= row) && (row <
numRows()));
1025 CAROM_ASSERT((0 <= col) && (col <
numColumns()));
1026 return d_mat[row*d_num_cols+col];
1047 CAROM_ASSERT((0 <= row) && (row <
numRows()));
1048 CAROM_ASSERT((0 <= col) && (col <
numColumns()));
1049 return d_mat[row*d_num_cols+col];
1064 return item(row, col);
1081 return item(row, col);
1090 void print(
const char * prefix)
const;
1098 void write(
const std::string& base_file_name)
const;
1106 void read(
const std::string& base_file_name);
1115 void local_read(
const std::string& base_file_name,
int rank);
1153 calculateNumDistributedRows();
1170 qrcp_pivots_transpose_serial(
int* row_pivot,
1171 int* row_pivot_owner,
1172 int pivots_requested)
const;
1189 qrcp_pivots_transpose_distributed(
int* row_pivot,
1190 int* row_pivot_owner,
1191 int pivots_requested)
const;
1194 qrcp_pivots_transpose_distributed_scalapack(
int* row_pivot,
1195 int* row_pivot_owner,
1196 int pivots_requested)
const;
1214 qrcp_pivots_transpose_distributed_elemental(
int* row_pivot,
1215 int* row_pivot_owner,
1216 int pivots_requested)
const;
1234 qrcp_pivots_transpose_distributed_elemental_balanced
1235 (
int* row_pivot,
int* row_pivot_owner,
int pivots_requested)
const;
1253 qrcp_pivots_transpose_distributed_elemental_unbalanced
1254 (
int* row_pivot,
int* row_pivot_owner,
int pivots_requested)
const;
1267 void qr_factorize(
bool computeR, std::vector<Matrix*> & QR)
const;
1282 int d_num_distributed_rows;
1412 std::vector<std::complex<double>>
eigs;
1486 const std::vector<double> *tscale=NULL,
1487 const bool At0at0=
false,
const bool Bt0at0=
false,
1488 const bool lagB=
false,
const bool skip0=
false);
void local_read(const std::string &base_file_name, int rank)
read a single rank of a distributed Matrix into (a) HDF file(s).
void setSize(int num_rows, int num_cols)
Sets the number of rows and columns of the matrix and reallocates storage if needed....
Matrix * getFirstNColumns(int n) const
Get the first N columns of a matrix.
double * getData() const
Get the matrix data as a pointer.
bool balanced() const
Returns true if rows of matrix are load-balanced.
Matrix * inverse() const
Computes and returns the inverse of this.
const double & operator()(int row, int col) const
Const Matrix member access.
void transposePseudoinverse()
Computes the transposePseudoinverse of this.
bool distributed() const
Returns true if the Matrix is distributed.
void rescale_rows_max()
Rescale every matrix row by its maximum absolute value.
Matrix & operator=(const Matrix &rhs)
Assignment operator.
void write(const std::string &base_file_name) const
write Matrix into (a) HDF file(s).
int numColumns() const
Returns the number of columns in the Matrix. This method will return the same value from each process...
void distribute(const int &local_num_rows)
Distribute this matrix rows among MPI processes, based on the specified local number of rows....
void gather()
Gather all the distributed rows among MPI processes. This becomes not distributed after this function...
Matrix * elementwise_mult(const Matrix &other) const
Multiplies two matrices element-wise.
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.
std::unique_ptr< Matrix > qr_factorize() const
Computes and returns the Q of the QR factorization of this.
Matrix & operator+=(const Matrix &rhs)
Addition operator.
Matrix * mult(const Matrix *other) const
Multiplies this Matrix with other and returns the product, pointer version.
int numDistributedRows() const
Returns the number of rows of the Matrix across all processors.
Vector * getColumn(int column) const
Returns a column of the matrix (not owned by Matrix).
Matrix & operator-=(const Matrix &rhs)
Subtraction operator.
void orthogonalize(bool double_pass=false, double zero_tol=1.0e-15)
Orthonormalizes the matrix.
void multPlus(Vector &a, const Vector &b, double c) const
Computes a += this*b*c.
Vector * transposeMult(const Vector &other) const
Multiplies the transpose of this Matrix with other and returns the product, reference version.
double & item(int row, int col)
Non-const Matrix member access. Matrix data is stored in row-major format.
void print(const char *prefix) const
print Matrix into (a) ascii file(s).
Matrix * transposeMult(const Matrix &other) const
Multiplies the transpose of this Matrix with other and returns the product, reference version.
void pointwise_mult(int this_row, const Vector &other, Vector &result) const
Multiplies a specified row of this Matrix with other pointwise.
Matrix * mult(const Matrix &other) const
Multiplies this Matrix with other and returns the product, reference version.
Vector * mult(const Vector *other) const
Multiplies this Matrix with other and returns the product, pointer version.
Matrix * elementwise_mult(const Matrix *other) const
Multiplies two matrices element-wise and returns the product, pointer version.
const double & item(int row, int col) const
Const Matrix member access. Matrix data is stored in row-major format.
Vector * mult(const Vector &other) const
Multiplies this Matrix with other and returns the product, reference version.
Matrix * transposeMult(const Matrix *other) const
Multiplies the transpose of this Matrix with other and returns the product, pointer version.
int numRows() const
Returns the number of rows of the Matrix on this processor.
void read(const std::string &base_file_name)
read Matrix into (a) HDF file(s).
Matrix * elementwise_square() const
Square every element in the matrix.
Vector * transposeMult(const Vector *other) const
Multiplies the transpose of this Matrix with other and returns the product, pointer version.
void transpose()
Replaces this Matrix with its transpose (in place), in the serial square case only.
void rescale_cols_max()
Rescale every matrix column by its maximum absolute value.
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...
Matrix IdentityMatrixFactory(const Vector &v)
Factory function to make an identity matrix with rows distributed in the same fashion as its Vector a...
struct ComplexEigenPair NonSymmetricRightEigenSolve(Matrix *A)
Computes the eigenvectors/eigenvalues of an NxN real nonsymmetric matrix.
struct EigenPair SymmetricRightEigenSolve(Matrix *A)
Computes the eigenvectors/eigenvalues of an NxN real symmetric matrix.
Matrix outerProduct(const Vector &v, const Vector &w)
Computes the outer product of two Vectors, v and w.
void SerialSVD(Matrix *A, Matrix *U, Vector *S, Matrix *V)
Computes the SVD of a undistributed matrix in serial.
Matrix DiagonalMatrixFactory(const Vector &v)
Factory function to make a diagonal matrix with nonzero entries as in its Vector argument....
Matrix * ev_real
The real component of the eigenvectors in matrix form.
std::vector< std::complex< double > > eigs
The corresponding complex eigenvalues.
Matrix * ev_imaginary
The imaginary component of the eigenvectors in matrix form.
std::vector< double > eigs
The corresponding eigenvalues.
Matrix * ev
The eigenvectors in matrix form.
Matrix * U
The left singular vectors.
Vector * S
The singular values.
Matrix * V
The right singular vectors.