imate
C++/CUDA Reference
cublas_interface.cu
Go to the documentation of this file.
1 /*
2  * SPDX-FileCopyrightText: Copyright 2021, Siavash Ameli <sameli@berkeley.edu>
3  * SPDX-License-Identifier: BSD-3-Clause
4  * SPDX-FileType: SOURCE
5  *
6  * This program is free software: you can redistribute it and/or modify it
7  * under the terms of the license found in the LICENSE.txt file in the root
8  * directory of this source tree.
9  */
10 
11 
12 // =======
13 // Headers
14 // =======
15 
16 #include "./cublas_interface.h"
17 
18 
19 // ================
20 // cublas interface
21 // ================
22 
28 
30 {
31 
32  // ===========
33  // cublasXgemv (float)
34  // ===========
35 
38 
39  template<>
40  cublasStatus_t cublasXgemv<float>(
41  cublasHandle_t handle,
42  cublasOperation_t trans,
43  int m,
44  int n,
45  const float* alpha,
46  const float* A,
47  int lda,
48  const float* x,
49  int incx,
50  const float* beta,
51  float* y,
52  int incy)
53  {
54  return cublasSgemv(handle, trans, m, n, alpha, A, lda, x, incx, beta,
55  y, incy);
56  }
57 
58 
59  // ===========
60  // cublasXgemv (double)
61  // ===========
62 
65 
66  template<>
67  cublasStatus_t cublasXgemv<double>(
68  cublasHandle_t handle,
69  cublasOperation_t trans,
70  int m,
71  int n,
72  const double* alpha,
73  const double* A,
74  int lda,
75  const double* x,
76  int incx,
77  const double* beta,
78  double* y,
79  int incy)
80  {
81  return cublasDgemv(handle, trans, m, n, alpha, A, lda, x, incx, beta,
82  y, incy);
83  }
84 
85 
86 
87  // ===========
88  // cublasXcopy (float)
89  // ===========
90 
93 
94  template<>
95  cublasStatus_t cublasXcopy<float>(
96  cublasHandle_t handle,
97  int n,
98  const float* x,
99  int incx,
100  float* y,
101  int incy)
102  {
103  return cublasScopy(handle, n, x, incx, y, incy);
104  }
105 
106 
107  // ===========
108  // cublasXcopy (double)
109  // ===========
110 
113 
114  template<>
115  cublasStatus_t cublasXcopy<double>(
116  cublasHandle_t handle,
117  int n,
118  const double* x,
119  int incx,
120  double* y,
121  int incy)
122  {
123  return cublasDcopy(handle, n, x, incx, y, incy);
124  }
125 
126 
127  // ===========
128  // cublasXaxpy (float)
129  // ===========
130 
133 
134  template<>
135  cublasStatus_t cublasXaxpy<float>(
136  cublasHandle_t handle,
137  int n,
138  const float *alpha,
139  const float *x,
140  int incx,
141  float *y,
142  int incy)
143  {
144  return cublasSaxpy(handle, n, alpha, x, incx, y, incy);
145  }
146 
147 
148  // ===========
149  // cublasXaxpy (double)
150  // ===========
151 
154 
155  template<>
156  cublasStatus_t cublasXaxpy<double>(
157  cublasHandle_t handle,
158  int n,
159  const double *alpha,
160  const double *x,
161  int incx,
162  double *y,
163  int incy)
164  {
165  return cublasDaxpy(handle, n, alpha, x, incx, y, incy);
166  }
167 
168 
169  // ==========
170  // cublasXdot (float)
171  // ==========
172 
175 
176  template<>
177  cublasStatus_t cublasXdot<float>(
178  cublasHandle_t handle,
179  int n,
180  const float *x,
181  int incx,
182  const float *y,
183  int incy,
184  float *result)
185  {
186  return cublasSdot(handle, n, x, incx, y, incy, result);
187  }
188 
189 
190  // ==========
191  // cublasXdot (double)
192  // ==========
193 
196 
197  template<>
198  cublasStatus_t cublasXdot<double>(
199  cublasHandle_t handle,
200  int n,
201  const double *x,
202  int incx,
203  const double *y,
204  int incy,
205  double *result)
206  {
207  return cublasDdot(handle, n, x, incx, y, incy, result);
208  }
209 
210 
211  // ===========
212  // cublasXnrm2 (float)
213  // ===========
214 
217 
218  template<>
219  cublasStatus_t cublasXnrm2<float>(
220  cublasHandle_t handle,
221  int n,
222  const float *x,
223  int incx,
224  float *result)
225  {
226  return cublasSnrm2(handle, n, x, incx, result);
227  }
228 
229 
230  // ===========
231  // cublasXnrm2 (double)
232  // ===========
233 
236 
237  template<>
238  cublasStatus_t cublasXnrm2<double>(
239  cublasHandle_t handle,
240  int n,
241  const double *x,
242  int incx,
243  double *result)
244  {
245  return cublasDnrm2(handle, n, x, incx, result);
246  }
247 
248 
249  // ===========
250  // cublasXscal (float)
251  // ===========
252 
255 
256  template<>
257  cublasStatus_t cublasXscal<float>(
258  cublasHandle_t handle,
259  int n,
260  const float *alpha,
261  float *x,
262  int incx)
263  {
264  return cublasSscal(handle, n, alpha, x, incx);
265  }
266 
267 
268  // ===========
269  // cublasXscal (double)
270  // ===========
271 
274 
275  template<>
276  cublasStatus_t cublasXscal<double>(
277  cublasHandle_t handle,
278  int n,
279  const double *alpha,
280  double *x,
281  int incx)
282  {
283  return cublasDscal(handle, n, alpha, x, incx);
284  }
285 } // namespace cublas_interface
cublasStatus_t cublasDcopy(cublasHandle_t handle, int n, const double *x, int incx, double *y, int incy)
Definition of CUDA's cublasDcopy function using dynamically loaded cublas library.
cublasStatus_t cublasSscal(cublasHandle_t handle, int n, const float *alpha, float *x, int incx)
Definition of CUDA's cublasSscal function using dynamically loaded cublas library.
cublasStatus_t cublasDscal(cublasHandle_t handle, int n, const double *alpha, double *x, int incx)
Definition of CUDA's cublasDscal function using dynamically loaded cublas library.
cublasStatus_t cublasSdot(cublasHandle_t handle, int n, const float *x, int incx, const float *y, int incy, float *result)
Definition of CUDA's cublasSdot function using dynamically loaded cublas library.
cublasStatus_t cublasSnrm2(cublasHandle_t handle, int n, const float *x, int incx, float *result)
Definition of CUDA's cublasSnrm2 function using dynamically loaded cublas library.
cublasStatus_t cublasSaxpy(cublasHandle_t handle, int n, const float *alpha, const float *x, int incx, float *y, int incy)
Definition of CUDA's cublasSaxpy function using dynamically loaded cublas library.
cublasStatus_t cublasDaxpy(cublasHandle_t handle, int n, const double *alpha, const double *x, int incx, double *y, int incy)
Definition of CUDA's cublasDaxpy function using dynamically loaded cublas library.
cublasStatus_t cublasScopy(cublasHandle_t handle, int n, const float *x, int incx, float *y, int incy)
Definition of CUDA's cublasScopy function using dynamically loaded cublas library.
cublasStatus_t cublasDnrm2(cublasHandle_t handle, int n, const double *x, int incx, double *result)
Definition of CUDA's cublasDnrm2 function using dynamically loaded cublas library.
cublasStatus_t cublasDdot(cublasHandle_t handle, int n, const double *x, int incx, const double *y, int incy, double *result)
Definition of CUDA's cublasDdot function using dynamically loaded cublas library.
A collection of templates to wrapper cublas functions.
cublasStatus_t cublasXgemv< float >(cublasHandle_t handle, cublasOperation_t trans, int m, int n, const float *alpha, const float *A, int lda, const float *x, int incx, const float *beta, float *y, int incy)
A template wrapper for cublasSgemv.
cublasStatus_t cublasXnrm2< float >(cublasHandle_t handle, int n, const float *x, int incx, float *result)
A template wrapper to cublasSnrm2.
cublasStatus_t cublasXaxpy< double >(cublasHandle_t handle, int n, const double *alpha, const double *x, int incx, double *y, int incy)
A template wrapper for cublasDaxpy.
cublasStatus_t cublasXaxpy< float >(cublasHandle_t handle, int n, const float *alpha, const float *x, int incx, float *y, int incy)
A template wrapper for cublasSaxpy.
cublasStatus_t cublasXscal< double >(cublasHandle_t handle, int n, const double *alpha, double *x, int incx)
A template wrapper for cublasDscal.
cublasStatus_t cublasXcopy< double >(cublasHandle_t handle, int n, const double *x, int incx, double *y, int incy)
A template wrapper for cublasDcopy.
cublasStatus_t cublasXscal< float >(cublasHandle_t handle, int n, const float *alpha, float *x, int incx)
A template wrapper for cublasSscal.
cublasStatus_t cublasXdot< float >(cublasHandle_t handle, int n, const float *x, int incx, const float *y, int incy, float *result)
A template wrapper for cublasSdot.
cublasStatus_t cublasXgemv< double >(cublasHandle_t handle, cublasOperation_t trans, int m, int n, const double *alpha, const double *A, int lda, const double *x, int incx, const double *beta, double *y, int incy)
A template wrapper for cublasDgemv.
cublasStatus_t cublasXcopy< float >(cublasHandle_t handle, int n, const float *x, int incx, float *y, int incy)
A template wrapper for cublasScopy.
cublasStatus_t cublasXdot< double >(cublasHandle_t handle, int n, const double *x, int incx, const double *y, int incy, double *result)
A template wrapper for cublasDdot.
cublasStatus_t cublasXnrm2< double >(cublasHandle_t handle, int n, const double *x, int incx, double *result)
A template wrapper to cublasDnrm2.