1.3.2
CCC
 
AMatrix Class Referenceabstract

#include <AMatrix.hpp>

Inheritance diagram for AMatrix:
AStringable ICloneable AMatrixDense MatrixSparse MatrixRectangular ProjMatrix AMatrixSquare Table MatrixSquareGeneral MatrixSquareSymmetric

Detailed Description

This class is the root of the Matrix organization in gstlearn A matrix is a 2-D organization: it is characterized by its number of rows and its number of columns. Although the user should not bother with this remark, the elements of a matrix processed in 'gstlearn' are stored in a Row-major format. This is to say that the internal rank of an element characterized by its row and column numbers is: (icol * getNRows() + irow)

Since gstlearn version v1.3:

  • Dense Matrices storage and algebra rely on Eigen3 library only
  • Sparse Matrices storage and algebra rely on Eigen3 or cs library (see MatrixSparse.hpp)

Public Member Functions

 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 void resetFromValue (int nrows, int ncols, double value)
 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)
 Reset the matrix from an array of double values. More...
 
virtual void resetFromVD (int nrows, int ncols, const VectorDouble &tab, bool byCol=true)
 Reset the matrix from a vector of double values. More...
 
virtual void resetFromVVD (const VectorVectorDouble &tab, bool byCol=true)
 Reset the matrix from an array of double values. More...
 
virtual String toString (const AStringFormat *strfmt=nullptr) const override
 Interface to AStringable. More...
 
virtual bool isDense () const =0
 Interface to AMatrix. More...
 
virtual bool isSparse () const =0
 
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 setColumn (int icol, const VectorDouble &tab)
 
virtual void setRow (int irow, const VectorDouble &tab)
 
virtual void setDiagonal (const VectorDouble &tab)
 
virtual void setDiagonalToConstant (double value=1.)
 
virtual void transposeInPlace ()
 
virtual AMatrixtranspose () const
 
virtual void addScalar (double v)
 
virtual void addScalarDiag (double v)
 
virtual void prodScalar (double v)
 
virtual void fill (double value)
 
virtual void multiplyRow (const VectorDouble &vec)
 
virtual void multiplyColumn (const VectorDouble &vec)
 
virtual void divideRow (const VectorDouble &vec)
 
virtual void divideColumn (const VectorDouble &vec)
 
virtual VectorDouble prodVecMat (const VectorDouble &x, bool transpose=false) const
 
virtual VectorDouble prodMatVec (const VectorDouble &x, bool transpose=false) const
 
virtual VectorDouble getRow (int irow) const
 
virtual VectorDouble getColumn (int icol) const
 
virtual void prodMatMatInPlace (const AMatrix *x, const AMatrix *y, bool transposeX=false, bool transposeY=false)
 
virtual NF_Triplet getMatrixToTriplet (int shiftRow=0, int shiftCol=0) 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...
 
virtual double getValue (int irow, int icol, bool flagCheck=true) const =0
 
virtual void setValue (int irow, int icol, double value, bool flagCheck=true)=0
 
virtual void updValue (int irow, int icol, const EOperator &oper, double value, bool flagCheck=true)=0
 
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
 

Constructor & Destructor Documentation

◆ AMatrix() [1/2]

AMatrix::AMatrix ( int  nrow = 0,
int  ncol = 0 
)

◆ AMatrix() [2/2]

AMatrix::AMatrix ( const AMatrix m)

◆ ~AMatrix()

AMatrix::~AMatrix ( )
virtual

Member Function Documentation

◆ addMatInPlace()

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

Add a matrix (multiplied by a constant)

Add the matrix 'y' to the current Matrix

Parameters
yMatrix to be added
cxMultiplicative parameter for this
cyMultiplicative parameter for y

◆ addScalar()

void AMatrix::addScalar ( double  v)
virtual

Add a value to each matrix component

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

Reimplemented in MatrixSparse, and AMatrixDense.

◆ addScalarDiag()

void AMatrix::addScalarDiag ( double  v)
virtual

Add value to matrix diagonal

Parameters
vAdd constant value to the diagonal of the current Matrix

Reimplemented in MatrixSparse, and AMatrixDense.

◆ addValue()

void AMatrix::addValue ( int  irow,
int  icol,
double  value 
)

Add a value to a matrix term

◆ compare()

double AMatrix::compare ( const AMatrix mat) const

Returns the sum of absolute difference between argument and this

◆ copyElements()

void AMatrix::copyElements ( const AMatrix m,
double  factor = 1. 
)

Copy the contents of matrix 'm' into 'this' Warning: matrices must have the same dimensions (not checked)

Parameters
mInput matrix
factorMultiplicative factor (applied to each element)

◆ copyReduce()

void AMatrix::copyReduce ( const AMatrix x,
const VectorInt activeRows,
const VectorInt activeCols 
)

◆ divideColumn()

void AMatrix::divideColumn ( const VectorDouble vec)
virtual

Divide a Matrix column-wise

Reimplemented in MatrixSparse, and AMatrixDense.

◆ divideRow()

void AMatrix::divideRow ( const VectorDouble vec)
virtual

Divide a Matrix row-wise

Reimplemented in MatrixSparse, and AMatrixDense.

◆ dumpElements()

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

Dump a specific range of samples from the internal storage

◆ empty()

bool AMatrix::empty ( ) const
inline

Returns if the current matrix is Empty

◆ fill()

void AMatrix::fill ( double  value)
virtual

Set all the values of the Matrix at once

Fill 'this' with the constant 'value'

Parameters
valueConstant value used for filling 'this'

Reimplemented in MatrixSparse, and AMatrixDense.

◆ fillRandom()

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

◆ getColumn()

VectorDouble AMatrix::getColumn ( int  icol) const
virtual

Extract a Column

Reimplemented in AMatrixDense.

◆ getDiagonal()

VectorDouble AMatrix::getDiagonal ( int  shift = 0) const

Extract a Diagonal (main or secondary) of this

◆ getMatrixToTriplet()

NF_Triplet AMatrix::getMatrixToTriplet ( int  shiftRow = 0,
int  shiftCol = 0 
) const
virtual

Extract the contents of the matrix

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

Reimplemented in MatrixSparse.

◆ getMaximum()

double AMatrix::getMaximum ( ) const

◆ getMeanByColumn()

double AMatrix::getMeanByColumn ( int  icol) const

◆ getMinimum()

double AMatrix::getMinimum ( ) const

◆ getNCols()

int AMatrix::getNCols ( ) const
inline

Returns the number of columns

◆ getNormInf()

double AMatrix::getNormInf ( ) const

◆ getNRows()

int AMatrix::getNRows ( ) const
inline

Returns the number of rows

◆ getNumberColumnDefined()

int AMatrix::getNumberColumnDefined ( ) const

Define the number of defined columns

◆ getNumberRowDefined()

int AMatrix::getNumberRowDefined ( ) const

Define the number of defined rows

◆ getRow()

VectorDouble AMatrix::getRow ( int  irow) const
virtual

Extract a Row

Reimplemented in AMatrixDense.

◆ getValue()

virtual double AMatrix::getValue ( int  irow,
int  icol,
bool  flagCheck = true 
) const
pure virtual

Gets the value at row 'irow' and column 'icol'

Implemented in MatrixSparse, and AMatrixDense.

◆ getValues()

VectorDouble AMatrix::getValues ( bool  byCol = true) const

Returns the contents of the whole matrix as a VectorDouble

◆ invert()

int AMatrix::invert ( )

Matrix inversion in place

◆ isColumnDefined()

bool AMatrix::isColumnDefined ( int  icol) const

Checks if a Column is valid (contains a non TEST value)

◆ isDense()

virtual bool AMatrix::isDense ( ) const
pure virtual

Interface to AMatrix.

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

Implemented in MatrixSparse, and AMatrixDense.

◆ isIdentity()

bool AMatrix::isIdentity ( bool  printWhyNot = false) const
virtual

Check if the matrix is square and Identity

Indicate if the current matrix is the Identity

Parameters
printWhyNotPrint the message is the answer is false
Returns
true if the current matrix is the Identity

◆ isNonNegative()

bool AMatrix::isNonNegative ( bool  verbose = false)

Check if the matrix does not contain any negative element

Check if all the elements of a matrix are non-negative

Returns
True if the matrix is non-negative; False otherwise
Parameters
[in]verboseTrue for the verbose option

◆ isRowDefined()

bool AMatrix::isRowDefined ( int  irow) const

Checks if a Row is valid (contains a non TEST value)

◆ isSame()

bool AMatrix::isSame ( const AMatrix m,
double  eps = EPSILON4,
bool  printWhyNot = false 
)

Check if a matrix is the same as me (norm L1)

Check that Matrix 'm' is equal to the current Matrix

Parameters
mMatrix to be compared to the current Matrix
epsEpsilon for double equality comparison
printWhyNotPrint the message is the answer is false
Returns
true if 'm' is equal to the current Matrix

◆ isSameSize()

bool AMatrix::isSameSize ( const AMatrix m) const

Check that both matrix have the same number of rows and columns

Check that Matrix 'm' share the same dimensions as current one

Parameters
mMatrix to be compared to the current Matrix
Returns
true if 'm' has same dimensions as the current Matrix

◆ isSparse()

virtual bool AMatrix::isSparse ( ) const
pure virtual

Returns if the current matrix is Sparse

Implemented in MatrixSparse, and AMatrixDense.

◆ isSquare()

bool AMatrix::isSquare ( bool  printWhyNot = false) const
virtual

Check if the matrix is (non empty) square

Indicate if the given matrce is a square matrix

Parameters
printWhyNotPrint the message is the answer is false
Returns
true if the matrix is square

Reimplemented in AMatrixSquare.

◆ isSymmetric()

bool AMatrix::isSymmetric ( bool  printWhyNot = false,
double  eps = EPSILON10 
) const
virtual

Check if the input matrix is (non empty and square) symmetric

Indicate if the current matrix is symmetric

Parameters
printWhyNotPrint the message is the answer is false
epsEpsilon for double equality comparison
Returns
true if the current matrix is symmetric

Reimplemented in MatrixSquareSymmetric.

◆ isValid()

bool AMatrix::isValid ( int  irow,
int  icol,
bool  printWhyNot = false 
) const
virtual

Indicate if the given indices are valid for the current matrix size

Indicate if the given indices are valid for the current matrix size

Parameters
irowRow index
icolColumn index
printWhyNotPrint the message is the answer is false
Returns
true if indices are valid for the current matrix size

◆ linearCombination()

void AMatrix::linearCombination ( double  val1,
const AMatrix mat1,
double  val2 = 1.,
const AMatrix mat2 = nullptr 
)

◆ makePositiveColumn()

void AMatrix::makePositiveColumn ( )

Modify the contents of the matrix so that each column has a positive sum of elements. If this is not the case, simply invert the sign of the column

◆ multiplyColumn()

void AMatrix::multiplyColumn ( const VectorDouble vec)
virtual

Multiply a Matrix column-wise

Reimplemented in MatrixSparse, and AMatrixDense.

◆ multiplyRow()

void AMatrix::multiplyRow ( const VectorDouble vec)
virtual

Multiply a Matrix row-wise

Reimplemented in MatrixSparse, and AMatrixDense.

◆ mustBeSymmetric()

virtual bool AMatrix::mustBeSymmetric ( ) const
inlinevirtual

Say if the matrix must be symmetric

Reimplemented in MatrixSquareGeneral, MatrixRectangular, and MatrixSquareSymmetric.

◆ operator()() [1/2]

double& AMatrix::operator() ( int  row,
int  col 
)
inline

Set value operator override

◆ operator()() [2/2]

double AMatrix::operator() ( int  row,
int  col 
) const
inline

Get value operator override

◆ operator=()

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

◆ prodMatInPlace()

void AMatrix::prodMatInPlace ( const AMatrix matY,
bool  transposeY = false 
)

Multiply 'this' by matrix 'y' and store in 'this'

◆ prodMatMatInPlace()

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

Multiply matrix 'x' by matrix 'y' and store the result in 'this'

Store the product of 'x'(or 't(x)') by 'y' (or 't(y') in this

Parameters
xFirst Matrix
ySecond matrix
transposeXTrue if first matrix must be transposed
transposeYTrue if second matrix must be transposed

Reimplemented in MatrixSparse, and AMatrixDense.

◆ prodMatVec()

VectorDouble AMatrix::prodMatVec ( const VectorDouble x,
bool  transpose = false 
) const
virtual

Perform 'this' * 'vec'

Reimplemented in MatrixSparse, and AMatrixDense.

◆ prodMatVecInPlace()

void AMatrix::prodMatVecInPlace ( const VectorDouble x,
VectorDouble y,
bool  transpose = false 
) const

Perform 'y' = 'this' * 'x'

Returns 'y' = 'this' %*% 'x'

Parameters
xInput vector
yOutput vector obtained by multiplying 'inv' by current Matrix
transposeTrue if the matrix 'this' must be transposed

◆ prodMatVecInPlacePtr()

void AMatrix::prodMatVecInPlacePtr ( const double *  x,
double *  y,
bool  transpose = false 
) const

◆ prodNormMatInPlace()

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

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

◆ prodNormMatMatInPlace()

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

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

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

Parameters
aMatrix A
mMatrix M
transposeTrue for first implementation, False for the second

◆ prodScalar()

void AMatrix::prodScalar ( double  v)
virtual

Multiply each matrix component by a value

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

Reimplemented in MatrixSparse, and AMatrixDense.

◆ prodVecMat()

VectorDouble AMatrix::prodVecMat ( const VectorDouble x,
bool  transpose = false 
) const
virtual

Perform 'vec' * 'this'

Reimplemented in MatrixSparse, and AMatrixDense.

◆ prodVecMatInPlace()

void AMatrix::prodVecMatInPlace ( const VectorDouble x,
VectorDouble y,
bool  transpose = false 
) const

Perform 'y' = 'x' * 'this'

◆ prodVecMatInPlacePtr()

void AMatrix::prodVecMatInPlacePtr ( const double *  x,
double *  y,
bool  transpose = false 
) const

◆ quadraticMatrix()

double AMatrix::quadraticMatrix ( const VectorDouble x,
const VectorDouble y 
)

Perform x %*% mat %*% y

◆ reset()

void AMatrix::reset ( int  nrows,
int  ncols 
)
virtual

Reimplemented in Table.

◆ resetFromArray()

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

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 in MatrixSparse.

◆ resetFromValue()

void AMatrix::resetFromValue ( int  nrows,
int  ncols,
double  value 
)
virtual

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 in MatrixSparse.

◆ resetFromVD()

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

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 in MatrixSparse.

◆ resetFromVVD()

void AMatrix::resetFromVVD ( const VectorVectorDouble tab,
bool  byCol = true 
)
virtual

Reset the matrix from an array of double values.

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

Reimplemented in MatrixSparse.

◆ resize()

void AMatrix::resize ( int  nrows,
int  ncols 
)

Resize the matrix to new dimensions (this method doesn't change the storage type)

Modify the dimension of the matrix (if needed)

Parameters
nrowsNew number of rows
ncolsNew number of columns

◆ setColumn()

void AMatrix::setColumn ( int  icol,
const VectorDouble tab 
)
virtual

Set the contents of a Column

Reimplemented in MatrixSparse, and AMatrixDense.

◆ setDiagonal()

void AMatrix::setDiagonal ( const VectorDouble tab)
virtual

Set the contents of the (main) Diagonal

Reset the contents of a matrix by setting all terms to 0 and update diagonal terms from the input argument 'tab'

Parameters
tabInput vector to be copied to the diagonal of the output matrix

Reimplemented in MatrixSparse, and AMatrixDense.

◆ setDiagonalToConstant()

void AMatrix::setDiagonalToConstant ( double  value = 1.)
virtual

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

Reimplemented in MatrixSparse, and AMatrixDense.

◆ setFlagCheckAddress()

void AMatrix::setFlagCheckAddress ( bool  flagCheckAddress)
inline

◆ setIdentity()

void AMatrix::setIdentity ( double  value = 1.)

Sets the matrix as Identity

◆ setRow()

void AMatrix::setRow ( int  irow,
const VectorDouble tab 
)
virtual

Set the contents of a Row

Reimplemented in MatrixSparse, and AMatrixDense.

◆ setValue()

virtual void AMatrix::setValue ( int  irow,
int  icol,
double  value,
bool  flagCheck = true 
)
pure virtual

Sets the value at row 'irow' and column 'icol'

Implemented in MatrixSparse, and AMatrixDense.

◆ setValues()

void AMatrix::setValues ( const VectorDouble values,
bool  byCol = true 
)

Filling the matrix with an array of values Note that this array is ALWAYS dimensioned to the total number of elements in the matrix. Kept for compatibility with old code where matrix contents was stored as a VectorDouble

Parameters
values
byColtrue for Column major; false for Row Major

◆ size()

int AMatrix::size ( ) const
inline

Get the total number of elements of the (full) matrix

◆ solve()

int AMatrix::solve ( const VectorDouble b,
VectorDouble x 
) const

Solving the Matrix Linear system

◆ toString()

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

Interface to AStringable.

Reimplemented from AStringable.

Reimplemented in Table, MatrixSparse, and ProjMatrix.

◆ transpose()

AMatrix * AMatrix::transpose ( ) const
virtual

Transpose the matrix and return it as a copy

Reimplemented in MatrixSparse.

◆ transposeInPlace()

void AMatrix::transposeInPlace ( )
virtual

Transpose the matrix in place

◆ updValue()

virtual void AMatrix::updValue ( int  irow,
int  icol,
const EOperator &  oper,
double  value,
bool  flagCheck = true 
)
pure virtual

Update the value at row 'irow' and column 'icol'

Implemented in MatrixSparse, and AMatrixDense.


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