TiledArray::math::linalg::non_distributed Namespace Reference

Functions

template<typename Tile , typename Policy >
auto rank_local_cholesky (const DistArray< Tile, Policy > &A)
 
template<typename Array , typename = std::enable_if_t<TiledArray::detail::is_array_v<Array>>>
auto cholesky (const Array &A, TiledRange l_trange=TiledRange())
 Compute the Cholesky factorization of a HPD rank-2 tensor. More...
 
template<typename ContiguousTensor , typename = std::enable_if_t< TiledArray::detail::is_contiguous_tensor_v<ContiguousTensor>>>
auto cholesky (const ContiguousTensor &A)
 Compute the Cholesky factorization of a HPD rank-2 tensor. More...
 
template<bool Both, typename Array , typename = std::enable_if_t<TiledArray::detail::is_array_v<Array>>>
auto cholesky_linv (const Array &A, TiledRange l_trange=TiledRange())
 Compute the inverse of the Cholesky factor of an HPD rank-2 tensor. Optionally return the Cholesky factor itself. More...
 
template<typename Array , typename = std::enable_if_t<TiledArray::detail::is_array_v<Array>>>
auto cholesky_solve (const Array &A, const Array &B, TiledRange x_trange=TiledRange())
 
template<typename Array , typename = std::enable_if_t<TiledArray::detail::is_array_v<Array>>>
auto cholesky_lsolve (Op transpose, const Array &A, const Array &B, TiledRange l_trange=TiledRange(), TiledRange x_trange=TiledRange())
 
template<typename Array >
auto heig (const Array &A, TiledRange evec_trange=TiledRange())
 Solve the standard eigenvalue problem with LAPACK. More...
 
template<typename ArrayA , typename ArrayB , typename EVecType = ArrayA>
auto heig (const ArrayA &A, const ArrayB &B, TiledRange evec_trange=TiledRange())
 Solve the generalized eigenvalue problem with LAPACK. More...
 
template<typename ArrayA , typename ArrayB >
auto lu_solve (const ArrayA &A, const ArrayB &B, TiledRange x_trange=TiledRange())
 Solve a linear system via LU factorization. More...
 
template<typename Array >
auto lu_inv (const Array &A, TiledRange ainv_trange=TiledRange())
 Invert a matrix via LU. More...
 
template<SVD::Vectors Vectors, typename Array >
auto svd (const Array &A, TiledRange u_trange=TiledRange(), TiledRange vt_trange=TiledRange())
 Compute the singular value decomposition (SVD) via ScaLAPACK. More...
 

Function Documentation

◆ cholesky() [1/2]

template<typename Array , typename = std::enable_if_t<TiledArray::detail::is_array_v<Array>>>
auto TiledArray::math::linalg::non_distributed::cholesky ( const Array A,
TiledRange  l_trange = TiledRange() 
)

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
Arraya DistArray type (i.e., is_array_v<Array> is true)
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()
Returns
The lower triangular Cholesky factor L in TA format
Note
this is a collective operation with respect to the world of A

Definition at line 72 of file cholesky.h.

Here is the call graph for this function:

◆ cholesky() [2/2]

template<typename ContiguousTensor , typename = std::enable_if_t< TiledArray::detail::is_contiguous_tensor_v<ContiguousTensor>>>
auto TiledArray::math::linalg::non_distributed::cholesky ( const ContiguousTensor &  A)

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
ContiguousTensora contiguous tensor type (i.e., is_contiguous_tensor_v<ContiguousTensor> is true)
Parameters
[in]AInput array to be diagonalized. Must be rank-2
Returns
The lower triangular Cholesky factor L as a ContiguousTensor
Note
this is a non-collective operation, only computes on the rank on which invoked

Definition at line 99 of file cholesky.h.

Here is the call graph for this function:

◆ cholesky_linv()

template<bool Both, typename Array , typename = std::enable_if_t<TiledArray::detail::is_array_v<Array>>>
auto TiledArray::math::linalg::non_distributed::cholesky_linv ( const Array A,
TiledRange  l_trange = TiledRange() 
)

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
Arraya DistArray type (i.e., is_array_v<Array> is true)
RetLWhether 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()
Returns
The inverse lower triangular Cholesky factor in TA format
Note
this is a collective operation with respect to the world of A

Definition at line 129 of file cholesky.h.

Here is the call graph for this function:

◆ cholesky_lsolve()

template<typename Array , typename = std::enable_if_t<TiledArray::detail::is_array_v<Array>>>
auto TiledArray::math::linalg::non_distributed::cholesky_lsolve ( Op  transpose,
const Array A,
const Array B,
TiledRange  l_trange = TiledRange(),
TiledRange  x_trange = TiledRange() 
)

Definition at line 176 of file cholesky.h.

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

◆ cholesky_solve()

template<typename Array , typename = std::enable_if_t<TiledArray::detail::is_array_v<Array>>>
auto TiledArray::math::linalg::non_distributed::cholesky_solve ( const Array A,
const Array B,
TiledRange  x_trange = TiledRange() 
)

Definition at line 156 of file cholesky.h.

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

◆ heig() [1/2]

template<typename Array >
auto TiledArray::math::linalg::non_distributed::heig ( const Array A,
TiledRange  evec_trange = TiledRange() 
)

Solve the standard eigenvalue problem with LAPACK.

A(i,k) X(k,j) = X(i,j) E(j)

Example Usage:

auto [E, X] = heig(A, ...)

Template Parameters
ArrayInput array type
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()
Returns
A tuple containing the eigenvalues and eigenvectors of input array as std::vector and in TA format, respectively.

Definition at line 54 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::non_distributed::heig ( const ArrayA &  A,
const ArrayB &  B,
TiledRange  evec_trange = TiledRange() 
)

Solve the generalized eigenvalue problem with LAPACK.

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
Parameters
[in]AInput array to be diagonalized. Must be rank-2
[in]BPositive-definite matrix
[in]evec_trangeTiledRange for resulting eigenvectors. If left empty, will default to array.trange()
Returns
A tuple containing the eigenvalues and eigenvectors of input array as std::vector and in TA format, respectively.

Definition at line 95 of file heig.h.

Here is the call graph for this function:

◆ lu_inv()

template<typename Array >
auto TiledArray::math::linalg::non_distributed::lu_inv ( const Array A,
TiledRange  ainv_trange = TiledRange() 
)

Invert a matrix via LU.

Definition at line 58 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::non_distributed::lu_solve ( const ArrayA &  A,
const ArrayB &  B,
TiledRange  x_trange = TiledRange() 
)

Solve a linear system via LU factorization.

Definition at line 40 of file lu.h.

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

◆ rank_local_cholesky()

template<typename Tile , typename Policy >
auto TiledArray::math::linalg::non_distributed::rank_local_cholesky ( const DistArray< Tile, Policy > &  A)

Definition at line 36 of file cholesky.h.

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

◆ svd()

template<SVD::Vectors Vectors, typename Array >
auto TiledArray::math::linalg::non_distributed::svd ( const Array A,
TiledRange  u_trange = TiledRange(),
TiledRange  vt_trange = TiledRange() 
)

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<SVDLeftstd::vectors> (A, ...) auto [S, VT] = svd<SVDRightstd::vectors>(A, ...) auto [S, U, VT] = svd<SVDAllstd::vectors> (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).
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 svd.h.

Here is the call graph for this function: