imate
C++/CUDA Reference
Diagonalization< DataType > Class Template 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...
 

Detailed Description

template<typename DataType>
class Diagonalization< DataType >

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.

See also
StochasticLanczosQuadrature

Definition at line 35 of file diagonalization.h.

Member Function Documentation

◆ eigh_tridiagonal()

template<typename DataType >
int Diagonalization< DataType >::eigh_tridiagonal ( DataType *  diagonals,
DataType *  subdiagonals,
DataType *  eigenvectors,
IndexType  matrix_size 
)
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.

Algorithm

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.

References

Note: Remove blank in the URLs below when opened in the browser.

Parameters
[in,out]diagonalsA 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]subdiagonalsA 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]eigenvectors1D 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_sizeThe size of square matrix.
Returns
The info result of the sstev and dstev subroutine. Zero indicates a successful computation.

Definition at line 89 of file diagonalization.cpp.

94 {
95  // 'V' computes both eigenvalues and eigenvectors
96  char jobz = 'V';
97 
98  // Workspace array
99  DataType* work = new DataType[2*matrix_size - 2];
100 
101  // Leading dimension of 2D eigenvectors array
102  int ldz = matrix_size;
103 
104  // matrix size
105  int n = static_cast<int>(matrix_size);
106 
107  // Error code output
108  int info;
109 
110  // Calling Fortran subroutine
111  lapack_xstev(&jobz, &n, diagonals, subdiagonals, eigenvectors, &ldz, work,
112  &info);
113 
114  delete[] work;
115  assert((info == 0, "?stev subroutine returned non-zero status."));
116 
117  return info;
118 }
void lapack_xstev(char *jobz, int *n, DataType *d, DataType *e, DataType *z, int *ldz, DataType *work, int *info)

References lapack_xstev().

Referenced by cTraceEstimator< DataType >::_c_stochastic_lanczos_quadrature(), and cuTraceEstimator< DataType >::_cu_stochastic_lanczos_quadrature().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ svd_bidiagonal()

template<typename DataType >
int Diagonalization< DataType >::svd_bidiagonal ( DataType *  diagonals,
DataType *  supdiagonals,
DataType *  U,
DataType *  Vt,
IndexType  matrix_size 
)
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.

Algorithm

This function uses LAPACK's dbdsdc Fortran subroutine.

References

Note: Remove blank in the URLs below when opened in the browser.

Parameters
[in,out]diagonalsA 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]supdiagonalsA 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]URight 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]VtTranspose 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_sizeThe size of square matrix.
Returns
Integer info indicates the status of the computation. Zero indicates a successful computation.

Definition at line 190 of file diagonalization.cpp.

196 {
197  // Code 'U' indicates the matrix is upper bi-diagonal
198  char UPLO = 'U';
199 
200  // Code 'I' computes both singular values and singular vectors
201  char COMPQ = 'I';
202 
203  // matrix size
204  int n = static_cast<int>(matrix_size);
205 
206  // Leading dimensions of arrays U and Vt
207  int LDU = matrix_size;
208  int LDVT = matrix_size;
209 
210  // There arrays are not referenced when COMPQ is set to 'I'
211  DataType* Q = NULL;
212  int* IQ = NULL;
213 
214  // Work arrays
215  DataType* work = new DataType[3*matrix_size*matrix_size + 4*matrix_size];
216  int* iwork = new int[8 * matrix_size];
217 
218  // Error code output
219  int info;
220 
221  // Calling Fortran subroutine
222  lapack_xbdsdc(&UPLO, &COMPQ, &n, diagonals, supdiagonals, U, &LDU, Vt,
223  &LDVT, Q, IQ, work, iwork, &info);
224 
225  delete[] work;
226  delete[] iwork;
227 
228  assert((info == 0, "?stev subroutine returned non-zero status."));
229 
230  return info;
231 }
void lapack_xbdsdc(char *uplo, char *compq, int *n, DataType *d, DataType *e, DataType *u, int *ldu, DataType *vt, int *ldvt, DataType *q, int *iq, DataType *work, int *iwork, int *info)

References lapack_xbdsdc().

Referenced by cTraceEstimator< DataType >::_c_stochastic_lanczos_quadrature(), and cuTraceEstimator< DataType >::_cu_stochastic_lanczos_quadrature().

Here is the call graph for this function:
Here is the caller graph for this function:

The documentation for this class was generated from the following files: