File kernel/src/utils/SiconosAlgebra/SiconosMatrix.hpp#
Go to the source code of this file
Interface for matrices handling.
Typedefs
-
typedef std::vector<int> VInt#
A STL vector of int.
-
union MATRIX_UBLAS_TYPE#
- #include <SiconosMatrix.hpp>
Union of DenseMat pointer, TriangMat pointer BandedMat, SparseMat, SymMat, Zero and Identity mat pointers.
-
class SiconosMatrix
- #include <SiconosMatrix.hpp>
Abstract class to provide interface for matrices handling.
Matrices can be either block or Simple. See Derived classes for details.
In Siconos, a “matrix” can be either a SimpleMatrix or a BlockMatrix, ie a container of several pointers to SiconosMatrix
You can find an overview on how to build and use vectors and matrices in siconos users guide .
Subclassed by BlockMatrix, SimpleMatrix
Public Functions
-
inline virtual ~SiconosMatrix()
Destructor.
-
inline bool isBlock(void) const
true if the matrix is block else false.
- Returns:
a bool
-
inline virtual bool isPLUInversed() const
determines if the matrix has been inversed in place
- Returns:
true if the matrix is inversed
-
inline bool isSymmetric() const
true if the matrix is symmetric (the flag is just returned)
- Returns:
true if the matrix is symmetric
-
inline void setIsSymmetric(bool b)
set the flag _isSymmetric
-
inline bool isPositiveDefinite() const
true if the matrix is definite positive (the flag is just returned)
- Returns:
true if the matrix is
-
inline void setIsPositiveDefinite(bool b)
set the flag _isPositiveDefinite
-
virtual bool checkSymmetry(double tol) const = 0
determines if the matrix is symmetric up to a given tolerance
- Returns:
true if the matrix is inversed
-
inline virtual bool isPLUFactorized() const
determines if the matrix has been PLU factorized
- Returns:
true if the matrix is factorized
-
inline virtual bool isPLUFactorizedInPlace() const
determines if the matrix has been PLU factorized in place
- Returns:
true if the matrix is factorized
-
inline virtual bool isCholeskyFactorized() const
determines if the matrix has been Cholesky factorized
- Returns:
true if the matrix is factorized
-
inline bool isQRFactorized() const
determines if the matrix has been QR factorized
- Returns:
true if the matrix is factorized
-
inline virtual bool isFactorized() const
determines if the matrix has been factorized
- Returns:
true if the matrix is factorized
-
virtual unsigned int size(unsigned int index) const = 0
get the number of rows or columns of the matrix
- Parameters:
index – 0 for rows, 1 for columns
- Returns:
an int
-
inline siconos::UBLAS_TYPE num() const
get the attribute num of current matrix
- Returns:
an unsigned int.
-
inline virtual unsigned int numberOfBlocks(unsigned int i) const
get the number of block (i=0, row, i=1 col)
- Parameters:
i – unsigned int(i=0, row, i=1 col)
- Returns:
an unsigned int. 1 as default for SimpleMatrix.
-
virtual const SP::Index tabRow() const
reserved to BlockMatrix - get the index tab for rows
- Returns:
a pointer to a standard vector of int
-
virtual const SP::Index tabCol() const
reserved to BlockMatrix - get the index tab of columns
- Returns:
a pointer to a standard vector of int
-
virtual const DenseMat getDense(unsigned int row = 0, unsigned int col = 0) const = 0
get DenseMat matrix
- Parameters:
row – an unsigned int position of the block (row) - Useless for SimpleMatrix
col – an unsigned int position of the block (column) - Useless for SimpleMatrix
- Returns:
a DenseMat
-
virtual const TriangMat getTriang(unsigned int row = 0, unsigned int col = 0) const = 0
get TriangMat matrix
- Parameters:
row – an unsigned int, position of the block (row) - Useless for SimpleMatrix
col – an unsigned int, position of the block (column) - Useless for SimpleMatrix
- Returns:
a TriangMat
-
virtual const SymMat getSym(unsigned int row = 0, unsigned int col = 0) const = 0
get SymMat matrix
- Parameters:
row – an unsigned int, position of the block (row) - Useless for SimpleMatrix
col – an unsigned int, position of the block (column) - Useless for SimpleMatrix
- Returns:
a SymMat
-
virtual const BandedMat getBanded(unsigned int row = 0, unsigned int col = 0) const = 0
get BandedMat matrix
- Parameters:
row – an unsigned int, position of the block (row) - Useless for SimpleMatrix
col – an unsigned int, position of the block (column) - Useless for SimpleMatrix
- Returns:
a BandedMat
-
virtual const SparseMat getSparse(unsigned int row = 0, unsigned int col = 0) const = 0
get SparseMat matrix
- Parameters:
row – an unsigned int, position of the block (row) - Useless for SimpleMatrix
col – an unsigned int, position of the block (column) - Useless for SimpleMatrix
- Returns:
a SparseMat
-
virtual const SparseCoordinateMat getSparseCoordinate(unsigned int row = 0, unsigned int col = 0) const = 0
get SparseCoordinateMat matrix
- Parameters:
row – an unsigned int, position of the block (row) - Useless for SimpleMatrix
col – an unsigned int, position of the block (column) - Useless for SimpleMatrix
- Returns:
a SparseCoordinateMat
-
virtual const ZeroMat getZero(unsigned int row = 0, unsigned int col = 0) const = 0
get ZeroMat matrix
- Parameters:
row – an unsigned int, position of the block (row) - Useless for SimpleMatrix
col – an unsigned int, position of the block (column) - Useless for SimpleMatrix
- Returns:
a ZeroMat
-
virtual const IdentityMat getIdentity(unsigned int row = 0, unsigned int col = 0) const = 0
get getIdentity matrix
- Parameters:
row – an unsigned int, position of the block (row) - Useless for SimpleMatrix
col – an unsigned int, position of the block (column) - Useless for SimpleMatrix
- Returns:
an IdentityMat
-
virtual DenseMat *dense(unsigned int row = 0, unsigned int col = 0) const = 0
get a pointer on DenseMat matrix
- Parameters:
row – an unsigned int, position of the block (row) - Useless for SimpleMatrix
col – an unsigned int, position of the block (column) - Useless for SimpleMatrix
- Returns:
a DenseMat*
-
virtual TriangMat *triang(unsigned int row = 0, unsigned int col = 0) const = 0
get a pointer on TriangMat matrix
- Parameters:
row – an unsigned int, position of the block (row) - Useless for SimpleMatrix
col – an unsigned int, position of the block (column) - Useless for SimpleMatrix
- Returns:
a TriangMat*
-
virtual SymMat *sym(unsigned int row = 0, unsigned int col = 0) const = 0
get a pointer on SymMat matrix
- Parameters:
row – an unsigned int, position of the block (row) - Useless for SimpleMatrix
col – an unsigned int, position of the block (column) - Useless for SimpleMatrix
- Returns:
a SymMat*
-
virtual BandedMat *banded(unsigned int row = 0, unsigned int col = 0) const = 0
get a pointer on BandedMat matrix
- Parameters:
row – an unsigned int, position of the block (row) - Useless for SimpleMatrix
col – an unsigned int, position of the block (column) - Useless for SimpleMatrix
- Returns:
a BandedMat*
-
virtual SparseMat *sparse(unsigned int row = 0, unsigned int col = 0) const = 0
get a pointer on SparseMat matrix
- Parameters:
row – an unsigned int, position of the block (row) - Useless for SimpleMatrix
col – an unsigned int, position of the block (column) - Useless for SimpleMatrix
- Returns:
a SparseMat*
-
virtual SparseCoordinateMat *sparseCoordinate(unsigned int row = 0, unsigned int col = 0) const = 0
get a pointer on SparseCoordinateMat matrix
- Parameters:
row – an unsigned int, position of the block (row) - Useless for SimpleMatrix
col – an unsigned int, position of the block (column) - Useless for SimpleMatrix
- Returns:
a SparseCoordinateMat*
-
virtual ZeroMat *zero_mat(unsigned int row = 0, unsigned int col = 0) const = 0
get a pointer on ZeroMat matrix
- Parameters:
row – an unsigned int, position of the block (row) - Useless for SimpleMatrix
col – an unsigned int, position of the block (column) - Useless for SimpleMatrix
- Returns:
a ZeroMat*
-
virtual IdentityMat *identity(unsigned int row = 0, unsigned int col = 0) const = 0
get a pointer on Identity matrix
- Parameters:
row – an unsigned int, position of the block (row) - Useless for SimpleMatrix
col – an unsigned int, position of the block (column) - Useless for SimpleMatrix
- Returns:
an IdentityMat*
-
virtual double *getArray(unsigned int row = 0, unsigned int col = 0) const = 0
return the address of the array of double values of the matrix ( for block(i,j) if this is a block matrix)
- Parameters:
row – position for the required block
col – position for the required block
- Returns:
double* : the pointer on the double array
-
virtual void zero() = 0
sets all the values of the matrix to 0.0
-
virtual void randomize() = 0
Initialize the matrix with random values.
-
virtual void eye() = 0
set an identity matrix
-
virtual void resize(unsigned int nbrow, unsigned int nbcol, unsigned int lower = 0, unsigned int upper = 0, bool preserve = true) = 0
resize the matrix with nbrow rows and nbcol columns, upper and lower are only useful for BandedMatrix .
The existing elements of the matrix are preseved when specified.
- Parameters:
nbrow –
nbcol –
lower, upper – for banded matrices
preserve –
-
virtual double normInf() const = 0
compute the infinite norm of the matrix
- Returns:
a double
-
virtual void display() const = 0
display data on standard output
-
virtual void displayExpert(bool brief = true) const = 0
display data on standard output
-
virtual std::string toString() const = 0
put data of the matrix into a std::string
- Returns:
std::string
-
virtual double &operator()(unsigned int i, unsigned int j) = 0
get or set the element matrix[i,j]
- Parameters:
i – an unsigned int i
j – an unsigned int j
- Returns:
the element matrix[i,j]
-
virtual double operator()(unsigned int i, unsigned int j) const = 0
get or set the element matrix[i,j]
- Parameters:
i – an unsigned int i
j – an unsigned int j
- Returns:
the element matrix[i,j]
-
virtual double getValue(unsigned int i, unsigned int j) const = 0
return the element matrix[i,j]
- Parameters:
i – an unsigned int i
j – an unsigned int j
- Returns:
a double
-
virtual void setValue(unsigned int i, unsigned int j, double value) = 0
set the element matrix[i,j]
- Parameters:
i – an unsigned int i
j – an unsigned int j
value –
-
inline virtual SP::SiconosMatrix block(unsigned int row = 0, unsigned int col = 0)
get block at position row-col if BlockMatrix, else if SimpleMatrix return this
- Parameters:
row – unsigned int row
col – unsigned int col
- Returns:
SP::SiconosMatrix
-
inline virtual SPC::SiconosMatrix block(unsigned int row = 0, unsigned int col = 0) const
get block at position row-col if BlockMatrix, else if SimpleMatrix return this
- Parameters:
row – unsigned int row
col – unsigned int col
- Returns:
SPC::SiconosMatrix
-
virtual void getRow(unsigned int index, SiconosVector &vOut) const = 0
get row index of current matrix and save it into vOut
- Parameters:
index – row we want to get
vOut – [out] SiconosVector that will contain the desired row
-
virtual void getCol(unsigned int index, SiconosVector &vOut) const = 0
get column index of current matrix and save it into vOut
- Parameters:
index – column we want to get
vOut – [out] SiconosVector that will contain the desired column
-
virtual void setRow(unsigned int index, const SiconosVector &vIn) = 0
set line row of the current matrix with vector v
- Parameters:
index – row we want to set
vIn – SiconosVector containing the new row
-
virtual void setCol(unsigned int index, const SiconosVector &vIn) = 0
set column col of the current matrix with vector v
- Parameters:
index – column we want to set
vIn – a SiconosVector containing the new column
-
virtual void trans() = 0
transpose in place: x->trans() is x = transpose of x.
-
virtual void trans(const SiconosMatrix &m) = 0
transpose a matrix: x->trans(m) is x = transpose of m.
- Parameters:
m – the matrix to be transposed.
-
virtual SiconosMatrix &operator=(const SiconosMatrix &m) = 0
operator =
- Parameters:
m – the matrix to be copied
- Returns:
-
virtual SiconosMatrix &operator=(const DenseMat &m) = 0
operator = to a DenseMat
- Parameters:
m – the DenseMat to be copied
- Returns:
-
virtual SiconosMatrix &operator+=(const SiconosMatrix &m) = 0
operator +=
- Parameters:
m – a matrix to add
- Returns:
-
virtual SiconosMatrix &operator-=(const SiconosMatrix &m) = 0
operator -=
- Parameters:
m – a matrix to subtract
- Returns:
-
virtual void updateNumericsMatrix() = 0#
-
inline virtual NumericsMatrix *numericsMatrix() const#
-
virtual void PLUFactorizationInPlace() = 0
computes a LU factorization of a general M-by-N matrix with partial pivoting and row interchanges.
The result is returned in this (InPlace). Based on Blas dgetrf function for dense matrix and ublas cholesky decomposition for sparse matrix (work only for a symmetric matrix and very slow because it uses matric accessor) use preferably PLUFactorize()
-
virtual void Factorize() = 0
computes a factorization of a general M-by-N matrix The implementation is based on an internal NumericsMatrix
-
virtual void PLUInverseInPlace() = 0
compute inverse of this thanks to LU factorization with partial pivoting.
This method inverts U and then computes inv(A) by solving the system inv(A)*L = inv(U) for inv(A). The result is returned in this (InPlace). Based on Blas dgetri function for dense function
-
virtual void PLUForwardBackwardInPlace(SiconosMatrix &B) = 0
solves a system of linear equations A * X = B (A=this) for a general N-by-N matrix A using the LU factorization computed by PLUFactorizationInPlace.
Based on Blas dgetrs function for dense matrix.
- Parameters:
B – [inout] on input the RHS matrix b; on output the result x
-
virtual void Solve(SiconosMatrix &B) = 0
solves a system of linear equations A * X = B (A=this) for a general N-by-N matrix A using the LU factorization computed by PLUFactorize.
- Parameters:
B – [inout] on input the RHS matrix b; on output the result x
-
virtual void PLUForwardBackwardInPlace(SiconosVector &B) = 0
solves a system of linear equations A * X = B (A=this) for a general N-by-N matrix A using the LU factorization computed by PLUFactorizationInPlace.
Based on Blas dgetrs function for dense matrix.
- Parameters:
B – [inout] on input the RHS matrix b; on output the result x
-
virtual void Solve(SiconosVector &B) = 0
solves a system of linear equations A * X = B (A=this) for a general N-by-N matrix A using the LU factorization computed by PLUFactorize.
- Parameters:
B – [inout] on input the RHS matrix b; on output the result x
-
inline virtual void resetLU()
set to false all LU indicators.
Useful in case of assignment for example.
-
inline virtual void resetFactorizationFlags()
set to false all factorization indicators.
Useful in case of assignment for example.
-
virtual size_t nnz(double tol = 1e-14)
return the number of non-zero in the matrix
- Parameters:
tol – the tolerance to consider a number zero (not used if the matrix is sparse)
- Returns:
the number of non-zeros
-
bool fillCSC(CSparseMatrix *csc, size_t row_off, size_t col_off, double tol = 1e-14)
Fill CSparseMatrix compresses column sparse matrix.
Warning
not clear that it works for an empty csr matrix with row_off =0 and col_off =0
- Parameters:
csc – the compressed column sparse matrix
row_off –
col_off –
tol – the tolerance under which a number is considered as equal to zero
- Returns:
true if function worked.
-
bool fillCSC(CSparseMatrix *csc, double tol = 1e-14)
Fill CSparseMatrix compresses column sparse matrix.
- Parameters:
csc – the compressed column sparse matrix
tol – the tolerance under which a number is considered as equal to zero
- Returns:
true if function worked.
-
bool fromCSC(CSparseMatrix *csc)#
-
bool fillTriplet(CSparseMatrix *csc, size_t row_off, size_t col_off, double tol = 1e-14)
return the number of non-zero in the matrix
- Parameters:
csc – the compressed column sparse matrix
row_off –
col_off –
tol – the tolerance to consider a number zero (not used if the matrix is sparse)
- Returns:
the number of non-zeros
-
VIRTUAL_ACCEPT_VISITORS(SiconosMatrix)#
Protected Functions
-
ACCEPT_SERIALIZATION(SiconosMatrix)#
-
inline SiconosMatrix()#
default constructor
-
SiconosMatrix(siconos::UBLAS_TYPE type)#
basic constructor
- Parameters:
type – unsigned int type-number of the vector
-
void private_prod(unsigned int startRow, const SiconosVector &x, SiconosVector &y, bool init) const#
computes y = subA*x (init =true) or += subA * x (init = false), subA being a submatrix of A (all columns, and rows between start and start+sizeY).
If x is a block vector, it call the present function for all blocks.
- Parameters:
A – a pointer to SiconosMatrix
startRow – an int, sub-block position
x – a pointer to a SiconosVector
y – a pointer to a SiconosVector
init – a bool
-
void private_addprod(unsigned int startRow, unsigned int startCol, const SiconosVector &x, SiconosVector &res) const#
computes res = subA*x +res, subA being a submatrix of A (rows from startRow to startRow+sizeY and columns between startCol and startCol+sizeX).
If x is a block vector, it call the present function for all blocks.
- Parameters:
A – a pointer to SiconosMatrix
startRow – an int, sub-block position
startCol – an int, sub-block position
x – a pointer to a SiconosVector
res – a DenseVect
Protected Attributes
-
siconos::UBLAS_TYPE _num#
A number to specify the type of the matrix: (block or ublas-type) 0-> BlockMatrix, 1 -> DenseMat, 2 -> TriangMat, 3 -> SymMat, 4->SparseMat, 5->BandedMat, 6->zeroMat, 7->IdentityMat.
-
bool _isSymmetric#
bool _isSymmetric; Boolean = true if the Matrix is symmetric
-
bool _isPositiveDefinite#
bool _isPositiveDefinite; Boolean = true if the Matrix is positive definite
Friends
-
friend void prod(const SiconosMatrix &A, const SiconosVector &x, BlockVector &y, bool init)
prod(A, x, y, init) computes y = A*x or y += A*x if init = false
- Parameters:
A – a SiconosMatrix
x – a SiconosVector
y – [inout] a SiconosVector
init – a bool (default = true)
-
friend void prod(const SiconosMatrix &A, const BlockVector &x, SiconosVector &y, bool init)#
-
friend std::ostream &operator<<(std::ostream &os, const SiconosMatrix &sm)
send data of the matrix to an ostream
- Parameters:
os – An output stream
sm – a SiconosMatrix
- Returns:
The same output stream
-
friend SiconosMatrix &operator*=(SiconosMatrix &m, const double &s)
multiply the current matrix with a scalar
- Parameters:
m – the matrix to operate on
s – the scalar
- Returns:
-
friend SiconosMatrix &operator/=(SiconosMatrix &m, const double &s)
divide the current matrix with a scalar
- Parameters:
m – the matrix to operate on
s – the scalar
- Returns:
-
inline virtual ~SiconosMatrix()