SkMatrix< T_ > Class Template Reference

To handle square matrices in skyline storage format. More...

#include <SkMatrix.h>

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

Public Member Functions

 SkMatrix ()
 Default constructor.
 
 SkMatrix (size_t size, int is_diagonal=false)
 Constructor that initializes a dense symmetric matrix.
 
 SkMatrix (Mesh &mesh, size_t dof=0, int is_diagonal=false)
 Constructor using mesh to initialize skyline structure of matrix.
 
 SkMatrix (const vector< size_t > &ColHt)
 Constructor that initializes skyline structure of matrix using vector of column heights.
 
 SkMatrix (const SkMatrix< T_ > &m)
 Copy Constructor.
 
 ~SkMatrix ()
 Destructor.
 
void setMesh (Mesh &mesh, size_t dof=0)
 Determine mesh graph and initialize matrix.
 
void setSkyline (Mesh &mesh)
 Determine matrix structure.
 
void setDiag ()
 Store diagonal entries in a separate internal vector.
 
void setDOF (size_t i)
 Choose DOF to activate.
 
void set (size_t i, size_t j, const T_ &val)
 Assign a value to an entry ofthe matrix.
 
void Axpy (T_ a, const SkMatrix< 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.
 
void MultAdd (const Vect< T_ > &x, Vect< T_ > &y) const
 Multiply matrix by vector x and add to y.
 
void TMultAdd (const Vect< T_ > &x, Vect< T_ > &y) const
 Multiply transpose of matrix by vector x and add to y.
 
void MultAdd (T_ a, const Vect< T_ > &x, Vect< T_ > &y) const
 Multiply matrix by a vector and add to another one.
 
void Mult (const Vect< T_ > &x, Vect< T_ > &y) const
 Multiply matrix by vector x and save in y.
 
void TMult (const Vect< T_ > &x, Vect< T_ > &y) const
 Multiply transpose of matrix by vector x and save in y.
 
void add (size_t i, size_t j, const T_ &val)
 Add a constant value to an entry ofthe matrix.
 
size_t getColHeight (size_t i) const
 Return column height.
 
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).
 
T_ & operator() (size_t i, size_t j)
 Operator () (Non constant version).
 
void DiagPrescribe (Mesh &mesh, Vect< T_ > &b, const Vect< T_ > &u, int flag=0)
 Impose an essential boundary condition.
 
void DiagPrescribe (Vect< T_ > &b, const Vect< T_ > &u, int flag=0)
 Impose an essential boundary condition using the Mesh instance provided by the constructor.
 
SkMatrix< T_ > & operator= (const SkMatrix< T_ > &m)
 Operator =.
 
SkMatrix< T_ > & operator= (const T_ &x)
 Operator =.
 
SkMatrix< T_ > & operator+= (const SkMatrix< T_ > &m)
 Operator +=.
 
SkMatrix< T_ > & operator+= (const T_ &x)
 Operator +=.
 
SkMatrix< T_ > & operator*= (const T_ &x)
 Operator *=.
 
int setLU ()
 Factorize the matrix (LU factorization)
 
int solve (Vect< T_ > &b, bool fact=true)
 Solve linear system.
 
int solve (const Vect< T_ > &b, Vect< T_ > &x, bool fact=true)
 Solve linear system.
 
T_ * get () const
 Return C-Array.
 
T_ get (size_t i, size_t j) const
 Return entry (i,j) of matrix if this one is stored, 0 else.
 
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.
 
virtual void reset ()
 Set matrix to 0 and reset factorization parameter.
 
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::SkMatrix< T_ >

To handle square matrices in skyline storage format.

This template class allows storing and manipulating a matrix in skyline storage format.

The matrix entries are stored in 2 vectors column by column as in the following example:

/                    \        /                       \ 
| l0            .    |        | u0   u1    0   0   u7 |
| l1  l2        .    |        |      u2   u3   0   u8 |
|  0  l3  l4    .    |        | ...       u4  u5   u9 |
|  0   0  l5  l6     |        |               u6  u10 |
| l7  l8  l9 l10 l11 |        |                   u11 |
\                    /        \                       /
Template Parameters
T_Data type (double, float, complex<double>, ...)
Author
Rachid Touzani

Constructor & Destructor Documentation

◆ SkMatrix() [1/4]

template<class T_ >
SkMatrix ( )

Default constructor.

Initializes a zero-dimension matrix

◆ SkMatrix() [2/4]

template<class T_ >
SkMatrix ( size_t  size,
int  is_diagonal = false 
)

Constructor that initializes a dense symmetric matrix.

Normally, for a dense matrix this is not the right class.

Parameters
[in]sizeNumber of matrix rows (and columns).
[in]is_diagonalBoolean to select if the matrix is diagonal or not [Default: false]

◆ SkMatrix() [3/4]

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

Constructor using mesh to initialize skyline 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.

◆ SkMatrix() [4/4]

template<class T_ >
SkMatrix ( const vector< size_t > &  ColHt)

Constructor that initializes skyline structure of matrix using vector of column heights.

Parameters
[in]ColHtVect instance that contains rows lengths of matrix.

Member Function Documentation

◆ add()

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

Add a constant value to an entry ofthe matrix.

Parameters
[in]iRow index
[in]jColumn index
[in]valConstant value to add to a(i,j)

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 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_ >.

◆ Axpy() [2/2]

template<class T_ >
void Axpy ( T_  a,
const SkMatrix< 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

◆ DiagPrescribe() [1/2]

template<class T_ >
void DiagPrescribe ( Mesh mesh,
Vect< T_ > &  b,
const Vect< T_ > &  u,
int  flag = 0 
)

Impose an essential boundary condition.

This member function modifies diagonal terms in matrix and terms in vector that correspond to degrees of freedom with nonzero code in order to impose a boundary condition. It can be modified by member function setPenal(..).

Parameters
[in]meshMesh instance from which information is extracted.
[in]bVect instance that contains right-hand side.
[in]uVect instance that conatins imposed valued at DOFs where they are to be imposed.
[in]flagParameter to determine whether only the right-hand side is to be modified (dof>0)
or both matrix and right-hand side (dof=0, default value).

◆ DiagPrescribe() [2/2]

template<class T_ >
void DiagPrescribe ( Vect< T_ > &  b,
const Vect< T_ > &  u,
int  flag = 0 
)

Impose an essential boundary condition using the Mesh instance provided by the constructor.

This member function modifies diagonal terms in matrix and terms in vector that correspond to degrees of freedom with nonzero code in order to impose a boundary condition. It can be modified by member function setPenal(..).

Parameters
[in]bVect instance that contains right-hand side.
[in]uVect instance that conatins imposed valued at DOFs where they are to be imposed.
[in]flagParameter to determine whether only the right-hand side is to be modified (dof>0)
or both matrix and right-hand side (dof=0, default value).

◆ get()

template<class T_ >
T_ * get ( ) const

Return C-Array.

Skyline of matrix is stored row by row.

◆ getColHeight()

template<class T_ >
size_t getColHeight ( size_t  i) const

Return column height.

Column height at entry i is returned.

◆ Mult()

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

Multiply matrix by vector x and save in y.

Parameters
[in]xVector to multiply by matrix
[out]yVector that contains on output the result.

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 to y.

Parameters
[in]xVector to multiply by matrix
[in,out]yVector to add to the result. y contains on output the 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 a vector and add to another one.

Parameters
[in]aConstant to multiply by matrix
[in]xVector to multiply by matrix
[in,out]yVector to add to the result. y contains on output the result.

Implements Matrix< T_ >.

◆ operator()() [1/2]

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

Operator () (Non constant version).

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).

Parameters
[in]iRow index
[in]jColumn index

Implements Matrix< T_ >.

◆ operator*=()

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

Operator *=.

Premultiply matrix entries by constant value x.

◆ operator+=() [1/2]

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

Operator +=.

Add matrix m to current matrix instance.

◆ operator+=() [2/2]

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

Operator +=.

Add constant value x to matrix entries.

◆ operator=() [1/2]

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

Operator =.

Copy matrix m to current matrix instance.

◆ operator=() [2/2]

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

Operator =.

define the matrix as a diagonal one with all diagonal entries equal to x.

◆ set()

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

Assign a value to an entry ofthe matrix.

Parameters
[in]iRow index (starting at i=1)
[in]jColumn index (starting at i=1)
[in]valValue to assign to entry a(i,j)

Implements Matrix< T_ >.

◆ setDOF()

template<class T_ >
void setDOF ( size_t  i)

Choose DOF to activate.

This function is available only if variable dof is equal to 1 in the constructor

Parameters
[in]iIndex of the DOF

◆ 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.

◆ setMesh()

template<class T_ >
void setMesh ( Mesh mesh,
size_t  dof = 0 
)

Determine mesh graph and initialize matrix.

This member function is called by constructor with the same arguments

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.

◆ setSkyline()

template<class T_ >
void setSkyline ( Mesh mesh)

Determine matrix structure.

This member function calculates matrix structure using a Mesh instance.

Parameters
[in]meshMesh instance

◆ 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_ >.

◆ TMult()

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

Multiply transpose of matrix by vector x and save in y.

Parameters
[in]xVector to multiply by matrix
[out]yVector that contains on output the result.

Implements Matrix< T_ >.

◆ TMultAdd()

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

Multiply transpose of matrix by vector x and add to y.

Parameters
[in]xVector to multiply by matrix
[in,out]yVector to add to the result. y contains on output the result.