// ================================================================================================= // 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 // // 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 #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 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 = CommandQueue(*queue); auto event_cpp = Event(*event); auto routine = Xaxpy(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(x_buffer), x_offset, x_inc, Buffer(y_buffer), y_offset, y_inc); } template StatusCode Axpy(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(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(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(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 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 = CommandQueue(*queue); auto event_cpp = Event(*event); auto routine = Xgemv(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(a_buffer), a_offset, a_ld, Buffer(x_buffer), x_offset, x_inc, beta, Buffer(y_buffer), y_offset, y_inc); } template StatusCode Gemv(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(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(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(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 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 = CommandQueue(*queue); auto event_cpp = Event(*event); auto routine = Xgemm(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(a_buffer), a_offset, a_ld, Buffer(b_buffer), b_offset, b_ld, beta, Buffer(c_buffer), c_offset, c_ld); } template StatusCode Gemm(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(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(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(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 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 = CommandQueue(*queue); auto event_cpp = Event(*event); auto routine = Xsymm(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(a_buffer), a_offset, a_ld, Buffer(b_buffer), b_offset, b_ld, beta, Buffer(c_buffer), c_offset, c_ld); } template StatusCode Symm(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(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(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(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 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 = CommandQueue(*queue); auto event_cpp = Event(*event); auto routine = Xhemm(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(a_buffer), a_offset, a_ld, Buffer(b_buffer), b_offset, b_ld, beta, Buffer(c_buffer), c_offset, c_ld); } template StatusCode Hemm(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(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 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 = CommandQueue(*queue); auto event_cpp = Event(*event); auto routine = Xsyrk(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(a_buffer), a_offset, a_ld, beta, Buffer(c_buffer), c_offset, c_ld); } template StatusCode Syrk(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(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(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(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 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 = CommandQueue(*queue); auto event_cpp = Event(*event); auto routine = Xherk,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(a_buffer), a_offset, a_ld, beta, Buffer(c_buffer), c_offset, c_ld); } template StatusCode Herk(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(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 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 = CommandQueue(*queue); auto event_cpp = Event(*event); auto routine = Xsyr2k(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(a_buffer), a_offset, a_ld, Buffer(b_buffer), b_offset, b_ld, beta, Buffer(c_buffer), c_offset, c_ld); } template StatusCode Syr2k(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(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(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(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 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 = CommandQueue(*queue); auto event_cpp = Event(*event); auto routine = Xher2k(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(a_buffer), a_offset, a_ld, Buffer(b_buffer), b_offset, b_ld, beta, Buffer(c_buffer), c_offset, c_ld); } template StatusCode Her2k(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(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 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 = CommandQueue(*queue); auto event_cpp = Event(*event); auto routine = Xtrmm(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(a_buffer), a_offset, a_ld, Buffer(b_buffer), b_offset, b_ld); } template StatusCode Trmm(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(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(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(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 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 = CommandQueue(*queue); auto event_cpp = Event(*event); auto routine = Xtrsm(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(a_buffer), a_offset, a_ld, Buffer(b_buffer), b_offset, b_ld); } template StatusCode Trsm(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(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(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(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