imate
C++/CUDA Reference
|
A static class to find eigenvalues and eigenvectors (diagonalize) tridiagonal and bidiagonal matrices. This class acts as a templated namespace, where all member methods are public and static. More...
#include <diagonalization.h>
Static Public Member Functions | |
static int | eigh_tridiagonal (DataType *diagonals, DataType *subdiagonals, DataType *eigenvectors, IndexType matrix_size) |
Computes all eigenvalues and eigenvectors of a real and symmetric tridiagonal matrix. More... | |
static int | svd_bidiagonal (DataType *diagonals, DataType *subdiagonals, DataType *U, DataType *Vt, IndexType matrix_size) |
Computes all singular-values and left and right eigenvectors of a real and symmetric upper bi-diagonal matrix. More... | |
A static class to find eigenvalues and eigenvectors (diagonalize) tridiagonal and bidiagonal matrices. This class acts as a templated namespace, where all member methods are public and static.
Definition at line 35 of file diagonalization.h.
|
static |
Computes all eigenvalues and eigenvectors of a real and symmetric tridiagonal matrix.
The symmetric tridiagonal matrix \( \mathbf{A} \) is decomposed in to:
\[ \mathbf{A} = \mathbf{V} \boldsymbol{\Lambda} \mathbf{V}^{\intercal} \]
where \( \mathbf{V} \) is unitary and \( \boldsymbol{\Lambda} \) is diagonal.
This function is equivalent to scipy.linalg.eigh_tridigonal
which wraps around LAPACK's sstev
and dstev
subroutine. Except, this function does not acquire python's GIL, whereas the scipy's function does.
Note: Remove blank in the URLs below when opened in the browser.
[in,out] | diagonals | A 1D array of the length matrix_size containing the diagonal elements of the matrix. This array will be written in-place with the computed eigenvalues. This array is both input and output of this function. |
[in] | subdiagonals | A 1D array of the length matrix_size-1 containing the sub-diagonal elements of the matrix. This array will be written in-place with intermediate computation values, but is not an output of this function. |
[out] | eigenvectors | 1D array of size matrix_size*matrix_size and represents a 2D array of the shape (matrix_size,matrix_size). The second index of the equivalent 2D array iterates over the column vectors. The array has Fortran ordering, meaning that the first index is contiguous. Thus the i-th element of the j-th vector should be accessed by eigenvectors [j*matrix_size+i]. This array is written in-place and is the output of this function. |
[in] | matrix_size | The size of square matrix. |
info
result of the sstev
and dstev
subroutine. Zero indicates a successful computation. Definition at line 89 of file diagonalization.cpp.
References lapack_xstev().
Referenced by cTraceEstimator< DataType >::_c_stochastic_lanczos_quadrature(), and cuTraceEstimator< DataType >::_cu_stochastic_lanczos_quadrature().
|
static |
Computes all singular-values and left and right eigenvectors of a real and symmetric upper bi-diagonal matrix.
The symmetric upper bi-diagonal matrix \( \mathbf{A} \) is decomposed in to
\[ \mathbf{A} = \mathbf{U} \boldsymbol{\Sigma} \mathbf{V}^{\intercal} \]
where \( \mathbf{U} \) and \( \mathbf{V} \) are unitary and \( \boldsymbol{\Sigma} \) is diagonal with positive entries.
This function uses LAPACK's dbdsdc
Fortran subroutine.
Note: Remove blank in the URLs below when opened in the browser.
[in,out] | diagonals | A 1D array of the length matrix_size containing the diagonal elements of the matrix. This array will be written in-place with the computed eigenvalues. This array is both input and output of this function. |
[in] | supdiagonals | A 1D array of the length matrix_size-1 containing the sub-diagonal elements of the matrix. This array will be written in-place with intermediate computation values, but is not an output of this function. |
[out] | U | Right eigenvectors represented by 1D array of the length matrix_size*matrix_size which denotes a 2D array of the shape (matrix_size, matrix_size). The second index of the matrix iterates over the column vectors. The array has Fortran ordering, meaning that the first index is contiguous. Thus, the i-th element of the j-th vector should be accessed by U [j*matrix_size+i]. This array is written in place and is the output of this function. |
[out] | Vt | Transpose of left eigenvectors represented by 1D array of the length matrix_size*matrix_size which denotes a 2D array of the shape (matrix_size, matrix_size). The second index of the matrix iterates over the column vectors. The array has Fortran ordering, meaning that the first index is contiguous. Thus, the i-th element of the j-th vector should be accessed by Vt [j*matrix_size+i]. This array is written in place and is the output of this function. |
[in] | matrix_size | The size of square matrix. |
info
indicates the status of the computation. Zero indicates a successful computation. Definition at line 190 of file diagonalization.cpp.
References lapack_xbdsdc().
Referenced by cTraceEstimator< DataType >::_c_stochastic_lanczos_quadrature(), and cuTraceEstimator< DataType >::_cu_stochastic_lanczos_quadrature().