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_(clblasdgemmex_low, CLBLASDGEMMEX_LOW)(int * order,
78 int * transA,
79 int * transB,
80 cl_long * M,
81 cl_long * N,
82 cl_long * K,
83 double * alpha,
84 const cl_mem * A,
85 cl_long * offA,
86 cl_long * lda,
87 const cl_mem * B,
88 cl_long * offB,
89 cl_long * ldb,
90 double * beta,
91 cl_mem * C,
92 cl_long * offC,
93 cl_long * ldc,
94 cl_command_queue * CommandQueue,
95 int * status){
96
97 *status = clblasDgemm((clblasOrder) *order, (clblasTranspose) *transA, (clblasTranspose) *transB,
98 (size_t) *M, (size_t) *N, (size_t) *K, *alpha,
99 *A, (size_t) *offA, (size_t) *lda,
100 *B, (size_t) *offB, (size_t) *ldb, *beta,
101 *C, (size_t) *offC, (size_t) *ldc,
102 1, CommandQueue, 0, NULL, NULL);
103}
104
105void FC_FUNC_(clblaszgemmex_low, CLBLASZGEMMEX_LOW)(int * order,
106 int * transA,
107 int * transB,
108 cl_long * M,
109 cl_long * N,
110 cl_long * K,
111 DoubleComplex * alpha,
112 const cl_mem * A,
113 cl_long * offA,
114 cl_long * lda,
115 const cl_mem * B,
116 cl_long * offB,
117 cl_long * ldb,
118 DoubleComplex * beta,
119 cl_mem * C,
120 cl_long * offC,
121 cl_long * ldc,
122 cl_command_queue * CommandQueue,
123 int * status){
124
125 *status = clblasZgemm((clblasOrder) *order, (clblasTranspose) *transA, (clblasTranspose) *transB,
126 (size_t) *M, (size_t) *N, (size_t) *K, *alpha,
127 *A, (size_t) *offA, (size_t) *lda,
128 *B, (size_t) *offB, (size_t) *ldb, *beta,
129 *C, (size_t) *offC, (size_t) *ldc,
130 1, CommandQueue, 0, NULL, NULL);
131}
132
133void FC_FUNC_(clblasdsyrkex_low, CLBLASDSYRKEX_LOW)(int * order,
134 int * uplo,
135 int * transA,
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 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 = clblasDsyrk((clblasOrder) *order, (clblasUplo) *uplo, (clblasTranspose) *transA,
150 (size_t) *N, (size_t) *K, *alpha,
151 *A, (size_t) *offA, (size_t) *lda, *beta,
152 *C, (size_t) *offC, (size_t) *ldc,
153 1, CommandQueue, 0, NULL, NULL);
154}
155
156void FC_FUNC_(clblaszherkex_low, CLBLASZHERKEX_LOW)(int * order,
157 int * uplo,
158 int * transA,
159 cl_long * N,
160 cl_long * K,
161 double * alpha,
162 const cl_mem * A,
163 cl_long * offA,
164 cl_long * lda,
165 double * beta,
166 cl_mem * C,
167 cl_long * offC,
168 cl_long * ldc,
169 cl_command_queue * CommandQueue,
170 int * status){
171
172 *status = clblasZherk((clblasOrder) *order, (clblasUplo) *uplo, (clblasTranspose) *transA,
173 (size_t) *N, (size_t) *K, *alpha,
174 *A, (size_t) *offA, (size_t) *lda, *beta,
175 *C, (size_t) *offC, (size_t) *ldc,
176 1, CommandQueue, 0, NULL, NULL);
177}
178
179
180clblasStatus FC_FUNC_(clblasddot_low, CLBLASDDOT_LOW)(cl_long * N,
181 cl_mem * dotProduct,
182 cl_long * offDP,
183 const cl_mem *X,
184 cl_long * offx,
185 int *incx,
186 const cl_mem * Y,
187 cl_long * offy,
188 int * incy,
189 cl_mem * scratchBuff,
190 cl_command_queue * CommandQueue,
191 int * status){
192
193
194 *status = clblasDdot((size_t) *N, *dotProduct, (size_t) *offDP, *X, (size_t) *offx, *incx,
195 *Y, (size_t) *offy, *incy, *scratchBuff, 1, CommandQueue, 0, NULL, NULL);
196
197}
198
199
200clblasStatus FC_FUNC_(clblaszdotc_low, CLBLASZDOTC_LOW)(cl_long * N,
201 cl_mem * dotProduct,
202 cl_long * offDP,
203 const cl_mem *X,
204 cl_long * offx,
205 int *incx,
206 const cl_mem * Y,
207 cl_long * offy,
208 int * incy,
209 cl_mem * scratchBuff,
210 cl_command_queue * CommandQueue,
211 int * status){
212
213
214 *status = clblasZdotc((size_t) *N, *dotProduct, (size_t) *offDP, *X, (size_t) *offx, *incx,
215 *Y, (size_t) *offy, *incy, *scratchBuff, 1, CommandQueue, 0, NULL, NULL);
216
217}
218
219
220clblasStatus FC_FUNC_(clblasdnrm2_low, CLBLASDNRM2_LOW)(cl_long * N,
221 cl_mem * NRM2,
222 cl_long * offNRM2,
223 const cl_mem *X,
224 cl_long * offx,
225 int *incx,
226 cl_mem * scratchBuff,
227 cl_command_queue * CommandQueue,
228 int * status){
229
230
231 *status = clblasDnrm2((size_t) *N, *NRM2, (size_t) *offNRM2, *X, (size_t) *offx, *incx,
232 *scratchBuff, 1, CommandQueue, 0, NULL, NULL);
233
234}
235
236clblasStatus FC_FUNC_(clblasdznrm2_low, CLBLASDZNRM2_LOW)(cl_long * N,
237 cl_mem * NRM2,
238 cl_long * offNRM2,
239 const cl_mem *X,
240 cl_long * offx,
241 int *incx,
242 cl_mem * scratchBuff,
243 cl_command_queue * CommandQueue,
244 int * status){
245
246
247 *status = clblasDznrm2((size_t) *N, *NRM2, (size_t) *offNRM2, *X, (size_t) *offx, *incx,
248 *scratchBuff, 1, CommandQueue, 0, NULL, NULL);
249
250}
251
252#endif