14#if GUANAQO_WITH_HL_BLAS_TRACING
15#define GUANAQO_TRACE_HL_BLAS(name, gflops) GUANAQO_TRACE_LINALG(name, gflops)
17#define GUANAQO_TRACE_HL_BLAS(name, gflops) GUANAQO_NOOP()
28template <
class T,
class I>
37 xgemv<T, I>(CblasColMajor, CblasNoTrans, A.rows, A.cols, alpha, A.data,
38 A.outer_stride, x.data, I{1}, beta, y.
data, I{1});
41template <
class T,
class I>
50 xgemv<T, I>(CblasColMajor, CblasTrans, A.rows, A.cols, alpha, A.data,
51 A.outer_stride, x.data, I{1}, beta, y.
data, I{1});
54template <
class T,
class I>
63 xgemv<T, I>(CblasRowMajor, CblasNoTrans, A.rows, A.cols, alpha, A.data,
64 A.outer_stride, x.data, I{1}, beta, y.
data, I{1});
67template <
class T,
class I>
76 xgemv<T, I>(CblasRowMajor, CblasTrans, A.rows, A.cols, alpha, A.data,
77 A.outer_stride, x.data, I{1}, beta, y.
data, I{1});
80template <
class T,
class I, StorageOrder O>
90 CblasNoTrans, CblasNoTrans, C.rows, C.cols, A.cols, alpha,
91 A.data, A.outer_stride, B.data, B.outer_stride, beta, C.data,
95template <
class T,
class I, StorageOrder O>
105 CblasTrans, CblasNoTrans, C.rows, C.cols, A.rows, alpha, A.data,
106 A.outer_stride, B.data, B.outer_stride, beta, C.data,
110template <
class T,
class I, StorageOrder O>
120 CblasTrans, CblasTrans, C.rows, C.cols, A.rows, alpha, A.data,
121 A.outer_stride, B.data, B.outer_stride, beta, C.data,
125template <
class T,
class I, StorageOrder O>
135 CblasNoTrans, CblasTrans, C.rows, C.cols, A.cols, alpha, A.data,
136 A.outer_stride, B.data, B.outer_stride, beta, C.data,
140template <
class T,
class I>
150 A.cols, alpha, A.data, A.outer_stride, B.data, B.outer_stride,
154template <
class T,
class I>
164 A.rows, alpha, A.data, A.outer_stride, B.data, B.outer_stride,
168template <
class T,
class I>
178 A.rows, alpha, A.data, A.outer_stride, B.data, B.outer_stride,
182template <
class T,
class I>
192 A.cols, alpha, A.data, A.outer_stride, B.data, B.outer_stride,
196template <
class T,
class I, StorageOrder O>
205 CblasLower, A.rows, alpha, A.data, A.outer_stride, x.data, I{1},
209template <
class T,
class I>
216 xtrmv<T, I>(CblasColMajor, CblasLower, CblasNoTrans, CblasNonUnit, A.rows,
217 A.data, A.outer_stride, x.
data, I{1});
220template <
class T,
class I>
227 xtrmv<T, I>(CblasColMajor, CblasLower, CblasTrans, CblasNonUnit, A.rows,
228 A.data, A.outer_stride, x.
data, I{1});
231template <
class T,
class I>
238 xtrsv<T, I>(CblasColMajor, CblasLower, CblasNoTrans, CblasNonUnit, A.rows,
239 A.data, A.outer_stride, x.
data, I{1});
242template <
class T,
class I>
249 xtrsv<T, I>(CblasColMajor, CblasLower, CblasTrans, CblasNonUnit, A.rows,
250 A.data, A.outer_stride, x.
data, I{1});
253template <
class T,
class I>
259 xtrmm<T, I>(CblasColMajor, CblasLeft, CblasLower, CblasNoTrans,
260 CblasNonUnit, B.
rows, B.
cols, alpha, A.data, A.outer_stride,
264template <
class T,
class I>
270 xtrmm<T, I>(CblasColMajor, CblasLeft, CblasLower, CblasTrans, CblasNonUnit,
275template <
class T,
class I>
281 xtrmm<T, I>(CblasColMajor, CblasRight, CblasLower, CblasNoTrans,
282 CblasNonUnit, B.
rows, B.
cols, alpha, A.data, A.outer_stride,
286template <
class T,
class I>
292 xtrmm<T, I>(CblasColMajor, CblasRight, CblasLower, CblasTrans, CblasNonUnit,
297template <
class T,
class I>
303 xsyrk<T, I>(CblasColMajor, CblasLower, CblasNoTrans, C.
rows, A.cols, alpha,
307template <
class T,
class I>
313 xsyrk<T, I>(CblasColMajor, CblasLower, CblasTrans, C.
rows, A.rows, alpha,
317template <
class T,
class I>
323 xtrsm<T, I>(CblasColMajor, CblasLeft, CblasLower, CblasNoTrans,
324 CblasNonUnit, B.
rows, B.
cols, alpha, A.data, A.outer_stride,
328template <
class T,
class I>
334 xtrsm<T, I>(CblasColMajor, CblasLeft, CblasLower, CblasTrans, CblasNonUnit,
339template <
class T,
class I>
345 xtrsm<T, I>(CblasColMajor, CblasRight, CblasLower, CblasNoTrans,
346 CblasNonUnit, B.
rows, B.
cols, alpha, A.data, A.outer_stride,
350template <
class T,
class I>
356 xtrsm<T, I>(CblasColMajor, CblasRight, CblasLower, CblasTrans, CblasNonUnit,
362template <
class T,
class I>
363void xsytrf_rk(
const char *uplo,
const I *n, T *a,
const I *lda, T *e, I *ipiv,
364 T *work,
const I *lwork, I *info);
367template <
class T,
class I>
368void xtrtrs(
const char *uplo,
const char *trans,
const char *diag,
const I *n,
369 const I *nrhs,
const T *A,
const I *ldA, T *B,
const I *ldB,
373template <
class T,
class I>
384template <
class T,
class I>
396template <
class T,
class I>
Assertion and assumption macros with debug/release semantics.
This file provides simple overloaded wrappers around standard BLAS functions.
std::integral_constant< I, 1 > UnitStride
void xgemm_NN(T alpha, std::type_identity_t< MatrixView< const T, I, UnitStride< I >, O > > A, std::type_identity_t< MatrixView< const T, I, UnitStride< I >, O > > B, T beta, MatrixView< T, I, UnitStride< I >, O > C)
void xgemmt_LNN(T alpha, std::type_identity_t< MatrixView< const T, I > > A, std::type_identity_t< MatrixView< const T, I > > B, T beta, MatrixView< T, I > C)
void xtrsm_LLNN(T alpha, std::type_identity_t< MatrixView< const T, I > > A, MatrixView< T, I > B)
void xgemm_TT(T alpha, std::type_identity_t< MatrixView< const T, I, UnitStride< I >, O > > A, std::type_identity_t< MatrixView< const T, I, UnitStride< I >, O > > B, T beta, MatrixView< T, I, UnitStride< I >, O > C)
void xtrsm_RLNN(T alpha, std::type_identity_t< MatrixView< const T, I > > A, MatrixView< T, I > B)
void xgemmt_LTT(T alpha, std::type_identity_t< MatrixView< const T, I > > A, std::type_identity_t< MatrixView< const T, I > > B, T beta, MatrixView< T, I > C)
void xtrmv_LNN(std::type_identity_t< MatrixView< const T, I > > A, MatrixView< T, I > x)
void xtrsv_LNN(std::type_identity_t< MatrixView< const T, I > > A, MatrixView< T, I > x)
void xtrmm_LLTN(T alpha, std::type_identity_t< MatrixView< const T, I > > A, MatrixView< T, I > B)
void xtrsv_LTN(std::type_identity_t< MatrixView< const T, I > > A, MatrixView< T, I > x)
void xtrtri_LN(MatrixView< T, I > A)
void xpotrf_L(MatrixView< T, I > A)
void xgemmt_LNT(T alpha, std::type_identity_t< MatrixView< const T, I > > A, std::type_identity_t< MatrixView< const T, I > > B, T beta, MatrixView< T, I > C)
void xtrsm_LLTN(T alpha, std::type_identity_t< MatrixView< const T, I > > A, MatrixView< T, I > B)
void xsyrk_LN(T alpha, std::type_identity_t< MatrixView< const T, I > > A, T beta, MatrixView< T, I > C)
void xtrmm_LLNN(T alpha, std::type_identity_t< MatrixView< const T, I > > A, MatrixView< T, I > B)
void xgemm_NT(T alpha, std::type_identity_t< MatrixView< const T, I, UnitStride< I >, O > > A, std::type_identity_t< MatrixView< const T, I, UnitStride< I >, O > > B, T beta, MatrixView< T, I, UnitStride< I >, O > C)
void xgemm_TN(T alpha, std::type_identity_t< MatrixView< const T, I, UnitStride< I >, O > > A, std::type_identity_t< MatrixView< const T, I, UnitStride< I >, O > > B, T beta, MatrixView< T, I, UnitStride< I >, O > C)
void xsymv_L(T alpha, MatrixView< const T, I, UnitStride< I >, O > A, std::type_identity_t< MatrixView< const T, I > > x, T beta, MatrixView< T, I > y)
void xsyrk_LT(T alpha, std::type_identity_t< MatrixView< const T, I > > A, T beta, MatrixView< T, I > C)
void xgemv_N(T alpha, std::type_identity_t< MatrixView< const T, I > > A, std::type_identity_t< MatrixView< const T, I > > x, T beta, MatrixView< T, I > y)
void xtrsm_RLTN(T alpha, std::type_identity_t< MatrixView< const T, I > > A, MatrixView< T, I > B)
void xgemmt_LTN(T alpha, std::type_identity_t< MatrixView< const T, I > > A, std::type_identity_t< MatrixView< const T, I > > B, T beta, MatrixView< T, I > C)
void xtrmm_RLNN(T alpha, std::type_identity_t< MatrixView< const T, I > > A, MatrixView< T, I > B)
void xtrmv_LTN(std::type_identity_t< MatrixView< const T, I > > A, MatrixView< T, I > x)
void xgemv_T(T alpha, std::type_identity_t< MatrixView< const T, I > > A, std::type_identity_t< MatrixView< const T, I > > x, T beta, MatrixView< T, I > y)
void xtrmm_RLTN(T alpha, std::type_identity_t< MatrixView< const T, I > > A, MatrixView< T, I > B)
void xlauum_L(MatrixView< T, I > A)
void xlauum(const char *uplo, I n, T *a, I lda, I *info)
void xtrsm(CBLAS_LAYOUT Layout, CBLAS_SIDE Side, CBLAS_UPLO Uplo, CBLAS_TRANSPOSE TransA, CBLAS_DIAG Diag, I M, I N, std::type_identity_t< T > alpha, std::type_identity_t< const T * > A, I lda, T *B, I ldb)
void xtrmm(CBLAS_LAYOUT Layout, CBLAS_SIDE Side, CBLAS_UPLO Uplo, CBLAS_TRANSPOSE TransA, CBLAS_DIAG Diag, I M, I N, std::type_identity_t< T > alpha, std::type_identity_t< const T * > A, I lda, T *B, I ldb)
void xpotrf(const char *uplo, I n, T *a, I lda, I *info)
void xsymv(CBLAS_LAYOUT Layout, CBLAS_UPLO uplo, I N, std::type_identity_t< T > alpha, std::type_identity_t< const T * > A, I lda, std::type_identity_t< const T * > X, I incX, std::type_identity_t< T > beta, T *Y, I incY)
void xtrmv(CBLAS_LAYOUT Layout, CBLAS_UPLO Uplo, CBLAS_TRANSPOSE TransA, CBLAS_DIAG Diag, I N, std::type_identity_t< const T * > A, I lda, T *X, I incX)
void xtrtrs(const char *uplo, const char *trans, const char *diag, const I *n, const I *nrhs, std::type_identity_t< const T * > A, const I *ldA, T *B, const I *ldB, I *info)
void lapack_throw_on_err(Name &&name, index_t info)
void xtrsv(CBLAS_LAYOUT Layout, CBLAS_UPLO Uplo, CBLAS_TRANSPOSE TransA, CBLAS_DIAG Diag, I N, std::type_identity_t< const T * > A, I lda, T *X, I incX)
void xgemmt(CBLAS_LAYOUT Layout, CBLAS_UPLO uplo, CBLAS_TRANSPOSE TransA, CBLAS_TRANSPOSE TransB, I N, I K, std::type_identity_t< T > alpha, std::type_identity_t< const T * > A, I lda, std::type_identity_t< const T * > B, I ldb, std::type_identity_t< T > beta, T *C, I ldc)
void xsyrk(CBLAS_LAYOUT Layout, CBLAS_UPLO Uplo, CBLAS_TRANSPOSE Trans, I N, I K, std::type_identity_t< T > alpha, std::type_identity_t< const T * > A, I lda, std::type_identity_t< T > beta, T *C, I ldc)
void xsytrf_rk(const char *uplo, const I *n, T *a, const I *lda, T *e, I *ipiv, T *work, const I *lwork, I *info)
void xgemm(CBLAS_LAYOUT Layout, CBLAS_TRANSPOSE TransA, CBLAS_TRANSPOSE TransB, I M, I N, I K, std::type_identity_t< T > alpha, std::type_identity_t< const T * > A, I lda, std::type_identity_t< const T * > B, I ldb, std::type_identity_t< T > beta, T *C, I ldc)
void xgemv(CBLAS_LAYOUT Layout, CBLAS_TRANSPOSE TransA, I M, I N, std::type_identity_t< T > alpha, std::type_identity_t< const T * > A, I lda, std::type_identity_t< const T * > X, I incX, std::type_identity_t< T > beta, T *Y, I incY)
void xtrtri(const char *uplo, const char *diag, I n, T *a, I lda, I *info)
MatrixView< T, I, S, StorageOrder::RowMajor > MatrixViewRM
Convenience alias for a row-major MatrixView.
@ RowMajor
Row-major storage order (C-style, layout right, unit column stride).
#define GUANAQO_ASSUME(x)
Invokes undefined behavior if the expression x does not evaluate to true.
#define GUANAQO_TRACE_HL_BLAS(name, gflops)
A lightweight view of a 2D matrix.
Tracing logger and macros (ITT, Perfetto, or fallback).