summaryrefslogtreecommitdiff
path: root/src/clblast.cc
blob: eddb8022068907a0a241a22878f598a46c3b552f (plain)
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
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
// =================================================================================================
// This file is part of the CLBlast project. The project is licensed under Apache Version 2.0. This
// project loosely follows the Google C++ styleguide and uses a tab-size of two spaces and a max-
// width of 100 characters per line.
//
// Author(s):
//   Cedric Nugteren <www.cedricnugteren.nl>
//
// This file implements all the BLAS API calls. In all cases, it does not much more than creating
// a new object of the appropriate type, and calling the main routine on that object. It forwards
// all status codes to the caller.
//
// =================================================================================================

#include <string>

#include "clblast.h"

// BLAS level-1 includes
#include "internal/routines/level1/xaxpy.h"

// BLAS level-2 includes
#include "internal/routines/level2/xgemv.h"

// BLAS level-3 includes
#include "internal/routines/level3/xgemm.h"
#include "internal/routines/level3/xsymm.h"
#include "internal/routines/level3/xhemm.h"
#include "internal/routines/level3/xsyrk.h"
#include "internal/routines/level3/xherk.h"
#include "internal/routines/level3/xsyr2k.h"
#include "internal/routines/level3/xher2k.h"
#include "internal/routines/level3/xtrmm.h"

namespace clblast {
// =================================================================================================
// BLAS level-1 (vector-vector) routines

// AXPY
template <typename T>
StatusCode Axpy(const size_t n, const T alpha,
                const cl_mem x_buffer, const size_t x_offset, const size_t x_inc,
                cl_mem y_buffer, const size_t y_offset, const size_t y_inc,
                cl_command_queue* queue, cl_event* event) {
  auto queue_cpp = Queue(*queue);
  auto event_cpp = Event(*event);
  auto routine = Xaxpy<T>(queue_cpp, event_cpp);

  // Compiles the routine's device kernels
  auto status = routine.SetUp();
  if (status != StatusCode::kSuccess) { return status; }

  // Runs the routine
  return routine.DoAxpy(n, alpha,
                        Buffer<T>(x_buffer), x_offset, x_inc,
                        Buffer<T>(y_buffer), y_offset, y_inc);
}
template StatusCode Axpy<float>(const size_t, const float,
                                const cl_mem, const size_t, const size_t,
                                cl_mem, const size_t, const size_t,
                                cl_command_queue*, cl_event*);
template StatusCode Axpy<double>(const size_t, const double,
                                 const cl_mem, const size_t, const size_t,
                                 cl_mem, const size_t, const size_t,
                                 cl_command_queue*, cl_event*);
template StatusCode Axpy<float2>(const size_t, const float2,
                                 const cl_mem, const size_t, const size_t,
                                 cl_mem, const size_t, const size_t,
                                 cl_command_queue*, cl_event*);
template StatusCode Axpy<double2>(const size_t, const double2,
                                  const cl_mem, const size_t, const size_t,
                                  cl_mem, const size_t, const size_t,
                                  cl_command_queue*, cl_event*);

// =================================================================================================
// BLAS level-2 (matrix-vector) routines

// GEMV
template <typename T>
StatusCode Gemv(const Layout layout, const Transpose a_transpose,
                const size_t m, const size_t n, const T alpha,
                const cl_mem a_buffer, const size_t a_offset, const size_t a_ld,
                const cl_mem x_buffer, const size_t x_offset, const size_t x_inc, const T beta,
                cl_mem y_buffer, const size_t y_offset, const size_t y_inc,
                cl_command_queue* queue, cl_event* event) {
  
  auto queue_cpp = Queue(*queue);
  auto event_cpp = Event(*event);
  auto routine = Xgemv<T>(queue_cpp, event_cpp);

  // Compiles the routine's device kernels
  auto status = routine.SetUp();
  if (status != StatusCode::kSuccess) { return status; }

  // Runs the routine
  return routine.DoGemv(layout, a_transpose, m, n, alpha,
                        Buffer<T>(a_buffer), a_offset, a_ld,
                        Buffer<T>(x_buffer), x_offset, x_inc, beta,
                        Buffer<T>(y_buffer), y_offset, y_inc);
}
template StatusCode Gemv<float>(const Layout, const Transpose,
                                const size_t, const size_t, const float,
                                const cl_mem, const size_t, const size_t,
                                const cl_mem, const size_t, const size_t, const float,
                                cl_mem, const size_t, const size_t,
                                cl_command_queue*, cl_event*);
template StatusCode Gemv<double>(const Layout, const Transpose,
                                 const size_t, const size_t, const double,
                                 const cl_mem, const size_t, const size_t,
                                 const cl_mem, const size_t, const size_t, const double,
                                 cl_mem, const size_t, const size_t,
                                 cl_command_queue*, cl_event*);
template StatusCode Gemv<float2>(const Layout, const Transpose,
                                 const size_t, const size_t, const float2,
                                 const cl_mem, const size_t, const size_t,
                                 const cl_mem, const size_t, const size_t, const float2,
                                 cl_mem, const size_t, const size_t,
                                 cl_command_queue*, cl_event*);
template StatusCode Gemv<double2>(const Layout, const Transpose,
                                  const size_t, const size_t, const double2,
                                  const cl_mem, const size_t, const size_t,
                                  const cl_mem, const size_t, const size_t, const double2,
                                  cl_mem, const size_t, const size_t,
                                  cl_command_queue*, cl_event*);

// =================================================================================================
// BLAS level-3 (matrix-matrix) routines

// GEMM
template <typename T>
StatusCode Gemm(const Layout layout, const Transpose a_transpose, const Transpose b_transpose,
                const size_t m, const size_t n, const size_t k, const T alpha,
                const cl_mem a_buffer, const size_t a_offset, const size_t a_ld,
                const cl_mem b_buffer, const size_t b_offset, const size_t b_ld, const T beta,
                cl_mem c_buffer, const size_t c_offset, const size_t c_ld,
                cl_command_queue* queue, cl_event* event) {
  auto queue_cpp = Queue(*queue);
  auto event_cpp = Event(*event);
  auto routine = Xgemm<T>(queue_cpp, event_cpp);

  // Compiles the routine's device kernels
  auto status = routine.SetUp();
  if (status != StatusCode::kSuccess) { return status; }

  // Runs the routine
  return routine.DoGemm(layout, a_transpose, b_transpose, m, n, k, alpha,
                        Buffer<T>(a_buffer), a_offset, a_ld,
                        Buffer<T>(b_buffer), b_offset, b_ld, beta,
                        Buffer<T>(c_buffer), c_offset, c_ld);
}
template StatusCode Gemm<float>(const Layout, const Transpose, const Transpose,
                                const size_t, const size_t, const size_t, const float,
                                const cl_mem, const size_t, const size_t,
                                const cl_mem, const size_t, const size_t, const float,
                                cl_mem, const size_t, const size_t,
                                cl_command_queue*, cl_event*);
template StatusCode Gemm<double>(const Layout, const Transpose, const Transpose,
                                 const size_t, const size_t, const size_t, const double,
                                 const cl_mem, const size_t, const size_t,
                                 const cl_mem, const size_t, const size_t, const double,
                                 cl_mem, const size_t, const size_t,
                                 cl_command_queue*, cl_event*);
template StatusCode Gemm<float2>(const Layout, const Transpose, const Transpose,
                                 const size_t, const size_t, const size_t, const float2,
                                 const cl_mem, const size_t, const size_t,
                                 const cl_mem, const size_t, const size_t, const float2,
                                 cl_mem, const size_t, const size_t,
                                 cl_command_queue*, cl_event*);
template StatusCode Gemm<double2>(const Layout, const Transpose, const Transpose,
                                  const size_t, const size_t, const size_t, const double2,
                                  const cl_mem, const size_t, const size_t,
                                  const cl_mem, const size_t, const size_t, const double2,
                                  cl_mem, const size_t, const size_t,
                                  cl_command_queue*, cl_event*);

// =================================================================================================

// SYMM
template <typename T>
StatusCode Symm(const Layout layout, const Side side, const Triangle triangle,
                const size_t m, const size_t n, const T alpha,
                const cl_mem a_buffer, const size_t a_offset, const size_t a_ld,
                const cl_mem b_buffer, const size_t b_offset, const size_t b_ld, const T beta,
                cl_mem c_buffer, const size_t c_offset, const size_t c_ld,
                cl_command_queue* queue, cl_event* event) {
  auto queue_cpp = Queue(*queue);
  auto event_cpp = Event(*event);
  auto routine = Xsymm<T>(queue_cpp, event_cpp);

  // Compiles the routine's device kernels
  auto status = routine.SetUp();
  if (status != StatusCode::kSuccess) { return status; }

  // Runs the routine
  return routine.DoSymm(layout, side, triangle, m, n, alpha,
                        Buffer<T>(a_buffer), a_offset, a_ld,
                        Buffer<T>(b_buffer), b_offset, b_ld, beta,
                        Buffer<T>(c_buffer), c_offset, c_ld);
}
template StatusCode Symm<float>(const Layout, const Side, const Triangle,
                                const size_t, const size_t, const float,
                                const cl_mem, const size_t, const size_t,
                                const cl_mem, const size_t, const size_t, const float,
                                cl_mem, const size_t, const size_t,
                                cl_command_queue*, cl_event*);
template StatusCode Symm<double>(const Layout, const Side, const Triangle,
                                 const size_t, const size_t, const double,
                                 const cl_mem, const size_t, const size_t,
                                 const cl_mem, const size_t, const size_t, const double,
                                 cl_mem, const size_t, const size_t,
                                 cl_command_queue*, cl_event*);
template StatusCode Symm<float2>(const Layout, const Side, const Triangle,
                                 const size_t, const size_t, const float2,
                                 const cl_mem, const size_t, const size_t,
                                 const cl_mem, const size_t, const size_t, const float2,
                                 cl_mem, const size_t, const size_t,
                                 cl_command_queue*, cl_event*);
template StatusCode Symm<double2>(const Layout, const Side, const Triangle,
                                  const size_t, const size_t, const double2,
                                  const cl_mem, const size_t, const size_t,
                                  const cl_mem, const size_t, const size_t, const double2,
                                  cl_mem, const size_t, const size_t,
                                  cl_command_queue*, cl_event*);

// =================================================================================================

// HEMM
template <typename T>
StatusCode Hemm(const Layout layout, const Side side, const Triangle triangle,
                const size_t m, const size_t n, const T alpha,
                const cl_mem a_buffer, const size_t a_offset, const size_t a_ld,
                const cl_mem b_buffer, const size_t b_offset, const size_t b_ld, const T beta,
                cl_mem c_buffer, const size_t c_offset, const size_t c_ld,
                cl_command_queue* queue, cl_event* event) {
  auto queue_cpp = Queue(*queue);
  auto event_cpp = Event(*event);
  auto routine = Xhemm<T>(queue_cpp, event_cpp);

  // Compiles the routine's device kernels
  auto status = routine.SetUp();
  if (status != StatusCode::kSuccess) { return status; }

  // Runs the routine
  return routine.DoHemm(layout, side, triangle, m, n, alpha,
                        Buffer<T>(a_buffer), a_offset, a_ld,
                        Buffer<T>(b_buffer), b_offset, b_ld, beta,
                        Buffer<T>(c_buffer), c_offset, c_ld);
}
template StatusCode Hemm<float2>(const Layout, const Side, const Triangle,
                                 const size_t, const size_t, const float2,
                                 const cl_mem, const size_t, const size_t,
                                 const cl_mem, const size_t, const size_t, const float2,
                                 cl_mem, const size_t, const size_t,
                                 cl_command_queue*, cl_event*);
template StatusCode Hemm<double2>(const Layout, const Side, const Triangle,
                                  const size_t, const size_t, const double2,
                                  const cl_mem, const size_t, const size_t,
                                  const cl_mem, const size_t, const size_t, const double2,
                                  cl_mem, const size_t, const size_t,
                                  cl_command_queue*, cl_event*);

// =================================================================================================

// SYRK
template <typename T>
StatusCode Syrk(const Layout layout, const Triangle triangle, const Transpose a_transpose,
                const size_t n, const size_t k, const T alpha,
                const cl_mem a_buffer, const size_t a_offset, const size_t a_ld, const T beta,
                cl_mem c_buffer, const size_t c_offset, const size_t c_ld,
                cl_command_queue* queue, cl_event* event) {
  auto queue_cpp = Queue(*queue);
  auto event_cpp = Event(*event);
  auto routine = Xsyrk<T>(queue_cpp, event_cpp);

  // Compiles the routine's device kernels
  auto status = routine.SetUp();
  if (status != StatusCode::kSuccess) { return status; }

  // Runs the routine
  return routine.DoSyrk(layout, triangle, a_transpose, n, k, alpha,
                        Buffer<T>(a_buffer), a_offset, a_ld, beta,
                        Buffer<T>(c_buffer), c_offset, c_ld);
}
template StatusCode Syrk<float>(const Layout, const Triangle, const Transpose,
                                const size_t, const size_t, const float,
                                const cl_mem, const size_t, const size_t, const float,
                                cl_mem, const size_t, const size_t,
                                cl_command_queue*, cl_event*);
template StatusCode Syrk<double>(const Layout, const Triangle, const Transpose,
                                 const size_t, const size_t, const double,
                                 const cl_mem, const size_t, const size_t, const double,
                                 cl_mem, const size_t, const size_t,
                                 cl_command_queue*, cl_event*);
template StatusCode Syrk<float2>(const Layout, const Triangle, const Transpose,
                                 const size_t, const size_t, const float2,
                                 const cl_mem, const size_t, const size_t, const float2,
                                 cl_mem, const size_t, const size_t,
                                 cl_command_queue*, cl_event*);
template StatusCode Syrk<double2>(const Layout, const Triangle, const Transpose,
                                  const size_t, const size_t, const double2,
                                  const cl_mem, const size_t, const size_t, const double2,
                                  cl_mem, const size_t, const size_t,
                                  cl_command_queue*, cl_event*);

// =================================================================================================

// HERK
template <typename T>
StatusCode Herk(const Layout layout, const Triangle triangle, const Transpose a_transpose,
                const size_t n, const size_t k, const T alpha,
                const cl_mem a_buffer, const size_t a_offset, const size_t a_ld, const T beta,
                cl_mem c_buffer, const size_t c_offset, const size_t c_ld,
                cl_command_queue* queue, cl_event* event) {
  auto queue_cpp = Queue(*queue);
  auto event_cpp = Event(*event);
  auto routine = Xherk<std::complex<T>,T>(queue_cpp, event_cpp);

  // Compiles the routine's device kernels
  auto status = routine.SetUp();
  if (status != StatusCode::kSuccess) { return status; }

  // Runs the routine
  return routine.DoHerk(layout, triangle, a_transpose, n, k, alpha,
                        Buffer<std::complex<T>>(a_buffer), a_offset, a_ld, beta,
                        Buffer<std::complex<T>>(c_buffer), c_offset, c_ld);
}
template StatusCode Herk<float>(const Layout, const Triangle, const Transpose,
                                const size_t, const size_t, const float,
                                const cl_mem, const size_t, const size_t, const float,
                                cl_mem, const size_t, const size_t,
                                cl_command_queue*, cl_event*);
template StatusCode Herk<double>(const Layout, const Triangle, const Transpose,
                                 const size_t, const size_t, const double,
                                 const cl_mem, const size_t, const size_t, const double,
                                 cl_mem, const size_t, const size_t,
                                 cl_command_queue*, cl_event*);

// =================================================================================================

// SYR2K
template <typename T>
StatusCode Syr2k(const Layout layout, const Triangle triangle, const Transpose ab_transpose,
                 const size_t n, const size_t k, const T alpha,
                 const cl_mem a_buffer, const size_t a_offset, const size_t a_ld,
                 const cl_mem b_buffer, const size_t b_offset, const size_t b_ld, const T beta,
                 cl_mem c_buffer, const size_t c_offset, const size_t c_ld,
                 cl_command_queue* queue, cl_event* event) {
  auto queue_cpp = Queue(*queue);
  auto event_cpp = Event(*event);
  auto routine = Xsyr2k<T>(queue_cpp, event_cpp);

  // Compiles the routine's device kernels
  auto status = routine.SetUp();
  if (status != StatusCode::kSuccess) { return status; }

  // Runs the routine
  return routine.DoSyr2k(layout, triangle, ab_transpose, n, k, alpha,
                         Buffer<T>(a_buffer), a_offset, a_ld,
                         Buffer<T>(b_buffer), b_offset, b_ld, beta,
                         Buffer<T>(c_buffer), c_offset, c_ld);
}
template StatusCode Syr2k<float>(const Layout, const Triangle, const Transpose,
                                 const size_t, const size_t, const float,
                                 const cl_mem, const size_t, const size_t,
                                 const cl_mem, const size_t, const size_t, const float,
                                 cl_mem, const size_t, const size_t,
                                 cl_command_queue*, cl_event*);
template StatusCode Syr2k<double>(const Layout, const Triangle, const Transpose,
                                  const size_t, const size_t, const double,
                                  const cl_mem, const size_t, const size_t,
                                  const cl_mem, const size_t, const size_t, const double,
                                  cl_mem, const size_t, const size_t,
                                  cl_command_queue*, cl_event*);
template StatusCode Syr2k<float2>(const Layout, const Triangle, const Transpose,
                                  const size_t, const size_t, const float2,
                                  const cl_mem, const size_t, const size_t,
                                  const cl_mem, const size_t, const size_t, const float2,
                                  cl_mem, const size_t, const size_t,
                                  cl_command_queue*, cl_event*);
template StatusCode Syr2k<double2>(const Layout, const Triangle, const Transpose,
                                   const size_t, const size_t, const double2,
                                   const cl_mem, const size_t, const size_t,
                                   const cl_mem, const size_t, const size_t, const double2,
                                   cl_mem, const size_t, const size_t,
                                   cl_command_queue*, cl_event*);

// =================================================================================================

// SYR2K
template <typename T, typename U>
StatusCode Her2k(const Layout layout, const Triangle triangle, const Transpose ab_transpose,
                 const size_t n, const size_t k, const T alpha,
                 const cl_mem a_buffer, const size_t a_offset, const size_t a_ld,
                 const cl_mem b_buffer, const size_t b_offset, const size_t b_ld, const U beta,
                 cl_mem c_buffer, const size_t c_offset, const size_t c_ld,
                 cl_command_queue* queue, cl_event* event) {
  auto queue_cpp = Queue(*queue);
  auto event_cpp = Event(*event);
  auto routine = Xher2k<T,U>(queue_cpp, event_cpp);

  // Compiles the routine's device kernels
  auto status = routine.SetUp();
  if (status != StatusCode::kSuccess) { return status; }

  // Runs the routine
  return routine.DoHer2k(layout, triangle, ab_transpose, n, k, alpha,
                         Buffer<T>(a_buffer), a_offset, a_ld,
                         Buffer<T>(b_buffer), b_offset, b_ld, beta,
                         Buffer<T>(c_buffer), c_offset, c_ld);
}
template StatusCode Her2k<float2,float>(const Layout, const Triangle, const Transpose,
                                        const size_t, const size_t, const float2,
                                        const cl_mem, const size_t, const size_t,
                                        const cl_mem, const size_t, const size_t, const float,
                                        cl_mem, const size_t, const size_t,
                                        cl_command_queue*, cl_event*);
template StatusCode Her2k<double2,double>(const Layout, const Triangle, const Transpose,
                                          const size_t, const size_t, const double2,
                                          const cl_mem, const size_t, const size_t,
                                          const cl_mem, const size_t, const size_t, const double,
                                          cl_mem, const size_t, const size_t,
                                          cl_command_queue*, cl_event*);

// =================================================================================================

// TRMM
template <typename T>
StatusCode Trmm(const Layout layout, const Side side, const Triangle triangle,
                const Transpose a_transpose, const Diagonal diagonal,
                const size_t m, const size_t n,
                const T alpha,
                const cl_mem a_buffer, const size_t a_offset, const size_t a_ld,
                cl_mem b_buffer, const size_t b_offset, const size_t b_ld,
                cl_command_queue* queue, cl_event* event) {
  auto queue_cpp = Queue(*queue);
  auto event_cpp = Event(*event);
  auto routine = Xtrmm<T>(queue_cpp, event_cpp);

  // Compiles the routine's device kernels
  auto status = routine.SetUp();
  if (status != StatusCode::kSuccess) { return status; }

  // Runs the routine
  return routine.DoTrmm(layout, side, triangle, a_transpose, diagonal, m, n, alpha,
                        Buffer<T>(a_buffer), a_offset, a_ld,
                        Buffer<T>(b_buffer), b_offset, b_ld);
}
template StatusCode Trmm<float>(const Layout, const Side, const Triangle,
                                const Transpose, const Diagonal,
                                const size_t, const size_t, const float,
                                const cl_mem, const size_t, const size_t,
                                cl_mem, const size_t, const size_t,
                                cl_command_queue*, cl_event*);
template StatusCode Trmm<double>(const Layout, const Side, const Triangle,
                                 const Transpose, const Diagonal,
                                 const size_t, const size_t, const double,
                                 const cl_mem, const size_t, const size_t,
                                 cl_mem, const size_t, const size_t,
                                 cl_command_queue*, cl_event*);
template StatusCode Trmm<float2>(const Layout, const Side, const Triangle,
                                 const Transpose, const Diagonal,
                                 const size_t, const size_t, const float2,
                                 const cl_mem, const size_t, const size_t,
                                 cl_mem, const size_t, const size_t,
                                 cl_command_queue*, cl_event*);
template StatusCode Trmm<double2>(const Layout, const Side, const Triangle,
                                  const Transpose, const Diagonal,
                                  const size_t, const size_t, const double2,
                                  const cl_mem, const size_t, const size_t,
                                  cl_mem, const size_t, const size_t,
                                  cl_command_queue*, cl_event*);

// =================================================================================================

// TRSM
/*
template <typename T>
StatusCode Trsm(const Layout layout, const Side side, const Triangle triangle,
                const Transpose a_transpose, const Diagonal diagonal,
                const size_t m, const size_t n,
                const T alpha,
                const cl_mem a_buffer, const size_t a_offset, const size_t a_ld,
                cl_mem b_buffer, const size_t b_offset, const size_t b_ld,
                cl_command_queue* queue, cl_event* event) {
  auto queue_cpp = Queue(*queue);
  auto event_cpp = Event(*event);
  auto routine = Xtrsm<T>(queue_cpp, event_cpp);

  // Compiles the routine's device kernels
  auto status = routine.SetUp();
  if (status != StatusCode::kSuccess) { return status; }

  // Runs the routine
  return routine.DoTrsm(layout, side, triangle, a_transpose, diagonal, m, n, alpha,
                        Buffer<T>(a_buffer), a_offset, a_ld,
                        Buffer<T>(b_buffer), b_offset, b_ld);
}
template StatusCode Trsm<float>(const Layout, const Side, const Triangle,
                                const Transpose, const Diagonal,
                                const size_t, const size_t, const float,
                                const cl_mem, const size_t, const size_t,
                                cl_mem, const size_t, const size_t,
                                cl_command_queue*, cl_event*);
template StatusCode Trsm<double>(const Layout, const Side, const Triangle,
                                 const Transpose, const Diagonal,
                                 const size_t, const size_t, const double,
                                 const cl_mem, const size_t, const size_t,
                                 cl_mem, const size_t, const size_t,
                                 cl_command_queue*, cl_event*);
template StatusCode Trsm<float2>(const Layout, const Side, const Triangle,
                                 const Transpose, const Diagonal,
                                 const size_t, const size_t, const float2,
                                 const cl_mem, const size_t, const size_t,
                                 cl_mem, const size_t, const size_t,
                                 cl_command_queue*, cl_event*);
template StatusCode Trsm<double2>(const Layout, const Side, const Triangle,
                                  const Transpose, const Diagonal,
                                  const size_t, const size_t, const double2,
                                  const cl_mem, const size_t, const size_t,
                                  cl_mem, const size_t, const size_t,
                                  cl_command_queue*, cl_event*);
*/
// =================================================================================================
} // namespace clblast