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