Class xlifepp::MatrixEigenDense#
-
template<typename K>
class MatrixEigenDense : public xlifepp::Matrix<K>#
-
Inheritence diagram for xlifepp::MatrixEigenDense:
Collaboration diagram for xlifepp::MatrixEigenDense:
Class (row) dense matrix to deal with direct eigen problem solver.
This class provides some basic methods to support calculation of eigen problem dense matrix.
- Template Parameters:
-
Type – of scalar (real_t, complex_t)
Public Functions
-
inline MatrixEigenDense()#
-
default constructor
-
inline MatrixEigenDense(const dimen_t r, const dimen_t c, const K v)#
-
ocnstructor from dimensions and value
-
inline MatrixEigenDense(const MatrixEigenDense<K> &mat, const int_t rIdx, const int_t cIdx, const dimen_t rSize, const dimen_t cSize)#
-
Construct matrix from an existing matrix (Extract a sub-matrix from an existing one) This is a very simple version which uses copy operation, and needs optimizing.
- Parameters:
-
mat – [in] Existing matrix from which we extract sub-matrix
rIdx – [in] Row position (top left corner) of matrix at which sub-matrix is extracted
cIdx – [in] Column position (top left corner) of matrix at which sub-matrix is extracted
rSize – [in] Number of row of sub-matrix
cSize – [in] Number of column sub-matrix
-
inline void addAssignCol(const dimen_t c, const VectorEigenDense<K> &col, K scale)#
-
Computes the addition assignment of column c of matrix A and (column) vector X column(c) += scale * X Content of the column will be replaced with a new one!!!!
- Parameters:
-
c – [in] column index
col – [in] (column) vector X
scale – [in] scale factor
-
inline void addAssignRow(const dimen_t r, const VectorEigenDense<K> &row, K scale)#
-
Computes the addition assignment a row r of matrix A and the (row) vector X row(r) += scale * X Content of the row will be replaced with a new one!!!!
- Parameters:
-
r – [in] row index
row – [in] (row) vector X
scale – [in] scale factor
-
inline void applyHouseholderOnTheLeft(const VectorEigenDense<K> &essential, const K &tau)#
-
Apply the elementary reflector H given by \( H = I - tau v v^*\) with \( v^T = [1 essential^T] \) from the left to a matrix.
On input:
- Parameters:
-
essential – the essential part of the vector
v
tau – the scaling factor of the Householder transformation
-
inline void applyHouseholderOnTheRight(const VectorEigenDense<K> &essential, const K &tau)#
-
Apply the elementary reflector H given by \( H = I - tau v v^*\) with \( v^T = [1 essential^T] \) from the right to a matrix.
On input:
- Parameters:
-
essential – the essential part of the vector
v
tau – the scaling factor of the Householder transformation
-
template<typename OtherScalar>
inline void applyOnTheLeft(Index p, Index q, const JacobiRotation<OtherScalar> &j)#
-
Applies the rotation in the plane j to the rows p and q of
*this
, i.e., it computes B = J * B, with \( B = \left ( \begin{array}{cc} \mbox{*this.row(p)} \\ \mbox{*this.row(q)} \end{array} \right ) \).
-
template<typename OtherScalar>
inline void applyOnTheRight(Index p, Index q, const JacobiRotation<OtherScalar> &j)#
-
Applies the rotation in the plane j to the columns p and q of
*this
, i.e., it computes B = B * J with \( B = \left ( \begin{array}{cc} \mbox{*this.col}(p) & \mbox{*this.col}(q) \end{array} \right ) \).
-
inline VectorEigenDense<K> blockCol(int_t rowIdx, int_t colIdx, int_t size) const#
-
Extract a block of a column of a matrix.
- Parameters:
-
rowIdx – row index (beginning position of block)
colIdx – column index (beginning position of block)
size – size of (column) block
- Returns:
-
(column) vector
-
inline VectorEigenDense<K> blockRow(int_t rowIdx, int_t colIdx, int_t size) const#
-
Extract a block of a row of a matrix.
- Parameters:
-
rowIdx – row index (beginning position of block)
colIdx – column index (beginning position of block)
size – size of (row) block
- Returns:
-
(row) vector
-
inline MatrixEigenDense bottomRightCorner(const Index rowSize, const Index colSize)#
-
Extract the bottom right corner sub-matrix.
- Parameters:
-
rowSize – [in] number of row of sub-matrix
colSize – [in] number of row of sub-matrix
- Returns:
-
sub-matrix
-
inline void bottomRightCorner(const Index rowSize, const Index colSize, MatrixEigenDense<K> &subMat)#
-
Replace the bottom right corner of a matrix with a matrix.
- Parameters:
-
rowSize – [in] number of row of sub-matrix
colSize – [in] number of row of sub-matrix
subMat – [in] matrix to replace
-
inline MatrixEigenDense<K> cholesky() const#
-
Standard Cholesky decomposition (LL^T) of a matrix and associated features.
This function performs a LL^T Cholesky decomposition of a symmetric, positive definite matrix A such that A = LL^* = U^*U, where L is lower triangular.
While the Cholesky decomposition is particularly useful to solve selfadjoint problems like D^*D x = b, Nevertheless, this standard Cholesky decomposition remains useful in many other situations like generalised eigen problems with hermitian matrices.
Remember that Cholesky decompositions are not rank-revealing. This LLT decomposition is only stable on positive definite matrices, Also, do not use a Cholesky decomposition to determine whether a system of equations has a solution.
- Returns:
-
Lower triangular matrix
-
inline K coeff(const number_t rowIdx, const number_t colIdx) const#
-
Return coefficient of a matrix.
- Parameters:
-
rowIdx – row index of coefficient
colIdx – column index of coefficient
- Returns:
-
coefficient
-
inline K &coeffRef(const number_t rowIdx, const number_t colIdx)#
-
Return reference of a coefficient of matrix.
- Parameters:
-
rowIdx – row index of coefficient
colIdx – column index of coefficient
- Returns:
-
reference of coefficient
-
inline VectorEigenDense<K> columnVector(const dimen_t c) const#
-
Return a column of matrix (column index begins from 0)
- Parameters:
-
c – [in] column to return
- Returns:
-
c-th column of matrix
-
inline void columnVector(const dimen_t c, VectorEigenDense<K> &v)#
-
Set a column of matrix with a column vector.
- Parameters:
-
c – [in] column to be set
v – [in] column vector
-
inline VectorEigenDense<K> diagonal() const#
-
Return diagonal of a matrix in a column vector.
-
inline real_t maxCoeff() const#
-
Return maximum absolute coefficient of a matrix.
-
inline real_t minCoeff() const#
-
Return minimum absolute coefficient of a matrix.
-
template<typename T>
inline void multSubMatVecVec(const Indexing &idx, const VectorEigenDense<T> &vec, VectorEigenDense<T> &res)#
-
Computes the product of a sub-matrix and a vector with: Matrix A | B C | D Vector X, Y Product: D x X = Y.
- Parameters:
-
idx – [in] index specify block of sub-matrix inside the main matrix. idx[0]: row position of top-left point idx[1]: col position of top-left point idx[2]: number of row of sub-matrix idx[3]: number of column of sub-matrix
vec – [in] vector X
res – [out] vector Y: the result of D * X
-
template<typename T>
inline void multVecSubMatVec(const Indexing &idx, const VectorEigenDense<T> &vec, VectorEigenDense<T> &res)#
-
Computes the product of a vector and sub-matrix with: Matrix A | B C | D Vector X, Y Product: X * D = Y.
- Parameters:
-
idx – [in] index specify block of sub-matrix inside the main matrix. idx[0]: row position of top-left point idx[1]: column position of top-left point idx[2]: number of row of sub-matrix idx[3]: number of column of sub-matrix
vec – [in] vector X
res – [out] vector Y: the result of X * D
-
inline void multVecVecSubMat(const Indexing &idx, const VectorEigenDense<K> &vec1, const VectorEigenDense<K> &vec2)#
-
Computes the product of a vector and a vector, result is returned into a sub-matrix with: Matrix A | B C | D Vector X, Y Product: D= X * Y * scale The content of sub-matrix will be replaced with the new one!!!!
- Parameters:
-
idx – [in] index specify block of sub-matrix inside the main matrix. idx[0]: row position of top-left point idx[1]: column position of top-left point idx[2]: number of row of sub-matrix idx[3]: number of column of sub-matrix
vec1 – [in] vector X
vec2 – [in] vector Y
-
inline void multVecVecSubMatAdditionAssign(const Indexing &idx, const VectorEigenDense<K> &vec1, const VectorEigenDense<K> &vec2, const K &scale)#
-
Computes the product of a vector and a vector, then addition assignment and the result is returned into a sub-matrix with: Matrix A | B C | D Vector X, Y Product: D += X * Y * scale The content of sub-matrix will be replaced with a new one!!!!
- Parameters:
-
idx – [in] index specify block of sub-matrix inside the main matrix. idx[0]: row position of top-left point idx[1]: column position of top-left point idx[2]: number of row of sub-matrix idx[3]: number of column of sub-matrix
vec1 – [in] vector X
vec2 – [in] vector Y
scale – [in] scale factor
-
inline MatrixEigenDense<K> &operator=(const MatrixEigenDense<K> &mat)#
-
assignment operator
-
inline void rankUpdate(const Indexing &idx, const VectorEigenDense<K> &vec1, const VectorEigenDense<K> &vec2, const K &scale)#
-
Computes the submatrix += scale * vec1*vec2^T + conj(scale)*vec2*vec1^T Matrix A | B C | D Vector X, Y Product: D += scaleX * X * Y^T + conjugate(scale) * Y * X^T The content of sub-matrix will be replaced with a new one!!!! This implementation can be optimized!!!
- Parameters:
-
idx – [in] index specify block of sub-matrix inside the main matrix. idx[0]: row position of top-left point idx[1]: column position of top-left point idx[2]: number of row of sub-matrix idx[3]: number of column of sub-matrix
vec1 – [in] colum vector X
vec2 – [in] colum vector Y
scale – [in] scale factor
-
inline void replace(const MatrixEigenDense<K> &mat, const int_t rIdx, const int_t cIdx, const int_t rSize, const int_t cSize)#
-
Replace content of a block of a matrix with another matrix.
- Parameters:
-
mat – [in] replacing matrix
rIdx – [in] row top left corner index
cIdx – [in] column top left corner index
rSize – [in] row of replaced block (and replacing matrix)
cSize – [in] column of replaced block (and replacing matrix)
-
inline VectorEigenDense<K> rowVector(const dimen_t r) const#
-
Return a row of matrix (row index begins from 0)
- Parameters:
-
r – [in] row to return
- Returns:
-
r-th row of matrix
-
inline void rowVector(const dimen_t r, VectorEigenDense<K> &v)#
-
Set a row of matrix (beginning from 0) with a row vector.
- Parameters:
-
r – [in] row to be set
v – [in] row vector
-
template<typename T>
inline MatrixEigenDense<K> &scale(const T s) const#
-
Scale a matrix.
- Parameters:
-
s – scaling factor
- Returns:
-
scaled matrix TODO need to do with round-off error. Limit using this function!!!
-
inline void setZeroLowerTriangle()#
-
Zeroing the lower part of a matrix TODO Change to large matrix to have better performance?
-
inline void setZeroStrictUpper()#
-
Zeroing the strict upper of a matrix TODO Change to large matrix to have better performance?
-
inline void sizeMisMatch(const string_t &s, const size_t r, const size_t c) const#
-
error: in matrix operation s dimensions disagree dim(op1)= … dim(op2)=…
-
template<typename T>
inline void solveCholeskyInplaceLower(const MatrixEigenDense<T> &mat, MatrixEigenDense<K> &res)#
-
Linear solver with Cholesky decomposed matrix Solve linear equation: Lx = b with L: lower triangular matrix !!!This function should only be called after the function cholesky()
- Parameters:
-
mat – [in] right-hand matrix (b)
res – [out] result matrix (x)
-
template<typename T>
inline void solveCholeskyInplaceUpper(const MatrixEigenDense<T> &mat, MatrixEigenDense<K> &res)#
-
Linear solver with Cholesky decomposed matrix Solve linear equation: XU = b with U: upper triangular matrix !!!This function should only be called after the function cholesky()
- Parameters:
-
mat – [in] right-hand matrix (b)
res – [out] result matrix (x)
-
inline VectorEigenDense<K> subDiagonal() const#
-
Return (lower) sub-diagonal of a matrix in a column vector.