ViennaCL - The Vienna Computing Library  1.7.0
Free open-source GPU-accelerated linear algebra and solver library.
blas3_opencl.cpp
Go to the documentation of this file.
1 /* =========================================================================
2  Copyright (c) 2010-2014, Institute for Microelectronics,
3  Institute for Analysis and Scientific Computing,
4  TU Wien.
5  Portions of this software are copyright by UChicago Argonne, LLC.
6 
7  -----------------
8  ViennaCL - The Vienna Computing Library
9  -----------------
10 
11  Project Head: Karl Rupp rupp@iue.tuwien.ac.at
12 
13  (A list of authors and contributors can be found in the PDF manual)
14 
15  License: MIT (X11), see file LICENSE in the base directory
16 ============================================================================= */
17 
18 // include necessary system headers
19 #include <iostream>
20 
21 #include "viennacl.hpp"
22 #include "viennacl_private.hpp"
23 
24 #include "blas3.hpp"
25 
26 //include basic scalar and vector types of ViennaCL
27 #include "viennacl/scalar.hpp"
28 #include "viennacl/vector.hpp"
29 #include "viennacl/matrix.hpp"
31 #include "viennacl/linalg/prod.hpp"
32 
33 
34 #ifdef VIENNACL_WITH_OPENCL
35 
36 
37 
38 //
39 // xGEMV
40 //
41 
42 namespace detail
43 {
44  template <typename NumericT>
45  ViennaCLStatus ViennaCLOpenCLgemm_impl(ViennaCLBackend backend,
46  ViennaCLOrder orderA, ViennaCLTranspose transA,
47  ViennaCLOrder orderB, ViennaCLTranspose transB,
48  ViennaCLOrder orderC,
50  NumericT alpha,
51  cl_mem A, ViennaCLInt offA_row, ViennaCLInt offA_col, ViennaCLInt incA_row, ViennaCLInt incA_col, ViennaCLInt lda,
52  cl_mem B, ViennaCLInt offB_row, ViennaCLInt offB_col, ViennaCLInt incB_row, ViennaCLInt incB_col, ViennaCLInt ldb,
53  NumericT beta,
54  cl_mem C, ViennaCLInt offC_row, ViennaCLInt offC_col, ViennaCLInt incC_row, ViennaCLInt incC_col, ViennaCLInt ldc)
55  {
56  typedef typename viennacl::matrix_base<NumericT>::size_type size_type;
57  typedef typename viennacl::matrix_base<NumericT>::size_type difference_type;
58 
59  size_type A_size1 = static_cast<size_type>((transA == ViennaCLTrans) ? k : m);
60  size_type A_size2 = static_cast<size_type>((transA == ViennaCLTrans) ? m : k);
61 
62  size_type B_size1 = static_cast<size_type>((transB == ViennaCLTrans) ? n : k);
63  size_type B_size2 = static_cast<size_type>((transB == ViennaCLTrans) ? k : n);
64 
65  bool A_row_major = (orderA == ViennaCLRowMajor);
66  bool B_row_major = (orderB == ViennaCLRowMajor);
67  bool C_row_major = (orderC == ViennaCLRowMajor);
68 
70  A_size1, size_type(offA_row), difference_type(incA_row), size_type(A_row_major ? m : lda),
71  A_size2, size_type(offA_col), difference_type(incA_col), size_type(A_row_major ? lda : k), A_row_major);
72 
74  B_size1, size_type(offB_row), difference_type(incB_row), size_type(B_row_major ? k : ldb),
75  B_size2, size_type(offB_col), difference_type(incB_col), size_type(B_row_major ? ldb : n), B_row_major);
76 
78  size_type(m), size_type(offC_row), difference_type(incC_row), size_type(C_row_major ? m : ldc),
79  size_type(n), size_type(offC_col), difference_type(incC_col), size_type(C_row_major ? ldc : n), C_row_major);
80 
81  detail::gemm_dispatch(alpha, matA, transA, matB, transB, beta, matC);
82 
83  return ViennaCLSuccess;
84  }
85 
86 }
87 
88 
90  ViennaCLOrder orderA, ViennaCLTranspose transA,
91  ViennaCLOrder orderB, ViennaCLTranspose transB,
92  ViennaCLOrder orderC,
94  float alpha,
95  cl_mem A, ViennaCLInt offA_row, ViennaCLInt offA_col, ViennaCLInt incA_row, ViennaCLInt incA_col, ViennaCLInt lda,
96  cl_mem B, ViennaCLInt offB_row, ViennaCLInt offB_col, ViennaCLInt incB_row, ViennaCLInt incB_col, ViennaCLInt ldb,
97  float beta,
98  cl_mem C, ViennaCLInt offC_row, ViennaCLInt offC_col, ViennaCLInt incC_row, ViennaCLInt incC_col, ViennaCLInt ldc)
99 {
100  return detail::ViennaCLOpenCLgemm_impl<float>(backend,
101  orderA, transA,
102  orderB, transB,
103  orderC,
104  m, n, k,
105  alpha,
106  A, offA_row, offA_col, incA_row, incA_col, lda,
107  B, offB_row, offB_col, incB_row, incB_col, ldb,
108  beta,
109  C, offC_row, offC_col, incC_row, incC_col, ldc);
110 }
111 
112 VIENNACL_EXPORTED_FUNCTION ViennaCLStatus ViennaCLOpenCLDgemm(ViennaCLBackend backend,
113  ViennaCLOrder orderA, ViennaCLTranspose transA,
114  ViennaCLOrder orderB, ViennaCLTranspose transB,
115  ViennaCLOrder orderC,
117  double alpha,
118  cl_mem A, ViennaCLInt offA_row, ViennaCLInt offA_col, ViennaCLInt incA_row, ViennaCLInt incA_col, ViennaCLInt lda,
119  cl_mem B, ViennaCLInt offB_row, ViennaCLInt offB_col, ViennaCLInt incB_row, ViennaCLInt incB_col, ViennaCLInt ldb,
120  double beta,
121  cl_mem C, ViennaCLInt offC_row, ViennaCLInt offC_col, ViennaCLInt incC_row, ViennaCLInt incC_col, ViennaCLInt ldc)
122 {
123  return detail::ViennaCLOpenCLgemm_impl<double>(backend,
124  orderA, transA,
125  orderB, transB,
126  orderC,
127  m, n, k,
128  alpha,
129  A, offA_row, offA_col, incA_row, incA_col, lda,
130  B, offB_row, offB_col, incB_row, incB_col, ldb,
131  beta,
132  C, offC_row, offC_col, incC_row, incC_col, ldc);
133 }
134 
135 
136 #endif
ViennaCLTranspose
Definition: viennacl.hpp:68
Generic backend for CUDA, OpenCL, host-based stuff.
void gemm_dispatch(ScalarType alpha, MatrixTypeA const &A, ViennaCLTranspose transA, MatrixTypeB const &B, ViennaCLTranspose transB, ScalarType beta, MatrixTypeC &C)
Definition: blas3.hpp:39
Generic interface for matrix-vector and matrix-matrix products. See viennacl/linalg/vector_operations...
Implementation of the dense matrix class.
ViennaCLStatus
Definition: viennacl.hpp:97
#define VIENNACL_EXPORTED_FUNCTION
Definition: viennacl.hpp:40
float NumericT
Definition: bisect.cpp:40
Definition: blas3.hpp:36
ViennaCLOpenCLBackend_impl opencl_backend
ViennaCLOrder
Definition: viennacl.hpp:61
int ViennaCLInt
Definition: viennacl.hpp:48
Implementations of dense direct solvers are found here.
The vector type with operator-overloads and proxy classes is defined here. Linear algebra operations ...
viennacl::ocl::context & get_context(long i)
Convenience function for returning the current context.
Definition: backend.hpp:225
Implementation of the ViennaCL scalar class.