Class BlockMatrix#
Defined in Program listing for file kernel/src/utils/SiconosAlgebra/BlockMatrix.hpp
-
class BlockMatrix : public SiconosMatrix#
“Block” matrix, ie container of matrices
A BlockMatrix is a boost::ublas::compressed_matrix of SP::SiconosMatrix.
The blocks positions are given by two Index objects, tabRow and tabCol.
If block 1 is n1xm1, block2 n2xm2, block3 n3xm3 …, then:
tabRow = [ n1 n1+n2 n1+n2+n3 …]
tabCol = [ m1 m1+m2 m1+m2+m3 …]Public Functions
-
BlockMatrix(const SiconosMatrix &m)#
copy constructor
- Parameters:
m – a SiconosMatrix
-
BlockMatrix(const BlockMatrix &m)#
copy constructor
- Parameters:
m – a BlockMatrix
-
BlockMatrix(const std::vector<SP::SiconosMatrix> &m, unsigned int row, unsigned int col)#
constructor with a list of pointer to SiconosMatrix (!links with pointer, no copy!)
- Parameters:
m – a vector of SiconosMatrix
row – number of blocks in a row
col – number of col in a row
-
BlockMatrix(SP::SiconosMatrix A, SP::SiconosMatrix B, SP::SiconosMatrix C, SP::SiconosMatrix D)#
contructor with a list of 4 pointer to SiconosMatrix (!links with pointer, no copy!)
- Parameters:
A – block (0,0)
B – block (0,1)
C – block (1,0)
D – block (1,1)
-
~BlockMatrix(void) noexcept#
destructor
-
inline virtual bool checkSymmetry(double tol) const override#
determines if the matrix is symmetric up to a given tolerance
- Returns:
true if the matrix is inversed
-
virtual unsigned int numberOfBlocks(unsigned int i) const override#
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
-
virtual const DenseMat getDense(unsigned int row = 0, unsigned int col = 0) const override#
get DenseMat matrix
- Parameters:
row – an unsigned int, position of the block (row)
col – an unsigned int, position of the block (column)
- Returns:
a DenseMat
-
virtual const TriangMat getTriang(unsigned int row = 0, unsigned int col = 0) const override#
get TriangMat matrix
- Parameters:
row – an unsigned int, position of the block (row)
col – an unsigned int, position of the block (column)
- Returns:
a TriangMat
-
virtual const SymMat getSym(unsigned int row = 0, unsigned int col = 0) const override#
get SymMat matrix
- Parameters:
row – an unsigned int, position of the block (row)
col – an unsigned int, position of the block (column)
- Returns:
a SymMat
-
virtual const BandedMat getBanded(unsigned int row = 0, unsigned int col = 0) const override#
get BandedMat matrix
- Parameters:
row – an unsigned int, position of the block (row)
col – an unsigned int, position of the block (column)
- Returns:
a BandedMat
-
virtual const SparseMat getSparse(unsigned int row = 0, unsigned int col = 0) const override#
get SparseMat matrix
- Parameters:
row – an unsigned int, position of the block (row)
col – an unsigned int, position of the block (column)
- Returns:
a SparseMat
-
virtual const SparseCoordinateMat getSparseCoordinate(unsigned int row = 0, unsigned int col = 0) const override#
get SparseCoordinateMat matrix
- Parameters:
row – an unsigned int, position of the block (row)
col – an unsigned int, position of the block (column)
- Returns:
a SparseCoordinateMat
-
virtual const ZeroMat getZero(unsigned int row = 0, unsigned int col = 0) const override#
get ZeroMat matrix
- Parameters:
row – an unsigned int, position of the block (row)
col – an unsigned int, position of the block (column)
- Returns:
a ZeroMat
-
virtual const IdentityMat getIdentity(unsigned int row = 0, unsigned int col = 0) const override#
get getIdentity matrix
- Parameters:
row – an unsigned int, position of the block (row)
col – an unsigned int, position of the block (column)
- Returns:
an IdentityMat
-
virtual DenseMat *dense(unsigned int row = 0, unsigned int col = 0) const override#
get a pointer on DenseMat matrix
- Parameters:
row – an unsigned int, position of the block (row)
col – an unsigned int, position of the block (column)
- Returns:
a DenseMat*
-
virtual TriangMat *triang(unsigned int row = 0, unsigned int col = 0) const override#
get a pointer on TriangMat matrix
- Parameters:
row – an unsigned int, position of the block (row)
col – an unsigned int, position of the block (column)
- Returns:
a TriangMat*
-
virtual SymMat *sym(unsigned int row = 0, unsigned int col = 0) const override#
get a pointer on SymMat matrix
- Parameters:
row – an unsigned int, position of the block (row)
col – `an unsigned int, position of the block (column)
- Returns:
a SymMat*
-
virtual BandedMat *banded(unsigned int row = 0, unsigned int col = 0) const override#
get a pointer on BandedMat matrix
- Parameters:
row – an unsigned int, position of the block (row)
col – an unsigned int, position of the block (column)
- Returns:
a BandedMat*
-
virtual SparseMat *sparse(unsigned int row = 0, unsigned int col = 0) const override#
get a pointer on SparseMat matrix
- Parameters:
row – an unsigned int, position of the block (row)
col – an unsigned int, position of the block (column)
- Returns:
a SparseMat*
-
virtual SparseCoordinateMat *sparseCoordinate(unsigned int row = 0, unsigned int col = 0) const override#
get a pointer on SparseCoordinateMat matrix
- Parameters:
row – an unsigned int, position of the block (row)
col – an unsigned int, position of the block (column)
- Returns:
a SparseCoordinateMat*
-
virtual ZeroMat *zero_mat(unsigned int row = 0, unsigned int col = 0) const override#
get a pointer on ZeroMat matrix
- Parameters:
row – an unsigned int, position of the block (row)
col – an unsigned int, position of the block (column)
- Returns:
a ZeroMat*
-
virtual IdentityMat *identity(unsigned int row = 0, unsigned int col = 0) const override#
get a pointer on Identity matrix
- Parameters:
row – an unsigned int, position of the block (row)
col – an unsigned int, position of the block (column)
- Returns:
an IdentityMat*
-
virtual double *getArray(unsigned int row = 0, unsigned int col = 0) const override#
return the address of the array of double values of the matrix
- Parameters:
row – position for the required block ->useless for SimpleMatrix
col – position for the required block ->useless for SimpleMatrix
- Returns:
double* : the pointer on the double array
-
virtual void zero() override#
sets all the values of the matrix to 0.0
-
virtual void randomize() override#
Initialize the matrix with random values.
-
virtual void eye() override#
set an identity matrix
-
virtual unsigned int size(unsigned int index) const override#
get the number of rows or columns of the matrix
- Parameters:
index – 0 for rows, 1 for columns
- Returns:
an int
-
virtual void resize(unsigned int nbrow, unsigned int nbcol, unsigned int lower = 0, unsigned int upper = 0, bool b = true) override#
resize the matrix with nbrow rows and nbcol columns, lower and upper are useful only for SparseMat.The existing elements of the Block matrix are preseved when specified.
- Parameters:
nbrow –
nbcol –
lower –
upper –
b –
-
virtual double normInf() const override#
compute the infinite norm of the Block matrix
- Returns:
a double
-
virtual void display() const override#
display data on standard output
-
virtual void displayExpert(bool brief = true) const override#
display data on standard output
-
virtual std::string toString() const override#
put data of the matrix into a std::string
- Returns:
std::string
-
virtual double &operator()(unsigned int i, unsigned int j) override#
get or set the element matrix[i,j]
- Parameters:
i – an unsigned int
j – an unsigned int
- Returns:
the element matrix[i,j]
-
virtual double operator()(unsigned int i, unsigned int j) const override#
get or set the element matrix[i,j]
- Parameters:
i – an unsigned int
j – an unsigned int
- Returns:
the element matrix[i,j]
-
virtual double getValue(unsigned int i, unsigned int j) const override#
return the element matrix[i,j]
- Parameters:
i – an unsigned int
j – an unsigned int
- Returns:
a double
-
virtual void setValue(unsigned int i, unsigned int j, double value) override#
set the element matrix[i,j]
- Parameters:
i – an unsigned int i
j – an unsigned int j
value –
-
virtual void trans(const SiconosMatrix &m) override#
transpose a matrix: x->trans(m) is x = transpose of m.
- Parameters:
m – the matrix to be transposed.
-
inline virtual const SP::Index tabRow() const override#
get the vector tabRow
- Returns:
a pointer to vector of int
-
inline virtual const SP::Index tabCol() const override#
get the vector tabCol
- Returns:
a pointer to vector of int
-
virtual SP::SiconosMatrix block(unsigned int row = 0, unsigned int col = 0) override#
get block at position row-col
- Parameters:
row – unsigned int
col – unsigned int
- Returns:
SP::SiconosMatrix the requested block
-
virtual SPC::SiconosMatrix block(unsigned int row = 0, unsigned int col = 0) const override#
get block at position row-col
- Parameters:
row – unsigned int
col – unsigned int
- Returns:
SP::SiconosMatrix the requested block
-
virtual void getRow(unsigned int r, SiconosVector &v) const override#
get row index of current matrix and save it in v
- Parameters:
r – index of required line
v – [out] a vector
-
virtual void setRow(unsigned int r, const SiconosVector &v) override#
set line row of the current matrix with vector v
- Parameters:
r – index of required line
v – a vector
-
virtual void getCol(unsigned int c, SiconosVector &v) const override#
get column index of current matrix and save it into vOut
- Parameters:
c – index of required column
v – [out] a vector
-
virtual void setCol(unsigned int c, const SiconosVector &v) override#
set column col of the current matrix with vector
- Parameters:
c – index of required column
v – a vector
-
void addSimple(unsigned int &i, unsigned int &j, const SiconosMatrix &m)#
add a part of the input matrix (starting from (i,j) pos) to the current matrix
- Parameters:
i – an unsigned int i (in-out)
j – an unsigned int j (in-out)
m – a SiconosMatrix (in-out)
-
void subSimple(unsigned int &i, unsigned int &j, const SiconosMatrix &m)#
subtract a part of the input matrix (starting from (i,j) pos) to the current matrix
- Parameters:
i – an unsigned int i (in-out)
j – an unsigned int j (in-out)
m – a SiconosMatrix (in-out)
-
virtual BlockMatrix &operator=(const SiconosMatrix &m) override#
assignment
- Parameters:
m – the matrix to be copied
- Returns:
-
BlockMatrix &operator=(const BlockMatrix &m)#
assignment
- Parameters:
m – the matrix to be copied
- Returns:
-
virtual BlockMatrix &operator=(const DenseMat &m) override#
assignment
- Parameters:
m – the matrix to be copied
- Returns:
-
virtual BlockMatrix &operator+=(const SiconosMatrix &m) override#
operator +=
- Parameters:
m – the matrix to add
- Returns:
-
virtual BlockMatrix &operator-=(const SiconosMatrix &m) override#
operator -=
- Parameters:
m – the matrix to subtract
- Returns:
-
virtual void PLUFactorizationInPlace() override#
computes an LU factorization of a general M-by-N matrix using partial pivoting with row interchanges.
The result is returned in this (InPlace). Based on Blas dgetrf function.
-
virtual void Factorize() override#
computes a factorization of a general M-by-N matrix The implementation is based on an internal NumericsMatrix
-
virtual void PLUInverseInPlace() override#
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.
-
virtual void PLUForwardBackwardInPlace(SiconosMatrix &B) override#
solves a system of linear equations A * X = B (A=this) with a general N-by-N matrix A using the LU factorization computed by PLUFactorizationInPlace.
Based on Blas dgetrs function.
- Parameters:
B – [inout] on input the RHS matrix b; on output: the result x
-
virtual void Solve(SiconosMatrix &B) override#
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) override#
solves a system of linear equations A * X = B (A=this) with a general N-by-N matrix A using the LU factorization computed by PLUFactorizationInPlace.
Based on Blas dgetrs function.
- Parameters:
B – [inout] on input the RHS matrix b; on output: the result x
-
virtual void Solve(SiconosVector &B) override#
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 size_t nnz(double tol = 1.e-14) override#
number of non-zero in the matrix
- Parameters:
tol – the tolerance under which a number is considered zero
Friends
-
friend std::ostream &operator<<(std::ostream &os, const BlockMatrix &bm)#
send data of the matrix to an ostream
- Parameters:
os – An output stream
bm – a BlockMatrix
- Returns:
The same output stream
-
friend void scal(double, const SiconosMatrix&, SiconosMatrix&, bool)#
multiplication of a matrix by a scalar, B = a*A (init = true) or B += a*A (init = false)
- Parameters:
a – a double
A – a SiconosMatrix
B – [inout] a SiconosMatrix
init – a bool
-
BlockMatrix(const SiconosMatrix &m)#