Namespaces | |
detail | |
Classes | |
class | BlockCyclicMatrix |
Functions | |
template<typename Array > | |
BlockCyclicMatrix< typename std::remove_cv_t< Array >::element_type > | array_to_block_cyclic (const Array &array, const blacspp::Grid &grid, size_t MB, size_t NB) |
Convert a dense DistArray to block-cyclic storage format. More... | |
template<typename Array > | |
std::remove_cv_t< Array > | block_cyclic_to_array (const BlockCyclicMatrix< typename std::remove_cv_t< Array >::element_type > &matrix, const TiledRange &trange) |
Convert a block-cyclic matrix to DistArray. More... | |
template<typename Array > | |
auto | cholesky (const Array &A, TiledRange l_trange=TiledRange(), size_t NB=default_block_size()) |
Compute the Cholesky factorization of a HPD rank-2 tensor. More... | |
template<bool Both, typename Array > | |
auto | cholesky_linv (const Array &A, TiledRange l_trange=TiledRange(), size_t NB=default_block_size()) |
Compute the inverse of the Cholesky factor of an HPD rank-2 tensor. Optionally return the Cholesky factor itself. More... | |
template<typename Array > | |
auto | cholesky_solve (const Array &A, const Array &B, TiledRange x_trange=TiledRange(), size_t NB=default_block_size()) |
template<typename Array > | |
auto | cholesky_lsolve (Op trans, const Array &A, const Array &B, TiledRange l_trange=TiledRange(), TiledRange x_trange=TiledRange(), size_t NB=default_block_size()) |
template<typename Array > | |
auto | heig (const Array &A, TiledRange evec_trange=TiledRange(), size_t NB=default_block_size()) |
Solve the standard eigenvalue problem with ScaLAPACK. More... | |
template<typename ArrayA , typename ArrayB , typename EVecType = ArrayA> | |
auto | heig (const ArrayA &A, const ArrayB &B, TiledRange evec_trange=TiledRange(), size_t NB=default_block_size()) |
Solve the generalized eigenvalue problem with ScaLAPACK. More... | |
template<typename ArrayA , typename ArrayB > | |
auto | lu_solve (const ArrayA &A, const ArrayB &B, TiledRange x_trange=TiledRange(), size_t NB=default_block_size(), size_t MB=default_block_size()) |
Solve a linear system via LU factorization. More... | |
template<typename Array > | |
auto | lu_inv (const Array &A, TiledRange ainv_trange=TiledRange(), size_t NB=default_block_size(), size_t MB=default_block_size()) |
Invert a matrix via LU. More... | |
template<SVD::Vectors Vectors, typename Array > | |
auto | svd (const Array &A, TiledRange u_trange, TiledRange vt_trange, size_t MB=default_block_size(), size_t NB=default_block_size()) |
Compute the singular value decomposition (SVD) via ScaLAPACK. More... | |
scalapackpp::TransposeFlag | to_scalapackpp_transposeflag (Op t) |
template<typename T > | |
void | zero_triangle (blacspp::Triangle tri, scalapack::BlockCyclicMatrix< T > &A, bool zero_diag=false) |
std::size_t | default_block_size () |
void | set_default_block_size (std::size_t NB) |
Function Documentation
◆ array_to_block_cyclic()
BlockCyclicMatrix<typename std::remove_cv_t<Array>::element_type> TiledArray::math::linalg::scalapack::array_to_block_cyclic | ( | const Array & | array, |
const blacspp::Grid & | grid, | ||
size_t | MB, | ||
size_t | NB | ||
) |
Convert a dense DistArray to block-cyclic storage format.
- Template Parameters
-
T Datatype of underlying tile
- Parameters
-
[in] array DistArray to be converted to block-cyclic format [in] grid BLACS grid context for block-cyclic matrix [in] MB Row blocking factor of resulting block-cyclic matrix [in] NB Column blocking factor of resulting block-cyclic matrix
- Returns
- Block-cyclic conversion of input DistArray
Definition at line 308 of file block_cyclic.h.
◆ block_cyclic_to_array()
std::remove_cv_t<Array> TiledArray::math::linalg::scalapack::block_cyclic_to_array | ( | const BlockCyclicMatrix< typename std::remove_cv_t< Array >::element_type > & | matrix, |
const TiledRange & | trange | ||
) |
Convert a block-cyclic matrix to DistArray.
- Template Parameters
-
Datatype of underlying tile
- Parameters
-
[in] matrix Block-cyclic matrix to convert to DistArray [in] trange Tiled ranges for the resulting DistArray
- Returns
- DistArray conversion of input block-cyclic matrix
Definition at line 325 of file block_cyclic.h.
◆ cholesky()
auto TiledArray::math::linalg::scalapack::cholesky | ( | const Array & | A, |
TiledRange | l_trange = TiledRange() , |
||
size_t | NB = default_block_size() |
||
) |
Compute the Cholesky factorization of a HPD rank-2 tensor.
A(i,j) = L(i,k) * conj(L(j,k))
Example Usage:
auto L = cholesky(A, ...)
- Template Parameters
-
Array Input array type, must be convertible to BlockCyclicMatrix
- Parameters
-
[in] A Input array to be diagonalized. Must be rank-2 [in] l_trange TiledRange for resulting Cholesky factor. If left empty, will default to array.trange() [in] NB ScaLAPACK block size. Defaults to 128
- Returns
- The lower triangular Cholesky factor L in TA format
Definition at line 60 of file cholesky.h.
◆ cholesky_linv()
auto TiledArray::math::linalg::scalapack::cholesky_linv | ( | const Array & | A, |
TiledRange | l_trange = TiledRange() , |
||
size_t | NB = default_block_size() |
||
) |
Compute the inverse of the Cholesky factor of an HPD rank-2 tensor. Optionally return the Cholesky factor itself.
A(i,j) = L(i,k) * conj(L(j,k)) -> compute Linv
Example Usage:
auto Linv = cholesky_Linv(A, ...) auto [L,Linv] = cholesky_Linv<decltype(A),true>(A, ...)
- Template Parameters
-
Array Input array type, must be convertible to BlockCyclicMatrix Both Whether or not to return the cholesky factor
- Parameters
-
[in] A Input array to be diagonalized. Must be rank-2 [in] l_trange TiledRange for resulting inverse Cholesky factor. If left empty, will default to array.trange() [in] NB ScaLAPACK block size. Defaults to 128
- Returns
- The inverse lower triangular Cholesky factor in TA format
Definition at line 114 of file cholesky.h.
◆ cholesky_lsolve()
auto TiledArray::math::linalg::scalapack::cholesky_lsolve | ( | Op | trans, |
const Array & | A, | ||
const Array & | B, | ||
TiledRange | l_trange = TiledRange() , |
||
TiledRange | x_trange = TiledRange() , |
||
size_t | NB = default_block_size() |
||
) |
◆ cholesky_solve()
auto TiledArray::math::linalg::scalapack::cholesky_solve | ( | const Array & | A, |
const Array & | B, | ||
TiledRange | x_trange = TiledRange() , |
||
size_t | NB = default_block_size() |
||
) |
◆ default_block_size()
|
inline |
◆ heig() [1/2]
auto TiledArray::math::linalg::scalapack::heig | ( | const Array & | A, |
TiledRange | evec_trange = TiledRange() , |
||
size_t | NB = default_block_size() |
||
) |
Solve the standard eigenvalue problem with ScaLAPACK.
A(i,k) X(k,j) = X(i,j) E(j)
Example Usage:
auto [E, X] = heig(A, ...)
- Template Parameters
-
Array Input array type, must be convertible to BlockCyclicMatrix
- Parameters
-
[in] A Input array to be diagonalized. Must be rank-2 [in] evec_trange TiledRange for resulting eigenvectors. If left empty, will default to array.trange() [in] NB ScaLAPACK block size. Defaults to 128
- Returns
- A tuple containing the eigenvalues and eigenvectors of input array as std::vector and in TA format, respectively.
Definition at line 59 of file heig.h.
◆ heig() [2/2]
auto TiledArray::math::linalg::scalapack::heig | ( | const ArrayA & | A, |
const ArrayB & | B, | ||
TiledRange | evec_trange = TiledRange() , |
||
size_t | NB = default_block_size() |
||
) |
Solve the generalized eigenvalue problem with ScaLAPACK.
A(i,k) X(k,j) = B(i,k) X(k,j) E(j)
with
X(k,i) B(k,l) X(l,j) = I(i,j)
Example Usage:
auto [E, X] = heig(A, B, ...)
- Template Parameters
-
Array Input array type, must be convertible to BlockCyclicMatrix
- Parameters
-
[in] A Input array to be diagonalized. Must be rank-2 [in] B Metric [in] evec_trange TiledRange for resulting eigenvectors. If left empty, will default to array.trange() [in] NB ScaLAPACK block size. Defaults to 128
- Returns
- A tuple containing the eigenvalues and eigenvectors of input array as std::vector and in TA format, respectively.
Definition at line 122 of file heig.h.
◆ lu_inv()
auto TiledArray::math::linalg::scalapack::lu_inv | ( | const Array & | A, |
TiledRange | ainv_trange = TiledRange() , |
||
size_t | NB = default_block_size() , |
||
size_t | MB = default_block_size() |
||
) |
◆ lu_solve()
auto TiledArray::math::linalg::scalapack::lu_solve | ( | const ArrayA & | A, |
const ArrayB & | B, | ||
TiledRange | x_trange = TiledRange() , |
||
size_t | NB = default_block_size() , |
||
size_t | MB = default_block_size() |
||
) |
◆ set_default_block_size()
|
inline |
◆ svd()
auto TiledArray::math::linalg::scalapack::svd | ( | const Array & | A, |
TiledRange | u_trange, | ||
TiledRange | vt_trange, | ||
size_t | MB = default_block_size() , |
||
size_t | NB = default_block_size() |
||
) |
Compute the singular value decomposition (SVD) via ScaLAPACK.
A(i,j) = S(k) U(i,k) conj(V(j,k))
Example Usage:
auto S = svd<SVDValuesOnly> (A, ...) auto [S, U] = svd<SVDLeftVectors> (A, ...) auto [S, VT] = svd<SVDRightVectors>(A, ...) auto [S, U, VT] = svd<SVDAllVectors> (A, ...)
- Template Parameters
-
Array Input array type, must be convertible to BlockCyclicMatrix
- Parameters
-
[in] A Input array to be decomposed. Must be rank-2 [in] u_trange TiledRange for resulting left singular vectors. [in] vt_trange TiledRange for resulting right singular vectors (transposed). [in] MB ScaLAPACK row block size. Defaults to 128 [in] NB ScaLAPACK column block size. Defaults to 128
- Returns
- A tuple containing the eigenvalues and eigenvectors of input array as std::vector and in TA format, respectively.
Definition at line 63 of file svd.h.
◆ to_scalapackpp_transposeflag()
|
inline |
◆ zero_triangle()
void TiledArray::math::linalg::scalapack::zero_triangle | ( | blacspp::Triangle | tri, |
scalapack::BlockCyclicMatrix< T > & | A, | ||
bool | zero_diag = false |
||
) |