Functions | |
| template<typename LhsScalar, typename RhsScalar, int KcFactor, typename SizeType> | |
| void | computeProductBlockingSizes (SizeType &k, SizeType &m, SizeType &n) |
| Computes the blocking parameters for a m x k times k x n matrix product. | |
| template<typename Scalar, int Flags, typename Index> | |
| MappedSparseMatrix< Scalar, Flags, Index > | map_superlu (SluMatrix &sluMat) |
| template<typename MatrixType, typename DiagonalType, typename SubDiagonalType> | |
| void | tridiagonalization_inplace (MatrixType &mat, DiagonalType &diag, SubDiagonalType &subdiag, bool extractQ) |
| Performs a full tridiagonalization in place. | |
This is defined in the Jacobi module.
Applies the clock wise 2D rotation j to the set of 2D vectors of cordinates x and y: 
| void computeProductBlockingSizes | ( | SizeType & | k, |
| SizeType & | m, | ||
| SizeType & | n ) |
Computes the blocking parameters for a m x k times k x n matrix product.
| [in,out] | k | Input: the third dimension of the product. Output: the blocking size along the same dimension. |
| [in,out] | m | Input: the number of rows of the left hand side. Output: the blocking size along the same dimension. |
| [in,out] | n | Input: the number of columns of the right hand side. Output: the blocking size along the same dimension. |
Given a m x k times k x n matrix product of scalar types LhsScalar and RhsScalar, this function computes the blocking size parameters along the respective dimensions for matrix products and related algorithms. The blocking sizes depends on various parameters:
| MappedSparseMatrix< Scalar, Flags, Index > map_superlu | ( | SluMatrix & | sluMat | ) |
View a Super LU matrix as an Eigen expression
References Eigen::ColMajor, and Eigen::RowMajor.
| void tridiagonalization_inplace | ( | MatrixType & | mat, |
| DiagonalType & | diag, | ||
| SubDiagonalType & | subdiag, | ||
| bool | extractQ ) |
Performs a full tridiagonalization in place.
| [in,out] | mat | On input, the selfadjoint matrix whose tridiagonal decomposition is to be computed. Only the lower triangular part referenced. The rest is left unchanged. On output, the orthogonal matrix Q in the decomposition if extractQ is true. |
| [out] | diag | The diagonal of the tridiagonal matrix T in the decomposition. |
| [out] | subdiag | The subdiagonal of the tridiagonal matrix T in the decomposition. |
| [in] | extractQ | If true, the orthogonal matrix Q in the decomposition is computed and stored in mat. |
Computes the tridiagonal decomposition of the selfadjoint matrix mat in place such that 


The tridiagonal matrix T is passed to the output parameters diag and subdiag. If extractQ is true, then the orthogonal matrix Q is passed to mat. Otherwise the lower part of the matrix mat is destroyed.
The vectors diag and subdiag are not resized. The function assumes that they are already of the correct size. The length of the vector diag should equal the number of rows in mat, and the length of the vector subdiag should be one left.
This implementation contains an optimized path for 3-by-3 matrices which is especially useful for plane fitting.
Example (this uses the same matrix as the example in Tridiagonalization::Tridiagonalization(const MatrixType&)):
Output:
Here is a random symmetric 5x5 matrix:
1.36 -0.816 0.521 1.43 -0.144
-0.816 -0.659 0.794 -0.173 -0.406
0.521 0.794 -0.541 0.461 0.179
1.43 -0.173 0.461 -1.43 0.822
-0.144 -0.406 0.179 0.822 -1.37
The orthogonal matrix Q is:
1 0 0 0 0
0 -0.471 0.127 -0.671 -0.558
0 0.301 -0.195 0.437 -0.825
0 0.825 0.0459 -0.563 -0.00872
0 -0.0832 -0.971 -0.202 0.0922
The diagonal of the tridiagonal matrix T is:
1.36
-1.2
-1.28
-1.69
0.164
The subdiagonal of the tridiagonal matrix T is:
1.73
-0.966
0.214
0.345