Class xlifepp::LargeMatrix#

template<typename T>
class LargeMatrix#

Inheritence diagram for xlifepp::LargeMatrix:

digraph { graph [bgcolor="#00000000"] node [shape=rectangle style=filled fillcolor="#FFFFFF" font=Helvetica padding=2] edge [color="#1414CE"] "2" [label="xlifepp::LargeMatrix< TA_ >" tooltip="xlifepp::LargeMatrix< TA_ >"] "3" [label="xlifepp::LargeMatrix< TB_ >" tooltip="xlifepp::LargeMatrix< TB_ >"] "4" [label="xlifepp::LargeMatrix< complex_t >" tooltip="xlifepp::LargeMatrix< complex_t >"] "5" [label="xlifepp::LargeMatrix< real_t >" tooltip="xlifepp::LargeMatrix< real_t >"] "7" [label="xlifepp::LargeMatrix< xlifepp::Matrix< complex_t > >" tooltip="xlifepp::LargeMatrix< xlifepp::Matrix< complex_t > >"] "6" [label="xlifepp::LargeMatrix< xlifepp::Matrix< real_t > >" tooltip="xlifepp::LargeMatrix< xlifepp::Matrix< real_t > >"] "1" [label="xlifepp::LargeMatrix< T >" tooltip="xlifepp::LargeMatrix< T >" fillcolor="#BFBFBF"] "2" -> "1" [dir=forward tooltip="template-instance"] "3" -> "1" [dir=forward tooltip="template-instance"] "4" -> "1" [dir=forward tooltip="template-instance"] "5" -> "1" [dir=forward tooltip="template-instance"] "7" -> "1" [dir=forward tooltip="template-instance"] "6" -> "1" [dir=forward tooltip="template-instance"] }

Collaboration diagram for xlifepp::LargeMatrix:

digraph { graph [bgcolor="#00000000"] node [shape=rectangle style=filled fillcolor="#FFFFFF" font=Helvetica padding=2] edge [color="#1414CE"] "3" [label="std::vector< T >" tooltip="std::vector< T >"] "2" [label="std::vector< number_t >" tooltip="std::vector< number_t >"] "5" [label="std::vector< xlifepp::MatrixStorage * >" tooltip="std::vector< xlifepp::MatrixStorage * >"] "1" [label="xlifepp::LargeMatrix< T >" tooltip="xlifepp::LargeMatrix< T >" fillcolor="#BFBFBF"] "4" [label="xlifepp::MatrixStorage" tooltip="xlifepp::MatrixStorage"] "2" -> "3" [dir=forward tooltip="template-instance"] "5" -> "3" [dir=forward tooltip="template-instance"] "1" -> "2" [dir=forward tooltip="usage"] "1" -> "4" [dir=forward tooltip="usage"] "4" -> "5" [dir=forward tooltip="usage"] }

templated class to deal with large matrices, it handles a vector of values (non zero coefficients) and a pointer to a MatrixStorage

Public Types

typedef LargeMatrixTrait<T>::ScalarType ScalarType#

useful typedef for LargeMatrix class

Public Functions

LargeMatrix()#

default constructor

template<typename K>
LargeMatrix(const LargeMatrix<K>&, bool storageCopy = false)#

generalized copy constructor (shallow copy), with option to copy or not the storage

LargeMatrix(const LargeMatrix<T>&, bool storageCopy = false)#

copy constructor (shadow copy), with option to copy or not the storage

LargeMatrix(const string_t &fn, StorageType fst, number_t nr, number_t nc, StorageType sto = _dense, AccessType at = _row, number_t nsr = 1, number_t nsc = 1)#

construct matrix from file

LargeMatrix(const string_t &fn, StorageType fst, number_t nr, StorageType sto = _dense, SymType sy = _symmetric, number_t nsr = 1)#

construct “symmetric” matrix from file

LargeMatrix(MatrixStorage*, const T &val = T(0), SymType sy = _noSymmetry)#

construct matrix from storage and constant value

LargeMatrix(MatrixStorage*, SymType)#

construct matrix from storage and symmetry type (scalar values assumed)

LargeMatrix(MatrixStorage *ms, dimPair dims, SymType sy = _noSymmetry)#

construct matrix from storage and symmetry type (matrix values)

LargeMatrix(number_t nr, number_t nc, StorageType sto = _dense, SymType = _symmetric, const T &val = T())#

construct “symmetric” matrix from dimensions, storage type, symmetry type a constant value

LargeMatrix(SpecialMatrix sp, StorageType st, AccessType at, number_t nbr, number_t nbc, const T v)#

construct Special matrix (diagonal)

LargeMatrix(ValueType vt, StrucType st, number_t nbr, number_t nbc, SymType sy, dimen_t nbrs, dimen_t nbcs, const T &val, const string_t &na, MatrixStorage *sp = nullptr)#

explicit constructor

~LargeMatrix()#

destructor construct matrix from dimensions, storage type and a constant value

inline AccessType accessType() const#

return the access type

template<typename K, typename C>
LargeMatrix<T> &add(const LargeMatrix<K>&, const std::vector<number_t>&, const std::vector<number_t>&, C)#

add scaled LargeMatrix at ranks (rows, cols) of current LargeMatrix

add to current LargeMatrix a scaled LargeMatrix mat (submatrix) given by its ranks (rows, cols) in current matrix and a scaling factor (scale) a void rows or cols means that the numbering of mat is the same as current matrix in that case, if the storage of mat and current matrix are the same, direct summation is used else access to position is used (expansive method) no storage update !!!

LargeMatrix<T> &add(const std::vector<T>&, const std::vector<number_t>&, const std::vector<number_t>&)#

add a matrix at ranks (rows, cols)

LargeMatrix<T> &add(const T &c, const std::vector<number_t> &rows, const std::vector<number_t> &cols)#

add a constant to submatrix given by its ranks (rows, cols)

template<typename iterator>
LargeMatrix<T> &add(iterator&, iterator&)#

add values given by iterator from beginning

void addColToCol(number_t c1, number_t c2, complex_t a, bool updateStorage = false)#

combine cols of matrix: c2 = c2 + a *c1 (indices start from 1) with update storage option

void addRowToRow(number_t r1, number_t r2, complex_t a, bool updateStorage = false)#

combine rows of matrix: r2 = r2 + a *r1 (indices start from 1) with update storage option

template<typename K>
LargeMatrix<T> &assign(const LargeMatrix<K>&, const std::vector<int_t>&, const std::vector<int_t>&)#

assign LargeMatrix at ranks (rows, cols) to current LargeMatrix

assign to current LargeMatrix a LargeMatrix mat (submatrix) given by its ranks (rows, cols) in current matrix a void rows or cols means that the numbering of mat is the same as current matrix in that case, if the storage of mat and current matrix are the same, direct summation is used if rows (resp cols) = [-n], it means a trivial numbering 1,2,…,n else more expansive method is used: use positions function !!! target storage has to be update and T = K has to be available NO CHECK!!!

inline void clear()#

deallocate memory used by matrix values, storage is not deallocated!

void clearall()#

clear values and storage pointer (not the attributes)

std::vector<T> col(number_t c) const#

get col c (>=1) as vector

void copyVal(const LargeMatrix<T>&, const std::vector<number_t>&, const std::vector<number_t>&)#

copy in current matrix values of mat located at mat_adrs at position cur_adrs

copy in current matrix values of mat located at mat_adrs at position cur_adrs it is assumed that mat_adrs and cur_adrs have the same size (address map!)

void deleteCols(number_t, number_t)#

delete cols c1,…, c2 (may be expansive for some storage)

void deleteRows(number_t, number_t)#

delete rows r1,…, r2 (may be expansive for some storage)

inline dimPair dimValues() const#

return dimensions of values, (1,1) when scalar

string_t encodeFileName(const string_t&, StorageType) const#

encode file name with additional matrix informations

LargeMatrix<T> *extract(const std::vector<number_t> &rowIndex, const std::vector<number_t> &colIndex) const#

extract LargeMatrix from current LargeMatrix (return a pointer to LargeMatrix)

extract sub matrix of current LargeMatrix and return it as a LargeMatrix rowIndex: row indices to be extracted (starts from 1) colIndex: col indices to be extracted (starts from 1) return a pointer to extracted LargeMatrix in the same storage type of current LargeMatrix

LargeMatrix<T> *extract(number_t r1, number_t r2, number_t c1, number_t c2) const#

extract LargeMatrix from current LargeMatrix (return a pointer to LargeMatrix)

extract block matrix [r1,r2]x[c1,c2] of current LargeMatrix and return it as a LargeMatrix return a pointer to extracted LargeMatrix in the same storage type of current LargeMatrix r1<=r2,c1<=c2 indices starting from 1

void extract2UmfPack(std::vector<int_t> &colPointer, std::vector<int_t> &rowIndex, std::vector<T> &matA) const#

extract and convert largeMatrix (and its storage) to Umfpack

inline T get(number_t i, number_t j) const#

access to entry (i,j) (do not use in heavy computation)

inline std::vector<std::pair<number_t, number_t>> getCol(number_t c, number_t r1 = 1, number_t r2 = 0) const#

get (row indices, adresses) of col c in set [r1,r2] (r2=0 means nbOfRows)

bool getCsColUmfPack(std::vector<int_t> &colPointer, std::vector<int_t> &rowIndex, const T *&matA) const#

get largeMatrix and its storage to Umfpack (only for cs col)

inline std::vector<std::pair<number_t, number_t>> getRow(number_t r, number_t c1 = 1, number_t c2 = 0) const#

get (col indices, adresses) of row r in set [c1,c2] (c2=0 means nbOfCols)

void ildlstarFactorize()#

Incomplete Factorization LDL*.

void ildltFactorize()#

Incomplete Factorization LDLt.

void illstarFactorize()#

Incomplete Factorization LL*.

void illtFactorize()#

Incomplete Factorization LLt.

void iluFactorize()#

Incomplete Factorization LU.

void init(MatrixStorage *ms, const T &val, SymType sy)#

initialize matrix from storage

bool isDiagonal() const#

true if matrix is in fact diagonal

bool isId(const double &tol = 0.) const#

true if matrix is in fact id with a tolerance

void ldlstarFactorize()#

Factorization LDL*.

template<typename S>
void ldlstarSolve(std::vector<S> &vec, std::vector<T> &res) const#

templated factorization LDL* solver

void ldltFactorize()#

return the condition number if it is a factorized matrix

Factorization LDLt

template<typename S1, typename S2>
void ldltSolve(std::vector<S1> &vec, std::vector<S2> &res) const#

templated factorization LDLt solver

template<typename S1, typename S2>
void lltSolve(std::vector<S1> &vec, std::vector<S2> &res) const#

templated factorization LLt solver

void loadFromFile(const string_t&, StorageType)#

load matrix from file

LargeMatrix<T> &lowerD1Solve(const LargeMatrix<T> &M, const LargeMatrix<T> &X) const#

lower and upper solvers on LargeMatrix: LD1^{-1}*M, U^{-1}*M, M*LD1^{-1}, M*U^{-1}, result is stored as col dense

std::vector<T> &lowerD1Solve(std::vector<T> &b, std::vector<T> &x) const#

lower and upper solvers on vector: LD1^{-1}*b, U^{-1}*b, b*LD1^{-1}, b*U^{-1}, b may be modified

void luFactorize(bool withPermutation = true)#

Factorization LU.

do LU factorization of skyline or dense matrix LU with row permutation is available only for dense matrices if matrix is stored as symmetric it is extended to store upper part (toUnsymmetric function)

template<typename S1, typename S2>
void luLeftSolve(std::vector<S1> &vec, std::vector<S2> &res) const#

templated factorization LU left solver

template<typename S1, typename S2>
void luSolve(std::vector<S1> &vec, std::vector<S2> &res) const#

templated factorization LU solver

void multLeftMatrixCol(T *pM, T *pR, number_t p) const#

product dense col matrix * LargeMatrix passed as pointer

do the product R = M * A where A is any largeMatrix pM pointer to the first element of the dense col matrix pR pointer to the first element of the dense col result matrix p the number of rows of M assuming the number of rows of M is equal to the number of cols of A

void multLeftMatrixRow(T *pM, T *pR, number_t p) const#

product dense row matrix * LargeMatrix passed as pointer

do the product R = M * A where A is any largeMatrix pM pointer to the first element of the dense row matrix pR pointer to the first element of the dense row result matrix p the number of rows of M assuming the number of rows of M is equal to the number of cols of A

void multMatrixCol(T *pM, T *pR, number_t p) const#

product LargeMatrix * dense col matrix passed as pointer

do the product R = A * M where A is any largeMatrix pM pointer to the first element of the dense col matrix pR pointer to the first element of the dense col result matrix p the number of cols of M assuming the number of rows of M is equal to the number of cols of A

void multMatrixRow(T *pM, T *pR, number_t p) const#

product LargeMatrix * dense row matrix passed as pointer

do the product R = A * M where A is any largeMatrix pM pointer to the first element of the dense row matrix pR pointer to the first element of the dense row result matrix p the number of cols of M assuming the number of rows of M is equal to the number of cols of A

inline number_t nbNonZero() const#

number of non zero coefficients stored (given by storage)

real_t norm2() const#

Frobenius norm of a LargeMatrix.

real_t norminfty() const#

infinite norm of a LargeMatrix

inline number_t numberOfCols() const#

return number of columns counted in T

inline number_t numberOfRows() const#

return number of rows counted in T

inline T &operator()(number_t adr)#

access to entry at adress adr

inline T operator()(number_t adr) const#

access to entry at adress adr

inline T operator()(number_t i, number_t j) const#

access to entry (i,j) (do not use in heavy computation)

inline T &operator()(number_t i, number_t j, bool errorOn = false)#

access to entry (i,j) with check (do not use in heavy computation)

template<typename K>
LargeMatrix<T> &operator*=(const K&)#

LargeMatrix<T> *= v.

LargeMatrix<T> &operator+=(const LargeMatrix<T>&)#

LargeMatrix<T> += LargeMatrix<T>

LargeMatrix<T> &operator-=(const LargeMatrix<T>&)#

LargeMatrix<T> -= LargeMatrix<T>

template<typename K>
LargeMatrix<T> &operator/=(const K&)#

LargeMatrix<T> /= v.

template<typename K>
LargeMatrix<T> &operator=(const LargeMatrix<K>&)#

generalized assign operator

LargeMatrix<T> &operator=(const LargeMatrix<T>&)#

assign operator

real_t partialNormOfCol(number_t c, number_t r1, number_t r2) const#

partial norm of column c, only components in range [r1,r2]

inline void positions(const std::vector<number_t> &rows, const std::vector<number_t> &cols, std::vector<number_t> &adrs, bool errorOn = false) const#

positions in storage of a submatrix

void print(std::ostream&) const#

print on output stream

inline void resetZero()#

reset first value to 0

void roundToZero(real_t aszero = 10 * theEpsilon)#

round to 0 all coefficients smaller (in norm) than a aszero, storage is not modified

std::vector<T> row(number_t r) const#

get row r (>=1) as vector

inline void saveToFile(const char *fn, StorageType st, bool encode, real_t tol = theTolerance) const#

save matrix to a file along different format

void saveToFile(const char *fn, StorageType st, number_t prec = fullPrec, bool encode = false, real_t tol = theTolerance) const#

save matrix to a file along different format

inline void saveToFile(const string_t &fn, StorageType st, bool encode, real_t tol = theTolerance) const#

save matrix to a file along different format

void saveToFile(const string_t &fn, StorageType st, number_t prec = fullPrec, bool encode = false, real_t tol = theTolerance) const#

save matrix to a file along different format

void setCol(const T &val, number_t c1 = 0, number_t c2 = 0)#

set to val cols c1,…, c2 (index start from 1), no col deletion

set to val cols/rows cr1,…, cr2 (index start from 1) no col/row deletion, no storage modification (zeros stay zeros !) generic algorithm use getCol/getRow of storage, may be slightly improved by using storage child versions, not yet done

void setColToZero(number_t c1 = 0, number_t c2 = 0)#

set to 0 cols c1,…, c2 (index start from 1), no col deletion

void setRow(const T &val, number_t r1 = 0, number_t r2 = 0)#

set to val rows r1,…, r2 (index start from 1), no col deletion

void setRowToZero(number_t r1 = 0, number_t r2 = 0)#

set to 0 rows r1,…, r2 (index start from 1), no row deletion

template<typename S1, typename S2>
void sorDiagonalMatrixVector(const std::vector<S1> &vec, std::vector<S2> &res, const real_t w) const#

templated factorization SOR diagonal solver

template<typename S1, typename S2>
void sorDiagonalSolve(const std::vector<S1> &vec, std::vector<S2> &res, const real_t w) const#

templated factorization SOR diagonal solver

template<typename S1, typename S2>
void sorLowerMatrixVector(const std::vector<S1> &vec, std::vector<S2> &res, const real_t w) const#

templated factorization SOR lower part solver

template<typename S1, typename S2>
void sorLowerSolve(const std::vector<S1> &vec, std::vector<S2> &res, const real_t w) const#

templated factorization SOR lower part solver

template<typename S1, typename S2>
void sorUpperMatrixVector(const std::vector<S1> &vec, std::vector<S2> &res, const real_t w) const#

templated factorization SOR upper part solver

template<typename S1, typename S2>
void sorUpperSolve(const std::vector<S1> &vec, std::vector<S2> &res, const real_t w) const#

templated factorization SOR upper part solver

inline real_t squaredNorm() const#

square of Frobenius norm of a LargeMatrix

inline MatrixStorage *&storagep()#

return the storage pointer

inline MatrixStorage *storagep() const#

return the storage pointer (const)

inline StorageType storageType() const#

return the storage type

inline SymType symmetry() const#

return kind of symmetry

LargeMatrix<T> &toConj()#

conjugate matrix

template<typename K>
LargeMatrix<K> *toScalar(K)#

create new scalar matrix entry from non scalar matrix entry

template<typename K>
void toScalarEntries(const std::vector<Matrix<K>>&, std::vector<K>&, const MatrixStorage&)#

copy matrix values to scalar values

copy matrix values to scalar values, generic algorithm using getrow or getcol K is the valueType of both matrices, sto is the matrixstorage of the scalar matrix

void toSkyline()#

convert current matrix storage to skyline storage

void toStorage(MatrixStorage*)#

convert current matrix storage to an other matrix storage

void toUnsymmetric(AccessType at = _sym)#

convert current matrix to unsymmetric representation if it has a symmetric one

convert current matrix to unsymmetric representation if it has a symmetric one convert only matrix using a symmetric access storage, if at = _sym symmetric storage is not changed (default behaviour) if at = _dual we moved to a dual storage values vector is extended to store upper part

template<typename S1, typename S2>
void umfluLeftSolve(std::vector<S1> &vec, std::vector<S2> &res) const#

templated factorization LU left solver (UmfPack)

template<typename S1, typename S2>
void umfluSolve(std::vector<S1> &vec, std::vector<S2> &res) const#

templated factorization LU solver (UmfPack)

void umfpackFactorize()#

Factorization using umfpack.

template<typename S>
real_t umfpackSolve(const std::vector<int_t> &colPointer, const std::vector<int_t> &rowIndex, const std::vector<T> &values, const std::vector<S> &vec, std::vector<typename Conditional<NumTraits<ScalarType>::IsComplex, ScalarType, S>::type> &res, bool soa = true)#

solve linear system using umfpack

template<typename S>
real_t umfpackSolve(const std::vector<S> &vec, std::vector<typename Conditional<NumTraits<ScalarType>::IsComplex, ScalarType, S>::type> &res, bool soa = true)#

templated UMFPACK solver

solve linear system using umfpack

template<typename S>
real_t umfpackSolve(const std::vector<std::vector<S>*> &bs, std::vector<std::vector<typename Conditional<NumTraits<ScalarType>::IsComplex, ScalarType, S>::type>*> &xs, bool soa = true)#

solve linear system using umfpack, multiple rhs

inline std::vector<T> &values()#

access to values

inline const std::vector<T> &values() const#

access to values (const)

inline ValueType valueType() const#

return type of values (real, complex)

void viewStorage(std::ostream&) const#

visualize storage on output stream

Public Members

std::vector<number_t> colPermutation_#

col permutation given by some factorizations (LU in col dense storage for instance)

FactorizationType factorization_#

one of _noFactorization, _lu, _ldlt, _ldlstar, _llt, _llstar, _qr, _umfpack

string_t name#

optionnal name, useful for documentation

number_t nbCols#

number of columns counted in T

dimen_t nbColsSub#

number of columns cof submatrices (1 if scalar values)

number_t nbRows#

number of rows counted in T

dimen_t nbRowsSub#

number of rows of submatrices (1 if scalar values)

std::vector<number_t> rowPermutation_#

row permutation given by some factorizations (LU in row dense storage for instance)

StrucType strucType#

struc of values (scalar, vector, matrix)

SymType sym#

type of symmetry

ValueType valueType_#

type of values (real, complex)

Friends

template<typename SA, typename SB, typename SR>
friend void multMatrixMatrix(const LargeMatrix<SA>&, const LargeMatrix<SB>&, LargeMatrix<SR>&)#

extern template product matrix*matrix the storage of the resulting matrix is ALWAYS dense there is no chance to get a sparse matrix except in case of product with a diagonal matrix where the storage is unchanged