DMatrix< T_ > Class Template Reference

To handle dense matrices. More...

#include <DMatrix.h>

Inheritance diagram for DMatrix< T_ >:
Matrix< T_ >

Public Member Functions

 DMatrix ()
 Default constructor.
 
 DMatrix (size_t nr)
 Constructor for a matrix with nr rows and nr columns.
 
 DMatrix (size_t nr, size_t nc)
 Constructor for a matrix with nr rows and nc columns.
 
 DMatrix (Vect< T_ > &v)
 Constructor that uses a Vect instance. The class uses the memory space occupied by this vector.
 
 DMatrix (const DMatrix< T_ > &m)
 Copy Constructor.
 
 DMatrix (Mesh &mesh, size_t dof=0, int is_diagonal=false)
 Constructor using mesh to initialize structure of matrix.
 
 ~DMatrix ()
 Destructor.
 
void setDiag ()
 Store diagonal entries in a separate internal vector.
 
void setDiag (const T_ &a)
 Set matrix as diagonal and assign its diagonal entries as a constant.
 
void setDiag (const vector< T_ > &d)
 Set matrix as diagonal and assign its diagonal entries.
 
void setSize (size_t size)
 Set size (number of rows) of matrix.
 
void setSize (size_t nr, size_t nc)
 Set size (number of rows and columns) of matrix.
 
void getColumn (size_t j, Vect< T_ > &v) const
 Get j-th column vector.
 
Vect< T_ > getColumn (size_t j) const
 Get j-th column vector.
 
void getRow (size_t i, Vect< T_ > &v) const
 Get i-th row vector.
 
Vect< T_ > getRow (size_t i) const
 Get i-th row vector.
 
void set (size_t i, size_t j, const T_ &val)
 Assign a constant value to an entry of the matrix.
 
void reset ()
 Set matrix to 0 and reset factorization parameter.
 
void setRow (size_t i, const Vect< T_ > &v)
 Copy a given vector to a prescribed row in the matrix.
 
void setColumn (size_t j, const Vect< T_ > &v)
 Copy a given vector to a prescribed column in the matrix.
 
void MultAdd (T_ a, const Vect< T_ > &x, Vect< T_ > &y) const
 Multiply matrix by vector a*x and add result to y.
 
void MultAdd (const Vect< T_ > &x, Vect< T_ > &y) const
 Multiply matrix by vector x and add result to y.
 
void Mult (const Vect< T_ > &x, Vect< T_ > &y) const
 Multiply matrix by vector x and save result in y.
 
void TMult (const Vect< T_ > &x, Vect< T_ > &y) const
 Multiply transpose of matrix by vector x and add result in y.
 
void add (size_t i, size_t j, const T_ &val)
 Add constant val to entry (i,j) of the matrix.
 
void Axpy (T_ a, const DMatrix< T_ > &m)
 Add to matrix the product of a matrix by a scalar.
 
void Axpy (T_ a, const Matrix< T_ > *m)
 Add to matrix the product of a matrix by a scalar.
 
int setQR ()
 Construct a QR factorization of the matrix.
 
int setTransQR ()
 Construct a QR factorization of the transpose of the matrix.
 
int solveQR (const Vect< T_ > &b, Vect< T_ > &x)
 Solve a linear system by QR decomposition.
 
int solveTransQR (const Vect< T_ > &b, Vect< T_ > &x)
 Solve a transpose linear system by QR decomposition.
 
T_ at (size_t i, size_t j)
 Return a value of a matrix entry.
 
T_ operator() (size_t i, size_t j) const
 Operator () (Constant version). Return a(i,j)
 
T_ & operator() (size_t i, size_t j)
 Operator () (Non constant version). Return a(i,j)
 
int setLU ()
 Factorize the matrix (LU factorization)
 
int setTransLU ()
 Factorize the transpose of the matrix (LU factorization)
 
int solve (Vect< T_ > &b, bool fact=true)
 Solve linear system.
 
int solveTrans (Vect< T_ > &b, bool fact=true)
 Solve the transpose linear system.
 
int solve (const Vect< T_ > &b, Vect< T_ > &x, bool fact=true)
 Solve linear system.
 
int solveTrans (const Vect< T_ > &b, Vect< T_ > &x, bool fact=true)
 Solve the transpose linear system.
 
DMatrixoperator= (DMatrix< T_ > &m)
 Operator =
 
DMatrixoperator+= (const DMatrix< T_ > &m)
 Operator +=.
 
DMatrixoperator-= (const DMatrix< T_ > &m)
 Operator -=.
 
DMatrixoperator= (const T_ &x)
 Operator =
 
DMatrixoperator*= (const T_ &x)
 Operator *=
 
DMatrixoperator+= (const T_ &x)
 Operator +=
 
DMatrixoperator-= (const T_ &x)
 Operator -=
 
const T_ * getArray ()
 Return matrix as C-Array.
 
T_ get (size_t i, size_t j) const
 Return entry (i,j) of matrix.
 
void add (size_t i, const T_ &val)
 Add val to entry i.
 
- Public Member Functions inherited from Matrix< T_ >
 Matrix ()
 Default constructor.
 
 Matrix (const Matrix< T_ > &m)
 Copy Constructor.
 
virtual ~Matrix ()
 Destructor.
 
size_t getNbRows () const
 Return number of rows.
 
size_t getNbColumns () const
 Return number of columns.
 
string getName () const
 Return name of matrix.
 
MatrixSize getMatrixSize () const
 Return storage type.
 
void setPenal (real_t p)
 Set Penalty Parameter (For boundary condition prescription).
 
void setDiagonal ()
 Set the matrix as diagonal.
 
T_ getDiag (size_t k) const
 Return k-th diagonal entry of matrix.
 
size_t size () const
 Return matrix dimension (Number of rows and columns).
 
void setDiagonal (Mesh &mesh)
 Initialize matrix storage in the case where only diagonal terms are stored.
 
virtual void clear ()
 brief Set all matrix entries to zero
 
void Assembly (const Element &el, T_ *a)
 Assembly of element matrix into global matrix.
 
void Assembly (const Side &sd, T_ *a)
 Assembly of side matrix into global matrix.
 
void Prescribe (Vect< T_ > &b, const Vect< T_ > &u, int flag=0)
 Impose by a penalty method an essential boundary condition, using the Mesh instance provided by the constructor.
 
void Prescribe (int dof, int code, Vect< T_ > &b, const Vect< T_ > &u, int flag=0)
 Impose by a penalty method an essential boundary condition to a given degree of freedom for a given code.
 
void Prescribe (Vect< T_ > &b, int flag=0)
 Impose by a penalty method a homegeneous (=0) essential boundary condition.
 
void Prescribe (size_t dof, Vect< T_ > &b, const Vect< T_ > &u, int flag=0)
 Impose by a penalty method an essential boundary condition when only one DOF is treated.
 
void PrescribeSide ()
 Impose by a penalty method an essential boundary condition when DOFs are supported by sides.
 
virtual int Factor ()=0
 Factorize matrix. Available only if the storage class enables it.
 
int FactorAndSolve (Vect< T_ > &b)
 Factorize matrix and solve the linear system.
 
int FactorAndSolve (const Vect< T_ > &b, Vect< T_ > &x)
 Factorize matrix and solve the linear system.
 
size_t getLength () const
 Return number of stored terms in matrix.
 
int isDiagonal () const
 Say if matrix is diagonal or not.
 
int isFactorized () const
 Say if matrix is factorized or not.
 
virtual size_t getColInd (size_t i) const
 Return Column index for column i (See the description for class SpMatrix).
 
virtual size_t getRowPtr (size_t i) const
 Return Row pointer for row i (See the description for class SpMatrix).
 
T_ operator() (size_t i) const
 Operator () with one argument (Constant version).
 
T_ & operator() (size_t i)
 Operator () with one argument (Non Constant version).
 
T_ & operator[] (size_t k)
 Operator [] (Non constant version).
 
T_ operator[] (size_t k) const
 Operator [] (Constant version).
 
Matrixoperator= (Matrix< T_ > &m)
 Operator =.
 
Matrixoperator+= (const Matrix< T_ > &m)
 Operator +=.
 
Matrixoperator-= (const Matrix< T_ > &m)
 Operator -=.
 
Matrixoperator= (const T_ &x)
 Operator =.
 
Matrixoperator*= (const T_ &x)
 Operator *=.
 
Matrixoperator+= (const T_ &x)
 Operator +=.
 
Matrixoperator-= (const T_ &x)
 Operator -=.
 

Detailed Description

template<class T_>
class OFELI::DMatrix< T_ >

To handle dense matrices.

This class enables storing and manipulating general dense matrices. Matrices can be square or rectangle ones.

Template Parameters
T_Data type (double, float, complex<double>, ...)
Author
Rachid Touzani

Constructor & Destructor Documentation

◆ DMatrix() [1/6]

template<class T_ >
DMatrix ( )

Default constructor.

Initializes a zero-dimension matrix.

◆ DMatrix() [2/6]

template<class T_ >
DMatrix ( size_t  nr)

Constructor for a matrix with nr rows and nr columns.

Matrix entries are set to 0.

◆ DMatrix() [3/6]

template<class T_ >
DMatrix ( size_t  nr,
size_t  nc 
)

Constructor for a matrix with nr rows and nc columns.

Matrix entries are set to 0.

◆ DMatrix() [4/6]

template<class T_ >
DMatrix ( Vect< T_ > &  v)

Constructor that uses a Vect instance. The class uses the memory space occupied by this vector.

Parameters
[in]vVector to copy

◆ DMatrix() [5/6]

template<class T_ >
DMatrix ( const DMatrix< T_ > &  m)

Copy Constructor.

Parameters
[in]mMatrix to copy

◆ DMatrix() [6/6]

template<class T_ >
DMatrix ( Mesh mesh,
size_t  dof = 0,
int  is_diagonal = false 
)

Constructor using mesh to initialize structure of matrix.

Parameters
[in]meshMesh instance for which matrix graph is determined.
[in]dofOption parameter, with default value 0.
dof=1 means that only one degree of freedom for each node (or element or side) is taken to determine matrix structure. The value dof=0 means that matrix structure is determined using all DOFs.
[in]is_diagonalBoolean argument to say is the matrix is actually a diagonal matrix or not.

Member Function Documentation

◆ add()

template<class T_ >
void add ( size_t  i,
size_t  j,
const T_ &  val 
)
virtual

Add constant val to entry (i,j) of the matrix.

Parameters
[in]irow index
[in]jcolumn index
[in]valConstant to add

Implements Matrix< T_ >.

◆ at()

template<class T_ >
T_ at ( size_t  i,
size_t  j 
)
virtual

Return a value of a matrix entry.

Parameters
[in]iRow index (starts at 1)
[in]jColumn index (starts at 1)

Implements Matrix< T_ >.

◆ Axpy() [1/2]

template<class T_ >
void Axpy ( T_  a,
const DMatrix< T_ > &  m 
)

Add to matrix the product of a matrix by a scalar.

Parameters
[in]aScalar to premultiply
[in]mMatrix by which a is multiplied. The result is added to current instance

◆ Axpy() [2/2]

template<class T_ >
void Axpy ( T_  a,
const Matrix< T_ > *  m 
)
virtual

Add to matrix the product of a matrix by a scalar.

Parameters
[in]aScalar to premultiply
[in]mMatrix by which a is multiplied. The result is added to current instance

Implements Matrix< T_ >.

◆ getArray()

template<class T_ >
const T_ * getArray ( )

Return matrix as C-Array.

Matrix is stored row by row.

◆ getColumn() [1/2]

template<class T_ >
Vect< T_ > getColumn ( size_t  j) const

Get j-th column vector.

Parameters
[in]jIndex of column to extract
Returns
Vect instance where the column is stored
Remarks
Vector v does not need to be sized before. It is resized in the function

◆ getColumn() [2/2]

template<class T_ >
void getColumn ( size_t  j,
Vect< T_ > &  v 
) const

Get j-th column vector.

Parameters
[in]jIndex of column to extract
[out]vReference to Vect instance where the column is stored
Remarks
Vector v does not need to be sized before. It is resized in the function

◆ getRow() [1/2]

template<class T_ >
Vect< T_ > getRow ( size_t  i) const

Get i-th row vector.

Parameters
[in]iIndex of row to extract
Returns
Vect instance where the row is stored
Remarks
Vector v does not need to be sized before. It is resized in the function

◆ getRow() [2/2]

template<class T_ >
void getRow ( size_t  i,
Vect< T_ > &  v 
) const

Get i-th row vector.

Parameters
[in]iIndex of row to extract
[out]vReference to Vect instance where the row is stored
Remarks
Vector v does not need to be sized before. It is resized in the function

◆ Mult()

template<class T_ >
void Mult ( const Vect< T_ > &  x,
Vect< T_ > &  y 
) const
virtual

Multiply matrix by vector x and save result in y.

Parameters
[in]xVector to add to y
[out]yResult.

Implements Matrix< T_ >.

◆ MultAdd() [1/2]

template<class T_ >
void MultAdd ( const Vect< T_ > &  x,
Vect< T_ > &  y 
) const
virtual

Multiply matrix by vector x and add result to y.

Parameters
[in]xVector to add to y
[in,out]yon input, vector to add to. On output, result.

Implements Matrix< T_ >.

◆ MultAdd() [2/2]

template<class T_ >
void MultAdd ( T_  a,
const Vect< T_ > &  x,
Vect< T_ > &  y 
) const
virtual

Multiply matrix by vector a*x and add result to y.

Parameters
[in]aconstant to multiply by
[in]xVector to multiply by a
[in,out]yon input, vector to add to. On output, result.

Implements Matrix< T_ >.

◆ operator()() [1/2]

template<class T_ >
T_ & operator() ( size_t  i,
size_t  j 
)
virtual

Operator () (Non constant version). Return a(i,j)

Parameters
[in]irow index
[in]jcolumn index

Implements Matrix< T_ >.

◆ operator()() [2/2]

template<class T_ >
T_ operator() ( size_t  i,
size_t  j 
) const
virtual

Operator () (Constant version). Return a(i,j)

Parameters
[in]irow index
[in]jcolumn index

Implements Matrix< T_ >.

◆ operator*=()

template<class T_ >
DMatrix & operator*= ( const T_ &  x)

Operator *=

Premultiply matrix entries by constant value x.

◆ operator+=() [1/2]

template<class T_ >
DMatrix & operator+= ( const DMatrix< T_ > &  m)

Operator +=.

Add matrix m to current matrix instance.

◆ operator+=() [2/2]

template<class T_ >
DMatrix & operator+= ( const T_ &  x)

Operator +=

Add constant value x to matrix entries

◆ operator-=() [1/2]

template<class T_ >
DMatrix & operator-= ( const DMatrix< T_ > &  m)

Operator -=.

Subtract matrix m from current matrix instance.

◆ operator-=() [2/2]

template<class T_ >
DMatrix & operator-= ( const T_ &  x)

Operator -=

Subtract constant value x from matrix entries.

◆ operator=() [1/2]

template<class T_ >
DMatrix & operator= ( const T_ &  x)

Operator =

Assign matrix to identity times x

◆ operator=() [2/2]

template<class T_ >
DMatrix & operator= ( DMatrix< T_ > &  m)

Operator =

Copy matrix m to current matrix instance.

◆ reset()

template<class T_ >
void reset ( )
virtual

Set matrix to 0 and reset factorization parameter.

Warning
This function must be used if after a factorization, the matrix has modified

Reimplemented from Matrix< T_ >.

◆ set()

template<class T_ >
void set ( size_t  i,
size_t  j,
const T_ &  val 
)
virtual

Assign a constant value to an entry of the matrix.

Parameters
[in]irow index of matrix
[in]jcolumn index of matrix
[in]valValue to assign to a(i,j).

Implements Matrix< T_ >.

◆ setColumn()

template<class T_ >
void setColumn ( size_t  j,
const Vect< T_ > &  v 
)

Copy a given vector to a prescribed column in the matrix.

Parameters
[in]jcolumn index to be assigned
[in]vVect instance to copy

◆ setDiag() [1/2]

template<class T_ >
void setDiag ( const T_ &  a)

Set matrix as diagonal and assign its diagonal entries as a constant.

Parameters
[in]aValue to assign to all diagonal entries

◆ setDiag() [2/2]

template<class T_ >
void setDiag ( const vector< T_ > &  d)

Set matrix as diagonal and assign its diagonal entries.

Parameters
[in]dVector entries to assign to matrix diagonal entries

◆ setLU()

template<class T_ >
int setLU ( )

Factorize the matrix (LU factorization)

LU factorization of the matrix is realized. Note that since this is an in place factorization, the contents of the matrix are modified.

Returns
  • 0 if factorization was normally performed,
  • n if the n-th pivot is null.
Remarks
A flag in this class indicates after factorization that this one has been realized, so that, if the member function solve is called after this no further factorization is done.

◆ setQR()

template<class T_ >
int setQR ( )

Construct a QR factorization of the matrix.

This function constructs the QR decomposition using the Householder method. The upper triangular matrix R is returned in the upper triangle of the current matrix, except for the diagonal elements of R which are stored in an internal vector. The orthogonal matrix Q is represented as a product of n-1 Householder matrices Q1 . . . Qn-1, where Qj = 1 - uj.uj /cj . The i-th component of uj is zero for i = 1, ..., j-1 while the nonzero components are returned in a[i][j] for i = j, ..., n.

Returns
0 if the decomposition was successful, k is the k-th row is singular
Remarks
The matrix can be square or rectangle

◆ setRow()

template<class T_ >
void setRow ( size_t  i,
const Vect< T_ > &  v 
)

Copy a given vector to a prescribed row in the matrix.

Parameters
[in]irow index to be assigned
[in]vVect instance to copy

◆ setSize() [1/2]

template<class T_ >
void setSize ( size_t  nr,
size_t  nc 
)

Set size (number of rows and columns) of matrix.

Parameters
[in]nrNumber of rows.
[in]ncNumber of columns.

◆ setSize() [2/2]

template<class T_ >
void setSize ( size_t  size)

Set size (number of rows) of matrix.

Parameters
[in]sizeNumber of rows and columns.

◆ setTransLU()

template<class T_ >
int setTransLU ( )

Factorize the transpose of the matrix (LU factorization)

LU factorization of the transpose of the matrix is realized. Note that since this is an in place factorization, the contents of the matrix are modified.

Returns
  • 0 if factorization was normally performed,
  • n if the n-th pivot is null.
Remarks
A flag in this class indicates after factorization that this one has been realized, so that, if the member function solve is called after this no further factorization is done.

◆ setTransQR()

template<class T_ >
int setTransQR ( )

Construct a QR factorization of the transpose of the matrix.

This function constructs the QR decomposition using the Householder method. The upper triangular matrix R is returned in the upper triangle of the current matrix, except for the diagonal elements of R which are stored in an internal vector. The orthogonal matrix Q is represented as a product of n-1 Householder matrices Q1 . . . Qn-1, where Qj = 1 - uj.uj /cj . The i-th component of uj is zero for i = 1, ..., j-1 while the nonzero components are returned in a[i][j] for i = j, ..., n.

Returns
0 if the decomposition was successful, k is the k-th row is singular
Remarks
The matrix can be square or rectangle

◆ solve() [1/2]

template<class T_ >
int solve ( const Vect< T_ > &  b,
Vect< T_ > &  x,
bool  fact = true 
)
virtual

Solve linear system.

The linear system having the current instance as a matrix is solved by using the LU decomposition. Solution is thus realized after a factorization step and a forward/backward substitution step. The factorization step is realized only if this was not already done.
Note that this function modifies the matrix contents is a factorization is performed. Naturally, if the the matrix has been modified after using this function, the user has to refactorize it using the function setLU. This is because the class has no non-expensive way to detect if the matrix has been modified. The function setLU realizes the factorization step only.

Parameters
[in]bVect instance that contains right-hand side.
[out]xVect instance that contains solution
[in]factSet true if matrix is to be factorized (Default value), false if not
Returns
  • 0 if solution was normally performed,
  • n if the n-th pivot is null.

Implements Matrix< T_ >.

◆ solve() [2/2]

template<class T_ >
int solve ( Vect< T_ > &  b,
bool  fact = true 
)
virtual

Solve linear system.

The linear system having the current instance as a matrix is solved by using the LU decomposition. Solution is thus realized after a factorization step and a forward/backward substitution step. The factorization step is realized only if this was not already done.
Note that this function modifies the matrix contents is a factorization is performed. Naturally, if the the matrix has been modified after using this function, the user has to refactorize it using the function setLU. This is because the class has no non-expensive way to detect if the matrix has been modified. The function setLU realizes the factorization step only.

Parameters
[in,out]bVect instance that contains right-hand side on input and solution on output.
[in]factSet true if matrix is to be factorized (Default value), false if not
Returns
  • 0 if solution was normally performed,
  • n if the n-th pivot is null.

Implements Matrix< T_ >.

◆ solveQR()

template<class T_ >
int solveQR ( const Vect< T_ > &  b,
Vect< T_ > &  x 
)

Solve a linear system by QR decomposition.

This function constructs the QR decomposition, if this was not already done by using the member function QR and solves the linear system

Parameters
[in]bRight-hand side vector
[out]xSolution vector. Must have been sized before using this function.
Returns
The same value as returned by the function QR

◆ solveTrans() [1/2]

template<class T_ >
int solveTrans ( const Vect< T_ > &  b,
Vect< T_ > &  x,
bool  fact = true 
)

Solve the transpose linear system.

The linear system having the current instance as a transpose matrix is solved by using the LU decomposition. Solution is thus realized after a factorization step and a forward/backward substitution step. The factorization step is realized only if this was not already done.
Note that this function modifies the matrix contents is a factorization is performed. Naturally, if the the matrix has been modified after using this function, the user has to refactorize it using the function setLU. This is because the class has no non-expensive way to detect if the matrix has been modified. The function setLU realizes the factorization step only.

Parameters
[in]bVect instance that contains right-hand side.
[out]xVect instance that contains solution
[in]factSet true if matrix is to be factorized (Default value), false if not
Returns
  • 0 if solution was normally performed,
  • n if the n-th pivot is null.

◆ solveTrans() [2/2]

template<class T_ >
int solveTrans ( Vect< T_ > &  b,
bool  fact = true 
)

Solve the transpose linear system.

The linear system having the current instance as a transpose matrix is solved by using the LU decomposition. Solution is thus realized after a factorization step and a forward/backward substitution step. The factorization step is realized only if this was not already done.
Note that this function modifies the matrix contents is a factorization is performed. Naturally, if the the matrix has been modified after using this function, the user has to refactorize it using the function setLU. This is because the class has no non-expensive way to detect if the matrix has been modified. The function setLU realizes the factorization step only.

Parameters
[in,out]bVect instance that contains right-hand side on input and solution on output.
[in]factSet true if matrix is to be factorized (Default value), false if not
Returns
  • 0 if solution was normally performed,
  • n if the n-th pivot is null.

◆ solveTransQR()

template<class T_ >
int solveTransQR ( const Vect< T_ > &  b,
Vect< T_ > &  x 
)

Solve a transpose linear system by QR decomposition.

This function constructs the QR decomposition, if this was not already done by using the member function QR and solves the linear system

Parameters
[in]bRight-hand side vector
[out]xSolution vector. Must have been sized before using this function.
Returns
The same value as returned by the function QR

◆ TMult()

template<class T_ >
void TMult ( const Vect< T_ > &  x,
Vect< T_ > &  y 
) const
virtual

Multiply transpose of matrix by vector x and add result in y.

Parameters
[in]xVector to add to y
[in,out]yon input, vector to add to. On output, result.

Implements Matrix< T_ >.