35 void scal<float>(
const int N,
40 cblas_sscal(N, alpha, X, incX);
44 void scal<double>(
const int N,
50 cblas_dscal(N, alpha, X, incX);
54 void scal<std::complex<float> >(
const int N,
55 const std::complex<float> & alpha,
56 std::complex<float> * X,
60 cblas_cscal(N, &alpha, X, incX);
65 void scal<std::complex<double> >(
const int N,
66 const std::complex<double> & alpha,
67 std::complex<double> * X,
71 cblas_zscal(N, &alpha, X, incX);
75 void gemm<float>(
const CBLAS_ORDER Order,
const CBLAS_TRANSPOSE TransA,
76 const CBLAS_TRANSPOSE TransB,
const int M,
const int N,
77 const int K,
const float & alpha,
const float *A,
78 const int lda,
const float *B,
const int ldb,
79 const float & beta,
float *C,
const int ldc)
81 cblas_sgemm(Order, TransA, TransB, M, N, K, alpha, A, lda, B, ldb, beta, C, ldc);
85 void gemm<double>(
const CBLAS_ORDER Order,
const CBLAS_TRANSPOSE TransA,
86 const CBLAS_TRANSPOSE TransB,
const int M,
const int N,
87 const int K,
const double & alpha,
const double *A,
88 const int lda,
const double *B,
const int ldb,
89 const double & beta,
double *C,
const int ldc)
91 cblas_dgemm(Order, TransA, TransB, M, N, K, alpha, A, lda, B, ldb, beta, C, ldc);
96 void gemm<std::complex<float> >(
const CBLAS_ORDER Order,
const CBLAS_TRANSPOSE TransA,
97 const CBLAS_TRANSPOSE TransB,
const int M,
const int N,
98 const int K,
const std::complex<float> & alpha,
const std::complex<float> *A,
99 const int lda,
const std::complex<float> *B,
const int ldb,
100 const std::complex<float> & beta, std::complex<float> *C,
const int ldc)
102 cblas_cgemm(Order, TransA, TransB, M, N, K, &alpha, A, lda, B, ldb, &beta, C, ldc);
106 void gemm<std::complex<double> >(
const CBLAS_ORDER Order,
const CBLAS_TRANSPOSE TransA,
107 const CBLAS_TRANSPOSE TransB,
const int M,
const int N,
108 const int K,
const std::complex<double> & alpha,
const std::complex<double> *A,
109 const int lda,
const std::complex<double> *B,
const int ldb,
110 const std::complex<double> & beta, std::complex<double> *C,
const int ldc)
112 cblas_zgemm(Order, TransA, TransB, M, N, K, &alpha, A, lda, B, ldb, &beta, C, ldc);
116 void syrk<float>(
const CBLAS_ORDER Order,
const CBLAS_UPLO Uplo,
117 const CBLAS_TRANSPOSE Trans,
const int N,
const int K,
118 const float & alpha,
const float *A,
const int lda,
119 const float & beta,
float *C,
const int ldc)
121 cblas_ssyrk( Order, Uplo, Trans, N, K, alpha, A, lda, beta, C, ldc);
125 void syrk<double>(
const CBLAS_ORDER Order,
const CBLAS_UPLO Uplo,
126 const CBLAS_TRANSPOSE Trans,
const int N,
const int K,
127 const double & alpha,
const double *A,
const int lda,
128 const double & beta,
double *C,
const int ldc)
130 cblas_dsyrk( Order, Uplo, Trans, N, K, alpha, A, lda, beta, C, ldc);
134 void syrk<std::complex<float> >(
const CBLAS_ORDER Order,
const CBLAS_UPLO Uplo,
135 const CBLAS_TRANSPOSE Trans,
const int N,
const int K,
136 const std::complex<float> & alpha,
const std::complex<float> *A,
const int lda,
137 const std::complex<float> & beta, std::complex<float> *C,
const int ldc)
139 cblas_csyrk( Order, Uplo, Trans, N, K, &alpha, A, lda, &beta, C, ldc);
143 void syrk<std::complex<double> >(
const CBLAS_ORDER Order,
const CBLAS_UPLO Uplo,
144 const CBLAS_TRANSPOSE Trans,
const int N,
const int K,
145 const std::complex<double> & alpha,
const std::complex<double> *A,
const int lda,
146 const std::complex<double> & beta, std::complex<double> *C,
const int ldc)
148 cblas_zsyrk( Order, Uplo, Trans, N, K, &alpha, A, lda, &beta, C, ldc);
Declares and defines templatized wrappers for the BLAS.