Loading [MathJax]/extensions/tex2jax.js
1.7.3
Geostatistics & Machine Learning toolbox | https://gstlearn.org
 
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
MatrixDense Class Reference

#include <MatrixDense.hpp>

Inheritance diagram for MatrixDense:
AMatrix AStringable ICloneable MatrixSquare Table MatrixSymmetric

Detailed Description

Dense Matrix This class provides all the functions that can be performed using a Matrix stored in "Dense" format (in opposition to the "Sparse" format). This class can be derived in the case the matrix is Square, and even more if it is Square and Symmetric.

Public Member Functions

 MatrixDense (int nrow=0, int ncol=0)
 
 MatrixDense (const MatrixDense &r)
 
 MatrixDense (const AMatrix &r)
 
MatrixDenseoperator= (const MatrixDense &r)
 
virtual ~MatrixDense ()
 
bool isDense () const override
 Cloneable interface.
 
bool isSparse () const override
 
bool mustBeSymmetric () const override
 
void setValue (int irow, int icol, double value, bool flagCheck=false) override
 
virtual double getValue (int irow, int icol, bool flagCheck=false) const override
 
void updValue (int irow, int icol, const EOperator &oper, double value, bool flagCheck=false) override
 
virtual void setColumn (int icol, const VectorDouble &tab, bool flagCheck=false) override
 
virtual void setRow (int irow, const VectorDouble &tab, bool flagCheck=false) override
 
virtual void setDiagonal (const VectorDouble &tab, bool flagCheck=false) override
 
virtual void setDiagonalToConstant (double value=1.) 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 VectorDouble getRow (int irow) const override
 
virtual VectorDouble getColumn (int icol) const override
 
constvect getColumnPtr (int icol) const
 
virtual void prodMatMatInPlace (const AMatrix *x, const AMatrix *y, bool transposeX=false, bool transposeY=false) override
 
void addMatInPlace (const MatrixDense &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 not generic of AMatrix anymore. WARNING: output matrix should not match any of input matrices (speed up).
 
virtual void prodNormMatMatInPlace (const MatrixDense *a, const MatrixDense *m, bool transpose=false)
 
virtual void prodNormMatVecInPlace (const MatrixDense &a, const VectorDouble &vec=VectorDouble(), bool transpose=false)
 
VectorDouble getEigenValues () const
 
const MatrixSquaregetEigenVectors () const
 
int invert2 (MatrixDense &res) const
 
void unsample (const AMatrix *A, const VectorInt &rowFetch, const VectorInt &colFetch, bool flagInvertRow=false, bool flagInvertCol=false)
 Set the values contained in 'A' into the current matrix.
 
MatrixDense compressMatLC (const MatrixDense &matLC, bool transpose=false)
 Perform the compressing product, according to 'transpose'.
 
void addRow (int nrow_added=1)
 
void addColumn (int ncolumn_added=1)
 
constvect getViewOnColumn (int icol) const
 
vect getViewOnColumnModify (int icol)
 
Eigen::Map< const Eigen::MatrixXd > getEigenMat () const
 
Eigen::Map< Eigen::MatrixXd > getEigenMat ()
 
- 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 void resetFromValue (int nrows, int ncols, double value)
 Reset the matrix to new dimensions and fill with a new value.
 
virtual void resetFromArray (int nrows, int ncols, const double *tab, bool byCol=true)
 Reset the matrix from an array of double values.
 
virtual void resetFromVD (int nrows, int ncols, const VectorDouble &tab, bool byCol=true)
 Reset the matrix from a vector of double values.
 
virtual void resetFromVVD (const VectorVectorDouble &tab, bool byCol=true)
 Reset the matrix from an array of double values.
 
virtual String toString (const AStringFormat *strfmt=nullptr) const override
 Interface to AStringable.
 
void clear ()
 
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 (double eps=EPSILON10, bool printWhyNot=false) const
 
virtual void transposeInPlace ()
 
virtual AMatrixtranspose () const
 
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 prodNormMatVecInPlace (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)
 
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 getNColDefined () const
 
int getNRowDefined () const
 
VectorDouble getColumnByRowRange (int icol, int rowFrom, int rowTo) const
 
bool isNonNegative (bool verbose=false) const
 
void prodMatVecInPlace (const VectorDouble &x, VectorDouble &y, bool transpose=false) const
 
int prodMatVecInPlace (const constvect x, vect 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 dumpStatistics (const String &title) const
 
void setIdentity (double value=1.)
 
void fillRandom (int seed=432432, double zeroPercent=0)
 
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 &validRows, const VectorInt &validCols)
 
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 val3=1., const AMatrix *mat3=nullptr)
 Perfom the algebraic equation this = val1 * mat1 + val2 * mat2 + val3 * mat3.
 
virtual int addProdMatVecInPlace (const constvect x, vect y, bool transpose=false) const
 
double operator() (int row, int col) const
 
double & operator() (int row, int col)
 
virtual bool needToReset (int nrows, int ncols)
 
- 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 MatrixDensecreate (const MatrixDense *mat)
 
static MatrixDensecreate (int nrow, int ncol)
 
static MatrixDensecreateFromVVD (const VectorVectorDouble &X)
 
static MatrixDensecreateFromVD (const VectorDouble &X, int nrow, int ncol, bool byCol=false, bool invertColumnOrder=false)
 
static MatrixDenseglue (const AMatrix *A1, const AMatrix *A2, bool flagShiftRow, bool flagShiftCol)
 
static MatrixDensesample (const AMatrix *A, const VectorInt &rowKeep=VectorInt(), const VectorInt &colKeep=VectorInt(), bool flagInvertRow=false, bool flagInvertCol=false)
 Create an output Rectangular Matrix by selecting some rows and columns of the Input matrix 'A'.
 
static void sum (const MatrixDense *mat1, const MatrixDense *mat2, MatrixDense *mat3)
 

Public Attributes

 DECLARE_TOTL
 Has a specific implementation in the Target language.
 

Constructor & Destructor Documentation

◆ MatrixDense() [1/3]

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

◆ MatrixDense() [2/3]

MatrixDense::MatrixDense ( const MatrixDense r)

◆ MatrixDense() [3/3]

MatrixDense::MatrixDense ( const AMatrix r)

◆ ~MatrixDense()

MatrixDense::~MatrixDense ( )
virtual

Member Function Documentation

◆ addColumn()

void MatrixDense::addColumn ( int  ncolumn_added = 1)

◆ addMatInPlace()

void MatrixDense::addMatInPlace ( const MatrixDense 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 not generic of AMatrix anymore. WARNING: output matrix should not match any of input matrices (speed up).

Add a matrix (multiplied by a constant)

◆ addRow()

void MatrixDense::addRow ( int  nrow_added = 1)

◆ addScalar()

void MatrixDense::addScalar ( double  v)
overridevirtual

Add a value to each matrix component

Reimplemented from AMatrix.

◆ addScalarDiag()

void MatrixDense::addScalarDiag ( double  v)
overridevirtual

Add value to matrix diagonal

Reimplemented from AMatrix.

◆ compressMatLC()

MatrixDense MatrixDense::compressMatLC ( const MatrixDense matLC,
bool  transpose = false 
)

Perform the compressing product, according to 'transpose'.

  • False: 'this'[nrows,ncols] %*% t('matLC')[ncolsCL,nrowsCL] -> mat[nrows,ncolsCL]
  • True: t('matLC')[ncolsCL,nrowsCL] %*% 'this'[nrows,ncols] -> mat[ncolsCL,ncols]
Parameters
matLC
transpose
Returns
MatrixDense

◆ create() [1/2]

MatrixDense * MatrixDense::create ( const MatrixDense mat)
static

◆ create() [2/2]

MatrixDense * MatrixDense::create ( int  nrow,
int  ncol 
)
static

◆ createFromVD()

MatrixDense * MatrixDense::createFromVD ( const VectorDouble X,
int  nrow,
int  ncol,
bool  byCol = false,
bool  invertColumnOrder = false 
)
static

◆ createFromVVD()

MatrixDense * MatrixDense::createFromVVD ( const VectorVectorDouble X)
static

Converts a VectorVectorDouble into a Matrix Note: the input argument is stored by row (if coming from [] specification)

Parameters
XInput VectorVectorDouble argument
Returns
The returned rectangular matrix
Remarks
: the matrix is transposed implicitly while reading

◆ divideColumn()

void MatrixDense::divideColumn ( const VectorDouble vec)
overridevirtual

Divide a Matrix column-wise

Reimplemented from AMatrix.

◆ divideRow()

void MatrixDense::divideRow ( const VectorDouble vec)
overridevirtual

Divide a Matrix row-wise

Reimplemented from AMatrix.

◆ fill()

void MatrixDense::fill ( double  value)
overridevirtual

Set all the values of the Matrix at once

Reimplemented from AMatrix.

◆ getColumn()

VectorDouble MatrixDense::getColumn ( int  icol) const
overridevirtual

Extract a Column

Reimplemented from AMatrix.

◆ getColumnPtr()

constvect MatrixDense::getColumnPtr ( int  icol) const

◆ getEigenMat() [1/2]

Eigen::Map< Eigen::MatrixXd > MatrixDense::getEigenMat ( )
inline

◆ getEigenMat() [2/2]

Eigen::Map< const Eigen::MatrixXd > MatrixDense::getEigenMat ( ) const
inline

◆ getEigenValues()

VectorDouble MatrixDense::getEigenValues ( ) const
inline

◆ getEigenVectors()

const MatrixSquare * MatrixDense::getEigenVectors ( ) const
inline

◆ getRow()

VectorDouble MatrixDense::getRow ( int  irow) const
overridevirtual

Extract a Row

Reimplemented from AMatrix.

◆ getValue()

double MatrixDense::getValue ( int  irow,
int  icol,
bool  flagCheck = false 
) const
overridevirtual

Get the value from a matrix cell

Implements AMatrix.

◆ getViewOnColumn()

constvect MatrixDense::getViewOnColumn ( int  icol) const

◆ getViewOnColumnModify()

vect MatrixDense::getViewOnColumnModify ( int  icol)

◆ glue()

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

◆ invert2()

int MatrixDense::invert2 ( MatrixDense res) const

TODO : check beforehand if matrix is invertible ?

◆ isDense()

bool MatrixDense::isDense ( ) const
inlineoverridevirtual

Cloneable interface.

Interface for AMatrix

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

Implements AMatrix.

◆ isSparse()

bool MatrixDense::isSparse ( ) const
inlineoverridevirtual

Returns if the current matrix is Sparse

Implements AMatrix.

◆ multiplyColumn()

void MatrixDense::multiplyColumn ( const VectorDouble vec)
overridevirtual

Multiply a Matrix column-wise

Reimplemented from AMatrix.

◆ multiplyRow()

void MatrixDense::multiplyRow ( const VectorDouble vec)
overridevirtual

Multiply a Matrix row-wise

Reimplemented from AMatrix.

◆ mustBeSymmetric()

bool MatrixDense::mustBeSymmetric ( ) const
inlineoverridevirtual

Say if the matrix must be symmetric

Reimplemented from AMatrix.

Reimplemented in MatrixSymmetric, and MatrixSquare.

◆ operator=()

MatrixDense & MatrixDense::operator= ( const MatrixDense r)

◆ prodMatMatInPlace()

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

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

Reimplemented from AMatrix.

◆ prodMatVec()

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

Perform 'this' * 'vec'

Reimplemented from AMatrix.

◆ prodNormMatMatInPlace()

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

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

Product of matrices, stored in 'this'

  • transpose = true: t('a') * 'm' * 'a'
  • transpose = false: 'a' * 'm' * t('a')
Parameters
aFirst input matrix
mSecond input matrix
transposeTrue if 'a' should be transposed beforehand
Note
: 'a' and 'm' may NOT coincide with 'this'

◆ prodNormMatVecInPlace()

void MatrixDense::prodNormMatVecInPlace ( const MatrixDense a,
const VectorDouble vec = VectorDouble(),
bool  transpose = false 
)
virtual

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

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

Parameters
aInput matrix
vecInput vector
transposeWhen True, the input Matrix is transposed

◆ prodScalar()

void MatrixDense::prodScalar ( double  v)
overridevirtual

Multiply each matrix component by a value

Reimplemented from AMatrix.

◆ prodVecMat()

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

Perform 'vec' * 'this'

Reimplemented from AMatrix.

◆ sample()

MatrixDense * MatrixDense::sample ( const AMatrix A,
const VectorInt rowKeep = VectorInt(),
const VectorInt colKeep = VectorInt(),
bool  flagInvertRow = false,
bool  flagInvertCol = false 
)
static

Create an output Rectangular Matrix by selecting some rows and columns of the Input matrix 'A'.

Parameters
AInput Rectangular Matrix
rowKeepSet of Rows to be kept (all if not defined)
colKeepSet of Columns to be kept (all if not defined)
flagInvertRowwhen True, transform 'rowKeep' into 'rowDrop'
flagInvertColwhen True, transform 'colKeep' into 'colDrop'
Returns
Newly created Rectangular Matrix

◆ setColumn()

void MatrixDense::setColumn ( int  icol,
const VectorDouble tab,
bool  flagCheck = false 
)
overridevirtual

Set the contents of a Column

Reimplemented from AMatrix.

◆ setDiagonal()

void MatrixDense::setDiagonal ( const VectorDouble tab,
bool  flagCheck = false 
)
overridevirtual

Set the contents of the (main) Diagonal

Reimplemented from AMatrix.

◆ setDiagonalToConstant()

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

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

Reimplemented from AMatrix.

◆ setRow()

void MatrixDense::setRow ( int  irow,
const VectorDouble tab,
bool  flagCheck = false 
)
overridevirtual

Set the contents of a Row

Reimplemented from AMatrix.

◆ setValue()

void MatrixDense::setValue ( int  irow,
int  icol,
double  value,
bool  flagCheck = false 
)
overridevirtual

Set the value for in a matrix cell

Implements AMatrix.

◆ sum()

void MatrixDense::sum ( const MatrixDense mat1,
const MatrixDense mat2,
MatrixDense mat3 
)
static

◆ unsample()

void MatrixDense::unsample ( const AMatrix A,
const VectorInt rowFetch,
const VectorInt colFetch,
bool  flagInvertRow = false,
bool  flagInvertCol = false 
)

Set the values contained in 'A' into the current matrix.

Parameters
AInput Matrix
rowFetchSet of row indices of 'this' where values of 'A' should be stored
colFetchSet of column indices of'this' where values of 'A' should be stored
flagInvertRowwhen True, transform 'rowFetch' into 'rowAvoid'
flagInvertColwhen True, transform 'colFetch' into 'colAvoid'

◆ updValue()

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

Update the contents of a matrix cell

Implements AMatrix.

Member Data Documentation

◆ DECLARE_TOTL

MatrixDense::DECLARE_TOTL

Has a specific implementation in the Target language.


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