strumpack::DenseMatrixWrapper< scalar_t > Class Template Reference

Like DenseMatrix, this class represents a matrix, stored in column major format, to allow direct use of BLAS/LAPACK routines. However, objects of the DenseMatrixWrapper class do not allocate, own or free any data. More...

#include <DenseMatrix.hpp>

Inheritance diagram for strumpack::DenseMatrixWrapper< scalar_t >:
Collaboration diagram for strumpack::DenseMatrixWrapper< scalar_t >:

Public Member Functions

 DenseMatrixWrapper ()
 
 DenseMatrixWrapper (std::size_t m, std::size_t n, scalar_t *D, std::size_t ld)
 
 DenseMatrixWrapper (std::size_t m, std::size_t n, DenseMatrix< scalar_t > &D, std::size_t i, std::size_t j)
 
virtual ~DenseMatrixWrapper ()
 
void clear () override
 
std::size_t memory () const override
 
std::size_t nonzeros () const override
 
 DenseMatrixWrapper (const DenseMatrixWrapper< scalar_t > &)=default
 
 DenseMatrixWrapper (const DenseMatrix< scalar_t > &)=delete
 
 DenseMatrixWrapper (DenseMatrixWrapper< scalar_t > &&)=default
 
 DenseMatrixWrapper (DenseMatrix< scalar_t > &&)=delete
 
DenseMatrixWrapper< scalar_t > & operator= (DenseMatrixWrapper< scalar_t > &&D)
 
DenseMatrix< scalar_t > & operator= (const DenseMatrix< scalar_t > &a) override
 
- Public Member Functions inherited from strumpack::DenseMatrix< scalar_t >
 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
 

Additional Inherited Members

- Static Public Member Functions inherited from strumpack::DenseMatrix< scalar_t >
static DenseMatrix< scalar_t > read (const std::string &fname)
 

Detailed Description

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

Like DenseMatrix, this class represents a matrix, stored in column major format, to allow direct use of BLAS/LAPACK routines. However, objects of the DenseMatrixWrapper class do not allocate, own or free any data.

The user has to make sure that the memory that was used to create this matrix stays valid for as long as this matrix wrapper, or any other wrappers derived from this wrapper (submatrices), is required. The DenseMatrixWrapper class is a subclass of DenseMatrix, so a DenseMatrixWrapper can be used as a DenseMatrix. A DenseMatrixWrapper can be used to wrap already allocated memory, or to create a submatrix of a DenseMatrix.

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

Constructor & Destructor Documentation

◆ DenseMatrixWrapper() [1/7]

template<typename scalar_t >
strumpack::DenseMatrixWrapper< scalar_t >::DenseMatrixWrapper ( )
inline

Default constructor. Creates an empty, 0x0 matrix.

◆ DenseMatrixWrapper() [2/7]

template<typename scalar_t >
strumpack::DenseMatrixWrapper< scalar_t >::DenseMatrixWrapper ( std::size_t  m,
std::size_t  n,
scalar_t *  D,
std::size_t  ld 
)
inline

Constructor. Create an m x n matrix wrapper using already allocated memory, pointed to by D, with leading dimension ld.

Parameters
mnumber of rows of the new (sub) matrix
nnumber of columns of the new matrix
Dpointer to memory representing matrix, this should point to at least ld*n bytes of allocated memory
ldleading dimension of matrix allocated at D. ld >= m

◆ DenseMatrixWrapper() [3/7]

template<typename scalar_t >
strumpack::DenseMatrixWrapper< scalar_t >::DenseMatrixWrapper ( std::size_t  m,
std::size_t  n,
DenseMatrix< scalar_t > &  D,
std::size_t  i,
std::size_t  j 
)
inline

Constructor. Create a DenseMatrixWrapper as a submatrix of size m x n, of a DenseMatrix (or DenseMatrixWrapper) D, at position i,j in D. The constructed DenseMatrixWrapper will be the submatrix D(i:i+m,j:j+n).

Parameters
mnumber of rows of the new (sub) matrix
nnumber of columns of the new matrix
Dmatrix from which to take a submatrix
irow offset in D of the top left corner of the submatrix
jcolumns offset in D of the top left corner of the submatrix

◆ ~DenseMatrixWrapper()

template<typename scalar_t >
virtual strumpack::DenseMatrixWrapper< scalar_t >::~DenseMatrixWrapper ( )
inlinevirtual

Virtual destructor. Since a DenseMatrixWrapper does not actually own it's memory, put just keeps a pointer, this will not free any memory.

◆ DenseMatrixWrapper() [4/7]

template<typename scalar_t >
strumpack::DenseMatrixWrapper< scalar_t >::DenseMatrixWrapper ( const DenseMatrixWrapper< scalar_t > &  )
default

Default copy constructor, from another DenseMatrixWrapper.

◆ DenseMatrixWrapper() [5/7]

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

Constructing a DenseMatrixWrapper from a DenseMatrixWrapper is not allowed. TODO Why not??!! just delegate to DenseMatrixWrapper(m, n, D, i, j)??

◆ DenseMatrixWrapper() [6/7]

template<typename scalar_t >
strumpack::DenseMatrixWrapper< scalar_t >::DenseMatrixWrapper ( DenseMatrixWrapper< scalar_t > &&  )
default

Default move constructor.

◆ DenseMatrixWrapper() [7/7]

template<typename scalar_t >
strumpack::DenseMatrixWrapper< scalar_t >::DenseMatrixWrapper ( DenseMatrix< scalar_t > &&  )
delete

Moving from a DenseMatrix is not allowed.

Member Function Documentation

◆ clear()

template<typename scalar_t >
void strumpack::DenseMatrixWrapper< scalar_t >::clear ( )
inlineoverridevirtual

Clear the DenseMatrixWrapper. Ie, set to an empty matrix. This will not affect the original matrix, to which this is a wrapper, only the wrapper itself is reset. No memory is released.

Reimplemented from strumpack::DenseMatrix< scalar_t >.

◆ memory()

template<typename scalar_t >
std::size_t strumpack::DenseMatrixWrapper< scalar_t >::memory ( ) const
inlineoverridevirtual

Return the amount of memory taken by this wrapper, ie, 0. (since the wrapper itself does not own the memory). The memory will likely be owned by a DenseMatrix, while this DenseMatrixWrapper is just a submatrix of that existing matrix. Returning 0 here avoids counting memory double.

See also
nonzeros

Reimplemented from strumpack::DenseMatrix< scalar_t >.

◆ nonzeros()

template<typename scalar_t >
std::size_t strumpack::DenseMatrixWrapper< scalar_t >::nonzeros ( ) const
inlineoverridevirtual

Return the number of nonzeros taken by this wrapper, ie, 0. (since the wrapper itself does not own the memory). The memory will likely be owned by a DenseMatrix, while this DenseMatrixWrapper is just a submatrix of that existing matrix. Returning 0 here avoids counting nonzeros double.

See also
memory

Reimplemented from strumpack::DenseMatrix< scalar_t >.

◆ operator=() [1/2]

template<typename scalar_t >
DenseMatrix< scalar_t > & strumpack::DenseMatrixWrapper< scalar_t >::operator= ( const DenseMatrix< scalar_t > &  a)
inlineoverridevirtual

Assignment operator, from a DenseMatrix. Assign the memory of the DenseMatrix to the matrix wrapped by this DenseMatrixWrapper object.

Parameters
amatrix to copy from, should be a.rows() == this->rows() and a.cols() == this->cols()

Reimplemented from strumpack::DenseMatrix< scalar_t >.

◆ operator=() [2/2]

template<typename scalar_t >
DenseMatrixWrapper< scalar_t > & strumpack::DenseMatrixWrapper< scalar_t >::operator= ( DenseMatrixWrapper< scalar_t > &&  D)
inline

Move assignment. This moves only the wrapper.

Parameters
Dmatrix wrapper to move from. This will not be modified.

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