TiledArray::math::linalg::scalapack Namespace Reference

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< Arrayblock_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()

template<typename Array >
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
TDatatype of underlying tile
Parameters
[in]arrayDistArray to be converted to block-cyclic format
[in]gridBLACS grid context for block-cyclic matrix
[in]MBRow blocking factor of resulting block-cyclic matrix
[in]NBColumn blocking factor of resulting block-cyclic matrix
Returns
Block-cyclic conversion of input DistArray

Definition at line 308 of file block_cyclic.h.

Here is the caller graph for this function:

◆ block_cyclic_to_array()

template<typename 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
Datatypeof underlying tile
Parameters
[in]matrixBlock-cyclic matrix to convert to DistArray
[in]trangeTiled ranges for the resulting DistArray
Returns
DistArray conversion of input block-cyclic matrix

Definition at line 325 of file block_cyclic.h.

◆ cholesky()

template<typename Array >
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
ArrayInput array type, must be convertible to BlockCyclicMatrix
Parameters
[in]AInput array to be diagonalized. Must be rank-2
[in]l_trangeTiledRange for resulting Cholesky factor. If left empty, will default to array.trange()
[in]NBScaLAPACK block size. Defaults to 128
Returns
The lower triangular Cholesky factor L in TA format

Definition at line 60 of file cholesky.h.

Here is the call graph for this function:

◆ cholesky_linv()

template<bool Both, typename Array >
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
ArrayInput array type, must be convertible to BlockCyclicMatrix
BothWhether or not to return the cholesky factor
Parameters
[in]AInput array to be diagonalized. Must be rank-2
[in]l_trangeTiledRange for resulting inverse Cholesky factor. If left empty, will default to array.trange()
[in]NBScaLAPACK block size. Defaults to 128
Returns
The inverse lower triangular Cholesky factor in TA format

Definition at line 114 of file cholesky.h.

Here is the call graph for this function:

◆ cholesky_lsolve()

template<typename Array >
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() 
)

Definition at line 218 of file cholesky.h.

Here is the call graph for this function:

◆ cholesky_solve()

template<typename Array >
auto TiledArray::math::linalg::scalapack::cholesky_solve ( const Array A,
const Array B,
TiledRange  x_trange = TiledRange(),
size_t  NB = default_block_size() 
)

Definition at line 169 of file cholesky.h.

Here is the call graph for this function:

◆ default_block_size()

std::size_t TiledArray::math::linalg::scalapack::default_block_size ( )
inline

Definition at line 88 of file util.h.

Here is the call graph for this function:

◆ heig() [1/2]

template<typename Array >
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
ArrayInput array type, must be convertible to BlockCyclicMatrix
Parameters
[in]AInput array to be diagonalized. Must be rank-2
[in]evec_trangeTiledRange for resulting eigenvectors. If left empty, will default to array.trange()
[in]NBScaLAPACK 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.

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

◆ heig() [2/2]

template<typename ArrayA , typename ArrayB , typename EVecType = ArrayA>
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
ArrayInput array type, must be convertible to BlockCyclicMatrix
Parameters
[in]AInput array to be diagonalized. Must be rank-2
[in]BMetric
[in]evec_trangeTiledRange for resulting eigenvectors. If left empty, will default to array.trange()
[in]NBScaLAPACK 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.

Here is the call graph for this function:

◆ lu_inv()

template<typename Array >
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() 
)

Invert a matrix via LU.

Definition at line 90 of file lu.h.

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

◆ lu_solve()

template<typename ArrayA , typename ArrayB >
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() 
)

Solve a linear system via LU factorization.

Definition at line 43 of file lu.h.

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

◆ set_default_block_size()

void TiledArray::math::linalg::scalapack::set_default_block_size ( std::size_t  NB)
inline

Definition at line 92 of file util.h.

Here is the call graph for this function:

◆ svd()

template<SVD::Vectors Vectors, typename Array >
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
ArrayInput array type, must be convertible to BlockCyclicMatrix
Parameters
[in]AInput array to be decomposed. Must be rank-2
[in]u_trangeTiledRange for resulting left singular vectors.
[in]vt_trangeTiledRange for resulting right singular vectors (transposed).
[in]MBScaLAPACK row block size. Defaults to 128
[in]NBScaLAPACK 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.

Here is the call graph for this function:

◆ to_scalapackpp_transposeflag()

scalapackpp::TransposeFlag TiledArray::math::linalg::scalapack::to_scalapackpp_transposeflag ( Op  t)
inline

Definition at line 36 of file util.h.

Here is the caller graph for this function:

◆ zero_triangle()

template<typename T >
void TiledArray::math::linalg::scalapack::zero_triangle ( blacspp::Triangle  tri,
scalapack::BlockCyclicMatrix< T > &  A,
bool  zero_diag = false 
)

Definition at line 50 of file util.h.

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