strumpack::DenseMatrix< scalar_t > Class Template Reference

This class represents a matrix, stored in column major format, to allow direct use of BLAS/LAPACK routines. More...

#include <DenseMatrix.hpp>

Inheritance diagram for strumpack::DenseMatrix< scalar_t >:

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)
 

Detailed Description

template<typename scalar_t>
class strumpack::DenseMatrix< scalar_t >

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.

Template Parameters
scalar_tPossible values for the scalar_t template type are: float, double, std::complex<float> and std::complex<double>.

Several BLAS-like interfaces are provided:

See also
strumpack::gemm
strumpack::trsm
strumpack::trmm
strumpack::gemv

Constructor & Destructor Documentation

◆ DenseMatrix() [1/7]

template<typename scalar_t >
strumpack::DenseMatrix< scalar_t >::DenseMatrix ( )

Default constructor, constucts 0x0 empty matrix, with leading dimension 1.

◆ DenseMatrix() [2/7]

template<typename scalar_t >
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).

Parameters
mNumber of rows in the constructed matrix.
nNumber of columns in the constructed matrix.

◆ DenseMatrix() [3/7]

template<typename scalar_t >
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).

Parameters
mNumber of rows in the constructed matrix.
nNumber of columns in the constructed matrix.
Aroutine to compute each element

◆ DenseMatrix() [4/7]

template<typename scalar_t >
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).

Parameters
mNumber of rows in the constructed matrix.
nNumber of columns in the constructed matrix.
Dpointer to data to be copied in newly allocated DenseMatrix. Cannot be null.
ldLeading dimension of logically 2D matrix pointed to by D. Should be >= m.

◆ DenseMatrix() [5/7]

template<typename scalar_t >
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().

Parameters
mNumber of rows in the constructed matrix.
nNumber of columns in the constructed matrix.
DMatrix from which the newly constructed DenseMatrix is a submatrix.
iRow-offset in D, denoting top side of submatrix to copy to new matrix.
jColumn-offset in D, denoting left side of submatrix to copy to new matrix.

◆ DenseMatrix() [6/7]

template<typename scalar_t >
strumpack::DenseMatrix< scalar_t >::DenseMatrix ( const DenseMatrix< scalar_t > &  D)

Copy constructor

◆ DenseMatrix() [7/7]

template<typename scalar_t >
strumpack::DenseMatrix< scalar_t >::DenseMatrix ( DenseMatrix< scalar_t > &&  D)

Move constructor

◆ ~DenseMatrix()

template<typename scalar_t >
virtual strumpack::DenseMatrix< scalar_t >::~DenseMatrix ( )
virtual

Destructor

Member Function Documentation

◆ add()

template<typename scalar_t >
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.

Parameters
Bmatrix to add to this matrix.
depthcurrent OpenMP task recursion depth

◆ Cholesky()

template<typename scalar_t >
int strumpack::DenseMatrix< scalar_t >::Cholesky ( int  depth = 0)

Compute a Cholesky factorization of this matrix in-place. This calls the LAPACK routine DPOTRF. Only the lower triangle is written. Only the lower triangle is referenced/stored.

Parameters
depthcurrent OpenMP task recursion depth
Returns
info from xpotrf
See also
LU, LDLt

◆ clear()

template<typename scalar_t >
virtual void strumpack::DenseMatrix< scalar_t >::clear ( )
virtual

Clear the matrix. Resets the number of rows and columns to 0.

Reimplemented in strumpack::DenseMatrixWrapper< scalar_t >.

◆ cols()

template<typename scalar_t >
std::size_t strumpack::DenseMatrix< scalar_t >::cols ( ) const
inline

Number of columns of the matrix

◆ copy() [1/2]

template<typename scalar_t >
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().

Parameters
BMatrix from which to copy
iRow-offset denoting the top of the submatrix of B to copy
jColumn-offset denoting the top of the submatrix of B to copy

◆ copy() [2/2]

template<typename scalar_t >
void strumpack::DenseMatrix< scalar_t >::copy ( const scalar_t *  B,
std::size_t  ldb 
)

Copy a submatrix of size rows() x cols() from matrix B, with leading dimension ld, into this matrix.

Parameters
BPointer to the matrix data to copy
ldLeading dimension of matrix pointed to by B, should be at least cols()

◆ data() [1/2]

template<typename scalar_t >
scalar_t* strumpack::DenseMatrix< scalar_t >::data ( )
inline

Pointer to the raw data used to represent this matrix.

◆ data() [2/2]

template<typename scalar_t >
const scalar_t* strumpack::DenseMatrix< scalar_t >::data ( ) const
inline

Const pointer to the raw data used to represent this matrix.

◆ div_rows()

template<typename scalar_t >
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].

Parameters
Dscalar factors, D.size() == rows()
depthcurrent OpenMP task recursion depth

◆ end() [1/2]

template<typename scalar_t >
scalar_t* strumpack::DenseMatrix< scalar_t >::end ( )
inline

Pointer to the element after the last one in this matrix. end() == data() + size()

◆ end() [2/2]

template<typename scalar_t >
const scalar_t* strumpack::DenseMatrix< scalar_t >::end ( ) const
inline

Pointer to the element after the last one in this matrix. end() == data() + size()

◆ extract()

template<typename scalar_t >
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]).

Parameters
Irow indices of elements to extract, I[i] < rows()
Jcolumn indices of elements to extract, J[j] < cols()

◆ extract_cols() [1/2]

template<typename scalar_t >
DenseMatrix<scalar_t> strumpack::DenseMatrix< scalar_t >::extract_cols ( const std::vector< std::size_t > &  I) const

Return a matrix with columns I of this matrix.

Parameters
Iset of indices of columns to extract from this matrix. The elements of I should not be sorted, but they should be I[i] < cols(). \params B matrix to extraxt to. B should have the correct size, ie. B.cols() == I.size() and B.rows() == rows()

◆ extract_cols() [2/2]

template<typename scalar_t >
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.

Parameters
Iset of indices of columns to extract from this matrix. The elements of I should not be sorted, but they should be I[i] < cols(). \params B matrix to extraxt to. B should have the correct size, ie. B.cols() == I.size() and B.rows() == rows()

◆ extract_rows() [1/2]

template<typename scalar_t >
DenseMatrix<scalar_t> strumpack::DenseMatrix< scalar_t >::extract_rows ( const std::vector< std::size_t > &  I) const

Return a matrix with rows I of this matrix.

Parameters
Iset of indices of rows to extract from this matrix. The elements of I should not be sorted, but they should be I[i] < rows(). \params B matrix to extraxt to. B should have the correct size, ie. B.cols() == cols() and B.rows() == I.size()

◆ extract_rows() [2/2]

template<typename scalar_t >
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.

Parameters
Iset of indices of rows to extract from this matrix. The elements of I should not be sorted, but they should be I[i] < rows(). \params B matrix to extraxt to. B should have the correct size, ie. B.cols() == cols() and B.rows() == I.size()

◆ eye()

template<typename scalar_t >
void strumpack::DenseMatrix< scalar_t >::eye ( )

Set the matrix to the identity matrix. Also works for rectangular matrices.

◆ fill() [1/2]

template<typename scalar_t >
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

Parameters
Aroutine, can be lambda, functor or function pointer. Takes row i, column j and should return value at i, j.

◆ fill() [2/2]

template<typename scalar_t >
void strumpack::DenseMatrix< scalar_t >::fill ( scalar_t  v)

Fill matrix with a constant value

Parameters
vvalue to set

◆ hconcat()

template<typename scalar_t >
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.

Parameters
bMatrix to concatenate to this matrix. The b matrix should have to same number of rows, rows == b.rows().

◆ ID_column()

template<typename scalar_t >
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]

Parameters
Xoutput, see description above, will have rank rows ???
pivpivot vector resulting from the column pivoted QR
indset of columns selected from this matrix by RRQR
rel_tolrelative tolerance used in the RRQR (column pivoted QR)
abs_tolabsolute tolerance used in the RRQR (column pivoted QR)
max_rankmaximum rank for RRQR
depthcurrent OpenMP task recursion depth
See also
ID_row

◆ ID_row()

template<typename scalar_t >
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.

Parameters
Xoutput, see description in ID_column, will have rank columns ???
pivpivot vector resulting from the column pivoted QR
indset of columns selected from this matrix by RRQR
rel_tolrelative tolerance used in the RRQR (column pivoted QR)
abs_tolabsolute tolerance used in the RRQR (column pivoted QR)
max_rankmaximum rank for RRQR
depthcurrent OpenMP task recursion depth
See also
ID_column

◆ lapmr()

template<typename scalar_t >
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.

Parameters
Ppermutation vector, should be of size rows()
fwdapply permutation, or inverse permutation, see above

◆ lapmt()

template<typename scalar_t >
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.

Parameters
Ppermutation vector, should be of size cols()
fwdapply permutation, or inverse permutation, see above

◆ laswp() [1/2]

template<typename scalar_t >
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.

Parameters
Pvector with row interchanges, this is assumed to be of size rows()
fwdif fwd is false, the pivots are applied in reverse order

◆ laswp() [2/2]

template<typename scalar_t >
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.

Parameters
Pvector with row interchanges, this is assumed to be of size rows()
fwdif fwd is false, the pivots are applied in reverse order

◆ ld()

template<typename scalar_t >
std::size_t strumpack::DenseMatrix< scalar_t >::ld ( ) const
inline

Leading dimension used to store the matrix, typically set to max(1, rows())

◆ LDLt()

template<typename scalar_t >
std::vector<int> strumpack::DenseMatrix< scalar_t >::LDLt ( int  depth = 0)

Compute an LDLt factorization of this matrix in-place. This calls the LAPACK routine sytrf. Only the lower triangle is referenced/stored.

Parameters
depthcurrent OpenMP task recursion depth
Returns
info from xsytrf
See also
LU, Cholesky, solve_LDLt

◆ low_rank()

template<typename scalar_t >
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

Parameters
Umatrix U, low-rank factor. Will be of size this->rows() x rank. Does not need to be allocated beforehand.
Vmatrix V, low-rank factor. Will be of size rank x this->cols(). Does not need to be allcoated beforehand.
rel_tolrelative tolerance used in the RRQR (column pivoted QR)
abs_tolabsolute tolerance used in the RRQR (column pivoted QR)
max_rankmaximum rank for RRQR
depthcurrent OpenMP task recursion depth

◆ LQ()

template<typename scalar_t >
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.

Parameters
binput, right hand side vector/matrix. On output this will be the solution.
pivpivot vector returned by LU factorization
depthcurrent OpenMP task recursion depth
See also
LDLt_rook, LDLt, solve_LDLt_in_place, LU, solve_LU_in_place Compute an LQ (lower triangular, unitary) factorization of this matrix. This matrix is not modified, the L and Q factors are returned by reference in the L and Q arguments. L and Q do not have to be allocated, they will be resized accordingly.
Parameters
Llower triangular matrix, output argument. Does not have to be allocated.
Qunitary matrix, not necessarily square. Does not have to be allocated.
depthcurrent OpenMP task recursion depth

◆ LU() [1/2]

template<typename scalar_t >
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.

Parameters
depthcurrent OpenMP task recursion depth
Returns
the pivot vector
See also
laswp, solve

◆ LU() [2/2]

template<typename scalar_t >
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.

Parameters
pivpivot vector, will be resized if necessary
depthcurrent OpenMP task recursion depth
Returns
if nonzero, the pivot in this column was exactly zero
See also
laswp, solve

◆ memory()

template<typename scalar_t >
virtual std::size_t strumpack::DenseMatrix< scalar_t >::memory ( ) const
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 >.

◆ nonzeros()

template<typename scalar_t >
virtual std::size_t strumpack::DenseMatrix< scalar_t >::nonzeros ( ) const
inlinevirtual

Return the number of nonzeros in this matrix, ie, simply rows()*cols().

Reimplemented in strumpack::DenseMatrixWrapper< scalar_t >.

◆ norm()

template<typename 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.

◆ norm1()

template<typename scalar_t >
real_t strumpack::DenseMatrix< scalar_t >::norm1 ( ) const

Return the 1-norm of this matrix.

◆ normF()

template<typename scalar_t >
real_t strumpack::DenseMatrix< scalar_t >::normF ( ) const

Return the Frobenius norm of this matrix.

◆ normI()

template<typename scalar_t >
real_t strumpack::DenseMatrix< scalar_t >::normI ( ) const

Return the infinity norm of this matrix.

◆ operator()() [1/2]

template<typename scalar_t >
scalar_t& strumpack::DenseMatrix< scalar_t >::operator() ( std::size_t  i,
std::size_t  j 
)
inline

Reference to element (i,j) in the matrix. This will do a bounds check with assertions, which are enabled in Debug mode, disabled in Release mode.

Parameters
iRow index, i < rows()
jColumn index, j < cols()

◆ operator()() [2/2]

template<typename scalar_t >
const scalar_t& strumpack::DenseMatrix< scalar_t >::operator() ( std::size_t  i,
std::size_t  j 
) const
inline

Const reference to element (i,j) in the matrix. This will do a bounds check with assertions, which are enabled in Debug mode, disabled in Release mode.

Parameters
iRow index, i < rows()
jColumn index, j < cols()

◆ operator=() [1/2]

template<typename scalar_t >
virtual DenseMatrix<scalar_t>& strumpack::DenseMatrix< scalar_t >::operator= ( const DenseMatrix< scalar_t > &  D)
virtual

Copy operator, expensive operation for large matrices.

Reimplemented in strumpack::DenseMatrixWrapper< scalar_t >.

◆ operator=() [2/2]

template<typename scalar_t >
virtual DenseMatrix<scalar_t>& strumpack::DenseMatrix< scalar_t >::operator= ( DenseMatrix< scalar_t > &&  D)
virtual

Move operator

Parameters
DMatrix to be moved into this object, will be emptied.

◆ orthogonalize()

template<typename scalar_t >
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.

Parameters
r_maxmaximum value encountered on the diagonal of R (as returned from QR factorization)
r_minminimum value encountered on the diagonal of R (as returned from QR factorization)
depthcurrent OpenMP task recursion depth

◆ print() [1/2]

template<typename scalar_t >
void strumpack::DenseMatrix< scalar_t >::print ( ) const
inline

Print the matrix to std::cout, in a format interpretable by Matlab/Octave. The matrix will only be printed in full when not too big, i.e., when rows() <= 20 and cols() <= 32. Otherwise only its sizes and its norm are printed. Useful for debugging.

◆ print() [2/2]

template<typename scalar_t >
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.

Parameters
nameName to use when printing.
allIf true, print all values, if false, only print all values when not too big. Defaults to false.
paramwidth Specifies how many digits to use for printing floating point values, defaults to 8.

◆ print_to_file()

template<typename scalar_t >
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.

Parameters
nameName to use for the printed matrix.
filenameName of the file to write to
widthNumber of digits to use to represent floating point values, defaults to 8.

◆ ptr() [1/2]

template<typename scalar_t >
scalar_t* strumpack::DenseMatrix< scalar_t >::ptr ( std::size_t  i,
std::size_t  j 
)
inline

Pointer to element (i,j) in the matrix. This will do a bounds check with assertions, which are enabled in Debug mode, disabled in Release mode.

Parameters
iRow index, i < rows()
jColumn index, j < cols()

◆ ptr() [2/2]

template<typename scalar_t >
const scalar_t* strumpack::DenseMatrix< scalar_t >::ptr ( std::size_t  i,
std::size_t  j 
) const
inline

Const pointer to element (i,j) in the matrix. This will do a bounds check with assertions, which are enabled in Debug mode, disabled in Release mode.

Parameters
iRow index, i < rows()
jColumn index, j < cols()

◆ random() [1/2]

template<typename scalar_t >
void strumpack::DenseMatrix< scalar_t >::random ( )

Fill the matrix with random numbers, using random number generator/distribution random::make_default_random_generator<real_t>()

◆ random() [2/2]

template<typename scalar_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.

◆ read()

template<typename scalar_t >
static DenseMatrix<scalar_t> strumpack::DenseMatrix< scalar_t >::read ( const std::string &  fname)
static

Read a DenseMatrix<scalar_t> from a binary file, called fname.

See also
write

◆ resize()

template<typename scalar_t >
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.

Parameters
mNumber of rows after resizing.
mNumber of columns after resizing.

◆ rows()

template<typename scalar_t >
std::size_t strumpack::DenseMatrix< scalar_t >::rows ( ) const
inline

Number of rows of the matrix

◆ scale()

template<typename scalar_t >
DenseMatrix<scalar_t>& strumpack::DenseMatrix< scalar_t >::scale ( scalar_t  alpha,
int  depth = 0 
)

Scale this matrix by a constant factor.

Parameters
alphascaling factor

◆ scale_and_add()

template<typename scalar_t >
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.

Parameters
alphascalar factor
Bmatrix to add, should be the same size of this matrix
depthcurrent OpenMP task recursion depth

◆ scale_rows() [1/2]

template<typename scalar_t >
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].

Parameters
Dscaling vector
depthcurrent OpenMP task recursion depth
See also
scale_rows

◆ scale_rows() [2/2]

template<typename scalar_t >
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].

Parameters
Dscaling vector, D.size() == rows()
depthcurrent OpenMP task recursion depth

◆ scaled_add()

template<typename scalar_t >
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.

Parameters
alphascalar factor
Bmatrix to add, should be the same size of this matrix
depthcurrent OpenMP task recursion depth

◆ scatter_rows_add()

template<typename scalar_t >
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.

Parameters
Iindex set in this matrix, where to add the rows of B, I[i] < rows()
Bmatrix with rows to scatter in this matrix
depthcurrent OpenMP task recursion depth

◆ shift()

template<typename scalar_t >
void strumpack::DenseMatrix< scalar_t >::shift ( scalar_t  sigma)

Shift the diagonal with a value sigma, ie. add a scaled identity matrix.

Parameters
sigmascalar value to add to diagonal

◆ singular_values()

template<typename scalar_t >
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.

◆ solve()

template<typename scalar_t >
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.

Parameters
depthcurrent OpenMP task recursion depth
Returns
info from xsytrf_rook
See also
LU, Cholesky, solve_LDLt 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. The solution is returned by value.
Parameters
binput, right hand side vector/matrix
pivpivot vector returned by LU factorization
depthcurrent OpenMP task recursion depth
Returns
the solution x
See also
LU, solve_LU_in_place, solve_LDLt_in_place, solve_LDLt_rook_in_place

◆ solve_LDLt_in_place()

template<typename scalar_t >
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.

Parameters
binput, right hand side vector/matrix. On output this will be the solution.
pivpivot vector returned by LU factorization
depthcurrent OpenMP task recursion depth
See also
LDLt, LDLt_rook, solve_LDLt_rook_in_place, LU, solve_LU_in_place

◆ solve_LU_in_place() [1/2]

template<typename scalar_t >
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.

Parameters
binput, right hand side vector/matrix. On output this will be the solution.
pivpivot vector returned by LU factorization
depthcurrent OpenMP task recursion depth
See also
LU, solve_LU_in_place, solve_LDLt_in_place, solve_LDLt_rook_in_place

◆ solve_LU_in_place() [2/2]

template<typename scalar_t >
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.

Parameters
binput, right hand side vector/matrix. On output this will be the solution.
pivpivot vector returned by LU factorization
depthcurrent OpenMP task recursion depth
See also
LU, solve_LU_in_place, solve_LDLt_in_place, solve_LDLt_rook_in_place

◆ sub()

template<typename scalar_t >
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.

Parameters
Bmatrix to subtract from this matrix
depthcurrent OpenMP task recursion depth

◆ syev()

template<typename scalar_t >
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.

Parameters
jobN: Compute eigenvalues only, V: Compute eigenvalues and eigenvectors.
ulU: Upper triangle is stored, L: Lower triangle is stored.
lambdaon exit this will contain all eigenvalues of this matrix.

◆ transpose() [1/2]

template<typename scalar_t >
DenseMatrix<scalar_t> strumpack::DenseMatrix< scalar_t >::transpose ( ) const

Return the transpose of this matrix

◆ transpose() [2/2]

template<typename scalar_t >
void strumpack::DenseMatrix< scalar_t >::transpose ( DenseMatrix< scalar_t > &  X) const

Set X to the transpose of this matrix.

◆ write()

template<typename scalar_t >
void strumpack::DenseMatrix< scalar_t >::write ( const std::string &  fname) const

Write this DenseMatrix<scalar_t> to a binary file, called fname.

See also
read

◆ zero()

template<typename scalar_t >
void strumpack::DenseMatrix< scalar_t >::zero ( )

Set all matrix elements to zero


The documentation for this class was generated from the following files: