ABACUS develop
Atomic-orbital Based Ab-initio Computation at UStc
Loading...
Searching...
No Matches
blas_connector.h
Go to the documentation of this file.
1#ifndef BLAS_CONNECTOR_H
2#define BLAS_CONNECTOR_H
3
4#include <complex>
6#include "../macros.h"
7
8// These still need to be linked in the header file
9// Because quite a lot of code will directly use the original cblas kernels.
10
11extern "C"
12{
13 // level 1: std::vector-std::vector operations, O(n) data and O(n) work.
14
15 // Peize Lin add ?scal 2016-08-04, to compute x=a*x
16 void sscal_(const int *N, const float *alpha, float *X, const int *incX);
17 void dscal_(const int *N, const double *alpha, double *X, const int *incX);
18 void cscal_(const int *N, const std::complex<float> *alpha, std::complex<float> *X, const int *incX);
19 void zscal_(const int *N, const std::complex<double> *alpha, std::complex<double> *X, const int *incX);
20
21 // Peize Lin add ?axpy 2016-08-04, to compute y=a*x+y
22 void saxpy_(const int *N, const float *alpha, const float *X, const int *incX, float *Y, const int *incY);
23 void daxpy_(const int *N, const double *alpha, const double *X, const int *incX, double *Y, const int *incY);
24 void caxpy_(const int *N, const std::complex<float> *alpha, const std::complex<float> *X, const int *incX, std::complex<float> *Y, const int *incY);
25 void zaxpy_(const int *N, const std::complex<double> *alpha, const std::complex<double> *X, const int *incX, std::complex<double> *Y, const int *incY);
26
27 void dcopy_(long const *n, const double *a, int const *incx, double *b, int const *incy);
28 void zcopy_(long const *n, const std::complex<double> *a, int const *incx, std::complex<double> *b, int const *incy);
29
30 //reason for passing results as argument instead of returning it:
31 //see https://www.numbercrunch.de/blog/2014/07/lost-in-translation/
32 // void zdotc_(std::complex<double> *result, const int *n, const std::complex<double> *zx,
33 // const int *incx, const std::complex<double> *zy, const int *incy);
34 // Peize Lin add ?dot 2017-10-27, to compute d=x*y
35 float sdot_(const int *N, const float *X, const int *incX, const float *Y, const int *incY);
36 double ddot_(const int *N, const double *X, const int *incX, const double *Y, const int *incY);
37
38 // Peize Lin add ?nrm2 2018-06-12, to compute out = ||x||_2 = \sqrt{ \sum_i x_i**2 }
39 float snrm2_( const int *n, const float *X, const int *incX );
40 double dnrm2_( const int *n, const double *X, const int *incX );
41 double dznrm2_( const int *n, const std::complex<double> *X, const int *incX );
42
43 // symmetric rank-k update
44 void dsyrk_(
45 const char* uplo,
46 const char* trans,
47 const int* n,
48 const int* k,
49 const double* alpha,
50 const double* a,
51 const int* lda,
52 const double* beta,
53 double* c,
54 const int* ldc
55 );
56
57 // level 2: matrix-std::vector operations, O(n^2) data and O(n^2) work.
58 void sgemv_(const char*const transa, const int*const m, const int*const n,
59 const float*const alpha, const float*const a, const int*const lda, const float*const x, const int*const incx,
60 const float*const beta, float*const y, const int*const incy);
61 void dgemv_(const char*const transa, const int*const m, const int*const n,
62 const double*const alpha, const double*const a, const int*const lda, const double*const x, const int*const incx,
63 const double*const beta, double*const y, const int*const incy);
64
65 void cgemv_(const char *trans, const int *m, const int *n, const std::complex<float> *alpha,
66 const std::complex<float> *a, const int *lda, const std::complex<float> *x, const int *incx,
67 const std::complex<float> *beta, std::complex<float> *y, const int *incy);
68
69 void zgemv_(const char *trans, const int *m, const int *n, const std::complex<double> *alpha,
70 const std::complex<double> *a, const int *lda, const std::complex<double> *x, const int *incx,
71 const std::complex<double> *beta, std::complex<double> *y, const int *incy);
72
73 void dsymv_(const char *uplo, const int *n,
74 const double *alpha, const double *a, const int *lda,
75 const double *x, const int *incx,
76 const double *beta, double *y, const int *incy);
77
78 // A := alpha x * y.T + A
79 void dger_(const int* m,
80 const int* n,
81 const double* alpha,
82 const double* x,
83 const int* incx,
84 const double* y,
85 const int* incy,
86 double* a,
87 const int* lda);
88 void zgerc_(const int* m,
89 const int* n,
90 const std::complex<double>* alpha,
91 const std::complex<double>* x,
92 const int* incx,
93 const std::complex<double>* y,
94 const int* incy,
95 std::complex<double>* a,
96 const int* lda);
97
98 // level 3: matrix-matrix operations, O(n^2) data and O(n^3) work.
99
100 // Peize Lin add ?gemm 2017-10-27, to compute C = a * A.? * B.? + b * C
101 // A is general
102 void sgemm_(const char *transa, const char *transb, const int *m, const int *n, const int *k,
103 const float *alpha, const float *a, const int *lda, const float *b, const int *ldb,
104 const float *beta, float *c, const int *ldc);
105 void dgemm_(const char *transa, const char *transb, const int *m, const int *n, const int *k,
106 const double *alpha, const double *a, const int *lda, const double *b, const int *ldb,
107 const double *beta, double *c, const int *ldc);
108 void cgemm_(const char *transa, const char *transb, const int *m, const int *n, const int *k,
109 const std::complex<float> *alpha, const std::complex<float> *a, const int *lda, const std::complex<float> *b, const int *ldb,
110 const std::complex<float> *beta, std::complex<float> *c, const int *ldc);
111 void zgemm_(const char *transa, const char *transb, const int *m, const int *n, const int *k,
112 const std::complex<double> *alpha, const std::complex<double> *a, const int *lda, const std::complex<double> *b, const int *ldb,
113 const std::complex<double> *beta, std::complex<double> *c, const int *ldc);
114
115 // A is symmetric. C = a * A.? * B.? + b * C
116 void ssymm_(const char *side, const char *uplo, const int *m, const int *n,
117 const float *alpha, const float *a, const int *lda, const float *b, const int *ldb,
118 const float *beta, float *c, const int *ldc);
119 void dsymm_(const char *side, const char *uplo, const int *m, const int *n,
120 const double *alpha, const double *a, const int *lda, const double *b, const int *ldb,
121 const double *beta, double *c, const int *ldc);
122 void csymm_(const char *side, const char *uplo, const int *m, const int *n,
123 const std::complex<float> *alpha, const std::complex<float> *a, const int *lda, const std::complex<float> *b, const int *ldb,
124 const std::complex<float> *beta, std::complex<float> *c, const int *ldc);
125 void zsymm_(const char *side, const char *uplo, const int *m, const int *n,
126 const std::complex<double> *alpha, const std::complex<double> *a, const int *lda, const std::complex<double> *b, const int *ldb,
127 const std::complex<double> *beta, std::complex<double> *c, const int *ldc);
128
129 // A is hermitian. C = a * A.? * B.? + b * C
130 void chemm_(char *side, char *uplo, int *m, int *n,std::complex<float> *alpha,
131 std::complex<float> *a, int *lda, std::complex<float> *b, int *ldb, std::complex<float> *beta, std::complex<float> *c, int *ldc);
132 void zhemm_(char *side, char *uplo, int *m, int *n,std::complex<double> *alpha,
133 std::complex<double> *a, int *lda, std::complex<double> *b, int *ldb, std::complex<double> *beta, std::complex<double> *c, int *ldc);
134
135 //solving triangular matrix with multiple right hand sides
136 void dtrsm_(char *side, char* uplo, char *transa, char *diag, int *m, int *n,
137 double* alpha, double* a, int *lda, double*b, int *ldb);
138 void ztrsm_(char *side, char* uplo, char *transa, char *diag, int *m, int *n,
139 std::complex<double>* alpha, std::complex<double>* a, int *lda, std::complex<double>*b, int *ldb);
140
141}
142
143// Class BlasConnector provide the connector to fortran lapack routine.
144// The entire function in this class are static and inline function.
145// Usage example: BlasConnector::functionname(parameter list).
147{
148public:
149
150 // Peize Lin add 2016-08-04
151 // y=a*x+y
152 static
153 void axpy( const int n, const float alpha, const float *X, const int incX, float *Y, const int incY, base_device::AbacusDevice_t device_type = base_device::AbacusDevice_t::CpuDevice);
154
155 static
156 void axpy( const int n, const double alpha, const double *X, const int incX, double *Y, const int incY, base_device::AbacusDevice_t device_type = base_device::AbacusDevice_t::CpuDevice);
157
158 static
159 void axpy( const int n, const std::complex<float> alpha, const std::complex<float> *X, const int incX, std::complex<float> *Y, const int incY, base_device::AbacusDevice_t device_type = base_device::AbacusDevice_t::CpuDevice);
160
161 static
162 void axpy( const int n, const std::complex<double> alpha, const std::complex<double> *X, const int incX, std::complex<double> *Y, const int incY, base_device::AbacusDevice_t device_type = base_device::AbacusDevice_t::CpuDevice);
163
164
165 // Peize Lin add 2016-08-04
166 // x=a*x
167 static
168 void scal( const int n, const float alpha, float *X, const int incX, base_device::AbacusDevice_t device_type = base_device::AbacusDevice_t::CpuDevice);
169
170 static
171 void scal( const int n, const double alpha, double *X, const int incX, base_device::AbacusDevice_t device_type = base_device::AbacusDevice_t::CpuDevice);
172
173 static
174 void scal( const int n, const std::complex<float> alpha, std::complex<float> *X, const int incX, base_device::AbacusDevice_t device_type = base_device::AbacusDevice_t::CpuDevice);
175
176 static
177 void scal( const int n, const std::complex<double> alpha, std::complex<double> *X, const int incX, base_device::AbacusDevice_t device_type = base_device::AbacusDevice_t::CpuDevice);
178
179
180 // Peize Lin add 2017-10-27
181 // d=x*y
182 static
183 float dot( const int n, const float*const X, const int incX, const float*const Y, const int incY, base_device::AbacusDevice_t device_type = base_device::AbacusDevice_t::CpuDevice);
184
185 static
186 double dot( const int n, const double*const X, const int incX, const double*const Y, const int incY, base_device::AbacusDevice_t device_type = base_device::AbacusDevice_t::CpuDevice);
187
188 // d=x*y
189 static
190 float dotu( const int n, const float*const X, const int incX, const float*const Y, const int incY, base_device::AbacusDevice_t device_type = base_device::AbacusDevice_t::CpuDevice);
191
192 static
193 double dotu( const int n, const double*const X, const int incX, const double*const Y, const int incY, base_device::AbacusDevice_t device_type = base_device::AbacusDevice_t::CpuDevice);
194
195 static
196 std::complex<float> dotu( const int n, const std::complex<float>*const X, const int incX, const std::complex<float>*const Y, const int incY, base_device::AbacusDevice_t device_type = base_device::AbacusDevice_t::CpuDevice);
197
198 static
199 std::complex<double> dotu( const int n, const std::complex<double>*const X, const int incX, const std::complex<double>*const Y, const int incY, base_device::AbacusDevice_t device_type = base_device::AbacusDevice_t::CpuDevice);
200
201 // d=x.conj()*y
202 static
203 float dotc( const int n, const float*const X, const int incX, const float*const Y, const int incY, base_device::AbacusDevice_t device_type = base_device::AbacusDevice_t::CpuDevice);
204
205 static
206 double dotc( const int n, const double*const X, const int incX, const double*const Y, const int incY, base_device::AbacusDevice_t device_type = base_device::AbacusDevice_t::CpuDevice);
207
208 static
209 std::complex<float> dotc( const int n, const std::complex<float>*const X, const int incX, const std::complex<float>*const Y, const int incY, base_device::AbacusDevice_t device_type = base_device::AbacusDevice_t::CpuDevice);
210
211 static
212 std::complex<double> dotc( const int n, const std::complex<double>*const X, const int incX, const std::complex<double>*const Y, const int incY, base_device::AbacusDevice_t device_type = base_device::AbacusDevice_t::CpuDevice);
213
214 // Peize Lin add 2017-10-27, fix bug trans 2019-01-17
215 // C = a * A.? * B.? + b * C
216 // Row Major by default
217 static
218 void gemm(const char transa, const char transb, const int m, const int n, const int k,
219 const float alpha, const float *a, const int lda, const float *b, const int ldb,
220 const float beta, float *c, const int ldc, base_device::AbacusDevice_t device_type = base_device::AbacusDevice_t::CpuDevice);
221
222 static
223 void gemm(const char transa, const char transb, const int m, const int n, const int k,
224 const double alpha, const double *a, const int lda, const double *b, const int ldb,
225 const double beta, double *c, const int ldc, base_device::AbacusDevice_t device_type = base_device::AbacusDevice_t::CpuDevice);
226
227 static
228 void gemm(const char transa, const char transb, const int m, const int n, const int k,
229 const std::complex<float> alpha, const std::complex<float> *a, const int lda, const std::complex<float> *b, const int ldb,
230 const std::complex<float> beta, std::complex<float> *c, const int ldc, base_device::AbacusDevice_t device_type = base_device::AbacusDevice_t::CpuDevice);
231
232 static
233 void gemm(const char transa, const char transb, const int m, const int n, const int k,
234 const std::complex<double> alpha, const std::complex<double> *a, const int lda, const std::complex<double> *b, const int ldb,
235 const std::complex<double> beta, std::complex<double> *c, const int ldc, base_device::AbacusDevice_t device_type = base_device::AbacusDevice_t::CpuDevice);
236
237 // Col-Major if you need to use it
238
239 static
240 void gemm_cm(const char transa, const char transb, const int m, const int n, const int k,
241 const float alpha, const float *a, const int lda, const float *b, const int ldb,
242 const float beta, float *c, const int ldc, base_device::AbacusDevice_t device_type = base_device::AbacusDevice_t::CpuDevice);
243
244 static
245 void gemm_cm(const char transa, const char transb, const int m, const int n, const int k,
246 const double alpha, const double *a, const int lda, const double *b, const int ldb,
247 const double beta, double *c, const int ldc, base_device::AbacusDevice_t device_type = base_device::AbacusDevice_t::CpuDevice);
248
249 static
250 void gemm_cm(const char transa, const char transb, const int m, const int n, const int k,
251 const std::complex<float> alpha, const std::complex<float> *a, const int lda, const std::complex<float> *b, const int ldb,
252 const std::complex<float> beta, std::complex<float> *c, const int ldc, base_device::AbacusDevice_t device_type = base_device::AbacusDevice_t::CpuDevice);
253
254 static
255 void gemm_cm(const char transa, const char transb, const int m, const int n, const int k,
256 const std::complex<double> alpha, const std::complex<double> *a, const int lda, const std::complex<double> *b, const int ldb,
257 const std::complex<double> beta, std::complex<double> *c, const int ldc, base_device::AbacusDevice_t device_type = base_device::AbacusDevice_t::CpuDevice);
258
259 // side=='L': C = a * A * B + b * C.
260 // side=='R': C = a * B * A + b * C.
261 // A == A^T
262 // Because you cannot pack symm or hemm into a row-major kernel by exchanging parameters, so only col-major functions are provided.
263 static
264 void symm_cm(const char side, const char uplo, const int m, const int n,
265 const float alpha, const float *a, const int lda, const float *b, const int ldb,
266 const float beta, float *c, const int ldc, base_device::AbacusDevice_t device_type = base_device::AbacusDevice_t::CpuDevice);
267
268 static
269 void symm_cm(const char side, const char uplo, const int m, const int n,
270 const double alpha, const double *a, const int lda, const double *b, const int ldb,
271 const double beta, double *c, const int ldc, base_device::AbacusDevice_t device_type = base_device::AbacusDevice_t::CpuDevice);
272
273 static
274 void symm_cm(const char side, const char uplo, const int m, const int n,
275 const std::complex<float> alpha, const std::complex<float> *a, const int lda, const std::complex<float> *b, const int ldb,
276 const std::complex<float> beta, std::complex<float> *c, const int ldc, base_device::AbacusDevice_t device_type = base_device::AbacusDevice_t::CpuDevice);
277
278 static
279 void symm_cm(const char side, const char uplo, const int m, const int n,
280 const std::complex<double> alpha, const std::complex<double> *a, const int lda, const std::complex<double> *b, const int ldb,
281 const std::complex<double> beta, std::complex<double> *c, const int ldc, base_device::AbacusDevice_t device_type = base_device::AbacusDevice_t::CpuDevice);
282
283 // side=='L': C = a * A * B + b * C.
284 // side=='R': C = a * B * A + b * C.
285 // A == A^H
286 static
287 void hemm_cm(const char side, const char uplo, const int m, const int n,
288 const float alpha, const float *a, const int lda, const float *b, const int ldb,
289 const float beta, float *c, const int ldc, base_device::AbacusDevice_t device_type = base_device::AbacusDevice_t::CpuDevice);
290
291 static
292 void hemm_cm(const char side, const char uplo, const int m, const int n,
293 const double alpha, const double *a, const int lda, const double *b, const int ldb,
294 const double beta, double *c, const int ldc, base_device::AbacusDevice_t device_type = base_device::AbacusDevice_t::CpuDevice);
295
296 static
297 void hemm_cm(char side, char uplo, int m, int n,
298 std::complex<float> alpha, std::complex<float> *a, int lda, std::complex<float> *b, int ldb,
299 std::complex<float> beta, std::complex<float> *c, int ldc, base_device::AbacusDevice_t device_type = base_device::AbacusDevice_t::CpuDevice);
300
301 static
302 void hemm_cm(char side, char uplo, int m, int n,
303 std::complex<double> alpha, std::complex<double> *a, int lda, std::complex<double> *b, int ldb,
304 std::complex<double> beta, std::complex<double> *c, int ldc, base_device::AbacusDevice_t device_type = base_device::AbacusDevice_t::CpuDevice);
305
306 // y = A*x + beta*y
307 static
308 void gemv(const char trans, const int m, const int n,
309 const float alpha, const float* A, const int lda, const float* X, const int incx,
310 const float beta, float* Y, const int incy, base_device::AbacusDevice_t device_type = base_device::AbacusDevice_t::CpuDevice);
311
312 static
313 void gemv(const char trans, const int m, const int n,
314 const double alpha, const double* A, const int lda, const double* X, const int incx,
315 const double beta, double* Y, const int incy, base_device::AbacusDevice_t device_type = base_device::AbacusDevice_t::CpuDevice);
316
317 static
318 void gemv(const char trans, const int m, const int n,
319 const std::complex<float> alpha, const std::complex<float> *A, const int lda, const std::complex<float> *X, const int incx,
320 const std::complex<float> beta, std::complex<float> *Y, const int incy, base_device::AbacusDevice_t device_type = base_device::AbacusDevice_t::CpuDevice);
321
322 static
323 void gemv(const char trans, const int m, const int n,
324 const std::complex<double> alpha, const std::complex<double> *A, const int lda, const std::complex<double> *X, const int incx,
325 const std::complex<double> beta, std::complex<double> *Y, const int incy, base_device::AbacusDevice_t device_type = base_device::AbacusDevice_t::CpuDevice);
326
327 // Peize Lin add 2018-06-12
328 // out = ||x||_2
329 static
330 float nrm2( const int n, const float *X, const int, base_device::AbacusDevice_t device_type = base_device::AbacusDevice_t::CpuDevice );
331
332 static
333 double nrm2( const int n, const double *X, const int incX, base_device::AbacusDevice_t device_type = base_device::AbacusDevice_t::CpuDevice );
334
335 static
336 double nrm2( const int n, const std::complex<double> *X, const int incX, base_device::AbacusDevice_t device_type = base_device::AbacusDevice_t::CpuDevice );
337
338
339 // copies a into b
340 static
341 void copy(const long n, const double *a, const int incx, double *b, const int incy, base_device::AbacusDevice_t device_type = base_device::AbacusDevice_t::CpuDevice);
342
343 static
344 void copy(const long n, const std::complex<double> *a, const int incx, std::complex<double> *b, const int incy, base_device::AbacusDevice_t device_type = base_device::AbacusDevice_t::CpuDevice);
345
346 // There is some other operators needed, so implemented manually here
347 template <typename T>
348 static
349 void vector_mul_vector(const int& dim, T* result, const T* vector1, const T* vector2, base_device::AbacusDevice_t device_type = base_device::AbacusDevice_t::CpuDevice);
350
351 template <typename T>
352 static
353 void vector_div_vector(const int& dim, T* result, const T* vector1, const T* vector2, base_device::AbacusDevice_t device_type = base_device::AbacusDevice_t::CpuDevice);
354
355 // y = alpha * x + beta * y
356 static
357 void vector_add_vector(const int& dim, float *result, const float *vector1, const float constant1, const float *vector2, const float constant2, base_device::AbacusDevice_t device_type = base_device::AbacusDevice_t::CpuDevice);
358
359 static
360 void vector_add_vector(const int& dim, double *result, const double *vector1, const double constant1, const double *vector2, const double constant2, base_device::AbacusDevice_t device_type = base_device::AbacusDevice_t::CpuDevice);
361
362 static
363 void vector_add_vector(const int& dim, std::complex<float> *result, const std::complex<float> *vector1, const float constant1, const std::complex<float> *vector2, const float constant2, base_device::AbacusDevice_t device_type = base_device::AbacusDevice_t::CpuDevice);
364
365 static
366 void vector_add_vector(const int& dim, std::complex<double> *result, const std::complex<double> *vector1, const double constant1, const std::complex<double> *vector2, const double constant2, base_device::AbacusDevice_t device_type = base_device::AbacusDevice_t::CpuDevice);
367};
368
369#ifdef __CUDA
370
371#include <cuda_runtime.h>
372#include "cublas_v2.h"
373
374// If you want to use cublas, you need these functions to create and destroy the cublas/hipblas handle.
375// You also need to use these functions to translate the transpose parameter into cublas/hipblas datatype.
376
377namespace BlasUtils{
378
379 static cublasHandle_t cublas_handle = nullptr;
380
381 void createGpuBlasHandle(); // Create a cublas/hipblas handle.
382
383 void destoryBLAShandle(); // Destroy the cublas/hipblas handle. Do this when the software is about to end.
384
385 cublasOperation_t judge_trans(bool is_complex, const char& trans, const char* name); // Translate a normal transpose parameter to a cublas/hipblas type.
386
387 cublasSideMode_t judge_side(const char& trans); // Translate a normal side parameter to a cublas/hipblas type.
388
389 cublasFillMode_t judge_fill(const char& trans); // Translate a normal fill parameter to a cublas/hipblas type.
390
391}
392
393#endif
394
395// If GATHER_INFO is defined, the original function is replaced with a "i" suffix,
396// preventing changes on the original code.
397// The real function call is at gather_math_lib_info.cpp
398#ifdef GATHER_INFO
399
400#define zgemm_ zgemm_i
401void zgemm_i(const char *transa,
402 const char *transb,
403 const int *m,
404 const int *n,
405 const int *k,
406 const std::complex<double> *alpha,
407 const std::complex<double> *a,
408 const int *lda,
409 const std::complex<double> *b,
410 const int *ldb,
411 const std::complex<double> *beta,
412 std::complex<double> *c,
413 const int *ldc);
414
415#define zaxpy_ zaxpy_i
416void zaxpy_i(const int *N,
417 const std::complex<double> *alpha,
418 const std::complex<double> *X,
419 const int *incX,
420 std::complex<double> *Y,
421 const int *incY);
422
423/*
424#define zgemv_ zgemv_i
425
426void zgemv_i(const char *trans,
427 const int *m,
428 const int *n,
429 const std::complex<double> *alpha,
430 const std::complex<double> *a,
431 const int *lda,
432 const std::complex<double> *x,
433 const int *incx,
434 const std::complex<double> *beta,
435 std::complex<double> *y,
436 const int *incy);
437*/
438
439#endif // GATHER_INFO
440#endif // BLAS_CONNECTOR_H
void chemm_(char *side, char *uplo, int *m, int *n, std::complex< float > *alpha, std::complex< float > *a, int *lda, std::complex< float > *b, int *ldb, std::complex< float > *beta, std::complex< float > *c, int *ldc)
double ddot_(const int *N, const double *X, const int *incX, const double *Y, const int *incY)
void dscal_(const int *N, const double *alpha, double *X, const int *incX)
void cscal_(const int *N, const std::complex< float > *alpha, std::complex< float > *X, const int *incX)
void zgerc_(const int *m, const int *n, const std::complex< double > *alpha, const std::complex< double > *x, const int *incx, const std::complex< double > *y, const int *incy, std::complex< double > *a, const int *lda)
void sgemm_(const char *transa, const char *transb, const int *m, const int *n, const int *k, const float *alpha, const float *a, const int *lda, const float *b, const int *ldb, const float *beta, float *c, const int *ldc)
void dgemv_(const char *const transa, const int *const m, const int *const n, const double *const alpha, const double *const a, const int *const lda, const double *const x, const int *const incx, const double *const beta, double *const y, const int *const incy)
void sscal_(const int *N, const float *alpha, float *X, const int *incX)
void daxpy_(const int *N, const double *alpha, const double *X, const int *incX, double *Y, const int *incY)
void dsymm_(const char *side, const char *uplo, const int *m, const int *n, const double *alpha, const double *a, const int *lda, const double *b, const int *ldb, const double *beta, double *c, const int *ldc)
void csymm_(const char *side, const char *uplo, const int *m, const int *n, const std::complex< float > *alpha, const std::complex< float > *a, const int *lda, const std::complex< float > *b, const int *ldb, const std::complex< float > *beta, std::complex< float > *c, const int *ldc)
void zgemm_(const char *transa, const char *transb, const int *m, const int *n, const int *k, const std::complex< double > *alpha, const std::complex< double > *a, const int *lda, const std::complex< double > *b, const int *ldb, const std::complex< double > *beta, std::complex< double > *c, const int *ldc)
void dger_(const int *m, const int *n, const double *alpha, const double *x, const int *incx, const double *y, const int *incy, double *a, const int *lda)
void zscal_(const int *N, const std::complex< double > *alpha, std::complex< double > *X, const int *incX)
void dtrsm_(char *side, char *uplo, char *transa, char *diag, int *m, int *n, double *alpha, double *a, int *lda, double *b, int *ldb)
void zhemm_(char *side, char *uplo, int *m, int *n, std::complex< double > *alpha, std::complex< double > *a, int *lda, std::complex< double > *b, int *ldb, std::complex< double > *beta, std::complex< double > *c, int *ldc)
void zsymm_(const char *side, const char *uplo, const int *m, const int *n, const std::complex< double > *alpha, const std::complex< double > *a, const int *lda, const std::complex< double > *b, const int *ldb, const std::complex< double > *beta, std::complex< double > *c, const int *ldc)
void dsyrk_(const char *uplo, const char *trans, const int *n, const int *k, const double *alpha, const double *a, const int *lda, const double *beta, double *c, const int *ldc)
void zaxpy_(const int *N, const std::complex< double > *alpha, const std::complex< double > *X, const int *incX, std::complex< double > *Y, const int *incY)
void cgemm_(const char *transa, const char *transb, const int *m, const int *n, const int *k, const std::complex< float > *alpha, const std::complex< float > *a, const int *lda, const std::complex< float > *b, const int *ldb, const std::complex< float > *beta, std::complex< float > *c, const int *ldc)
void ztrsm_(char *side, char *uplo, char *transa, char *diag, int *m, int *n, std::complex< double > *alpha, std::complex< double > *a, int *lda, std::complex< double > *b, int *ldb)
void zcopy_(long const *n, const std::complex< double > *a, int const *incx, std::complex< double > *b, int const *incy)
void dcopy_(long const *n, const double *a, int const *incx, double *b, int const *incy)
float sdot_(const int *N, const float *X, const int *incX, const float *Y, const int *incY)
void dgemm_(const char *transa, const char *transb, const int *m, const int *n, const int *k, const double *alpha, const double *a, const int *lda, const double *b, const int *ldb, const double *beta, double *c, const int *ldc)
void dsymv_(const char *uplo, const int *n, const double *alpha, const double *a, const int *lda, const double *x, const int *incx, const double *beta, double *y, const int *incy)
double dznrm2_(const int *n, const std::complex< double > *X, const int *incX)
void caxpy_(const int *N, const std::complex< float > *alpha, const std::complex< float > *X, const int *incX, std::complex< float > *Y, const int *incY)
void ssymm_(const char *side, const char *uplo, const int *m, const int *n, const float *alpha, const float *a, const int *lda, const float *b, const int *ldb, const float *beta, float *c, const int *ldc)
float snrm2_(const int *n, const float *X, const int *incX)
void cgemv_(const char *trans, const int *m, const int *n, const std::complex< float > *alpha, const std::complex< float > *a, const int *lda, const std::complex< float > *x, const int *incx, const std::complex< float > *beta, std::complex< float > *y, const int *incy)
void sgemv_(const char *const transa, const int *const m, const int *const n, const float *const alpha, const float *const a, const int *const lda, const float *const x, const int *const incx, const float *const beta, float *const y, const int *const incy)
void saxpy_(const int *N, const float *alpha, const float *X, const int *incX, float *Y, const int *incY)
double dnrm2_(const int *n, const double *X, const int *incX)
void zgemv_(const char *trans, const int *m, const int *n, const std::complex< double > *alpha, const std::complex< double > *a, const int *lda, const std::complex< double > *x, const int *incx, const std::complex< double > *beta, std::complex< double > *y, const int *incy)
Definition blas_connector.h:147
static void copy(const long n, const double *a, const int incx, double *b, const int incy, base_device::AbacusDevice_t device_type=base_device::AbacusDevice_t::CpuDevice)
Definition blas_connector_vector.cpp:325
static void vector_add_vector(const int &dim, std::complex< double > *result, const std::complex< double > *vector1, const double constant1, const std::complex< double > *vector2, const double constant2, base_device::AbacusDevice_t device_type=base_device::AbacusDevice_t::CpuDevice)
static void vector_div_vector(const int &dim, T *result, const T *vector1, const T *vector2, base_device::AbacusDevice_t device_type=base_device::AbacusDevice_t::CpuDevice)
static void axpy(const int n, const float alpha, const float *X, const int incX, float *Y, const int incY, base_device::AbacusDevice_t device_type=base_device::AbacusDevice_t::CpuDevice)
Definition blas_connector_vector.cpp:18
static void gemm_cm(const char transa, const char transb, const int m, const int n, const int k, const float alpha, const float *a, const int lda, const float *b, const int ldb, const float beta, float *c, const int ldc, base_device::AbacusDevice_t device_type=base_device::AbacusDevice_t::CpuDevice)
Definition blas_connector_matrix.cpp:190
static void gemv(const char trans, const int m, const int n, const float alpha, const float *A, const int lda, const float *X, const int incx, const float beta, float *Y, const int incy, base_device::AbacusDevice_t device_type=base_device::AbacusDevice_t::CpuDevice)
Definition blas_connector_matrix.cpp:501
static void vector_add_vector(const int &dim, double *result, const double *vector1, const double constant1, const double *vector2, const double constant2, base_device::AbacusDevice_t device_type=base_device::AbacusDevice_t::CpuDevice)
static void vector_add_vector(const int &dim, std::complex< float > *result, const std::complex< float > *vector1, const float constant1, const std::complex< float > *vector2, const float constant2, base_device::AbacusDevice_t device_type=base_device::AbacusDevice_t::CpuDevice)
static float nrm2(const int n, const float *X, const int, base_device::AbacusDevice_t device_type=base_device::AbacusDevice_t::CpuDevice)
Definition blas_connector_vector.cpp:271
static float dotc(const int n, const float *const X, const int incX, const float *const Y, const int incY, base_device::AbacusDevice_t device_type=base_device::AbacusDevice_t::CpuDevice)
Definition blas_connector_vector.cpp:224
static void gemm(const char transa, const char transb, const int m, const int n, const int k, const float alpha, const float *a, const int lda, const float *b, const int ldb, const float beta, float *c, const int ldc, base_device::AbacusDevice_t device_type=base_device::AbacusDevice_t::CpuDevice)
Definition blas_connector_matrix.cpp:20
static void symm_cm(const char side, const char uplo, const int m, const int n, const float alpha, const float *a, const int lda, const float *b, const int ldb, const float beta, float *c, const int ldc, base_device::AbacusDevice_t device_type=base_device::AbacusDevice_t::CpuDevice)
Definition blas_connector_matrix.cpp:361
static void hemm_cm(const char side, const char uplo, const int m, const int n, const float alpha, const float *a, const int lda, const float *b, const int ldb, const float beta, float *c, const int ldc, base_device::AbacusDevice_t device_type=base_device::AbacusDevice_t::CpuDevice)
Definition blas_connector_matrix.cpp:445
static void vector_add_vector(const int &dim, float *result, const float *vector1, const float constant1, const float *vector2, const float constant2, base_device::AbacusDevice_t device_type=base_device::AbacusDevice_t::CpuDevice)
static void vector_mul_vector(const int &dim, T *result, const T *vector1, const T *vector2, base_device::AbacusDevice_t device_type=base_device::AbacusDevice_t::CpuDevice)
static void scal(const int n, const float alpha, float *X, const int incX, base_device::AbacusDevice_t device_type=base_device::AbacusDevice_t::CpuDevice)
Definition blas_connector_vector.cpp:80
static float dotu(const int n, const float *const X, const int incX, const float *const Y, const int incY, base_device::AbacusDevice_t device_type=base_device::AbacusDevice_t::CpuDevice)
Definition blas_connector_vector.cpp:177
static float dot(const int n, const float *const X, const int incX, const float *const Y, const int incY, base_device::AbacusDevice_t device_type=base_device::AbacusDevice_t::CpuDevice)
Definition blas_connector_vector.cpp:142
#define N
Definition exp.cpp:24
#define T
Definition exp.cpp:237
void zgemm_i(const char *transa, const char *transb, const int *m, const int *n, const int *k, const std::complex< double > *alpha, const std::complex< double > *a, const int *lda, const std::complex< double > *b, const int *ldb, const std::complex< double > *beta, std::complex< double > *c, const int *ldc)
Definition gather_math_lib_info.cpp:16
void zaxpy_i(const int *N, const std::complex< double > *alpha, const std::complex< double > *X, const int *incX, std::complex< double > *Y, const int *incY)
Definition gather_math_lib_info.cpp:36
AbacusDevice_t
Definition types.h:12
@ CpuDevice
Definition types.h:14