This class represents a matrix, stored in column major format, to allow direct use of BLAS/LAPACK routines. More...
#include <DenseMatrix.hpp>
Public Member Functions | |
DenseMatrix () | |
DenseMatrix (std::size_t m, std::size_t n) | |
DenseMatrix (std::size_t m, std::size_t n, const std::function< scalar_t(std::size_t, std::size_t)> &A) | |
DenseMatrix (std::size_t m, std::size_t n, const scalar_t *D, std::size_t ld) | |
DenseMatrix (std::size_t m, std::size_t n, const DenseMatrix< scalar_t > &D, std::size_t i, std::size_t j) | |
DenseMatrix (const DenseMatrix< scalar_t > &D) | |
DenseMatrix (DenseMatrix< scalar_t > &&D) | |
virtual | ~DenseMatrix () |
virtual DenseMatrix< scalar_t > & | operator= (const DenseMatrix< scalar_t > &D) |
virtual DenseMatrix< scalar_t > & | operator= (DenseMatrix< scalar_t > &&D) |
std::size_t | rows () const |
std::size_t | cols () const |
std::size_t | ld () const |
const scalar_t * | data () const |
scalar_t * | data () |
scalar_t * | end () |
const scalar_t * | end () const |
const scalar_t & | operator() (std::size_t i, std::size_t j) const |
const scalar_t * | ptr (std::size_t i, std::size_t j) const |
scalar_t & | operator() (std::size_t i, std::size_t j) |
scalar_t * | ptr (std::size_t i, std::size_t j) |
void | print () const |
void | print (std::string name, bool all=false, int width=8) const |
void | print_to_file (std::string name, std::string filename, int width=8) const |
void | random () |
void | random (random::RandomGeneratorBase< typename RealType< scalar_t >::value_type > &rgen) |
void | fill (scalar_t v) |
void | fill (const std::function< scalar_t(std::size_t, std::size_t)> &A) |
void | zero () |
void | eye () |
virtual void | clear () |
void | resize (std::size_t m, std::size_t n) |
void | hconcat (const DenseMatrix< scalar_t > &b) |
void | copy (const DenseMatrix< scalar_t > &B, std::size_t i=0, std::size_t j=0) |
void | copy (const scalar_t *B, std::size_t ldb) |
DenseMatrix< scalar_t > | transpose () const |
void | transpose (DenseMatrix< scalar_t > &X) const |
void | laswp (const std::vector< int > &P, bool fwd) |
void | laswp (const int *P, bool fwd) |
void | lapmr (const std::vector< int > &P, bool fwd) |
void | lapmt (const std::vector< int > &P, bool fwd) |
void | extract_rows (const std::vector< std::size_t > &I, DenseMatrix< scalar_t > &B) const |
DenseMatrix< scalar_t > | extract_rows (const std::vector< std::size_t > &I) const |
void | extract_cols (const std::vector< std::size_t > &I, DenseMatrix< scalar_t > &B) const |
DenseMatrix< scalar_t > | extract_cols (const std::vector< std::size_t > &I) const |
DenseMatrix< scalar_t > | extract (const std::vector< std::size_t > &I, const std::vector< std::size_t > &J) const |
DenseMatrix< scalar_t > & | scatter_rows_add (const std::vector< std::size_t > &I, const DenseMatrix< scalar_t > &B, int depth) |
DenseMatrix< scalar_t > & | add (const DenseMatrix< scalar_t > &B, int depth=0) |
DenseMatrix< scalar_t > & | sub (const DenseMatrix< scalar_t > &B, int depth=0) |
DenseMatrix< scalar_t > & | scale (scalar_t alpha, int depth=0) |
DenseMatrix< scalar_t > & | scaled_add (scalar_t alpha, const DenseMatrix< scalar_t > &B, int depth=0) |
DenseMatrix< scalar_t > & | scale_and_add (scalar_t alpha, const DenseMatrix< scalar_t > &B, int depth=0) |
DenseMatrix< scalar_t > & | scale_rows (const std::vector< scalar_t > &D, int depth=0) |
DenseMatrix< scalar_t > & | scale_rows_real (const std::vector< real_t > &D, int depth=0) |
DenseMatrix< scalar_t > & | scale_rows (const scalar_t *D, int depth=0) |
DenseMatrix< scalar_t > & | scale_rows_real (const real_t *D, int depth=0) |
DenseMatrix< scalar_t > & | div_rows (const std::vector< scalar_t > &D, int depth=0) |
real_t | norm () const |
real_t | normF () const |
real_t | norm1 () const |
real_t | normI () const |
virtual std::size_t | memory () const |
virtual std::size_t | nonzeros () const |
int | LU (std::vector< int > &piv, int depth=0) |
std::vector< int > | LU (int depth=0) |
int | Cholesky (int depth=0) |
std::vector< int > | LDLt (int depth=0) |
DenseMatrix< scalar_t > | solve (const DenseMatrix< scalar_t > &b, const std::vector< int > &piv, int depth=0) const |
void | solve_LU_in_place (DenseMatrix< scalar_t > &b, const std::vector< int > &piv, int depth=0) const |
void | solve_LU_in_place (DenseMatrix< scalar_t > &b, const int *piv, int depth=0) const |
void | solve_LDLt_in_place (DenseMatrix< scalar_t > &b, const std::vector< int > &piv, int depth=0) const |
void | LQ (DenseMatrix< scalar_t > &L, DenseMatrix< scalar_t > &Q, int depth) const |
void | orthogonalize (scalar_t &r_max, scalar_t &r_min, int depth) |
void | ID_column (DenseMatrix< scalar_t > &X, std::vector< int > &piv, std::vector< std::size_t > &ind, real_t rel_tol, real_t abs_tol, int max_rank, int depth) |
void | ID_row (DenseMatrix< scalar_t > &X, std::vector< int > &piv, std::vector< std::size_t > &ind, real_t rel_tol, real_t abs_tol, int max_rank, int depth) const |
void | low_rank (DenseMatrix< scalar_t > &U, DenseMatrix< scalar_t > &V, real_t rel_tol, real_t abs_tol, int max_rank, int depth) const |
std::vector< scalar_t > | singular_values () const |
void | shift (scalar_t sigma) |
int | syev (Jobz job, UpLo ul, std::vector< scalar_t > &lambda) |
void | write (const std::string &fname) const |
Static Public Member Functions | |
static DenseMatrix< scalar_t > | read (const std::string &fname) |
Friends | |
template<typename T > | |
class | DistributedMatrix |
template<typename T > | |
std::ofstream & | operator<< (std::ofstream &os, const DenseMatrix< T > &D) |
template<typename T > | |
std::ifstream & | operator>> (std::ifstream &is, DenseMatrix< T > &D) |
This class represents a matrix, stored in column major format, to allow direct use of BLAS/LAPACK routines.
This class represents a (2D) matrix, stored in column major format, to allow direct use of BLAS/LAPACK routines. A DenseMatrix allocates, owns and deallocates its memory. If you want to use pre-allocated memory to represent a dense matrix, use the DenseMatrixWrapper<scalar_t> class.
Several routines in this matrix perform some sort of bounds or size checking using assertions. These assertions can be removed by compiling with -DNDEBUG, which is added by default when using a CMake Release build.
Several routines in this class take a depth parameter. This refers to the depth of the nested OpenMP task spawning. No more tasks will be generated once the depth reaches a certain maximum level (params::task_recursion_cutoff_level), in order to limit the overhead of task spawning.
scalar_t | Possible values for the scalar_t template type are: float, double, std::complex<float> and std::complex<double>. |
Several BLAS-like interfaces are provided:
strumpack::DenseMatrix< scalar_t >::DenseMatrix | ( | ) |
Default constructor, constucts 0x0 empty matrix, with leading dimension 1.
strumpack::DenseMatrix< scalar_t >::DenseMatrix | ( | std::size_t | m, |
std::size_t | n | ||
) |
Constructs, and allocates, an m x n dense matrix, using column major storage. The leading dimension will be max(1, m).
m | Number of rows in the constructed matrix. |
n | Number of columns in the constructed matrix. |
strumpack::DenseMatrix< scalar_t >::DenseMatrix | ( | std::size_t | m, |
std::size_t | n, | ||
const std::function< scalar_t(std::size_t, std::size_t)> & | A | ||
) |
Constructs, and allocates, an m x n dense matrix, using column major storage. The leading dimension will be max(1, m).
m | Number of rows in the constructed matrix. |
n | Number of columns in the constructed matrix. |
A | routine to compute each element |
strumpack::DenseMatrix< scalar_t >::DenseMatrix | ( | std::size_t | m, |
std::size_t | n, | ||
const scalar_t * | D, | ||
std::size_t | ld | ||
) |
Construct/allocate a dense m x n matrix, and initialize it by copying the data pointed to by D (with leading dimension ld).
m | Number of rows in the constructed matrix. |
n | Number of columns in the constructed matrix. |
D | pointer to data to be copied in newly allocated DenseMatrix. Cannot be null. |
ld | Leading dimension of logically 2D matrix pointed to by D. Should be >= m. |
strumpack::DenseMatrix< scalar_t >::DenseMatrix | ( | std::size_t | m, |
std::size_t | n, | ||
const DenseMatrix< scalar_t > & | D, | ||
std::size_t | i, | ||
std::size_t | j | ||
) |
Construct a dense m x n matrix by copying a submatrix from matrix D. The copied submatrix has has top-left corner at position i,j in D, i.e. is D(i:i+m,j:j+n). The submatrix should be contained in D, i.e.: i+m <= D.rows() and j+n <= D.cols().
m | Number of rows in the constructed matrix. |
n | Number of columns in the constructed matrix. |
D | Matrix from which the newly constructed DenseMatrix is a submatrix. |
i | Row-offset in D, denoting top side of submatrix to copy to new matrix. |
j | Column-offset in D, denoting left side of submatrix to copy to new matrix. |
strumpack::DenseMatrix< scalar_t >::DenseMatrix | ( | const DenseMatrix< scalar_t > & | D | ) |
Copy constructor
strumpack::DenseMatrix< scalar_t >::DenseMatrix | ( | DenseMatrix< scalar_t > && | D | ) |
Move constructor
|
virtual |
Destructor
DenseMatrix<scalar_t>& strumpack::DenseMatrix< scalar_t >::add | ( | const DenseMatrix< scalar_t > & | B, |
int | depth = 0 |
||
) |
Add matrix B to this matrix. Return a reference to this matrix.
B | matrix to add to this matrix. |
depth | current OpenMP task recursion depth |
int strumpack::DenseMatrix< scalar_t >::Cholesky | ( | int | depth = 0 | ) |
|
virtual |
Clear the matrix. Resets the number of rows and columns to 0.
Reimplemented in strumpack::DenseMatrixWrapper< scalar_t >.
|
inline |
Number of columns of the matrix
void strumpack::DenseMatrix< scalar_t >::copy | ( | const DenseMatrix< scalar_t > & | B, |
std::size_t | i = 0 , |
||
std::size_t | j = 0 |
||
) |
Copy a submatrix of size rows() x cols() from B, at position (i,j) into this matrix. The following conditions should be satisfied: i+rows() <= B.rows() and j+cols() <= B.cols().
B | Matrix from which to copy |
i | Row-offset denoting the top of the submatrix of B to copy |
j | Column-offset denoting the top of the submatrix of B to copy |
void strumpack::DenseMatrix< scalar_t >::copy | ( | const scalar_t * | B, |
std::size_t | ldb | ||
) |
|
inline |
Pointer to the raw data used to represent this matrix.
|
inline |
Const pointer to the raw data used to represent this matrix.
DenseMatrix<scalar_t>& strumpack::DenseMatrix< scalar_t >::div_rows | ( | const std::vector< scalar_t > & | D, |
int | depth = 0 |
||
) |
Scale the rows of this matrix with the inverses of the scalar values in vector D, ie, this->operator()(i,j) /= D[i].
D | scalar factors, D.size() == rows() |
depth | current OpenMP task recursion depth |
|
inline |
|
inline |
DenseMatrix<scalar_t> strumpack::DenseMatrix< scalar_t >::extract | ( | const std::vector< std::size_t > & | I, |
const std::vector< std::size_t > & | J | ||
) | const |
Return a submatrix of this matrix defined by (I,J). The vectors I and J define the row and column indices of the submatrix. The extracted submatrix will be I.size() x J.size(). The extracted submatrix, lets call it B, satisfies B(i,j) = this->operator()(I[i],J[j]).
DenseMatrix<scalar_t> strumpack::DenseMatrix< scalar_t >::extract_cols | ( | const std::vector< std::size_t > & | I | ) | const |
void strumpack::DenseMatrix< scalar_t >::extract_cols | ( | const std::vector< std::size_t > & | I, |
DenseMatrix< scalar_t > & | B | ||
) | const |
Extract columns of this matrix to a specified matrix. Column I[i] in this matrix will become column i in matrix B.
DenseMatrix<scalar_t> strumpack::DenseMatrix< scalar_t >::extract_rows | ( | const std::vector< std::size_t > & | I | ) | const |
void strumpack::DenseMatrix< scalar_t >::extract_rows | ( | const std::vector< std::size_t > & | I, |
DenseMatrix< scalar_t > & | B | ||
) | const |
Extract rows of this matrix to a specified matrix. Row I[i] in this matrix will become row i in matrix B.
void strumpack::DenseMatrix< scalar_t >::eye | ( | ) |
Set the matrix to the identity matrix. Also works for rectangular matrices.
void strumpack::DenseMatrix< scalar_t >::fill | ( | const std::function< scalar_t(std::size_t, std::size_t)> & | A | ) |
Fill matrix using a routine to compute values
A | routine, can be lambda, functor or function pointer. Takes row i, column j and should return value at i, j. |
void strumpack::DenseMatrix< scalar_t >::fill | ( | scalar_t | v | ) |
Fill matrix with a constant value
v | value to set |
void strumpack::DenseMatrix< scalar_t >::hconcat | ( | const DenseMatrix< scalar_t > & | b | ) |
Horizontally concatenate a matrix to this matrix: [this b]. The resulting matrix will have the same number of rows as b and as this original matrix, and will have cols()+b.cols() columns.
b | Matrix to concatenate to this matrix. The b matrix should have to same number of rows, rows == b.rows(). |
void strumpack::DenseMatrix< scalar_t >::ID_column | ( | DenseMatrix< scalar_t > & | X, |
std::vector< int > & | piv, | ||
std::vector< std::size_t > & | ind, | ||
real_t | rel_tol, | ||
real_t | abs_tol, | ||
int | max_rank, | ||
int | depth | ||
) |
Compute an interpolative decomposition on the columns, ie, write this matrix as a linear combination of some of its columns. This is computed using QR with column pivoting (dgeqp3, modified to take tolerances), also refered to as a rank-revealing QR (RRQR), followed by a triangular solve.
TODO check this, also definition of ind and piv??!!
this ~= this(:,ind) * [eye(rank,rank) X] \ ~= (this(:,piv))(:,1:rank) * [eye(rank,rank) X]
X | output, see description above, will have rank rows ??? |
piv | pivot vector resulting from the column pivoted QR |
ind | set of columns selected from this matrix by RRQR |
rel_tol | relative tolerance used in the RRQR (column pivoted QR) |
abs_tol | absolute tolerance used in the RRQR (column pivoted QR) |
max_rank | maximum rank for RRQR |
depth | current OpenMP task recursion depth |
void strumpack::DenseMatrix< scalar_t >::ID_row | ( | DenseMatrix< scalar_t > & | X, |
std::vector< int > & | piv, | ||
std::vector< std::size_t > & | ind, | ||
real_t | rel_tol, | ||
real_t | abs_tol, | ||
int | max_rank, | ||
int | depth | ||
) | const |
Similar to ID_column, but transposed. This is implemented by calling ID_column on the transpose of this matrix, the resulting X is then again transposed.
X | output, see description in ID_column, will have rank columns ??? |
piv | pivot vector resulting from the column pivoted QR |
ind | set of columns selected from this matrix by RRQR |
rel_tol | relative tolerance used in the RRQR (column pivoted QR) |
abs_tol | absolute tolerance used in the RRQR (column pivoted QR) |
max_rank | maximum rank for RRQR |
depth | current OpenMP task recursion depth |
void strumpack::DenseMatrix< scalar_t >::lapmr | ( | const std::vector< int > & | P, |
bool | fwd | ||
) |
Apply the LAPACK routine xLAPMR to the matrix. xLAPMR rearranges the rows of the M by N matrix X as specified by the permutation K(1),K(2),...,K(M) of the integers 1,...,M. If fwd == true, forward permutation: X(K(I),*) is moved X(I,*) for I = 1,2,...,M. If fwd == false, backward permutation: X(I,*) is moved to X(K(I),*) for I = 1,2,...,M.
P | permutation vector, should be of size rows() |
fwd | apply permutation, or inverse permutation, see above |
void strumpack::DenseMatrix< scalar_t >::lapmt | ( | const std::vector< int > & | P, |
bool | fwd | ||
) |
Apply the LAPACK routines xLAPMT to the matrix. xLAPMT rearranges the columns of the M by N matrix X as specified by the permutation K(1),K(2),...,K(N) of the integers 1,...,N. If fwd == true, forward permutation: X(*,K(J)) is moved X(*,J) for J = 1,2,...,N. If fwd == false, backward permutation: X(*,J) is moved to X(*,K(J)) for J = 1,2,...,N.
P | permutation vector, should be of size cols() |
fwd | apply permutation, or inverse permutation, see above |
void strumpack::DenseMatrix< scalar_t >::laswp | ( | const int * | P, |
bool | fwd | ||
) |
Apply the LAPACK routine xLASWP to the matrix. xLASWP performs a series of row interchanges on the matrix. One row interchange is initiated for each of rows in the vector P.
P | vector with row interchanges, this is assumed to be of size rows() |
fwd | if fwd is false, the pivots are applied in reverse order |
void strumpack::DenseMatrix< scalar_t >::laswp | ( | const std::vector< int > & | P, |
bool | fwd | ||
) |
Apply the LAPACK routine xLASWP to the matrix. xLASWP performs a series of row interchanges on the matrix. One row interchange is initiated for each of rows in the vector P.
P | vector with row interchanges, this is assumed to be of size rows() |
fwd | if fwd is false, the pivots are applied in reverse order |
|
inline |
Leading dimension used to store the matrix, typically set to max(1, rows())
std::vector<int> strumpack::DenseMatrix< scalar_t >::LDLt | ( | int | depth = 0 | ) |
void strumpack::DenseMatrix< scalar_t >::low_rank | ( | DenseMatrix< scalar_t > & | U, |
DenseMatrix< scalar_t > & | V, | ||
real_t | rel_tol, | ||
real_t | abs_tol, | ||
int | max_rank, | ||
int | depth | ||
) | const |
Computes a low-rank factorization of this matrix, with specified relative and absolute tolerances, (used in column pivoted (rank-revealing) QR).
this ~= U * V
U | matrix U, low-rank factor. Will be of size this->rows() x rank. Does not need to be allocated beforehand. |
V | matrix V, low-rank factor. Will be of size rank x this->cols(). Does not need to be allcoated beforehand. |
rel_tol | relative tolerance used in the RRQR (column pivoted QR) |
abs_tol | absolute tolerance used in the RRQR (column pivoted QR) |
max_rank | maximum rank for RRQR |
depth | current OpenMP task recursion depth |
void strumpack::DenseMatrix< scalar_t >::LQ | ( | DenseMatrix< scalar_t > & | L, |
DenseMatrix< scalar_t > & | Q, | ||
int | depth | ||
) | const |
Solve a linear system Ax=b with this matrix, factored in its LDLt factors (in place), using LDLt_rook. There can be multiple right hand side vectors. The solution is returned by value.
Disabled for now because not supported on some systems with older LAPACK versions.
b | input, right hand side vector/matrix. On output this will be the solution. |
piv | pivot vector returned by LU factorization |
depth | current OpenMP task recursion depth |
L | lower triangular matrix, output argument. Does not have to be allocated. |
Q | unitary matrix, not necessarily square. Does not have to be allocated. |
depth | current OpenMP task recursion depth |
std::vector<int> strumpack::DenseMatrix< scalar_t >::LU | ( | int | depth = 0 | ) |
Compute an LU factorization of this matrix using partial pivoting with row interchanges. The factorization has the form A = P * L * U, where P is a permutation matrix, L is lower triangular with unit diagonal elements, and U is upper triangular. This calls the LAPACK routine DGETRF. The L and U factors are stored in place, the permutation is returned, and can be applied with the laswp() routine.
depth | current OpenMP task recursion depth |
int strumpack::DenseMatrix< scalar_t >::LU | ( | std::vector< int > & | piv, |
int | depth = 0 |
||
) |
Compute an LU factorization of this matrix using partial pivoting with row interchanges. The factorization has the form A = P * L * U, where P is a permutation matrix, L is lower triangular with unit diagonal elements, and U is upper triangular. This calls the LAPACK routine DGETRF. The L and U factors are stored in place, the permutation is returned, and can be applied with the laswp() routine.
piv | pivot vector, will be resized if necessary |
depth | current OpenMP task recursion depth |
|
inlinevirtual |
Return the (approximate) amount of memory taken by this matrix, in bytes. Simply nonzeros()*sizeof(scalar_t). The matrix metadata is not counted in this.
Reimplemented in strumpack::DenseMatrixWrapper< scalar_t >.
|
inlinevirtual |
Return the number of nonzeros in this matrix, ie, simply rows()*cols().
Reimplemented in strumpack::DenseMatrixWrapper< scalar_t >.
real_t strumpack::DenseMatrix< scalar_t >::norm | ( | ) | const |
Return default norm of this matrix. Currently the default is set to the Frobenius norm.
real_t strumpack::DenseMatrix< scalar_t >::norm1 | ( | ) | const |
Return the 1-norm of this matrix.
real_t strumpack::DenseMatrix< scalar_t >::normF | ( | ) | const |
Return the Frobenius norm of this matrix.
real_t strumpack::DenseMatrix< scalar_t >::normI | ( | ) | const |
Return the infinity norm of this matrix.
|
inline |
|
inline |
|
virtual |
Copy operator, expensive operation for large matrices.
Reimplemented in strumpack::DenseMatrixWrapper< scalar_t >.
|
virtual |
Move operator
D | Matrix to be moved into this object, will be emptied. |
void strumpack::DenseMatrix< scalar_t >::orthogonalize | ( | scalar_t & | r_max, |
scalar_t & | r_min, | ||
int | depth | ||
) |
Builds an orthonormal basis for the columns in this matrix, using QR factorization. It return the maximum and minimum elements on the diagonal of the upper triangular matrix in the QR factorization. These values can be used to check whether the matrix was rank deficient.
r_max | maximum value encountered on the diagonal of R (as returned from QR factorization) |
r_min | minimum value encountered on the diagonal of R (as returned from QR factorization) |
depth | current OpenMP task recursion depth |
|
inline |
void strumpack::DenseMatrix< scalar_t >::print | ( | std::string | name, |
bool | all = false , |
||
int | width = 8 |
||
) | const |
Print the matrix to std::cout, in a format interpretable by Matlab/Octave. The matrix is printed in full when not too big, i.e., when rows() <= 20 and cols() <= 32, or when all is set to true. Otherwise only its sizes and its norm are printed. Useful for debugging.
name | Name to use when printing. |
all | If true, print all values, if false, only print all values when not too big. Defaults to false. |
param | width Specifies how many digits to use for printing floating point values, defaults to 8. |
void strumpack::DenseMatrix< scalar_t >::print_to_file | ( | std::string | name, |
std::string | filename, | ||
int | width = 8 |
||
) | const |
Print the matrix to a file, in a format readable by Matlab/Octave.
name | Name to use for the printed matrix. |
filename | Name of the file to write to |
width | Number of digits to use to represent floating point values, defaults to 8. |
|
inline |
|
inline |
void strumpack::DenseMatrix< scalar_t >::random | ( | ) |
Fill the matrix with random numbers, using random number generator/distribution random::make_default_random_generator<real_t>()
void strumpack::DenseMatrix< scalar_t >::random | ( | random::RandomGeneratorBase< typename RealType< scalar_t >::value_type > & | rgen | ) |
Fill the matrix with random numbers, using the specified random number generator.
|
static |
Read a DenseMatrix<scalar_t> from a binary file, called fname.
void strumpack::DenseMatrix< scalar_t >::resize | ( | std::size_t | m, |
std::size_t | n | ||
) |
Resize the matrix. The relevant parts of the original matrix will be copied to the new matrix. The contents of new parts of the matrix are undefined.
m | Number of rows after resizing. |
m | Number of columns after resizing. |
|
inline |
Number of rows of the matrix
DenseMatrix<scalar_t>& strumpack::DenseMatrix< scalar_t >::scale | ( | scalar_t | alpha, |
int | depth = 0 |
||
) |
Scale this matrix by a constant factor.
alpha | scaling factor |
DenseMatrix<scalar_t>& strumpack::DenseMatrix< scalar_t >::scale_and_add | ( | scalar_t | alpha, |
const DenseMatrix< scalar_t > & | B, | ||
int | depth = 0 |
||
) |
Scale this matrix, and add a given matrix to this matrix, ie, this = alpha * this + B.
alpha | scalar factor |
B | matrix to add, should be the same size of this matrix |
depth | current OpenMP task recursion depth |
DenseMatrix<scalar_t>& strumpack::DenseMatrix< scalar_t >::scale_rows | ( | const scalar_t * | D, |
int | depth = 0 |
||
) |
Scale the rows of this matrix with the scalar values from the vector D. Row i in this matrix is scaled with D[i].
D | scaling vector |
depth | current OpenMP task recursion depth |
DenseMatrix<scalar_t>& strumpack::DenseMatrix< scalar_t >::scale_rows | ( | const std::vector< scalar_t > & | D, |
int | depth = 0 |
||
) |
Scale the rows of this matrix with the scalar values from the vector D. Row i in this matrix is scaled with D[i].
D | scaling vector, D.size() == rows() |
depth | current OpenMP task recursion depth |
DenseMatrix<scalar_t>& strumpack::DenseMatrix< scalar_t >::scaled_add | ( | scalar_t | alpha, |
const DenseMatrix< scalar_t > & | B, | ||
int | depth = 0 |
||
) |
Add a scalar multiple of a given matrix to this matrix, ie, this += alpha * B.
alpha | scalar factor |
B | matrix to add, should be the same size of this matrix |
depth | current OpenMP task recursion depth |
DenseMatrix<scalar_t>& strumpack::DenseMatrix< scalar_t >::scatter_rows_add | ( | const std::vector< std::size_t > & | I, |
const DenseMatrix< scalar_t > & | B, | ||
int | depth | ||
) |
Add the rows of matrix B into this matrix at the rows specified by vector I, ie, add row i of matrix B to row I[i] of this matrix. This is used in the sparse solver.
I | index set in this matrix, where to add the rows of B, I[i] < rows() |
B | matrix with rows to scatter in this matrix |
depth | current OpenMP task recursion depth |
void strumpack::DenseMatrix< scalar_t >::shift | ( | scalar_t | sigma | ) |
Shift the diagonal with a value sigma, ie. add a scaled identity matrix.
sigma | scalar value to add to diagonal |
std::vector<scalar_t> strumpack::DenseMatrix< scalar_t >::singular_values | ( | ) | const |
Return a vector with the singular values of this matrix. Used only for debugging puposes.
DenseMatrix<scalar_t> strumpack::DenseMatrix< scalar_t >::solve | ( | const DenseMatrix< scalar_t > & | b, |
const std::vector< int > & | piv, | ||
int | depth = 0 |
||
) | const |
Compute an LDLt factorization of this matrix in-place. This calls the LAPACK routine sytrf_rook. Only the lower triangle is referenced/stored.
Disabled for now because not supported on some systems with older LAPACK versions.
depth | current OpenMP task recursion depth |
b | input, right hand side vector/matrix |
piv | pivot vector returned by LU factorization |
depth | current OpenMP task recursion depth |
void strumpack::DenseMatrix< scalar_t >::solve_LDLt_in_place | ( | DenseMatrix< scalar_t > & | b, |
const std::vector< int > & | piv, | ||
int | depth = 0 |
||
) | const |
Solve a linear system Ax=b with this matrix, factored in its LDLt factors (in place). There can be multiple right hand side vectors. The solution is returned by value.
b | input, right hand side vector/matrix. On output this will be the solution. |
piv | pivot vector returned by LU factorization |
depth | current OpenMP task recursion depth |
void strumpack::DenseMatrix< scalar_t >::solve_LU_in_place | ( | DenseMatrix< scalar_t > & | b, |
const int * | piv, | ||
int | depth = 0 |
||
) | const |
Solve a linear system Ax=b with this matrix, factored in its LU factors (in place), using a call to this->LU. There can be multiple right hand side vectors.
b | input, right hand side vector/matrix. On output this will be the solution. |
piv | pivot vector returned by LU factorization |
depth | current OpenMP task recursion depth |
void strumpack::DenseMatrix< scalar_t >::solve_LU_in_place | ( | DenseMatrix< scalar_t > & | b, |
const std::vector< int > & | piv, | ||
int | depth = 0 |
||
) | const |
Solve a linear system Ax=b with this matrix, factored in its LU factors (in place), using a call to this->LU. There can be multiple right hand side vectors.
b | input, right hand side vector/matrix. On output this will be the solution. |
piv | pivot vector returned by LU factorization |
depth | current OpenMP task recursion depth |
DenseMatrix<scalar_t>& strumpack::DenseMatrix< scalar_t >::sub | ( | const DenseMatrix< scalar_t > & | B, |
int | depth = 0 |
||
) |
Subtract matrix B from this matrix. Return a reference to this matrix.
B | matrix to subtract from this matrix |
depth | current OpenMP task recursion depth |
int strumpack::DenseMatrix< scalar_t >::syev | ( | Jobz | job, |
UpLo | ul, | ||
std::vector< scalar_t > & | lambda | ||
) |
SYEV computes all eigenvalues and, optionally, eigenvectors of this matrix. If job is N, the matrix is destroyed. If job is V, on exit the matrix will contain all eigenvectors.
job | N: Compute eigenvalues only, V: Compute eigenvalues and eigenvectors. |
ul | U: Upper triangle is stored, L: Lower triangle is stored. |
lambda | on exit this will contain all eigenvalues of this matrix. |
DenseMatrix<scalar_t> strumpack::DenseMatrix< scalar_t >::transpose | ( | ) | const |
Return the transpose of this matrix
void strumpack::DenseMatrix< scalar_t >::transpose | ( | DenseMatrix< scalar_t > & | X | ) | const |
Set X to the transpose of this matrix.
void strumpack::DenseMatrix< scalar_t >::write | ( | const std::string & | fname | ) | const |
Write this DenseMatrix<scalar_t> to a binary file, called fname.
void strumpack::DenseMatrix< scalar_t >::zero | ( | ) |
Set all matrix elements to zero