23#include <cusp/detail/config.h>
25#include <cusp/lapack/detail/defs.h>
41template<
typename DerivedPolicy,
44void getrf(
const thrust::detail::execution_policy_base<DerivedPolicy> &exec,
93template<
typename Array2d,
99template<
typename DerivedPolicy,
101void potrf(
const thrust::detail::execution_policy_base<DerivedPolicy> &exec,
148template<
typename Array2d>
153template<
typename DerivedPolicy,
156void sytrf(
const thrust::detail::execution_policy_base<DerivedPolicy> &exec,
212template<
typename Array2d,
219template<
typename DerivedPolicy,
222void getrs(
const thrust::detail::execution_policy_base<DerivedPolicy> &exec,
283template<
typename Array2d,
291template<
typename DerivedPolicy,
293void potrs(
const thrust::detail::execution_policy_base<DerivedPolicy> &exec,
358template<
typename Array2d>
364template<
typename DerivedPolicy,
367void sytrs(
const thrust::detail::execution_policy_base<DerivedPolicy> &exec,
436template<
typename Array2d,
444template<
typename DerivedPolicy,
446void trtrs(
const thrust::detail::execution_policy_base<DerivedPolicy> &exec,
504template<
typename Array2d>
512template<
typename DerivedPolicy,
514void trtri(
const thrust::detail::execution_policy_base<DerivedPolicy> &exec,
553template<
typename Array2d>
559template<
typename DerivedPolicy,
562void syev(
const thrust::detail::execution_policy_base<DerivedPolicy> &exec,
613template<
typename Array2d,
621template<
typename DerivedPolicy,
626void stev(
const thrust::detail::execution_policy_base<DerivedPolicy> &exec,
627 const Array1d1& alphas,
628 const Array1d2& betas,
680template<
typename Array1d1,
684void stev(
const Array1d1& alphas,
685 const Array1d2& betas,
691template<
typename DerivedPolicy,
696void sygv(
const thrust::detail::execution_policy_base<DerivedPolicy> &exec,
753template<
typename Array2d1,
763template<
typename DerivedPolicy,
766void gesv(
const thrust::detail::execution_policy_base<DerivedPolicy> &exec,
822template<
typename Array2d,
834#include <cusp/lapack/detail/lapack.inl>
void potrf(Array2d &A, char uplo='U')
Computes the Cholesky factorization of a symmetric (Hermitian) positive-definite matrix.
void sytrs(const Array2d &A, const Array1d &piv, Array2d &B, char uplo='U')
Solves a system of linear equations with a UDU- or LDL-factored symmetric matrix.
void stev(const Array1d1 &alphas, const Array1d2 &betas, Array1d3 &eigvals, Array2d &eigvecs, char job='V')
Computes all eigenvalues and, optionally, eigenvectors of a real symmetric tridiagonal matrix.
void sygv(const Array2d1 &A, const Array2d2 &B, Array1d &eigvals, Array2d3 &eigvecs)
Computes all eigenvalues and, optionally, eigenvectors of a real generalized symmetric definite eigen...
void getrs(const Array2d &A, const Array1d &piv, Array2d &B, char trans='N')
Solves a system of linear equations with an LU-factored square matrix, with multiple right-hand sides...
void syev(const Array2d &A, Array1d &eigvals, Array2d &eigvecs, char uplo='U')
Computes all eigenvalues and eigenvectors of a real symmetric matrix.
void gesv(const Array2d &A, Array2d &B, Array1d &piv)
Computes the solution of a system of linear equations with a square matrix A and multiple right-hand ...
void trtrs(const Array2d &A, Array2d &B, char uplo='U', char trans='N', char diag='N')
Solves a system of linear equations with a triangular matrix, with multiple right-hand sides.
void potrs(const Array2d &A, Array2d &B, char uplo='U')
Solves a system of linear equations with a Cholesky-factored symmetric (Hermitian) positive-definite ...
void sytrf(Array2d &A, Array1d &piv, char uplo='U')
Computes the Bunch-Kaufman factorization of a symmetric matrix.
void getrf(Array2d &A, Array1d &piv)
Compute LU factorization of matrix.
void trtri(Array2d &A, char uplo='U', char diag='N')
This routine computes the inverse of a triangular matrix.