Class xlifepp::LowRankMatrix#
-
template<class T>
class LowRankMatrix : public xlifepp::ApproximateMatrix<T>#
-
Inheritence diagram for xlifepp::LowRankMatrix:
Collaboration diagram for xlifepp::LowRankMatrix:
class dealing with matrix represented by a small sets of vectors A = U D V* with U a (m,r) matrix, V a (n,r) matrix and D a (r,r) diagonal matrix stored as a vector when D=Id it is not stored
Public Functions
-
LowRankMatrix()#
-
default constructor
-
LowRankMatrix(const LargeMatrix<T>&, HMApproximationMethod = _r3svdCompression, number_t = 0, real_t = theTolerance)#
-
constructor from LargeMatrix to be compressed
-
LowRankMatrix(const Matrix<T>&, const Matrix<T>&, const string_t &na = "")#
-
constructor from matrices, no diag
-
LowRankMatrix(const Matrix<T>&, const Matrix<T>&, const Vector<T>&, const string_t &na = "")#
-
constructor from matrices
-
LowRankMatrix(number_t m, number_t n, number_t r, bool noDiag, const string_t &na = "")#
-
constructor from dimensions, no diag
-
LowRankMatrix(number_t m, number_t n, number_t r, const string_t &na = "")#
-
constructor from dimensions
-
template<typename ITM>
LowRankMatrix(number_t m, number_t n, number_t r, ITM itmu, ITM itmv, const string_t &na = "")#
-
constructor from matrice iterator, no diag
constructor from matrice iterator
-
template<typename ITM, typename ITV>
LowRankMatrix(number_t m, number_t n, number_t r, ITM itmu, ITM itmv, ITV itd, const string_t &na = "")#
-
constructor from matrice iterator
-
LowRankMatrix<T> &add(const LowRankMatrix<T>&, const T&, real_t = 0)#
-
add a scaled LowRank matrix to current one
-
virtual ApproximateMatrix<T> *clone() const#
-
create a clone
-
LowRankMatrix<T> &compress(HMApproximationMethod cm, number_t r = 0, real_t eps = theTolerance)#
-
compression of the LowRankMatrix
compression of LowRankMatrix cm: compression method one of _svdCompression, _rsvdCompression, _r3svdCompression, _acaFull, _acaPartial, _acaPlus r: prescribed rank eps: prescribed precision if r>0 use prescribed rank else use prescribed precision
-
virtual void extend(const std::vector<number_t> &rowIndex, const std::vector<number_t> &colIndex, number_t m = 0, number_t n = 0)#
-
extend LowRankMatrix to larger numbering
extend LowRankMatrix to row and col parent numberings let I, J the row and col numbering of current LowRankMatrix (UI*D*VJ) related to parent numberings the new LowRankMatrix is
|UI| |VJ|' |UI*D*VJ' 0 | | | |D| | | = | | |0 | |0 | | 0 0 |
rowIndex: indices of rows starting at 1 colIndex: indices of cols starting at 1 m: number of rows of extended matrix n: number of cols of extended matrix if m=0 it is deduced from the rowIndex vector: m = max rowIndex if n=0 it is deduced from the colIndex vector: n = max colIndex
-
virtual LowRankMatrix<T> *extract(const std::vector<number_t> &rowIndex, const std::vector<number_t> &colIndex) const#
-
extract part of ApproximateMatrix (virtual)
extract part of LowRankMatrix given by row indices I and col indices J if low rank matrix has the form |UI| |VJ|’ |UI*D*VJ’ UI*D*VL’| | | |D| | | = | | |UK| |VL| |UK*D*VJ’ UK*D*VL’| so extract (I, J) part consists in building the new low rank matrix UI*D*VJ’ in other words, extract rows I of U and rows J of V
-
void fromLargeMatrix(const LargeMatrix<T>&)#
-
build LowRankMatrix from LargeMatrix
build low rank matrix from LargeMatrix
-
virtual void luFactorize(bool withPermutation = false)#
-
LU factorization.
-
virtual void multLeftMatrixCol(T *M, T *R, number_t p) const#
-
left product by a col dense matrix (matCol*L)-> col matrix
do the product R = M * L where L is the current LowRankMatrix, L=U*D*Vt, U (m,r) matrix, V (n,r) matrix stored as dense row, D (r,r) diagonal matrix stored as vector M a dense col matrix given by its pointer to first value R a dense col matrix given by its pointer to first value p the number of rows of M
-
virtual void multLeftMatrixRow(T *M, T *R, number_t p) const#
-
left product by a row dense matrix (matRow*L)-> row matrix
do the product R = M * L where L is the current LowRankMatrix, L=U*D*Vt, U (m,r) matrix, V (n,r) matrix stored as dense row, D (r,r) diagonal matrix stored as vector M a dense col matrix given by its pointer to first value R a dense col matrix given by its pointer to first value p the number of rows of M
-
virtual void multMatrixCol(T *M, T *R, number_t p) const#
-
right product by a col dense matrix (L*matCol)-> col matrix
do the product R = L * M where L is the current LowRankMatrix, L=U*D*Vt, U (m,r) matrix, V (n,r) matrix stored as dense row, D (r,r) diagonal matrix stored as vector M a dense col matrix given by its pointer to first value R a dense col matrix given by its pointer to first value p the number of cols of M
-
virtual void multMatrixRow(T *M, T *R, number_t p) const#
-
right product by a row dense matrix (L*matRow)-> row matrix
do the product R = L * M where L is the current LowRankMatrix, L=U*D*Vt, U (m,r) matrix, V (n,r) matrix stored as dense row, D (r,r) diagonal matrix stored as vector M a dense row matrix given by its pointer to first value R a dense row matrix given by its pointer to first value p the number of cols of M
-
virtual std::vector<T> &multMatrixVector(const std::vector<T>&, std::vector<T>&) const#
-
virtual product matrix * vector
-
virtual std::vector<T> &multVectorMatrix(const std::vector<T>&, std::vector<T>&) const#
-
virtual product vector * matrix
-
virtual real_t norminfty() const#
-
infinite norm
-
LowRankMatrix<T> &operator*=(const T&)#
-
multiply by a scalar
-
LowRankMatrix<T> &operator+=(const LowRankMatrix<T>&)#
-
add a LowRank matrix to current one
-
LowRankMatrix<T> &operator-=(const LowRankMatrix<T>&)#
-
substract a LowRank matrix to current one
-
LowRankMatrix<T> &operator/=(const T&)#
-
divide by a scalar
-
virtual void print(std::ostream&) const#
-
print to stream
-
LowRankMatrix<T> &r3svd(real_t eps = theTolerance, number_t t = 0, number_t p = 0, number_t q = 0, number_t maxit = 0)#
-
r3svd of a low rank matrice
faster random svd for low rank matrix returns the approximate svd USV* as a low rank matrix performing a truncation of singular values using an energy criteria eps: prescribed precision t: number of sampling p: number of over sampling q: power number (q=0 gives the standard r3svd) maxit: maximum of iterations
-
virtual void restrict(const std::vector<number_t> &rowIndex, const std::vector<number_t> &colIndex)#
-
restrict LowRankMatrix to smaller numbering
restrict a LowRankMatrix to row indices I and col indices J if low rank matrix has the form |UI| |VJ|’ |UI*D*VJ’ UI*D*VL’| | | |D| | | = | | |UK| |VL| |UK*D*VJ’ UK*D*VL’| so extract (I, J) part consists in building the new low rank matrix UI*D*VJ’
-
LowRankMatrix<T> &rsvd(number_t r = 0, real_t eps = theTolerance)#
-
rsvd of a low rank matrice
random svd for low rank matrix returns the approximate svd USV* as a low rank matrix performing a truncation of singular values at a given rank or below eps r: prescribed rank eps: prescribed precision if r>0 use prescribed rank else use prescribed precision
-
virtual real_t squaredNorm() const#
-
squared Frobenius norm
-
LowRankMatrix<T> &svd(number_t r = 0, real_t eps = theTolerance)#
-
svd of a low rank matrice
specific svd for low rank matrix returns the svd L=USV* as a low rank matrix performing a truncation of singular values below eps, r: prescribed rank eps: prescribed precision if r>0 use prescribed rank else use prescribed precision
convert to and use svd of class
-
template<typename I>
void toDenseRow(I itm) const#
-
convert to dense row matrix passing iterator on the first matrix value
convert to a dense row matrix passed by iterator itm pointing to the first value of dense row matrix
-
LowRankMatrix<T> &toFirsts(number_t r)#
-
restrict U_, V_ to their r first cols leading to a LowRankMatrix of rank r
restrict LowRankMatrix vector sets to their r first cols leading to a LowRankMatrix of rank r
-
virtual LargeMatrix<T> toLargeMatrix(StorageType st = _dense, AccessType = _row) const#
-
convert to LargeMatrix
convert to LargeMatrix (dense row)
Public Members
-
HMApproximationMethod cm_#
-
compression method used when re-compress
-
real_t eps_#
-
prescribed precision
-
FactorizationType factorization_#
-
one of _noFactorization, _lu, _ldlt, _ldlstar, _llt, _llstar, _qr, _umfpack
-
LowRankMatrix()#