7void FC_FUNC_(clblasgetversion_low, CLBLASGETVERSION_LOW)(
int * major, 
int * minor, 
int * patch, 
int * status){
 
    8  cl_uint cl_major, cl_minor, cl_patch;
 
   10  *status = clblasGetVersion(&cl_major, &cl_minor, &cl_patch);
 
   16void FC_FUNC_(clblassetup_low, CLBLASSETUP_LOW)(
int * status){
 
   17  *status = clblasSetup();
 
   20void FC_FUNC_(clblasteardown_low, CLBLASTEARDOWN_LOW)(){
 
   24void FC_FUNC_(clblasdtrsmex_low, CLBLASDTRSMEX_LOW)(
int * order,
 
   38                cl_command_queue * CommandQueue, 
 
   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);
 
   51void FC_FUNC_(clblasztrsmex_low, CLBLASZTRSMEX_LOW)(
int * order,
 
   58                DoubleComplex * alpha,
 
   65                cl_command_queue * CommandQueue, 
 
   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);
 
   77void FC_FUNC_(clblasdgemmex_low, CLBLASDGEMMEX_LOW)(
int * order,
 
   94                cl_command_queue * CommandQueue,
 
   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);
 
  105void FC_FUNC_(clblaszgemmex_low, CLBLASZGEMMEX_LOW)(
int * order,
 
  111                DoubleComplex * alpha,
 
  118                DoubleComplex * beta, 
 
  122                cl_command_queue * CommandQueue,
 
  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);
 
  133void FC_FUNC_(clblasdsyrkex_low, CLBLASDSYRKEX_LOW)(
int * order,
 
  146                cl_command_queue * CommandQueue,
 
  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);
 
  156void FC_FUNC_(clblaszherkex_low, CLBLASZHERKEX_LOW)(
int * order,
 
  169                cl_command_queue * CommandQueue,
 
  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);
 
  180clblasStatus FC_FUNC_(clblasddot_low, CLBLASDDOT_LOW)(cl_long * N,
 
  189                  cl_mem * scratchBuff,
 
  190                  cl_command_queue * CommandQueue,
 
  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);
 
  200clblasStatus FC_FUNC_(clblaszdotc_low, CLBLASZDOTC_LOW)(cl_long * N,
 
  209              cl_mem * scratchBuff,
 
  210              cl_command_queue * CommandQueue,
 
  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);
 
  220clblasStatus FC_FUNC_(clblasdnrm2_low, CLBLASDNRM2_LOW)(cl_long * N,
 
  226              cl_mem * scratchBuff,
 
  227              cl_command_queue * CommandQueue,
 
  231  *status = clblasDnrm2((
size_t) *N, *NRM2, (
size_t) *offNRM2, *X, (
size_t) *offx, *incx,
 
  232      *scratchBuff, 1, CommandQueue, 0, NULL, NULL);
 
  236clblasStatus FC_FUNC_(clblasdznrm2_low, CLBLASDZNRM2_LOW)(cl_long * N,
 
  242                cl_mem * scratchBuff,
 
  243                cl_command_queue * CommandQueue,
 
  247  *status = clblasDznrm2((
size_t) *N, *NRM2, (
size_t) *offNRM2, *X, (
size_t) *offx, *incx,
 
  248       *scratchBuff, 1, CommandQueue, 0, NULL, NULL);