Einsums/hipBLAS.hpp#
Defined in header Einsums/hipBLAS.hpp.
See Overview for a list of names and headers that are part of the public Einsums API.
-
namespace einsums
-
-
namespace blas
-
namespace gpu#
Functions
-
template<typename T>
void gemm(char transa, char transb, int_t m, int_t n, int_t k, T alpha, T const *a, int_t lda, T const *b, int_t ldb, T beta, T *c, int_t ldc)# Perform a General Matrix Multiply (GEMM) operation.
This function computes the product of two matrices,
\[\mathbf{C} := \alpha \mathbf{A}\mathbf{B} + \beta\mathbf{C} \]where \(\mathbf{A}\), \(\mathbf{B}\), and \(\mathbf{C}\) are matrices, and \(\alpha\) and \(\beta\) are scalar values.Added in version 2.0.0.
Note
The function performs one of the following matrix operations:
If transA is ‘N’ or ‘n’ and transB is ‘N’ or ‘n’: \(\mathbf{C} = \alpha\mathbf{AB} + \beta\mathbf{C}\)
If transA is ‘N’ or ‘n’ and transB is ‘T’ or ‘t’: \(\mathbf{C} = \alpha\mathbf{A}\mathbf{B}^T + \beta\mathbf{C}\)
If transA is ‘T’ or ‘t’ and transB is ‘N’ or ‘n’: \(\mathbf{C} = \alpha\mathbf{A}^T\mathbf{B} + \beta\mathbf{C}\)
If transA is ‘T’ or ‘t’ and transB is ‘T’ or ‘t’: \(\mathbf{C} = \alpha\mathbf{A}^T\mathbf{B}^T + \beta\mathbf{C}\)
If transA is ‘C’ or ‘c’ and transB is ‘N’ or ‘n’: \(\mathbf{C} = \alpha\mathbf{A}^H\mathbf{B} + \beta\mathbf{C}\)
If transA is ‘C’ or ‘c’ and transB is ‘T’ or ‘t’: \(\mathbf{C} = \alpha\mathbf{A}^H\mathbf{B}^Y + \beta\mathbf{C}\)
etc.
- Template Parameters:
T – The datatype of the GEMM.
- Parameters:
transa – [in] Whether to transpose matrix a :
’N’ or ‘n’ for no transpose,
’T’ or ‘t’ for transpose,
’C’ or ‘c’ for conjugate transpose.
transb – [in] Whether to transpose matrix b .
m – [in] The number of rows in matrix A and C.
n – [in] The number of columns in matrix B and C.
k – [in] The number of columns in matrix A and rows in matrix B.
alpha – [in] The scalar alpha.
a – [in] A pointer to the matrix A with dimensions
(lda, k)when transa is ‘N’ or ‘n’, and(lda, m)otherwise.lda – [in] Leading dimension of A, specifying the distance between two consecutive columns.
b – [in] A pointer to the matrix B with dimensions
(ldb, n)when transB is ‘N’ or ‘n’, and(ldb, k)otherwise.ldb – [in] Leading dimension of B, specifying the distance between two consecutive columns.
beta – [in] The scalar beta.
c – [inout] A pointer to the matrix C with dimensions
(ldc, n).ldc – [in] Leading dimension of C, specifying the distance between two consecutive columns.
- Throws:
-
template<typename T>
void gemv(char transa, int_t m, int_t n, T alpha, T const *a, int_t lda, T const *x, int_t incx, T beta, T *y, int_t incy)# Computes a matrix-vector product using a general matrix.
The gemv routine performs a matrix-vector operation defined as:
\[\mathbf{y} := \alpha \mathbf{A} \mathbf{x} + \beta \mathbf{y} \]or\[\mathbf{y} := \alpha \mathbf{A}^T \mathbf{x} + \beta \mathbf{y} \]or\[\mathbf{y} := \alpha \mathbf{A}^H \mathbf{x} + \beta \mathbf{y} \]Added in version 2.0.0.
- Template Parameters:
T – the underlying data type of the matrix and vector
- Parameters:
transa – [in] what to do with
a- no trans, trans, conjgm – [in] specifies the number of rows of
an – [in] specifies the number of columns of
aalpha – [in] Specifies the scaler alpha
a – [in] Array, size lda * m
lda – [in] Specifies the leading dimension of
aas declared in the calling functionx – [in] array, vector x
incx – [in] Specifies the increment for the elements of
xbeta – [in] Specifies the scalar beta. When beta is set to zero, then
yneed not be set on input.y – [inout] array, vector y
incy – [in] Specifies the increment for the elements of
y.
- Throws:
-
template<typename T>
int_t syev(char job, char uplo, int_t n, T *a, int_t lda, T *w, T *work, int_t lwork)# Performs diagonalization of a symmetrix matrix.
The syev routine finds the matrices that satisfy the following equation.
\[\mathbf{A} = \mathbf{P} \mathbf{\Lambda} \mathbf{P}^T \]In the above equation, \( \mathbf{A} \) is a real symmetric matrix, \( \mathbf{P} \) is a real orthogonal matrix whose columns are the eigenvectors of \( \mathbf{A} \), and \( \mathbf{\Lambda} \) is a diagonal matrix, whose elements are the eigenvalues of \(\mathbf{A} \). The eigenvalues are stored in a vector form on exit.Added in version 2.0.0.
- Template Parameters:
T – The type the function will handle.
- Parameters:
job – [in] Whether to compute the eigenvectors. Can be either ‘n’ or ‘v’, case insensitive.
uplo – [in] Whether the matrix data is stored in the upper or lower triangle. Can be either ‘u’ or ‘l’, case insensitive.
n – [in] The number of rows/columns of the input matrix.
a – [inout] The input matrix. On output, it will be changed. If the eigenvectors are requested, then they will be placed in the columns of
aon exit.lda – [in] The leading dimension of the input matrix.
w – [out] The output vector for the eigenvalues.
work – [inout] A work array. If
lworkis -1, then no operations are performed and the first value in the work array is the optimal work buffer size.lwork – [in] The size of the work array. If
lworkis -1, then a workspace query is assumed. No operations will be performed and the optimal workspace size will be put into the first element ofwork.
- Returns:
0 on success. If positive, this means that the algorithm did not converge. The return value indicates the number of eigenvalues that were able to be computed. If negative, this means that one of the parameters was invalid. The absolute value indicates which parameter was bad.
-
template<typename T>
int_t heev(char job, char uplo, int_t n, std::complex<T> *a, int_t lda, T *w, std::complex<T> *work, int_t lwork)# Performs diagonalization of a Hermitian matrix.
The heev routine finds the matrices that satisfy the following equation.
\[\mathbf{A} = \mathbf{P} \mathbf{\Lambda} \mathbf{P}^H \]In the above equation, \(\mathbf{A}\) is a Hermitian matrix, \(\mathbf{P}\) is a unitary matrix whose columns are the eigenvectors of \(\mathbf{A}\), and \(\mathbf{\Lambda}\) is a diagonal matrix, whose elements are the eigenvalues of \(\mathbf{A}\). The eigenvalues are stored in a vector form on exit.Added in version 2.0.0.
- Template Parameters:
T – The type this function handles.
- Parameters:
job – [in] Whether to compute the eigenvectors. Can be either ‘n’ or ‘v’, case insensitive.
uplo – [in] Whether the matrix data is stored in the upper or lower triangle. Can be either ‘u’ or ‘l’, case insensitive.
n – [in] The number of rows/columns of the input matrix.
a – [inout] The input matrix. On output, it will be changed. If the eigenvectors are requested, then they will be placed in the columns of
aon exit.lda – [in] The leading dimension of the input matrix.
w – [out] The output vector for the eigenvalues.
work – [inout] A work array. If
lworkis -1, then no operations are performed and the first value in the work array is the optimal work buffer size.lwork – [in] The size of the work array. If
lworkis -1, then a workspace query is assumed. No operations will be performed and the optimal workspace size will be put into the first element ofwork.
- Returns:
0 on success. If positive, this means that the algorithm did not converge. The return value indicates the number of eigenvalues that were able to be computed. If negative, this means that one of the parameters was invalid. The absolute value indicates which parameter was bad.
-
template<typename T>
int_t gesv(int_t n, int_t nrhs, T *a, int_t lda, int *ipiv, T *b, int_t ldb, T *x, int_t ldx)# Solve a system of linear equations.
Solves equations of the following form.
\[\mathbf{A}\mathbf{x} = \mathbf{B} \]Added in version 2.0.0.
- Template Parameters:
T – The type this function handles.
- Parameters:
n – [in] The number of rows and columns of \(\mathbf{A}\) and rows \(\mathbf{B}\).
nrhs – [in] The number of columns of \(\mathbf{B}\)
a – [inout] The coefficient matrix. On exit, it contains the LU decomposition of
a, where the lower-triangle matrix has unit diagonal entries, which are not stored.lda – [in] The leading dimension of
a.ipiv – [out] A list of pivots used in the decomposition.
b – [in] The results matrix. Unlike normal LAPACK, this is not overwritten by the function.
ldb – [in] The leading dimension of
b.x – [out] The output matrix. On AMD cards, this is allowed — even recommended — to be the same as
b. On NVidia cards, it is not.ldx – [in] The leading dimension of the output matrix.
- Returns:
0 on success. If positive, then the matrix was singular. If negative, then a bad value was passed to the function. The absolute value indicates which parameter was bad.
-
template<typename T>
void scal(int_t n, T const alpha, T *vec, int_t inc)# Scales a vector by a value.
Added in version 2.0.0.
- Template Parameters:
T – The type this function handles.
- Parameters:
n – [in] The number of elements to scale the vector by.
alpha – [in] The scale factor.
vec – [inout] The vector to scale.
inc – [in] The spacing between elements of the vector.
-
template<Complex T>
void scal(int_t n, RemoveComplexT<T> const alpha, T *vec, int_t inc)# Scales a complex vector by a real value.
Added in version 2.0.0.
- Parameters:
n – [in] The number of elements to scale the vector by.
alpha – [in] The scale factor.
vec – [inout] The vector to scale.
inc – [in] The spacing between elements of the vector.
-
template<typename T>
void rscl(int_t n, T const alpha, T *vec, int_t inc)# Scales a vector by the reciprocal of a value.
Added in version 2.0.0.
- Template Parameters:
T – The type this function handles.
- Parameters:
n – [in] The number of elements in the vector.
alpha – [in] The value to divide all the elements in the vector by.
vec – [inout] The vector to scale.
inc – [in] The spacing between elements in the vector.
-
template<Complex T>
void rscl(int_t n, RemoveComplexT<T> const alpha, T *vec, int_t inc)# Scales a complex vector by the reciprocal of a real value.
Added in version 2.0.0.
- Template Parameters:
T – The type this function handles.
- Parameters:
n – [in] The number of elements in the vector.
alpha – [in] The value to divide all the elements in the vector by.
vec – [inout] The vector to scale.
inc – [in] The spacing between elements in the vector.
-
template<typename T>
T dot(int_t n, T const *x, int_t incx, T const *y, int_t incy)# Computes the dot product of two vectors. For complex vectors it is the non-conjugated dot product; (c|z)dotu in BLAS nomenclature.
Added in version 2.0.0.
- Template Parameters:
T – underlying data type
- Parameters:
n – [in] length of the vectors
x – [in] first vector
incx – [in] how many elements to skip in x
y – [in] second vector
incy – [in] how many elements to skip in yo
- Returns:
result of the dot product
-
template<typename T>
T dotc(int_t n, T const *x, int_t incx, T const *y, int_t incy)# Computes the dot product of two vectors. For complex vector it is the conjugated dot product; (c|z)dotc in BLAS nomenclature.
Added in version 2.0.0.
- Template Parameters:
T – underlying data type
- Parameters:
n – [in] length of the vectors
x – [in] first vector
incx – [in] how many elements to skip in x
y – [in] second vector
incy – [in] how many elements to skip in yo
- Returns:
result of the dot product
-
template<typename T>
void axpy(int_t n, T alpha_x, T const *x, int_t inc_x, T *y, int_t inc_y)# Adds two vectors together with a scale factor.
Computes the following.
\[\mathbf{y} := \alpha\mathbf{x} + \mathbf{y} \]Added in version 2.0.0.
- Template Parameters:
T – The type this function handles.
- Parameters:
n – [in] The number of elements in the vectors.
alpha_x – [in] The scale factor for the input vector.
x – [in] The input vector.
inc_x – [in] The skip value for the output vector. It can be negative to go in reverse, or zero to broadcast values to
y.y – [inout] The output vector.
inc_y – [in] The skip value for the output vector. It can be negative to go in reverse, or zero to sum over the elements of
x.
-
template<typename T>
void axpby(int_t n, T alpha_x, T const *x, int_t inc_x, T b, T *y, int_t inc_y)# Adds two vectors together with a scale factor.
Computes the following.
\[\mathbf{y} := \alpha\mathbf{x} + \beta\mathbf{y} \]Added in version 2.0.0.
- Template Parameters:
T – The type this function handles.
- Parameters:
n – [in] The number of elements in the vectors.
alpha_x – [in] The scale factor for the input vector.
x – [in] The input vector.
inc_x – [in] The skip value for the output vector. It can be negative to go in reverse, or zero to broadcast values to
y.b – [in] The scale factor for the output vector.
y – [inout] The output vector.
inc_y – [in] The skip value for the output vector. It can be negative to go in reverse, or zero to sum over the elements of
x.
-
template<typename T>
void ger(int_t m, int_t n, T alpha, T const *x, int_t inc_x, T const *y, int_t inc_y, T *a, int_t lda)# Performs a rank-1 update of a general matrix.
The ?ger routines perform a matrix-vector operator defined as
\[ \mathbf{A} := \alpha\mathbf{x}\mathbf{y}^T + \mathbf{A} \]Added in version 2.0.0.
- Template Parameters:
T – The type this function handles.
- Parameters:
m – [in] The number of entries in
x.n – [in] The number of entries in
y.alpha – [in] The scale factor for the outer product.
x – [in] The left input vector.
inc_x – [in] The skip value for the left input. May be negative to go in reverse.
y – [in] The right input vector.
inc_y – [in] The skip value for the right input. May be negative to go in reverse.
a – [inout] The output matrix.
lda – [in] The leading dimension of
a.
- Throws:
-
template<typename T>
void gerc(int_t m, int_t n, T alpha, T const *x, int_t inc_x, T const *y, int_t inc_y, T *a, int_t lda)# Performs a rank-1 update of a general matrix.
The ?gerc routines perform a matrix-vector operator defined as
\[ \mathbf{A} := \alpha\mathbf{x}\mathbf{y}^H + \mathbf{A} \]Added in version 2.0.0.
- Template Parameters:
T – The type this function handle.
- Parameters:
m – [in] The number of entries in
x.n – [in] The number of entries in
y.alpha – [in] The scale factor for the outer product.
x – [in] The left input vector.
inc_x – [in] The skip value for the left input. May be negative to go in reverse.
y – [in] The right input vector.
inc_y – [in] The skip value for the right input. May be negative to go in reverse.
a – [inout] The output matrix.
lda – [in] The leading dimension of
a.
- Throws:
-
template<typename T>
int_t getrf(int_t m, int_t n, T *a, int_t lda, int *ipiv)# Computes the LU factorization of a general M-by-N matrix A using partial pivoting with row int_terchanges.
The factorization has the form
\[ \mathbf{A} = \mathbf{PLU} \]where \(\mathbf{P}\) is a permutation matrix, \(\mathbf{L}\) is lower triangular with unit diagonal elements (lower trapezoidal if m > n) and \(\mathbf{U}\) is upper triangular (upper trapezoidal if m < n).Added in version 2.0.0.
- Template Parameters:
T – The type this function handles.
- Parameters:
m – [in] The number of rows in the input.
n – [in] The number of columns in the input.
a – [inout] The input matrix. On exit, it contains the upper and lower triangular matrices. The elemnts of the lower triangular matrix are not stored since they are all 1.
lda – [in] The leading dimension of the matrix.
ipiv – [out] The list of pivots.
- Returns:
0 on success. If positive, the matrix is singular and the result should not be used for solving systems of equations. The decomposition was performed, though. If negative, one of the inputs had an invalid value. The absolute value indicates which input it is.
-
template<typename T>
int_t getri(int_t n, T *a, int_t lda, int *ipiv, T *c, int_t ldc)# Computes the inverse of a matrix using the LU factorization computed by getrf.
Added in version 2.0.0.
- Template Parameters:
T – The type this function handles.
- Parameters:
n – [in] The number of rows and columns of the matrix.
a – [in] The input matrix after being processed by getrf.
lda – [in] The leading dimension of the matrix.
ipiv – [in] The pivots from getrf.
c – [out] The calculated inverse.
ldc – [in] The leading dimension of the output matrix.
- Returns:
0 on success. If positive, the matrix is singular and an inverse could not be computed. If negative, one of the inputs is invalid, and the absolute value indicates which input is bad.
-
template<typename T>
void lassq(int_t n, T const *x, int_t incx, RemoveComplexT<T> *scale, RemoveComplexT<T> *sumsq)# Compute the sum of the squares of the input vector without roundoff error.
\[scale^2 sumsq := \left|\mathbf{x}\right|^2 + scale^2 sumsq \]Added in version 2.0.0.
- Template Parameters:
T – The type this function handles.
- Parameters:
n – [in] The number of elements in the vector.
x – [in] The input vector.
incx – [in] The skip value for the vector.
scale – [inout] The scale value used to avoid overflow/underflow. It is also used as an input to continue a previous calculation.
sumsq – [inout] The result of the operation, scaled to avoid overflow/underflow. It is also used as an input to continue a previous calculation.
-
template<typename T>
RemoveComplexT<T> nrm2(int_t n, T const *x, int_t incx)# Compute the Euclidean norm of a vector.
Added in version 2.0.0.
- Template Parameters:
T – The type this function handles.
- Parameters:
n – [in] The number of elements in the vector.
x – [in] The input vector.
incx – [in] The skip value for the vector.
- Returns:
The Euclidean norm of the vector.
-
template<typename T>
auto gesvd(char jobu, char jobvt, int_t m, int_t n, T *a, int_t lda, RemoveComplexT<T> *s, T *u, int_t ldu, T *vt, int_t ldvt, T *superb)# Performs singular value decomposition for a matrix using the QR algorithm.
\[\mathbf{A} = \mathbf{U\Sigma V}^T \]Added in version 2.0.0.
- Template Parameters:
T – The type this function handles.
- Parameters:
jobu – [in] Whether to compute the U matrix. Case insensitive. Can be ‘a’, ‘s’, ‘o’, or ‘n’.
jobvt – [in] Whether to compute the transpose of the V matrix. Case insensitive. Can be ‘a’, ‘s’, ‘o’, or ‘n’.
m – [in] The number of rows of the input matrix.
n – [in] The number of columns of the input matrix.
a – [inout] The input matrix.
lda – [in] The leading dimension of the input matrix.
s – [out] The singular values output.
u – [out] The U matrix from the singular value decomposition.
ldu – [in] The leading dimension of U.
vt – [out] The transpose of the V matrix from the singular value decomposition.
ldvt – [in] The leading dimension of the transpose of the V matrix.
superb – [inout] Temporary storage area for intermediates in the computation.
- Returns:
0 on success. If positive, the algorithm did not converge. If negative, then one of the parameters had a bad value. The absolute value of the return gives the parameter.
-
template<typename T>
int_t geqrf(int_t m, int_t n, T *a, int_t lda, T *tau)# Set up for computing the QR decomposition of a matrix.
\[\mathbf{A} = \mathbf{QR} \]Here, \(\mathbf{Q}\) is an orthogonal matrix and \(\mathbf{R}\) is an upper triangular matrix.
Added in version 2.0.0.
- Template Parameters:
T – The type this function handles.
- Parameters:
m – [in] The number of rows in the input matrix.
n – [in] The number of columns in the input matrix.
a – [inout] The input matrix. On exit, contains the data needed to compute the Q and R matrices. The entries on and above the diagonal are the entries of the R matrix. The rest is needed to find the Q matrix.
lda – [in] The leading dimension of the input matrix.
tau – [out] On exit, holds the Householder reflector parameters for computing the Q matrix.
- Returns:
0 on success. If negative, one of the inputs had a bad value, and the absolute value of the return tells you which one it was.
-
template<typename T>
int_t orgqr(int_t m, int_t n, int_t k, T *a, int_t lda, T *tau)# Extract the Q matrix after a call to geqrf.
Added in version 2.0.0.
- Template Parameters:
T – The type this function handles.
- Parameters:
m – [in] The number of rows of the input matrix.
n – [in] The number of columns in the input matrix.
k – [in] The number of elementary reflectors used in the calculation.
a – [inout] The input matrix after being processed by geqrf.
lda – [in] The leading dimension of the input matrix.
tau – [in] The scales for the elementary reflectors from geqrf.
- Returns:
0 on success. If negative, then one of the inputs had an invalid value, and the absolute value indicates which parameter it is.
-
template<typename T>
void copy(int_t n, T const *x, int_t inc_x, T *y, int_t inc_y)# Copy data from one vector to another.
Added in version 2.0.0.
- Template Parameters:
T – The type this function handles.
- Parameters:
n – [in] The number of elements to copy.
x – [in] The input vector.
inc_x – [in] The skip value for the input vector. If negative, the vector is traversed backwards. If zero, the values are broadcast to the output vector.
y – [out] The output vector.
inc_y – [in] The skip value for the output vector. If negative, the vector is traversed backwards.
-
template<typename T>
int_t lascl(char type, int_t kl, int_t ku, T cfrom, T cto, int_t m, int_t n, T *A, int_t lda)# Scales a general matrix. The scale factor is
cto / cfrom, but the scale is performed without overflow/underflow.Added in version 2.0.0.
- Template Parameters:
T – The type this function handles.
- Parameters:
type – [in] The type of matrix. Case insensitive. ‘g’ is for general matrices, ‘l’ is for lower triangular matrices, ‘u’ if for upper triangular matrices, ‘h’ is for hessenberg matrices, ‘b’ is for symmetric band matrices with lower bandwidth
kland upper bandwidth ofkuand with only the lower half stored, ‘q’ is the same as ‘b’ but with the upper half stored instead, and ‘z’ is the same as ‘b’ but with a more complicated storage scheme.kl – [in] The lower bandwidth of the matrix. Only used if the type is ‘b’, ‘q’, or ‘z’.
ku – [in] The upper bandwidth of the matrix. Only used if the type is ‘b’, ‘q’, or ‘z’.
cfrom – [in] The denominator for the scale.
cto – [in] The numerator for the scale.
m – [in] The number of rows in the matrix.
n – [in] The number of columns in the matrix.
A – [inout] The matrix being scaled.
lda – [in] The leading dimension of the matrix.
- Returns:
0 on success. If negative, then one of the parameters had an invalid value. The absolute value of the return indicates which parameter it was.
-
template<typename T>
void dirprod(int_t n, T alpha, T const *x, int_t incx, T const *y, int_t incy, T *z, int_t incz)# Computes the direct product between two vectors.
\[z_i := z_i + \alpha x_i y_i \]Added in version 2.0.0.
- Template Parameters:
T – The type this function handles.
- Parameters:
n – [in] The number of elements in the vectors.
alpha – [in] The scale factor for the product.
x – [in] The first input vector.
incx – [in] The skip value for the first vector.
y – [in] The second input vector.
incy – [in] The skip value for the second vector.
z – [inout] The accumulation vector.
incz – [in] The skip value for the accumulation vector.
-
template<typename T>
RemoveComplexT<T> asum(int_t n, T const *x, int_t incx)# Computes the sum of the absolute values of the input vector. If the vector is complex, then it is the sum of the absolute values of the components, not the magnitudes.
Added in version 2.0.0.
- Template Parameters:
T – The type this function handles.
- Parameters:
n – [in] The number of elements.
x – [in] The vector to process.
incx – [in] The skip value for the vector.
- Returns:
The sum of the absolute values of the inputs as stated above.
-
template<typename T>
RemoveComplexT<T> sum1(int_t n, T const *x, int_t incx)# Computes the sum of the absolute values of the input vector. If the vector is complex, then it is the sum of the magnitudes.
Added in version 2.0.0.
- Template Parameters:
T – The type this function handles.
- Parameters:
n – [in] The number of elements.
x – [in] The vector to process.
incx – [in] The skip value for the vector.
- Returns:
The sum of the absolute values of the inputs as stated above.
-
template<typename T>
void lacgv(int_t n, T *x, int_t incx)# Take the conjugate of a vector. Does nothing if the vector is real.
Added in version 2.0.0.
- Template Parameters:
T – The type this function handles.
- Parameters:
n – [in] The number of elements in the vector.
x – [in] The input vector.
incx – [in] The skip value for the vector.
-
template<typename T>
-
namespace gpu#
-
namespace blas