1.3.2
CCC
 
MatrixSparse Class Reference

#include <MatrixSparse.hpp>

Inheritance diagram for MatrixSparse:
AMatrix AStringable ICloneable ProjMatrix

Detailed Description

Sparse Matrix

Handle a sparse matrix that can be symmetrical, square or not. Storage relies either on Eigen3 Library (see opt_eigen flag) or cs code. Default storage option can be set globally by using setGlobalFlagEigen

Public Member Functions

 MatrixSparse (int nrow=0, int ncol=0, int opt_eigen=-1)
 
 MatrixSparse (const cs *A)
 
 MatrixSparse (const MatrixSparse &m)
 
MatrixSparseoperator= (const MatrixSparse &m)
 
virtual ~MatrixSparse ()
 
bool isFlagEigen () const
 Cloneable interface. More...
 
bool isSparse () const override
 Interface for AMatrix. More...
 
bool isDense () const override
 
void setValue (int irow, int icol, double value, bool flagCheck=true) override
 
double getValue (int row, int col, bool flagCheck=true) const override
 
void updValue (int irow, int icol, const EOperator &oper, double value, bool flagCheck=true) override
 
virtual void setColumn (int icol, const VectorDouble &tab) override
 
virtual void setRow (int irow, const VectorDouble &tab) override
 
virtual void setDiagonal (const VectorDouble &tab) override
 
virtual void setDiagonalToConstant (double value=1.) override
 
virtual MatrixSparsetranspose () const override
 
virtual void addScalar (double v) override
 
virtual void addScalarDiag (double v) override
 
virtual void prodScalar (double v) override
 
virtual void fill (double value) override
 
virtual void multiplyRow (const VectorDouble &vec) override
 
virtual void multiplyColumn (const VectorDouble &vec) override
 
virtual void divideRow (const VectorDouble &vec) override
 
virtual void divideColumn (const VectorDouble &vec) override
 
virtual VectorDouble prodVecMat (const VectorDouble &x, bool transpose=false) const override
 
virtual VectorDouble prodMatVec (const VectorDouble &x, bool transpose=false) const override
 
virtual void prodMatMatInPlace (const AMatrix *x, const AMatrix *y, bool transposeX=false, bool transposeY=false) override
 
virtual NF_Triplet getMatrixToTriplet (int shiftRow=0, int shiftCol=0) const override
 
virtual String toString (const AStringFormat *strfmt=nullptr) const override
 Interface to AStringable. More...
 
virtual void addMatInPlace (const MatrixSparse &y, double cx=1., double cy=1.)
 The next functions use specific definition of matrix (to avoid dynamic_cast) rather than manipulating AMatrix. They are no more generic of AMatrix. More...
 
virtual void prodNormMatMatInPlace (const MatrixSparse &a, const MatrixSparse &m, bool transpose=false)
 
virtual void prodNormMatInPlace (const MatrixSparse &a, const VectorDouble &vec=VectorDouble(), bool transpose=false)
 
const cs * getCS () const
 
void setCS (cs *cs)
 
void freeCS ()
 
cs * getCSUnprotected () const
 
virtual void resetFromValue (int nrows, int ncols, double value) override
 Reset the matrix to new dimensions and fill with a new value. More...
 
virtual void resetFromArray (int nrows, int ncols, const double *tab, bool byCol=true) override
 Reset the matrix from an array of double values. More...
 
virtual void resetFromVD (int nrows, int ncols, const VectorDouble &tab, bool byCol=true) override
 Reset the matrix from a vector of double values. More...
 
virtual void resetFromVVD (const VectorVectorDouble &tab, bool byCol=true) override
 Reset the matrix from an array of double values. More...
 
void resetFromTriplet (const NF_Triplet &NF_T)
 
void dumpElements (const String &title, int ifrom, int ito) const
 
void fillRandom (int seed=432432, double zeroPercent=0.1)
 
int computeCholesky ()
 
int solveCholesky (const VectorDouble &b, VectorDouble &x)
 
int simulateCholesky (const VectorDouble &b, VectorDouble &x)
 
double computeCholeskyLogDeterminant ()
 
void addValue (int row, int col, double value)
 
double L1Norm () const
 
void getStats (int *nrows, int *ncols, int *count, double *percent) const
 
int scaleByDiag ()
 
int addVecInPlace (const VectorDouble &x, VectorDouble &y)
 
void setConstant (double value)
 
VectorDouble extractDiag (int oper_choice=1) const
 
void prodNormDiagVecInPlace (const VectorDouble &vec, int oper=1)
 
const Eigen::SparseMatrix< double > & getEigenMatrix () const
 
void setEigenMatrix (const Eigen::SparseMatrix< double > &eigenMatrix)
 
MatrixSparseextractSubmatrixByRanks (const VectorInt &rank_rows, const VectorInt &rank_cols)
 
MatrixSparseextractSubmatrixByColor (const VectorInt &colors, int ref_color, bool row_ok, bool col_ok)
 
VectorInt colorCoding ()
 
int getNonZeros () const
 
- Public Member Functions inherited from AMatrix
 AMatrix (int nrow=0, int ncol=0)
 
 AMatrix (const AMatrix &m)
 
AMatrixoperator= (const AMatrix &m)
 
virtual ~AMatrix ()
 
virtual void reset (int nrows, int ncols)
 
virtual bool isSquare (bool printWhyNot=false) const
 
virtual bool isValid (int irow, int icol, bool printWhyNot=false) const
 
virtual bool isIdentity (bool printWhyNot=false) const
 
virtual bool isSymmetric (bool printWhyNot=false, double eps=EPSILON10) const
 
virtual bool mustBeSymmetric () const
 
virtual void transposeInPlace ()
 
virtual VectorDouble getRow (int irow) const
 
virtual VectorDouble getColumn (int icol) const
 
void addMatInPlace (const AMatrix &y, double cx=1., double cy=1.)
 
void prodMatInPlace (const AMatrix *matY, bool transposeY=false)
 
void prodNormMatMatInPlace (const AMatrix &a, const AMatrix &m, bool transpose=false)
 
void prodNormMatInPlace (const AMatrix &a, const VectorDouble &vec=VectorDouble(), bool transpose=false)
 
void resize (int nrows, int ncols)
 Resize the matrix to new dimensions (this method doesn't change the storage type) More...
 
void addValue (int irow, int icol, double value)
 
bool isSame (const AMatrix &m, double eps=EPSILON4, bool printWhyNot=false)
 
bool isSameSize (const AMatrix &m) const
 
bool empty () const
 
double compare (const AMatrix &mat) const
 
int getNRows () const
 
int getNCols () const
 
int size () const
 
VectorDouble getValues (bool byCol=true) const
 
VectorDouble getDiagonal (int shift=0) const
 
bool isColumnDefined (int icol) const
 
bool isRowDefined (int irow) const
 
int getNumberColumnDefined () const
 
int getNumberRowDefined () const
 
bool isNonNegative (bool verbose=false)
 
void prodMatVecInPlace (const VectorDouble &x, VectorDouble &y, bool transpose=false) const
 
void prodMatVecInPlacePtr (const double *x, double *y, bool transpose=false) const
 
void prodVecMatInPlace (const VectorDouble &x, VectorDouble &y, bool transpose=false) const
 
void prodVecMatInPlacePtr (const double *x, double *y, bool transpose=false) const
 
double quadraticMatrix (const VectorDouble &x, const VectorDouble &y)
 
int invert ()
 
int solve (const VectorDouble &b, VectorDouble &x) const
 
void dumpElements (const String &title, int ifrom, int ito) const
 
void setIdentity (double value=1.)
 
void fillRandom (int seed=432432, double zeroPercent=0.1)
 
void setValues (const VectorDouble &values, bool byCol=true)
 
double getMeanByColumn (int icol) const
 
double getMinimum () const
 
double getMaximum () const
 
double getNormInf () const
 
void copyReduce (const AMatrix *x, const VectorInt &activeRows, const VectorInt &activeCols)
 
void copyElements (const AMatrix &m, double factor=1.)
 
void setFlagCheckAddress (bool flagCheckAddress)
 
void makePositiveColumn ()
 
void linearCombination (double val1, const AMatrix *mat1, double val2=1., const AMatrix *mat2=nullptr)
 
double operator() (int row, int col) const
 
double & operator() (int row, int col)
 
- Public Member Functions inherited from AStringable
 AStringable ()
 
 AStringable (const AStringable &r)
 
AStringableoperator= (const AStringable &r)
 
virtual ~AStringable ()
 
virtual void display (const AStringFormat *strfmt=nullptr) const final
 
virtual void display (int level) const final
 
- Public Member Functions inherited from ICloneable
 ICloneable ()
 
virtual ~ICloneable ()
 
virtual ICloneableclone () const =0
 

Static Public Member Functions

static MatrixSparsecreateFromTriplet (const NF_Triplet &NF_T, int nrow=0, int ncol=0, int opt_eigen=-1)
 
static MatrixSparseaddMatMat (const MatrixSparse *x, const MatrixSparse *y, double cx=1., double cy=1.)
 
static MatrixSparsediagVec (const VectorDouble &vec, int opt_eigen=-1)
 
static MatrixSparsediagConstant (int number, double value=1., int opt_eigen=-1)
 
static MatrixSparsediagMat (MatrixSparse *A, int oper_choice, int opt_eigen=-1)
 
static MatrixSparseglue (const MatrixSparse *A1, const MatrixSparse *A2, bool flagShiftRow, bool flagShiftCol)
 

Public Attributes

 DECLARE_TOTL
 Has a specific implementation in the Target language. More...
 

Constructor & Destructor Documentation

◆ MatrixSparse() [1/3]

MatrixSparse::MatrixSparse ( int  nrow = 0,
int  ncol = 0,
int  opt_eigen = -1 
)

◆ MatrixSparse() [2/3]

MatrixSparse::MatrixSparse ( const cs *  A)

◆ MatrixSparse() [3/3]

MatrixSparse::MatrixSparse ( const MatrixSparse m)

◆ ~MatrixSparse()

MatrixSparse::~MatrixSparse ( )
virtual

Member Function Documentation

◆ addMatInPlace()

void MatrixSparse::addMatInPlace ( const MatrixSparse y,
double  cx = 1.,
double  cy = 1. 
)
virtual

The next functions use specific definition of matrix (to avoid dynamic_cast) rather than manipulating AMatrix. They are no more generic of AMatrix.

Add a matrix (multiplied by a constant)

Updates the current Matrix as a linear combination of matrices as follows: this <- cx * this + cy * y

Parameters
cxCoefficient applied to the current Matrix
cyCoefficient applied to the Matrix 'y'
ySecond Matrix in the Linear combination

◆ addMatMat()

MatrixSparse * MatrixSparse::addMatMat ( const MatrixSparse x,
const MatrixSparse y,
double  cx = 1.,
double  cy = 1. 
)
static

◆ addScalar()

void MatrixSparse::addScalar ( double  v)
overridevirtual

Add a value to each matrix component

Parameters
vAdd a scalar value to all terms of the current matrix

Reimplemented from AMatrix.

◆ addScalarDiag()

void MatrixSparse::addScalarDiag ( double  v)
overridevirtual

Add value to matrix diagonal

Parameters
vAdd constant value to the diagonal of the current Matrix

Reimplemented from AMatrix.

◆ addValue()

void MatrixSparse::addValue ( int  row,
int  col,
double  value 
)

◆ addVecInPlace()

int MatrixSparse::addVecInPlace ( const VectorDouble x,
VectorDouble y 
)

◆ colorCoding()

VectorInt MatrixSparse::colorCoding ( )

◆ computeCholesky()

int MatrixSparse::computeCholesky ( )

◆ computeCholeskyLogDeterminant()

double MatrixSparse::computeCholeskyLogDeterminant ( )

◆ createFromTriplet()

MatrixSparse * MatrixSparse::createFromTriplet ( const NF_Triplet NF_T,
int  nrow = 0,
int  ncol = 0,
int  opt_eigen = -1 
)
static

◆ diagConstant()

MatrixSparse * MatrixSparse::diagConstant ( int  number,
double  value = 1.,
int  opt_eigen = -1 
)
static

◆ diagMat()

MatrixSparse * MatrixSparse::diagMat ( MatrixSparse A,
int  oper_choice,
int  opt_eigen = -1 
)
static

Construct a sparse matrix with the diagonal of 'A', where each element is transformed

Parameters
AInput sparse matrix
oper_choiceOperation on the diagonal term (see Utilities::operate_XXX)
opt_eigenOption for choosing Eigen Library or not
Returns

◆ diagVec()

MatrixSparse * MatrixSparse::diagVec ( const VectorDouble vec,
int  opt_eigen = -1 
)
static

◆ divideColumn()

void MatrixSparse::divideColumn ( const VectorDouble vec)
overridevirtual

Divide the matrix column-wise

Divide a Matrix column-wise

Reimplemented from AMatrix.

◆ divideRow()

void MatrixSparse::divideRow ( const VectorDouble vec)
overridevirtual

Divide the matrix row-wise

Divide a Matrix row-wise

Reimplemented from AMatrix.

◆ dumpElements()

void MatrixSparse::dumpElements ( const String title,
int  ifrom,
int  ito 
) const

Dump a specific range of samples from the internal storage

◆ extractDiag()

VectorDouble MatrixSparse::extractDiag ( int  oper_choice = 1) const

◆ extractSubmatrixByColor()

MatrixSparse * MatrixSparse::extractSubmatrixByColor ( const VectorInt colors,
int  ref_color,
bool  row_ok,
bool  col_ok 
)

◆ extractSubmatrixByRanks()

MatrixSparse * MatrixSparse::extractSubmatrixByRanks ( const VectorInt rank_rows,
const VectorInt rank_cols 
)

◆ fill()

void MatrixSparse::fill ( double  value)
overridevirtual

Set all the values of the matrix at once

Change any nonzero term to 'value'

Parameters
valueConstant value used for filling 'this'

Reimplemented from AMatrix.

◆ fillRandom()

void MatrixSparse::fillRandom ( int  seed = 432432,
double  zeroPercent = 0.1 
)

Set all the values of the Matrix with random values

◆ freeCS()

void MatrixSparse::freeCS ( )

◆ getCS()

const cs * MatrixSparse::getCS ( ) const

Returns a pointer to the Sparse storage

◆ getCSUnprotected()

cs * MatrixSparse::getCSUnprotected ( ) const

Temporary function to get the CS contents of Sparse Matrix

◆ getEigenMatrix()

const Eigen::SparseMatrix<double>& MatrixSparse::getEigenMatrix ( ) const
inline

◆ getMatrixToTriplet()

NF_Triplet MatrixSparse::getMatrixToTriplet ( int  shiftRow = 0,
int  shiftCol = 0 
) const
overridevirtual

Extract the contents of the matrix

From a matrix of any type, creates the triplet (specific format for creating efficiently a Sparse matrix) It only takes the only non-zero elements of the matrix

Reimplemented from AMatrix.

◆ getNonZeros()

int MatrixSparse::getNonZeros ( ) const
inline

◆ getStats()

void MatrixSparse::getStats ( int *  nrows,
int *  ncols,
int *  count,
double *  percent 
) const

◆ getValue()

double MatrixSparse::getValue ( int  row,
int  col,
bool  flagCheck = true 
) const
overridevirtual

Get the value from a matrix cell

Implements AMatrix.

◆ glue()

MatrixSparse * MatrixSparse::glue ( const MatrixSparse A1,
const MatrixSparse A2,
bool  flagShiftRow,
bool  flagShiftCol 
)
static

◆ isDense()

bool MatrixSparse::isDense ( ) const
inlineoverridevirtual

Returns if the matrix belongs to the MatrixSparse class (avoids dynamic_cast)

Implements AMatrix.

◆ isFlagEigen()

bool MatrixSparse::isFlagEigen ( ) const
inline

Cloneable interface.

◆ isSparse()

bool MatrixSparse::isSparse ( ) const
inlineoverridevirtual

Interface for AMatrix.

Returns if the current matrix is Sparse

Implements AMatrix.

◆ L1Norm()

double MatrixSparse::L1Norm ( ) const

◆ multiplyColumn()

void MatrixSparse::multiplyColumn ( const VectorDouble vec)
overridevirtual

Multiply the matrix column-wise

Multiply a Matrix column-wise

Reimplemented from AMatrix.

◆ multiplyRow()

void MatrixSparse::multiplyRow ( const VectorDouble vec)
overridevirtual

Multiply the matrix row-wise

Multiply a Matrix row-wise

Reimplemented from AMatrix.

◆ operator=()

MatrixSparse & MatrixSparse::operator= ( const MatrixSparse m)

◆ prodMatMatInPlace()

void MatrixSparse::prodMatMatInPlace ( const AMatrix x,
const AMatrix y,
bool  transposeX = false,
bool  transposeY = false 
)
overridevirtual

Multiply a matrix by another and stored in 'this'

Store the product of 'x' by 'y' in this

Parameters
xFirst Matrix
ySecond matrix
transposeXTrue if First matrix is transposed
transposeYTrue if Second matrix is transposed

Reimplemented from AMatrix.

◆ prodMatVec()

VectorDouble MatrixSparse::prodMatVec ( const VectorDouble x,
bool  transpose = false 
) const
overridevirtual

Perform y = 'this' %*% x

Reimplemented from AMatrix.

◆ prodNormDiagVecInPlace()

void MatrixSparse::prodNormDiagVecInPlace ( const VectorDouble vec,
int  oper_choice = 1 
)

Perform: 'this' = diag('vec') %*% 'A' %*% diag('vec')

Parameters
vecInput Vector
oper_choiceType of transformation

◆ prodNormMatInPlace()

void MatrixSparse::prodNormMatInPlace ( const MatrixSparse a,
const VectorDouble vec = VectorDouble(),
bool  transpose = false 
)
virtual

Product 't(A)' %*% ['vec'] %*% 'A' or 'A' %*% ['vec'] %*% 't(A)' stored in 'this'

◆ prodNormMatMatInPlace()

void MatrixSparse::prodNormMatMatInPlace ( const MatrixSparse a,
const MatrixSparse m,
bool  transpose = false 
)
virtual

Product 't(A)' %*% 'M' %*% 'A' or 'A' %*% 'M' %*% 't(A)' stored in 'this'

◆ prodScalar()

void MatrixSparse::prodScalar ( double  v)
overridevirtual

Multiply each matrix component by a value

Parameters
vMultiply all the terms of the matrix by the scalar 'v'

Reimplemented from AMatrix.

◆ prodVecMat()

VectorDouble MatrixSparse::prodVecMat ( const VectorDouble x,
bool  transpose = false 
) const
overridevirtual

Perform y = x %*% 'this'

Reimplemented from AMatrix.

◆ resetFromArray()

void MatrixSparse::resetFromArray ( int  nrows,
int  ncols,
const double *  tab,
bool  byCol = true 
)
overridevirtual

Reset the matrix from an array of double values.

Parameters
nrowsNew number of rows
ncolsNew number of columns
tabThe array of values
byColTrue if values are column-major in the array

Reimplemented from AMatrix.

◆ resetFromTriplet()

void MatrixSparse::resetFromTriplet ( const NF_Triplet NF_T)

◆ resetFromValue()

void MatrixSparse::resetFromValue ( int  nrows,
int  ncols,
double  value 
)
overridevirtual

Reset the matrix to new dimensions and fill with a new value.

Parameters
nrowsNew number of rows
ncolsNew number of columns
valueThe new value used to fill the matrix

Reimplemented from AMatrix.

◆ resetFromVD()

void MatrixSparse::resetFromVD ( int  nrows,
int  ncols,
const VectorDouble tab,
bool  byCol = true 
)
overridevirtual

Reset the matrix from a vector of double values.

Parameters
nrowsNew number of rows
ncolsNew number of columns
tabThe vector of values
byColTrue if values are column-major in the vector

Reimplemented from AMatrix.

◆ resetFromVVD()

void MatrixSparse::resetFromVVD ( const VectorVectorDouble tab,
bool  byCol = true 
)
overridevirtual

Reset the matrix from an array of double values.

Parameters
tabThe array of values
byColTrue if values are column-major in the array

Reimplemented from AMatrix.

◆ scaleByDiag()

int MatrixSparse::scaleByDiag ( )

◆ setColumn()

void MatrixSparse::setColumn ( int  icol,
const VectorDouble tab 
)
overridevirtual

Set the contents of a Column

Fill a column of an already existing Sparse matrix, using 'tab' as entry The input 'tab' corresponds to the whole column contents

Parameters
icolColumn rank
tabVector containing the information (Dimension: nrows)
Warning
: This method only copies the values at the non-zero existing entries

Reimplemented from AMatrix.

◆ setConstant()

void MatrixSparse::setConstant ( double  value)

◆ setCS()

void MatrixSparse::setCS ( cs *  cs)

◆ setDiagonal()

void MatrixSparse::setDiagonal ( const VectorDouble tab)
overridevirtual

Set the contents of the (main) Diagonal

Reimplemented from AMatrix.

◆ setDiagonalToConstant()

void MatrixSparse::setDiagonalToConstant ( double  value = 1.)
overridevirtual

Set the contents of the (main) Diagonal to a constant value

Reimplemented from AMatrix.

◆ setEigenMatrix()

void MatrixSparse::setEigenMatrix ( const Eigen::SparseMatrix< double > &  eigenMatrix)
inline

◆ setRow()

void MatrixSparse::setRow ( int  irow,
const VectorDouble tab 
)
overridevirtual

Set the contents of a Row

Fill a row of an already existing Sparse matrix, using 'tab' as entry The input 'tab' corresponds to the whole row contents

Parameters
irowRow rank
tabVector containing the information (Dimension: ncols)
Warning
: This method only copies the values at the non-zero existing entries

Reimplemented from AMatrix.

◆ setValue()

void MatrixSparse::setValue ( int  irow,
int  icol,
double  value,
bool  flagCheck = true 
)
overridevirtual

Set the value for a matrix cell

Implements AMatrix.

◆ simulateCholesky()

int MatrixSparse::simulateCholesky ( const VectorDouble b,
VectorDouble x 
)

◆ solveCholesky()

int MatrixSparse::solveCholesky ( const VectorDouble b,
VectorDouble x 
)

◆ toString()

String MatrixSparse::toString ( const AStringFormat strfmt = nullptr) const
overridevirtual

Interface to AStringable.

Reimplemented from AMatrix.

Reimplemented in ProjMatrix.

◆ transpose()

MatrixSparse * MatrixSparse::transpose ( ) const
overridevirtual

Transpose the matrix and return it as a copy

Reimplemented from AMatrix.

◆ updValue()

void MatrixSparse::updValue ( int  irow,
int  icol,
const EOperator &  oper,
double  value,
bool  flagCheck = true 
)
overridevirtual

Modifies the contents of a matrix cell

Implements AMatrix.

Member Data Documentation

◆ DECLARE_TOTL

MatrixSparse::DECLARE_TOTL

Has a specific implementation in the Target language.


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