Octopus
clblast_low.c
Go to the documentation of this file.
1/* Copyright (C) 2022 N. Tancogne-Dejean
2*
3* This program is free software; you can redistribute it and/or modify
4* it under the terms of the GNU General Public License as published by
5* the Free Software Foundation; either version 2, or (at your option)
6* any later version.
7*
8* This program is distributed in the hope that it will be useful,
9* but WITHOUT ANY WARRANTY; without even the implied warranty of
10* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
11* GNU General Public License for more details.
12*
13* You should have received a copy of the GNU General Public License
14* along with this program; if not, write to the Free Software
15* Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
16* 02110-1301, USA.
17*/
18
19#include <config.h>
20
21#ifdef HAVE_CLBLAST
22
23#include <clblast_c.h>
24
25void FC_FUNC_(clblasdtrsmex_low, CLBLASDTRSMEX_LOW)(int * order,
26 int * side,
27 int * uplo,
28 int * transA,
29 int * diag,
30 cl_long * M,
31 cl_long * N,
32 double * alpha,
33 const cl_mem * A,
34 size_t * offA,
35 size_t * lda,
36 cl_mem * B,
37 size_t * offB,
38 size_t * ldb,
39 cl_command_queue * CommandQueue,
40 int * status){
41
42
43 *status = CLBlastDtrsm((CLBlastLayout) *order, (CLBlastSide) *side, (CLBlastTriangle) *uplo,
44 (CLBlastTranspose) *transA, (CLBlastDiagonal) *diag,
45 (size_t) *M, (size_t) *N, *alpha,
46 *A, (size_t) *offA, (size_t) *lda,
47 *B, (size_t) *offB, (size_t) *ldb,
48 CommandQueue, NULL);
49}
50
51
52void FC_FUNC_(clblasztrsmex_low, CLBLASZTRSMEX_LOW)(int * order,
53 int * side,
54 int * uplo,
55 int * transA,
56 int * diag,
57 cl_long * M,
58 cl_long * N,
59 cl_double2 * alpha,
60 const cl_mem * A,
61 size_t * offA,
62 size_t * lda,
63 cl_mem * B,
64 size_t * offB,
65 size_t * ldb,
66 cl_command_queue * CommandQueue,
67 int * status){
68
69
70 *status = CLBlastZtrsm((CLBlastLayout) *order, (CLBlastSide) *side, (CLBlastTriangle) *uplo,
71 (CLBlastTranspose) *transA, (CLBlastDiagonal) *diag,
72 (size_t) *M, (size_t) *N, *alpha,
73 *A, (size_t) *offA, (size_t) *lda,
74 *B, (size_t) *offB, (size_t) *ldb,
75 CommandQueue, NULL);
76}
77
78void FC_FUNC_(clblasdgemvex_low, CLBLASDGEMVEX_LOW)(int * order,
79 int * transA,
80 cl_long * M,
81 cl_long * N,
82 double * alpha,
83 const cl_mem * A,
84 cl_long * offA,
85 cl_long * lda,
86 const cl_mem * X,
87 cl_long * offX,
88 cl_long * incx,
89 double * beta,
90 cl_mem * Y,
91 cl_long * offY,
92 cl_long * incy,
93 cl_command_queue * CommandQueue,
94 int * status){
95
96 *status = CLBlastDgemv((CLBlastLayout) *order, (CLBlastTranspose) *transA,
97 (size_t) *M, (size_t) *N, *alpha,
98 *A, (size_t) *offA, (size_t) *lda,
99 *X, (size_t) *offX, (size_t) *incx, *beta,
100 *Y, (size_t) *offY, (size_t) *incy,
101 CommandQueue, NULL);
102}
103
104void FC_FUNC_(clblaszgemvex_low, CLBLASZGEMVEX_LOW)(int * order,
105 int * transA,
106 cl_long * M,
107 cl_long * N,
108 cl_double2 * alpha,
109 const cl_mem * A,
110 cl_long * offA,
111 cl_long * lda,
112 const cl_mem * X,
113 cl_long * offX,
114 cl_long * incx,
115 cl_double2 * beta,
116 cl_mem * Y,
117 cl_long * offY,
118 cl_long * incy,
119 cl_command_queue * CommandQueue,
120 int * status){
121
122 *status = CLBlastZgemv((CLBlastLayout) *order, (CLBlastTranspose) *transA,
123 (size_t) *M, (size_t) *N, *alpha,
124 *A, (size_t) *offA, (size_t) *lda,
125 *X, (size_t) *offX, (size_t) *incx, *beta,
126 *Y, (size_t) *offY, (size_t) *incy,
127 CommandQueue, NULL);
128}
129
130
131
132void FC_FUNC_(clblasdgemmex_low, CLBLASDGEMMEX_LOW)(int * order,
133 int * transA,
134 int * transB,
135 cl_long * M,
136 cl_long * N,
137 cl_long * K,
138 double * alpha,
139 const cl_mem * A,
140 cl_long * offA,
141 cl_long * lda,
142 const cl_mem * B,
143 cl_long * offB,
144 cl_long * ldb,
145 double * beta,
146 cl_mem * C,
147 cl_long * offC,
148 cl_long * ldc,
149 cl_command_queue * CommandQueue,
150 int * status){
151
152 *status = CLBlastDgemm((CLBlastLayout) *order, (CLBlastTranspose) *transA, (CLBlastTranspose) *transB,
153 (size_t) *M, (size_t) *N, (size_t) *K, *alpha,
154 *A, (size_t) *offA, (size_t) *lda,
155 *B, (size_t) *offB, (size_t) *ldb, *beta,
156 *C, (size_t) *offC, (size_t) *ldc,
157 CommandQueue, NULL);
158}
159
160void FC_FUNC_(clblaszgemmex_low, CLBLASZGEMMEX_LOW)(int * order,
161 int * transA,
162 int * transB,
163 cl_long * M,
164 cl_long * N,
165 cl_long * K,
166 cl_double2 * alpha,
167 const cl_mem * A,
168 cl_long * offA,
169 cl_long * lda,
170 const cl_mem * B,
171 cl_long * offB,
172 cl_long * ldb,
173 cl_double2 * beta,
174 cl_mem * C,
175 cl_long * offC,
176 cl_long * ldc,
177 cl_command_queue * CommandQueue,
178 int * status){
179
180 *status = CLBlastZgemm((CLBlastLayout) *order, (CLBlastTranspose) *transA, (CLBlastTranspose) *transB,
181 (size_t) *M, (size_t) *N, (size_t) *K, *alpha,
182 *A, (size_t) *offA, (size_t) *lda,
183 *B, (size_t) *offB, (size_t) *ldb, *beta,
184 *C, (size_t) *offC, (size_t) *ldc,
185 CommandQueue, NULL);
186}
187
188void FC_FUNC_(clblasdsyrkex_low, CLBLASDSYRKEX_LOW)(int * order,
189 int * uplo,
190 int * transA,
191 cl_long * N,
192 cl_long * K,
193 double * alpha,
194 const cl_mem * A,
195 cl_long * offA,
196 cl_long * lda,
197 double * beta,
198 cl_mem * C,
199 cl_long * offC,
200 cl_long * ldc,
201 cl_command_queue * CommandQueue,
202 int * status){
203
204 *status = CLBlastDsyrk((CLBlastLayout) *order, (CLBlastTriangle) *uplo, (CLBlastTranspose) *transA,
205 (size_t) *N, (size_t) *K, *alpha,
206 *A, (size_t) *offA, (size_t) *lda, *beta,
207 *C, (size_t) *offC, (size_t) *ldc,
208 CommandQueue, NULL);
209}
210
211void FC_FUNC_(clblaszherkex_low, CLBLASZHERKEX_LOW)(int * order,
212 int * uplo,
213 int * transA,
214 cl_long * N,
215 cl_long * K,
216 double * alpha,
217 const cl_mem * A,
218 cl_long * offA,
219 cl_long * lda,
220 double * beta,
221 cl_mem * C,
222 cl_long * offC,
223 cl_long * ldc,
224 cl_command_queue * CommandQueue,
225 int * status){
226
227 *status = CLBlastZherk((CLBlastLayout) *order, (CLBlastTriangle) *uplo, (CLBlastTranspose) *transA,
228 (size_t) *N, (size_t) *K, *alpha,
229 *A, (size_t) *offA, (size_t) *lda, *beta,
230 *C, (size_t) *offC, (size_t) *ldc,
231 CommandQueue, NULL);
232}
233
234
235CLBlastStatusCode FC_FUNC_(clblasddot_low, CLBLASDDOT_LOW)(cl_long * N,
236 cl_mem * dotProduct,
237 cl_long * offDP,
238 const cl_mem *X,
239 cl_long * offx,
240 int *incx,
241 const cl_mem * Y,
242 cl_long * offy,
243 int * incy,
244 cl_mem * scratchBuff,
245 cl_command_queue * CommandQueue,
246 int * status){
247
248
249 *status = CLBlastDdot((size_t) *N, *dotProduct, (size_t) *offDP, *X, (size_t) *offx, *incx,
250 *Y, (size_t) *offy, *incy, CommandQueue, NULL);
251
252}
253
254
255CLBlastStatusCode FC_FUNC_(clblaszdotc_low, CLBLASZDOTc_LOW)(cl_long * N,
256 cl_mem * dotProduct,
257 cl_long * offDP,
258 const cl_mem *X,
259 cl_long * offx,
260 int *incx,
261 const cl_mem * Y,
262 cl_long * offy,
263 int * incy,
264 cl_mem * scratchBuff,
265 cl_command_queue * CommandQueue,
266 int * status){
267
268
269 *status = CLBlastZdotc((size_t) *N, *dotProduct, (size_t) *offDP, *X, (size_t) *offx, *incx,
270 *Y, (size_t) *offy, *incy, CommandQueue, NULL);
271
272}
273
274
275CLBlastStatusCode FC_FUNC_(clblaszdotu_low, CLBLASZDOTU_LOW)(cl_long * N,
276 cl_mem * dotProduct,
277 cl_long * offDP,
278 const cl_mem *X,
279 cl_long * offx,
280 int *incx,
281 const cl_mem * Y,
282 cl_long * offy,
283 int * incy,
284 cl_mem * scratchBuff,
285 cl_command_queue * CommandQueue,
286 int * status){
287
288
289 *status = CLBlastZdotu((size_t) *N, *dotProduct, (size_t) *offDP, *X, (size_t) *offx, *incx,
290 *Y, (size_t) *offy, *incy, CommandQueue, NULL);
291
292}
293
294
295CLBlastStatusCode FC_FUNC_(clblasdnrm2_low, CLBLASDNRM2_LOW)(cl_long * N,
296 cl_mem * NRM2,
297 cl_long * offNRM2,
298 const cl_mem *X,
299 cl_long * offx,
300 int *incx,
301 cl_mem * scratchBuff,
302 cl_command_queue * CommandQueue,
303 int * status){
304
305
306 *status = CLBlastDnrm2((size_t) *N, *NRM2, (size_t) *offNRM2, *X, (size_t) *offx, *incx,
307 CommandQueue, NULL);
308
309}
310
311CLBlastStatusCode FC_FUNC_(clblasdznrm2_low, CLBLASDZNRM2_LOW)(cl_long * N,
312 cl_mem * NRM2,
313 cl_long * offNRM2,
314 const cl_mem *X,
315 cl_long * offx,
316 int *incx,
317 cl_mem * scratchBuff,
318 cl_command_queue * CommandQueue,
319 int * status){
320
321
322 *status = CLBlastDznrm2((size_t) *N, *NRM2, (size_t) *offNRM2, *X, (size_t) *offx, *incx,
323 CommandQueue, NULL);
324
325}
326
327#endif
328