Class xlifepp::MatrixEntry#

class MatrixEntry#

Collaboration diagram for xlifepp::MatrixEntry:

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

small class handling current types of LargeMatrix shadowing template parameter of LargeMatrix

Public Functions

MatrixEntry()#

default constructor

MatrixEntry(const MatrixEntry&, bool storageCopy = false)#

copy constructor (full copy) with option to force the copy of storage

MatrixEntry(SpecialMatrix, StorageType, AccessType, number_t, number_t, const complex_t&)#

construct special matrix (constant diagonal)

MatrixEntry(SpecialMatrix, StorageType, AccessType, number_t, number_t, const Matrix<complex_t>&)#

construct special matrix (constant diagonal)

MatrixEntry(SpecialMatrix, StorageType, AccessType, number_t, number_t, const Matrix<real_t>&)#

construct special matrix (constant diagonal)

MatrixEntry(SpecialMatrix, StorageType, AccessType, number_t, number_t, const real_t&)#

construct special matrix (constant diagonal)

MatrixEntry(ValueType, StrucType, MatrixStorage*, dimPair valdim = dimPair(1, 1), SymType sy = _noSymmetry)#

constructor from a storage pointer

~MatrixEntry()#

destructor

AccessType accessType() const#

return the access type

void add(const MatrixEntry&, const std::vector<number_t>&, const std::vector<number_t>&, complex_t)#

add MatrixEntry to current MatrixEntry with scalar multiplication

add MatrixEntry m to current MatrixEntry with scalar multiplication num_r: row dof numbers of MatrixEntry m related to current MatrixEntry num_c: col dof numbers of MatrixEntry m related to current MatrixEntry if num_r (resp.

num_c) is void, it means that the row (resp. col) dof numbers of MatrixEntry m is the same as current MatrixEntry

template<typename iterator>
void add(iterator&, iterator&)#

add entries to current MatrixEntry

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)

combine cols of matrix: c2 = c2 + a *c1 (indices start from 1) if updatestorage is true, the storage is dynamically updated (very expansive operation !!!)

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)

combine rows of matrix: r2 = rc2 + a *r1 (indices start from 1) if updatestorage is true, the storage is dynamically updated (very expansive operation !!!)

std::vector<std::pair<complex_t, VectorEntry*>> arpackSolve(number_t nev, real_t tol, string_t which, EigenComputationalMode eigCompMode, complex_t sigma, const char mode, bool isRealShift)#

arpack solver

void assign(const MatrixEntry&, const std::vector<int_t>&, const std::vector<int_t>&)#

assign entries to current MatrixEntry

assign MatrixEntry m to current MatrixEntry num_r: row dof numbers of MatrixEntry m related to current MatrixEntry num_c: col dof numbers of MatrixEntry m related to current MatrixEntry if num_r (resp.

num_c) is void, it means that the row (resp. col) dof numbers of MatrixEntry m is the same as current MatrixEntry if num_r (resp num_c) = [-n], it means a trivial numbering 1,2,…,n

void clear()#

deallocate memory used by a MatrixEntry

void copyVal(const MatrixEntry&, 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 and current matrix are consistent in structure type scalar/scalar or matrix/matrix in value type real/real or complex/complex if not do conversion before it is also assumed that mat_adrs and cur_adrs have the same size (addresses map!)

void deleteCols(number_t, number_t)#

delete cols c1,…, c2 in matrix

void deleteRows(number_t, number_t)#

delete rows r1,…, r2 in matrix

std::vector<std::pair<complex_t, VectorEntry*>> eigenSolve(const MatrixEntry *pB, EigenComputationalMode, number_t nev, real_t tol, string_t which, complex_t sigma = complex_t(), bool isRealShift = true, bool isShift = false)#

Resolve an eigen problem with different method This function is a wrapper of different eigenSolvers of class LargeMatrix (e.g: eigenDavidsonSolver, eigenKrylovSchurSolver, …)

Parameters:
  • pB[in] MatrixEntry containing the matrix B in the eigen problem A*X = lambda*B*X

  • eigCompMode[in] Computation mode of eigensolver

  • nev[in] The number of eigen values and eigenSolver searched for

  • tol[in] Tolerance

  • which[in] Specification of which eigen values are returned. Largest-LM, Smallest-SM

  • sigma[in] shift value used in shift mode. In case of real shift, the real part of sigma will be used

  • isRealShift[in] Specify if real shift is used

  • isShift[in] Specify if real or complex shift is used

std::vector<std::pair<complex_t, VectorEntry*>> eigenSolve(EigenComputationalMode, number_t nev, real_t tol, string_t which, complex_t sigma = complex_t(), bool isRealShift = true, bool isShift = false)#

eigen solvers

Resolve an eigen problem with different method This function is a wrapper of different eigenSolvers of class LargeMatrix (e.g: eigenDavidsonSolver, eigenKrylovSchurSolver, …)

Parameters:
  • eigCompMode[in] enumeration of different eigenSolver(davidson, krylovSchur). By default, it should be KrylovSchur

  • nev[in] The number of eigen values and eigenSolver searched for

  • tol[in] Tolerance

  • which[in] Specification of which eigen values are returned. Largest-LM, Smallest-SM

  • sigma[in] shift value used in shift mode. In case of real shift, the real part of sigma will be used

  • isRealShift[in] Specify if real shift is used

  • isShift[in] Specify if real or complex shift is used

FactorizationType factorization() const#

return factorisation type

Value getEntry(number_t i, number_t j, dimen_t k = 0, dimen_t l = 0) const#

get (i,j) as a Value (indices start from 1)

return matrix coefficients (i,j) >=1 as a value when a matrix of matrices, if k!=0 returns the scalar coefficients (k, l)>=1 of the matrix coefficient (i,j) else return the matrix coefficient (i,j) when matrix is a scalar one, (k,l) have no effect

inline complex_t getEntry(number_t, number_t, complex_t&) const#

get (i,j) complex value (indices start from 1)

inline Matrix<complex_t> getEntry(number_t, number_t, Matrix<complex_t>&) const#

get (i,j) complex matrix value (indices start from 1)

inline Matrix<real_t> getEntry(number_t, number_t, Matrix<real_t>&) const#

get (i,j) real matrix value (indices start from 1)

inline real_t getEntry(number_t, number_t, real_t&) const#

get (i,j) real value (indices start from 1)

get (i,j) entries, these functions use pos(i,j) functions of storage

template<>
LargeMatrix<real_t> &getLargeMatrix() const#

return a reference to LargeMatrix object

template<typename T>
LargeMatrix<T> &getLargeMatrix() const#

return a reference to LargeMatrix object

VectorEntry getRowCol(number_t, AccessType) const#

return a row or a col as VectorEntry

void ildlstarFactorize()#

factorize matrix entry as ILDL*

void ildltFactorize()#

factorize matrix entry as ILDLt

void illstarFactorize()#

factorize matrix entry as ILL*

void illtFactorize()#

factorize matrix entry as ILLt

void iluFactorize()#

factorize matrix entry as ILU

bool isDiagonal() const#

true if matrix is a diagonal one

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

true if matrix is a diagonal one

inline bool isvoid() const#

true if matrix is void

void ldlstarFactorize()#

factorize matrix entry as LDL*

void ldlstarSolve(const VectorEntry&, VectorEntry&)#

solve linear system using LDL* factorization

void ldltFactorize()#

factorize matrix entry as LDLt

void ldltSolve(const VectorEntry&, VectorEntry&)#

solve linear system using LDLt factorization

void lltSolve(const VectorEntry&, VectorEntry&)#

solve linear system using LLt factorization

void luFactorize(bool withPermutation = true)#

factorize matrix entry as LU

void luLeftSolve(const VectorEntry&, VectorEntry&)#

solve transposed linear system using LU factorization

void luSolve(const VectorEntry&, VectorEntry&)#

solve linear system using LU factorization

number_t nbOfCols() const#

return the number of columns

number_t nbOfRows() const#

return the number of rows

real_t norm2() const#

Frobenius norm.

real_t norminfty() const#

infinite norm

MatrixEntry &operator*=(const complex_t&)#

MatrixEntry *= c.

MatrixEntry &operator*=(const real_t&)#

MatrixEntry *= r.

MatrixEntry &operator+=(const MatrixEntry&)#

MatrixEntry += MatrixEntry.

MatrixEntry &operator-=(const MatrixEntry&)#

MatrixEntry -= MatrixEntry.

MatrixEntry &operator/=(const complex_t&)#

MatrixEntry /= c.

MatrixEntry &operator/=(const real_t&)#

MatrixEntry /= r.

MatrixEntry &operator=(const MatrixEntry&)#

assign operator (full copy)

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

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

partial norm of column c, only components in range [r1,r2], say sqrt( sum_{r1 <= r <= r2} Mrc^2 )

void print(std::ostream&) const#

print entry on stream

void roundToZero(real_t aszero = 10 * theEpsilon)#

round to zero all coefficients close to 0 (|.

| < aszero)

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

save matrix entries to dense or sparse Matlab format (_coo)

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

save matrix entries to dense or sparse Matlab format (_coo)

inline number_t scalarSize() const#

size of matrixentry, counted in scalar

void setCol(const Value &val, number_t c1, number_t c2)#

set c1 to c2 rows to val (index start from 1), no storage modification!

set cols c1->c2 to val (index start from 1), NO STORAGE MODIFICATION

void setColToZero(number_t, number_t)#

set c1 to c2 cols to 0 (index start from 1)

set cols or rows to 0 (index start from 1)

inline void setEntry(number_t, number_t, const complex_t&, bool errorOn = true)#

set (i,j) complex value (indices start from 1)

inline void setEntry(number_t, number_t, const Matrix<complex_t>&, bool errorOn = true)#

set (i,j) complex matrix value (indices start from 1)

inline void setEntry(number_t, number_t, const Matrix<real_t>&, bool errorOn = true)#

set (i,j) real matrix value (indices start from 1)

inline void setEntry(number_t, number_t, const real_t&, bool errorOn = true)#

set (i,j) real value (indices start from 1)

set (i,j) entries, be care: unallocated pointers are not checked ! these functions use pos(i,j) functions of storage, it may be expansive if onError is true then stops on indices out of storage

void setEntry(number_t, number_t, const Value&, dimen_t k = 0, dimen_t l = 0, bool errorOn = true)#

set (i,j) value (indices start from 1)

set matrix coefficient i,j from a Value (index start from 1) when a matrix of matrices, if k!=0 set the scalar coefficients (k, l)>=1 of the matrix coefficient (i,j) else set the matrix coefficient (i,j) when matrix is a scalar one, (k,l) have no effect

void setNbOfCols(number_t n)#

set the number of columns

void setNbOfRows(number_t n)#

set the number of rows

void setRow(const Value &val, number_t r1, number_t r2)#

set r1 to r2 rows to val (index start from 1), no storage modification!

set rows r1->r2 to val (index start from 1), NO STORAGE MODIFICATION

void setRowToZero(number_t, number_t)#

set r1 to r2 rows to 0 (index start from 1)

inline number_t size() const#

size of matrixentry, counted in dofs

void sorDiagonalMatrixVector(const VectorEntry&, VectorEntry&, const real_t)#

solve linear system using SOR block diagonal solver

void sorDiagonalSolve(const VectorEntry&, VectorEntry&, const real_t)#

solve linear system using SOR diagonal solver

void sorLowerSolve(const VectorEntry&, VectorEntry&, const real_t)#

solve linear system using SOR lower part solver

void sorUpperSolve(const VectorEntry&, VectorEntry&, const real_t)#

solve linear system using SOR upper part solver

MatrixStorage *storagep()#

return storage pointer (non const)

MatrixStorage *storagep() const#

return storage pointer (const)

StorageType storageType() const#

return the storage type

inline StrucType &strucType()#

return strucType (non const)

inline StrucType strucType() const#

return strucType (const)

SymType &symmetry()#

return symmetry of matrix (r/w)

SymType symmetry() const#

return symmetry of matrix

MatrixEntry &toComplex()#

change type of entries rEntries-> cEntries or rmEntries -> cmEntries

MatrixEntry &toConj()#

modify current entries to conjugate entries if they are not real

void toMatrix(const std::vector<dimen_t>&, const std::vector<dimen_t>&)#

change rEntries-> rmEntries or cEntries -> rcEntries

change matrix entry representation: from current representation to matrix of nbr x nbc matrices rowPos: list of row indexes to be moved, nbr = rowPos.size colPos: list of col indexes to be moved, nbc = colPos.size for instance [0 0] rowPos=[ 0 2 0], colPos=[1 2] means to move from representation [x y] to [x y] [0 0] [0 0 x] rowPos=[ 1 0 0], colPos=[0 0 3] means to move from representation [x] to [0 0 0] [0 0 0] Note: rowPos and colPos has to be consistent with current representation (not checked here)

MatrixEntry &toReal(bool realPart = true)#

change type of entries cEntries-> rEntries or cmEntries -> rmEntries taking real or imag part

MatrixEntry *toScalar()#

create new scalar matrix entry from non scalar matrix entry

void toSkyline()#

convert current matrix storage to skyline storage (useful for factorization)

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

void umfluLeftSolve(const VectorEntry&, VectorEntry&)#

solve transposed linear system using umfpack LU factorization

void umfluSolve(const VectorEntry&, VectorEntry&)#

solve linear system using umfpack LU factorization

void umfpackFactorize()#

factorize matrix using umfpack

inline ValueType &valueType()#

return valueType (non const)

inline ValueType valueType() const#

return valueType (const)

void viewStorage(std::ostream&) const#

print storage on stream

Public Members

LargeMatrix<complex_t> *cEntries_p#

pointer to complex LargeMatrix

LargeMatrix<Matrix<complex_t>> *cmEntries_p#

pointer to LargeMatrix of complex matrices

dimPair nbOfComponents#

number of line and columns of matrix values, (1,1) in scalar case

LargeMatrix<real_t> *rEntries_p#

pointer to real LargeMatrix

LargeMatrix<Matrix<real_t>> *rmEntries_p#

pointer to LargeMatrix of real matrices

StrucType strucType_#

structure type of the entries (_scalar,_matrix)

ValueType valueType_#

value type of the entries (_real,_complex)

Friends

friend MatrixEntry operator*(const MatrixEntry&, const MatrixEntry&)#

product of MatrixEntry