BLAS.ecl 11 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285
  1. /*##############################################################################
  2. HPCC SYSTEMS software Copyright (C) 2016 HPCC Systems®.
  3. Licensed under the Apache License, Version 2.0 (the "License");
  4. you may not use this file except in compliance with the License.
  5. You may obtain a copy of the License at
  6. http://www.apache.org/licenses/LICENSE-2.0
  7. Unless required by applicable law or agreed to in writing, software
  8. distributed under the License is distributed on an "AS IS" BASIS,
  9. WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
  10. See the License for the specific language governing permissions and
  11. limitations under the License.
  12. ############################################################################## */
  13. IMPORT lib_eclblas AS LIB_ECLBLAS;
  14. EXPORT BLAS := MODULE
  15. // Types for the Block Basic Linear Algebra Sub-programs support
  16. EXPORT Types := MODULE
  17. EXPORT value_t := LIB_ECLBLAS.blas_value_t; //REAL8;
  18. EXPORT dimension_t := LIB_ECLBLAS.blas_dimension_t; //UNSIGNED4;
  19. EXPORT matrix_t := LIB_ECLBLAS.blas_matrix_t; //SET OF REAL8
  20. EXPORT Triangle := LIB_ECLBLAS.blas_Triangle; //ENUM(UNSIGNED1, Upper=1, Lower=2)
  21. EXPORT Diagonal := LIB_ECLBLAS.blas_Diagonal; //ENUM(UNSIGNED1, UnitTri=1, NotUnitTri=2)
  22. EXPORT Side := LIB_ECLBLAS.blas_Side; //ENUM(UNSIGNED1, Ax=1, xA=2)
  23. END;
  24. /**
  25. * Function prototype for Apply2Cell.
  26. * @param v the value
  27. * @param r the row ordinal
  28. * @param c the column ordinal
  29. * @return the updated value
  30. */
  31. EXPORT Types.value_t ICellFunc(Types.value_t v,
  32. Types.dimension_t r,
  33. Types.dimension_t c) := v;
  34. /**
  35. * Iterate matrix and apply function to each cell
  36. *@param m number of rows
  37. *@param n number of columns
  38. *@param x matrix
  39. *@param f function to apply
  40. *@return updated matrix
  41. */
  42. EXPORT Types.matrix_t Apply2Cells(Types.dimension_t m,
  43. Types.dimension_t n,
  44. Types.matrix_t x,
  45. ICellFunc f) := FUNCTION
  46. Cell := {Types.value_t v};
  47. Cell applyFunc(Cell v, UNSIGNED pos) := TRANSFORM
  48. r := ((pos-1) % m) + 1;
  49. c := ((pos-1) DIV m) + 1;
  50. SELF.v := f(v.v, r, c);
  51. END;
  52. dIn := DATASET(x, Cell);
  53. dOut := PROJECT(dIn, applyFunc(LEFT, COUNTER));
  54. RETURN SET(dOut, v);
  55. END;
  56. /**
  57. * Absolute sum, the 1 norm of a vector.
  58. *@param m the number of entries
  59. *@param x the column major matrix holding the vector
  60. *@param incx the increment for x, 1 in the case of an actual vector
  61. *@param skipped default is zero, the number of entries stepped over
  62. * to get to the first entry
  63. *@return the sum of the absolute values
  64. */
  65. EXPORT Types.value_t
  66. dasum(Types.dimension_t m, Types.matrix_t x,
  67. Types.dimension_t incx, Types.dimension_t skipped=0)
  68. := LIB_ECLBLAS.ECLBLAS.dasum(m, x, incx, skipped);
  69. /**
  70. * alpha*X + Y
  71. * @param N number of elements in vector
  72. * @param alpha the scalar multiplier
  73. * @param X the column major matrix holding the vector X
  74. * @param incX the increment or stride for the vector
  75. * @param Y the column major matrix holding the vector Y
  76. * @param incY the increment or stride of Y
  77. * @param x_skipped number of entries skipped to get to the first X
  78. * @param y_skipped number of entries skipped to get to the first Y
  79. * @return the updated matrix
  80. */
  81. EXPORT Types.matrix_t
  82. daxpy(Types.dimension_t N, Types.value_t alpha, Types.matrix_t X,
  83. Types.dimension_t incX, Types.matrix_t Y, Types.dimension_t incY,
  84. Types.dimension_t x_skipped=0, Types.dimension_t y_skipped=0)
  85. := LIB_ECLBLAS.ECLBLAS.daxpy(N, alpha, X, incX, Y, incY, x_skipped, y_skipped);
  86. /**
  87. * alpha*op(A) op(B) + beta*C where op() is transpose
  88. * @param transposeA true when transpose of A is used
  89. * @param transposeB true when transpose of B is used
  90. * @param M number of rows in product
  91. * @param N number of columns in product
  92. * @param K number of columns/rows for the multiplier/multiplicand
  93. * @param alpha scalar used on A
  94. * @param A matrix A
  95. * @param B matrix B
  96. * @param beta scalar for matrix C
  97. * @param C matrix C or empty
  98. */
  99. EXPORT Types.matrix_t
  100. dgemm(BOOLEAN transposeA, BOOLEAN transposeB,
  101. Types.dimension_t M, Types.dimension_t N, Types.dimension_t K,
  102. Types.value_t alpha, Types.matrix_t A, Types.matrix_t B,
  103. Types.value_t beta=0.0, Types.matrix_t C=[])
  104. := LIB_ECLBLAS.ECLBLAS.dgemm(transposeA, transposeB, M, N, K, alpha, A, B, beta, C);
  105. /**
  106. * Compute LU Factorization of matrix A.
  107. * @param m number of rows of A
  108. * @param n number of columns of A
  109. * @return composite matrix of factors, lower triangle has an
  110. * implied diagonal of ones. Upper triangle has the diagonal of the
  111. * composite.
  112. */
  113. EXPORT Types.matrix_t
  114. dgetf2(Types.dimension_t m, Types.dimension_t n, Types.matrix_t a)
  115. := LIB_ECLBLAS.ECLBLAS.dgetf2(m, n, a);
  116. /**
  117. * DPOTF2 computes the Cholesky factorization of a real symmetric
  118. * positive definite matrix A.
  119. *The factorization has the form
  120. * A = U**T * U , if UPLO = 'U', or
  121. * A = L * L**T, if UPLO = 'L',
  122. * where U is an upper triangular matrix and L is lower triangular.
  123. * This is the unblocked version of the algorithm, calling Level 2 BLAS.
  124. * @param tri indicate whether upper or lower triangle is used
  125. * @param r number of rows/columns in the square matrix
  126. * @param A the square matrix
  127. * @param clear clears the unused triangle
  128. * @return the triangular matrix requested.
  129. */
  130. EXPORT Types.matrix_t
  131. dpotf2(Types.Triangle tri, Types.dimension_t r, Types.matrix_t A,
  132. BOOLEAN clear=TRUE)
  133. := LIB_ECLBLAS.ECLBLAS.dpotf2(tri, r, A, clear);
  134. /**
  135. * Scale a vector alpha
  136. * @param N number of elements in the vector
  137. * @param alpha the scaling factor
  138. * @param X the column major matrix holding the vector
  139. * @param incX the stride to get to the next element in the vector
  140. * @param skipped the number of elements skipped to get to the first element
  141. * @return the updated matrix
  142. */
  143. EXPORT Types.matrix_t
  144. dscal(Types.dimension_t N, Types.value_t alpha, Types.matrix_t X,
  145. Types.dimension_t incX, Types.dimension_t skipped=0)
  146. := LIB_ECLBLAS.ECLBLAS.dscal(N, alpha, X, incX, skipped);
  147. /**
  148. * Implements symmetric rank update C <- alpha A**T A + beta C or
  149. * c <- alpha A A**T + beta C. C is N x N.
  150. * @param tri update upper or lower triangle
  151. * @param transposeA Transpose the A matrix to be NxK
  152. * @param N number of rows
  153. * @param K number of columns in the update matrix or transpose
  154. * @param alpha the alpha scalar
  155. * @param A the update matrix, either NxK or KxN
  156. * @param beta the beta scalar
  157. * @param C the matrix to update
  158. * @param clear clear the triangle that is not updated. BLAS assumes
  159. * that symmetric matrices have only one of the triangles and this
  160. * option lets you make that true.
  161. */
  162. EXPORT Types.matrix_t
  163. dsyrk(Types.Triangle tri, BOOLEAN transposeA,
  164. Types.dimension_t N, Types.dimension_t K,
  165. Types.value_t alpha, Types.matrix_t A,
  166. Types.value_t beta, Types.matrix_t C, BOOLEAN clear=FALSE)
  167. := LIB_ECLBLAS.ECLBLAS.dsyrk(tri, transposeA, N, K, alpha, A, beta, C, clear);
  168. /**
  169. * Triangular matrix solver. op(A) X = alpha B or X op(A) = alpha B
  170. * where op is Transpose, X and B is MxN
  171. * @param side side for A, Side.Ax is op(A) X = alpha B
  172. * @param tri Says whether A is Upper or Lower triangle
  173. * @param transposeA is op(A) the transpose of A
  174. * @param diag is the diagonal an implied unit diagonal or supplied
  175. * @param M number of rows
  176. * @param N number of columns
  177. * @param lda the leading dimension of the A matrix, either M or N
  178. * @param alpha the scalar multiplier for B
  179. * @param A a triangular matrix
  180. * @param B the matrix of values for the solve
  181. * @return the matrix of coefficients to get B.
  182. */
  183. EXPORT Types.matrix_t
  184. dtrsm(Types.Side side, Types.Triangle tri,
  185. BOOLEAN transposeA, Types.Diagonal diag,
  186. Types.dimension_t M, Types.dimension_t N, Types.dimension_t lda,
  187. Types.value_t alpha, Types.matrix_t A, Types.matrix_t B)
  188. := LIB_ECLBLAS.ECLBLAS.dtrsm(side, tri, transposeA, diag, M, N, lda, alpha, A, B);
  189. /**
  190. * Extract the diagonal of he matrix
  191. * @param m number of rows
  192. * @param n number of columns
  193. * @param x matrix from which to extract the diagonal
  194. * @return diagonal matrix
  195. */
  196. EXPORT Types.matrix_t extract_diag(Types.dimension_t m,
  197. Types.dimension_t n,
  198. Types.matrix_t x) := FUNCTION
  199. cell := {Types.value_t v};
  200. cell ext(cell v, UNSIGNED pos) := TRANSFORM
  201. r := ((pos-1) % m) + 1;
  202. c := ((pos-1) DIV m) + 1;
  203. SELF.v := IF(r=c AND r<=m AND c<=n, v.v, SKIP);
  204. END;
  205. diag := SET(PROJECT(DATASET(x, cell), ext(LEFT, COUNTER)), v);
  206. RETURN diag;
  207. END;
  208. /**
  209. * Extract the upper or lower triangle. Diagonal can be actual or implied
  210. * unit diagonal.
  211. * @param m number of rows
  212. * @param n number of columns
  213. * @param tri Upper or Lower specifier, Triangle.Lower or Triangle.Upper
  214. * @param dt Use Diagonal.NotUnitTri or Diagonal.UnitTri
  215. * @param a Matrix, usually a composite from factoring
  216. * @return the triangle
  217. */
  218. EXPORT Types.matrix_t
  219. extract_tri(Types.dimension_t m, Types.dimension_t n,
  220. Types.Triangle tri, Types.Diagonal dt,
  221. Types.matrix_t a)
  222. := LIB_ECLBLAS.ECLBLAS.extract_tri(m, n, tri, dt, a);
  223. /**
  224. * Generate a diagonal matrix.
  225. * @param m number of diagonal entries
  226. * @param v option value, defaults to 1
  227. * @param X optional input of diagonal values, multiplied by v.
  228. * @return a diagonal matrix
  229. */
  230. EXPORT Types.matrix_t
  231. make_diag(Types.dimension_t m, Types.value_t v=1.0,
  232. Types.matrix_t X=[]) := LIB_ECLBLAS.ECLBLAS.make_diag(m, v, X);
  233. // make_vec helpers
  234. Cell := RECORD
  235. Types.value_t v;
  236. END;
  237. Cell makeCell(Types.value_t v) := TRANSFORM
  238. SELF.v := v;
  239. END;
  240. vec_dataset(Types.dimension_t m, Types.value_t v) := DATASET(m, makeCell(v));
  241. /**
  242. * Make a vector of dimension m
  243. * @param m number of elements
  244. * @param v the values, defaults to 1
  245. * @return the vector
  246. */
  247. EXPORT Types.matrix_t
  248. make_vector(Types.dimension_t m,
  249. Types.value_t v=1.0) := SET(vec_dataset(m, v), v);
  250. /**
  251. * The trace of the input matrix
  252. * @param m number of rows
  253. * @param n number of columns
  254. * @param x the matrix
  255. * @return the trace (sum of the diagonal entries)
  256. */
  257. EXPORT Types.value_t
  258. trace(Types.dimension_t m,
  259. Types.dimension_t n,
  260. Types.matrix_t x)
  261. := SUM(extract_diag(m,n,x));
  262. END;