blas.c 7.8 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242
  1. /*
  2. * StarPU
  3. * Copyright (C) INRIA 2008-2009 (see AUTHORS file)
  4. *
  5. * This program is free software; you can redistribute it and/or modify
  6. * it under the terms of the GNU Lesser General Public License as published by
  7. * the Free Software Foundation; either version 2.1 of the License, or (at
  8. * your option) any later version.
  9. *
  10. * This program is distributed in the hope that it will be useful, but
  11. * WITHOUT ANY WARRANTY; without even the implied warranty of
  12. * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
  13. *
  14. * See the GNU Lesser General Public License in COPYING.LGPL for more details.
  15. */
  16. #include <ctype.h>
  17. #include <stdio.h>
  18. #include <starpu.h>
  19. #include "blas.h"
  20. /*
  21. This files contains BLAS wrappers for the different BLAS implementations
  22. (eg. REFBLAS, ATLAS, GOTOBLAS ...). We assume a Fortran orientation as most
  23. libraries do not supply C-based ordering.
  24. */
  25. #ifdef ATLAS
  26. inline void SGEMM(char *transa, char *transb, int M, int N, int K,
  27. float alpha, float *A, int lda, float *B, int ldb,
  28. float beta, float *C, int ldc)
  29. {
  30. enum CBLAS_TRANSPOSE ta = (toupper(transa[0]) == 'N')?CblasNoTrans:CblasTrans;
  31. enum CBLAS_TRANSPOSE tb = (toupper(transb[0]) == 'N')?CblasNoTrans:CblasTrans;
  32. cblas_sgemm(CblasColMajor, ta, tb,
  33. M, N, K, alpha, A, lda, B, ldb, beta, C, ldc);
  34. }
  35. inline float SASUM(int N, float *X, int incX)
  36. {
  37. return cblas_sasum(N, X, incX);
  38. }
  39. void SSCAL(int N, float alpha, float *X, int incX)
  40. {
  41. cblas_sscal(N, alpha, X, incX);
  42. }
  43. void STRSM (const char *side, const char *uplo, const char *transa,
  44. const char *diag, const int m, const int n,
  45. const float alpha, const float *A, const int lda,
  46. float *B, const int ldb)
  47. {
  48. enum CBLAS_SIDE side_ = (toupper(side[0]) == 'L')?CblasLeft:CblasRight;
  49. enum CBLAS_UPLO uplo_ = (toupper(uplo[0]) == 'U')?CblasUpper:CblasLower;
  50. enum CBLAS_TRANSPOSE transa_ = (toupper(transa[0]) == 'N')?CblasNoTrans:CblasTrans;
  51. enum CBLAS_DIAG diag_ = (toupper(diag[0]) == 'N')?CblasNonUnit:CblasUnit;
  52. cblas_strsm(CblasColMajor, side_, uplo_, transa_, diag_, m, n, alpha, A, lda, B, ldb);
  53. }
  54. void SSYR (const char *uplo, const int n, const float alpha,
  55. const float *x, const int incx, float *A, const int lda)
  56. {
  57. enum CBLAS_UPLO uplo_ = (toupper(uplo[0]) == 'U')?CblasUpper:CblasLower;
  58. cblas_ssyr(CblasColMajor, uplo_, n, alpha, x, incx, A, lda);
  59. }
  60. void SSYRK (const char *uplo, const char *trans, const int n,
  61. const int k, const float alpha, const float *A,
  62. const int lda, const float beta, float *C,
  63. const int ldc)
  64. {
  65. enum CBLAS_UPLO uplo_ = (toupper(uplo[0]) == 'U')?CblasUpper:CblasLower;
  66. enum CBLAS_TRANSPOSE trans_ = (toupper(trans[0]) == 'N')?CblasNoTrans:CblasTrans;
  67. cblas_ssyrk(CblasColMajor, uplo_, trans_, n, k, alpha, A, lda, beta, C, ldc);
  68. }
  69. void SGER (const int m, const int n, const float alpha,
  70. const float *x, const int incx, const float *y,
  71. const int incy, float *A, const int lda)
  72. {
  73. cblas_sger(CblasRowMajor, m, n, alpha, x, incx, y, incy, A, lda);
  74. }
  75. void STRSV (const char *uplo, const char *trans, const char *diag,
  76. const int n, const float *A, const int lda, float *x,
  77. const int incx)
  78. {
  79. enum CBLAS_UPLO uplo_ = (toupper(uplo[0]) == 'U')?CblasUpper:CblasLower;
  80. enum CBLAS_TRANSPOSE trans_ = (toupper(trans[0]) == 'N')?CblasNoTrans:CblasTrans;
  81. enum CBLAS_DIAG diag_ = (toupper(diag[0]) == 'N')?CblasNonUnit:CblasUnit;
  82. cblas_strsv(CblasColMajor, uplo_, trans_, diag_, n, A, lda, x, incx);
  83. }
  84. void STRMM(const char *side, const char *uplo, const char *transA,
  85. const char *diag, const int m, const int n,
  86. const float alpha, const float *A, const int lda,
  87. float *B, const int ldb)
  88. {
  89. enum CBLAS_SIDE side_ = (toupper(side[0]) == 'L')?CblasLeft:CblasRight;
  90. enum CBLAS_UPLO uplo_ = (toupper(uplo[0]) == 'U')?CblasUpper:CblasLower;
  91. enum CBLAS_TRANSPOSE transA_ = (toupper(transA[0]) == 'N')?CblasNoTrans:CblasTrans;
  92. enum CBLAS_DIAG diag_ = (toupper(diag[0]) == 'N')?CblasNonUnit:CblasUnit;
  93. cblas_strmm(CblasColMajor, side_, uplo_, transA_, diag_, m, n, alpha, A, lda, B, ldb);
  94. }
  95. void STRMV(const char *uplo, const char *transA, const char *diag,
  96. const int n, const float *A, const int lda, float *X,
  97. const int incX)
  98. {
  99. enum CBLAS_UPLO uplo_ = (toupper(uplo[0]) == 'U')?CblasUpper:CblasLower;
  100. enum CBLAS_TRANSPOSE transA_ = (toupper(transA[0]) == 'N')?CblasNoTrans:CblasTrans;
  101. enum CBLAS_DIAG diag_ = (toupper(diag[0]) == 'N')?CblasNonUnit:CblasUnit;
  102. cblas_strmv(CblasColMajor, uplo_, transA_, diag_, n, A, lda, X, incX);
  103. }
  104. void SAXPY(const int n, const float alpha, float *X, const int incX, float *Y, const int incY)
  105. {
  106. cblas_saxpy(n, alpha, X, incX, Y, incY);
  107. }
  108. int ISAMAX (const int n, float *X, const int incX)
  109. {
  110. int retVal;
  111. retVal = cblas_isamax(n, X, incX);
  112. return retVal;
  113. }
  114. float SDOT(const int n, const float *x, const int incx, const float *y, const int incy)
  115. {
  116. return cblas_sdot(n, x, incx, y, incy);
  117. }
  118. #elif defined(GOTO) || defined(SYSTEM_BLAS)
  119. inline void SGEMM(char *transa, char *transb, int M, int N, int K,
  120. float alpha, float *A, int lda, float *B, int ldb,
  121. float beta, float *C, int ldc)
  122. {
  123. sgemm_(transa, transb, &M, &N, &K, &alpha,
  124. A, &lda, B, &ldb,
  125. &beta, C, &ldc);
  126. }
  127. inline float SASUM(int N, float *X, int incX)
  128. {
  129. return sasum_(&N, X, &incX);
  130. }
  131. void SSCAL(int N, float alpha, float *X, int incX)
  132. {
  133. sscal_(&N, &alpha, X, &incX);
  134. }
  135. void STRSM (const char *side, const char *uplo, const char *transa,
  136. const char *diag, const int m, const int n,
  137. const float alpha, const float *A, const int lda,
  138. float *B, const int ldb)
  139. {
  140. strsm_(side, uplo, transa, diag, &m, &n, &alpha, A, &lda, B, &ldb);
  141. }
  142. void SSYR (const char *uplo, const int n, const float alpha,
  143. const float *x, const int incx, float *A, const int lda)
  144. {
  145. ssyr_(uplo, &n, &alpha, x, &incx, A, &lda);
  146. }
  147. void SSYRK (const char *uplo, const char *trans, const int n,
  148. const int k, const float alpha, const float *A,
  149. const int lda, const float beta, float *C,
  150. const int ldc)
  151. {
  152. ssyrk_(uplo, trans, &n, &k, &alpha, A, &lda, &beta, C, &ldc);
  153. }
  154. void SGER (const int m, const int n, const float alpha,
  155. const float *x, const int incx, const float *y,
  156. const int incy, float *A, const int lda)
  157. {
  158. sger_(&m, &n, &alpha, x, &incx, y, &incy, A, &lda);
  159. }
  160. void STRSV (const char *uplo, const char *trans, const char *diag,
  161. const int n, const float *A, const int lda, float *x,
  162. const int incx)
  163. {
  164. strsv_(uplo, trans, diag, &n, A, &lda, x, &incx);
  165. }
  166. void STRMM(const char *side, const char *uplo, const char *transA,
  167. const char *diag, const int m, const int n,
  168. const float alpha, const float *A, const int lda,
  169. float *B, const int ldb)
  170. {
  171. strmm_(side, uplo, transA, diag, &m, &n, &alpha, A, &lda, B, &ldb);
  172. }
  173. void STRMV(const char *uplo, const char *transA, const char *diag,
  174. const int n, const float *A, const int lda, float *X,
  175. const int incX)
  176. {
  177. strmv_(uplo, transA, diag, &n, A, &lda, X, &incX);
  178. }
  179. void SAXPY(const int n, const float alpha, float *X, const int incX, float *Y, const int incY)
  180. {
  181. saxpy_(&n, &alpha, X, &incX, Y, &incY);
  182. }
  183. int ISAMAX (const int n, float *X, const int incX)
  184. {
  185. int retVal;
  186. retVal = isamax_ (&n, X, &incX);
  187. return retVal;
  188. }
  189. float SDOT(const int n, const float *x, const int incx, const float *y, const int incy)
  190. {
  191. float retVal = 0;
  192. /* GOTOBLAS will return a FLOATRET which is a double, not a float */
  193. retVal = (float)sdot_(&n, x, &incx, y, &incy);
  194. return retVal;
  195. }
  196. #else
  197. #error "no BLAS lib available..."
  198. #endif