Octopus
clblas_low.c
Go to the documentation of this file.
1#include <config.h>
2
3#ifdef HAVE_CLBLAS
4
5#include <clBLAS.h>
6
7void FC_FUNC_(clblasgetversion_low, CLBLASGETVERSION_LOW)(int * major, int * minor, int * patch, int * status){
8 cl_uint cl_major, cl_minor, cl_patch;
9
10 *status = clblasGetVersion(&cl_major, &cl_minor, &cl_patch);
11 *major = cl_major;
12 *minor = cl_minor;
13 *patch = cl_patch;
14}
15
16void FC_FUNC_(clblassetup_low, CLBLASSETUP_LOW)(int * status){
17 *status = clblasSetup();
18}
19
20void FC_FUNC_(clblasteardown_low, CLBLASTEARDOWN_LOW)(){
21 clblasTeardown();
22}
23
24void FC_FUNC_(clblasdtrsmex_low, CLBLASDTRSMEX_LOW)(int * order,
25 int * side,
26 int * uplo,
27 int * transA,
28 int * diag,
29 cl_long * M,
30 cl_long * N,
31 double * alpha,
32 const cl_mem * A,
33 size_t * offA,
34 size_t * lda,
35 cl_mem * B,
36 size_t * offB,
37 size_t * ldb,
38 cl_command_queue * CommandQueue,
39 int * status){
40
41
42 *status = clblasDtrsm((clblasOrder) *order, (clblasSide) *side, (clblasUplo) *uplo,
43 (clblasTranspose) *transA, (clblasDiag) *diag,
44 (size_t) *M, (size_t) *N, *alpha,
45 *A, (size_t) *offA, (size_t) *lda,
46 *B, (size_t) *offB, (size_t) *ldb,
47 1, CommandQueue, 0, NULL, NULL);
48}
49
50
51void FC_FUNC_(clblasztrsmex_low, CLBLASZTRSMEX_LOW)(int * order,
52 int * side,
53 int * uplo,
54 int * transA,
55 int * diag,
56 cl_long * M,
57 cl_long * N,
58 DoubleComplex * alpha,
59 const cl_mem * A,
60 size_t * offA,
61 size_t * lda,
62 cl_mem * B,
63 size_t * offB,
64 size_t * ldb,
65 cl_command_queue * CommandQueue,
66 int * status){
67
68
69 *status = clblasZtrsm((clblasOrder) *order, (clblasSide) *side, (clblasUplo) *uplo,
70 (clblasTranspose) *transA, (clblasDiag) *diag,
71 (size_t) *M, (size_t) *N, *alpha,
72 *A, (size_t) *offA, (size_t) *lda,
73 *B, (size_t) *offB, (size_t) *ldb,
74 1, CommandQueue, 0, NULL, NULL);
75}
76
77void FC_FUNC_(clblasdgemvex_low, CLBLASDGEMVEX_LOW)(int * order,
78 int * transA,
79 cl_long * M,
80 cl_long * N,
81 double * alpha,
82 const cl_mem * A,
83 cl_long * offA,
84 cl_long * lda,
85 const cl_mem * X,
86 cl_long * offX,
87 cl_long * incx,
88 double * beta,
89 cl_mem * Y,
90 cl_long * offY,
91 cl_long * incy,
92 cl_command_queue * CommandQueue,
93 int * status){
94
95 *status = clblasDgemv((clblasOrder) *order, (clblasTranspose) *transA,
96 (size_t) *M, (size_t) *N, *alpha,
97 *A, (size_t) *offA, (size_t) *lda,
98 *X, (size_t) *offX, (size_t) *incx, *beta,
99 *Y, (size_t) *offY, (size_t) *incy,
100 1, CommandQueue, 0, NULL, NULL);
101}
102
103void FC_FUNC_(clblaszgemvex_low, CLBLASZGEMVEX_LOW)(int * order,
104 int * transA,
105 cl_long * M,
106 cl_long * N,
107 cl_double2 * alpha,
108 const cl_mem * A,
109 cl_long * offA,
110 cl_long * lda,
111 const cl_mem * X,
112 cl_long * offX,
113 cl_long * incx,
114 cl_double2 * beta,
115 cl_mem * Y,
116 cl_long * offY,
117 cl_long * incy,
118 cl_command_queue * CommandQueue,
119 int * status){
120
121 *status = clblasZgemv((clblasOrder) *order, (clblasTranspose) *transA,
122 (size_t) *M, (size_t) *N, *alpha,
123 *A, (size_t) *offA, (size_t) *lda,
124 *X, (size_t) *offX, (size_t) *incx, *beta,
125 *Y, (size_t) *offY, (size_t) *incy,
126 1, CommandQueue, 0, NULL, NULL);
127}
128
129void FC_FUNC_(clblasdgemmex_low, CLBLASDGEMMEX_LOW)(int * order,
130 int * transA,
131 int * transB,
132 cl_long * M,
133 cl_long * N,
134 cl_long * K,
135 double * alpha,
136 const cl_mem * A,
137 cl_long * offA,
138 cl_long * lda,
139 const cl_mem * B,
140 cl_long * offB,
141 cl_long * ldb,
142 double * beta,
143 cl_mem * C,
144 cl_long * offC,
145 cl_long * ldc,
146 cl_command_queue * CommandQueue,
147 int * status){
148
149 *status = clblasDgemm((clblasOrder) *order, (clblasTranspose) *transA, (clblasTranspose) *transB,
150 (size_t) *M, (size_t) *N, (size_t) *K, *alpha,
151 *A, (size_t) *offA, (size_t) *lda,
152 *B, (size_t) *offB, (size_t) *ldb, *beta,
153 *C, (size_t) *offC, (size_t) *ldc,
154 1, CommandQueue, 0, NULL, NULL);
155}
156
157void FC_FUNC_(clblaszgemmex_low, CLBLASZGEMMEX_LOW)(int * order,
158 int * transA,
159 int * transB,
160 cl_long * M,
161 cl_long * N,
162 cl_long * K,
163 DoubleComplex * alpha,
164 const cl_mem * A,
165 cl_long * offA,
166 cl_long * lda,
167 const cl_mem * B,
168 cl_long * offB,
169 cl_long * ldb,
170 DoubleComplex * beta,
171 cl_mem * C,
172 cl_long * offC,
173 cl_long * ldc,
174 cl_command_queue * CommandQueue,
175 int * status){
176
177 *status = clblasZgemm((clblasOrder) *order, (clblasTranspose) *transA, (clblasTranspose) *transB,
178 (size_t) *M, (size_t) *N, (size_t) *K, *alpha,
179 *A, (size_t) *offA, (size_t) *lda,
180 *B, (size_t) *offB, (size_t) *ldb, *beta,
181 *C, (size_t) *offC, (size_t) *ldc,
182 1, CommandQueue, 0, NULL, NULL);
183}
184
185void FC_FUNC_(clblasdsyrkex_low, CLBLASDSYRKEX_LOW)(int * order,
186 int * uplo,
187 int * transA,
188 cl_long * N,
189 cl_long * K,
190 double * alpha,
191 const cl_mem * A,
192 cl_long * offA,
193 cl_long * lda,
194 double * beta,
195 cl_mem * C,
196 cl_long * offC,
197 cl_long * ldc,
198 cl_command_queue * CommandQueue,
199 int * status){
200
201 *status = clblasDsyrk((clblasOrder) *order, (clblasUplo) *uplo, (clblasTranspose) *transA,
202 (size_t) *N, (size_t) *K, *alpha,
203 *A, (size_t) *offA, (size_t) *lda, *beta,
204 *C, (size_t) *offC, (size_t) *ldc,
205 1, CommandQueue, 0, NULL, NULL);
206}
207
208void FC_FUNC_(clblaszherkex_low, CLBLASZHERKEX_LOW)(int * order,
209 int * uplo,
210 int * transA,
211 cl_long * N,
212 cl_long * K,
213 double * alpha,
214 const cl_mem * A,
215 cl_long * offA,
216 cl_long * lda,
217 double * beta,
218 cl_mem * C,
219 cl_long * offC,
220 cl_long * ldc,
221 cl_command_queue * CommandQueue,
222 int * status){
223
224 *status = clblasZherk((clblasOrder) *order, (clblasUplo) *uplo, (clblasTranspose) *transA,
225 (size_t) *N, (size_t) *K, *alpha,
226 *A, (size_t) *offA, (size_t) *lda, *beta,
227 *C, (size_t) *offC, (size_t) *ldc,
228 1, CommandQueue, 0, NULL, NULL);
229}
230
231
232clblasStatus FC_FUNC_(clblasddot_low, CLBLASDDOT_LOW)(cl_long * N,
233 cl_mem * dotProduct,
234 cl_long * offDP,
235 const cl_mem *X,
236 cl_long * offx,
237 int *incx,
238 const cl_mem * Y,
239 cl_long * offy,
240 int * incy,
241 cl_mem * scratchBuff,
242 cl_command_queue * CommandQueue,
243 int * status){
244
245
246 *status = clblasDdot((size_t) *N, *dotProduct, (size_t) *offDP, *X, (size_t) *offx, *incx,
247 *Y, (size_t) *offy, *incy, *scratchBuff, 1, CommandQueue, 0, NULL, NULL);
248
249}
250
251
252clblasStatus FC_FUNC_(clblaszdotc_low, CLBLASZDOTC_LOW)(cl_long * N,
253 cl_mem * dotProduct,
254 cl_long * offDP,
255 const cl_mem *X,
256 cl_long * offx,
257 int *incx,
258 const cl_mem * Y,
259 cl_long * offy,
260 int * incy,
261 cl_mem * scratchBuff,
262 cl_command_queue * CommandQueue,
263 int * status){
264
265
266 *status = clblasZdotc((size_t) *N, *dotProduct, (size_t) *offDP, *X, (size_t) *offx, *incx,
267 *Y, (size_t) *offy, *incy, *scratchBuff, 1, CommandQueue, 0, NULL, NULL);
268
269}
270
271clblasStatus FC_FUNC_(clblaszdotu_low, CLBLASZDOTU_LOW)(cl_long * N,
272 cl_mem * dotProduct,
273 cl_long * offDP,
274 const cl_mem *X,
275 cl_long * offx,
276 int *incx,
277 const cl_mem * Y,
278 cl_long * offy,
279 int * incy,
280 cl_mem * scratchBuff,
281 cl_command_queue * CommandQueue,
282 int * status){
283
284
285 *status = clblasZdotu((size_t) *N, *dotProduct, (size_t) *offDP, *X, (size_t) *offx, *incx,
286 *Y, (size_t) *offy, *incy, *scratchBuff, 1, CommandQueue, 0, NULL, NULL);
287
288}
289
290clblasStatus FC_FUNC_(clblasdnrm2_low, CLBLASDNRM2_LOW)(cl_long * N,
291 cl_mem * NRM2,
292 cl_long * offNRM2,
293 const cl_mem *X,
294 cl_long * offx,
295 int *incx,
296 cl_mem * scratchBuff,
297 cl_command_queue * CommandQueue,
298 int * status){
299
300
301 *status = clblasDnrm2((size_t) *N, *NRM2, (size_t) *offNRM2, *X, (size_t) *offx, *incx,
302 *scratchBuff, 1, CommandQueue, 0, NULL, NULL);
303
304}
305
306clblasStatus FC_FUNC_(clblasdznrm2_low, CLBLASDZNRM2_LOW)(cl_long * N,
307 cl_mem * NRM2,
308 cl_long * offNRM2,
309 const cl_mem *X,
310 cl_long * offx,
311 int *incx,
312 cl_mem * scratchBuff,
313 cl_command_queue * CommandQueue,
314 int * status){
315
316
317 *status = clblasDznrm2((size_t) *N, *NRM2, (size_t) *offNRM2, *X, (size_t) *offx, *incx,
318 *scratchBuff, 1, CommandQueue, 0, NULL, NULL);
319
320}
321
322#endif