You can not select more than 25 topics Topics must start with a letter or number, can include dashes ('-') and can be up to 35 characters long.

652 lines
35 KiB

  1. // This file is part of Eigen, a lightweight C++ template library
  2. // for linear algebra.
  3. //
  4. // Copyright (C) 2009-2010 Gael Guennebaud <gael.guennebaud@inria.fr>
  5. //
  6. // This Source Code Form is subject to the terms of the Mozilla
  7. // Public License v. 2.0. If a copy of the MPL was not distributed
  8. // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
  9. #include <iostream>
  10. #include "common.h"
  11. int EIGEN_BLAS_FUNC(gemm)(char *opa, char *opb, int *m, int *n, int *k, RealScalar *palpha, RealScalar *pa, int *lda, RealScalar *pb, int *ldb, RealScalar *pbeta, RealScalar *pc, int *ldc)
  12. {
  13. // std::cerr << "in gemm " << *opa << " " << *opb << " " << *m << " " << *n << " " << *k << " " << *lda << " " << *ldb << " " << *ldc << " " << *palpha << " " << *pbeta << "\n";
  14. typedef void (*functype)(DenseIndex, DenseIndex, DenseIndex, const Scalar *, DenseIndex, const Scalar *, DenseIndex, Scalar *, DenseIndex, Scalar, internal::level3_blocking<Scalar,Scalar>&, Eigen::internal::GemmParallelInfo<DenseIndex>*);
  15. static functype func[12];
  16. static bool init = false;
  17. if(!init)
  18. {
  19. for(int i=0; i<12; ++i)
  20. func[i] = 0;
  21. func[NOTR | (NOTR << 2)] = (internal::general_matrix_matrix_product<DenseIndex,Scalar,ColMajor,false,Scalar,ColMajor,false,ColMajor>::run);
  22. func[TR | (NOTR << 2)] = (internal::general_matrix_matrix_product<DenseIndex,Scalar,RowMajor,false,Scalar,ColMajor,false,ColMajor>::run);
  23. func[ADJ | (NOTR << 2)] = (internal::general_matrix_matrix_product<DenseIndex,Scalar,RowMajor,Conj, Scalar,ColMajor,false,ColMajor>::run);
  24. func[NOTR | (TR << 2)] = (internal::general_matrix_matrix_product<DenseIndex,Scalar,ColMajor,false,Scalar,RowMajor,false,ColMajor>::run);
  25. func[TR | (TR << 2)] = (internal::general_matrix_matrix_product<DenseIndex,Scalar,RowMajor,false,Scalar,RowMajor,false,ColMajor>::run);
  26. func[ADJ | (TR << 2)] = (internal::general_matrix_matrix_product<DenseIndex,Scalar,RowMajor,Conj, Scalar,RowMajor,false,ColMajor>::run);
  27. func[NOTR | (ADJ << 2)] = (internal::general_matrix_matrix_product<DenseIndex,Scalar,ColMajor,false,Scalar,RowMajor,Conj, ColMajor>::run);
  28. func[TR | (ADJ << 2)] = (internal::general_matrix_matrix_product<DenseIndex,Scalar,RowMajor,false,Scalar,RowMajor,Conj, ColMajor>::run);
  29. func[ADJ | (ADJ << 2)] = (internal::general_matrix_matrix_product<DenseIndex,Scalar,RowMajor,Conj, Scalar,RowMajor,Conj, ColMajor>::run);
  30. init = true;
  31. }
  32. Scalar* a = reinterpret_cast<Scalar*>(pa);
  33. Scalar* b = reinterpret_cast<Scalar*>(pb);
  34. Scalar* c = reinterpret_cast<Scalar*>(pc);
  35. Scalar alpha = *reinterpret_cast<Scalar*>(palpha);
  36. Scalar beta = *reinterpret_cast<Scalar*>(pbeta);
  37. int info = 0;
  38. if(OP(*opa)==INVALID) info = 1;
  39. else if(OP(*opb)==INVALID) info = 2;
  40. else if(*m<0) info = 3;
  41. else if(*n<0) info = 4;
  42. else if(*k<0) info = 5;
  43. else if(*lda<std::max(1,(OP(*opa)==NOTR)?*m:*k)) info = 8;
  44. else if(*ldb<std::max(1,(OP(*opb)==NOTR)?*k:*n)) info = 10;
  45. else if(*ldc<std::max(1,*m)) info = 13;
  46. if(info)
  47. return xerbla_(SCALAR_SUFFIX_UP"GEMM ",&info,6);
  48. if (*m == 0 || *n == 0)
  49. return 0;
  50. if(beta!=Scalar(1))
  51. {
  52. if(beta==Scalar(0)) matrix(c, *m, *n, *ldc).setZero();
  53. else matrix(c, *m, *n, *ldc) *= beta;
  54. }
  55. if(*k == 0)
  56. return 0;
  57. internal::gemm_blocking_space<ColMajor,Scalar,Scalar,Dynamic,Dynamic,Dynamic> blocking(*m,*n,*k,1,true);
  58. int code = OP(*opa) | (OP(*opb) << 2);
  59. func[code](*m, *n, *k, a, *lda, b, *ldb, c, *ldc, alpha, blocking, 0);
  60. return 0;
  61. }
  62. int EIGEN_BLAS_FUNC(trsm)(char *side, char *uplo, char *opa, char *diag, int *m, int *n, RealScalar *palpha, RealScalar *pa, int *lda, RealScalar *pb, int *ldb)
  63. {
  64. // std::cerr << "in trsm " << *side << " " << *uplo << " " << *opa << " " << *diag << " " << *m << "," << *n << " " << *palpha << " " << *lda << " " << *ldb<< "\n";
  65. typedef void (*functype)(DenseIndex, DenseIndex, const Scalar *, DenseIndex, Scalar *, DenseIndex, internal::level3_blocking<Scalar,Scalar>&);
  66. static functype func[32];
  67. static bool init = false;
  68. if(!init)
  69. {
  70. for(int i=0; i<32; ++i)
  71. func[i] = 0;
  72. func[NOTR | (LEFT << 2) | (UP << 3) | (NUNIT << 4)] = (internal::triangular_solve_matrix<Scalar,DenseIndex,OnTheLeft, Upper|0, false,ColMajor,ColMajor>::run);
  73. func[TR | (LEFT << 2) | (UP << 3) | (NUNIT << 4)] = (internal::triangular_solve_matrix<Scalar,DenseIndex,OnTheLeft, Lower|0, false,RowMajor,ColMajor>::run);
  74. func[ADJ | (LEFT << 2) | (UP << 3) | (NUNIT << 4)] = (internal::triangular_solve_matrix<Scalar,DenseIndex,OnTheLeft, Lower|0, Conj, RowMajor,ColMajor>::run);
  75. func[NOTR | (RIGHT << 2) | (UP << 3) | (NUNIT << 4)] = (internal::triangular_solve_matrix<Scalar,DenseIndex,OnTheRight,Upper|0, false,ColMajor,ColMajor>::run);
  76. func[TR | (RIGHT << 2) | (UP << 3) | (NUNIT << 4)] = (internal::triangular_solve_matrix<Scalar,DenseIndex,OnTheRight,Lower|0, false,RowMajor,ColMajor>::run);
  77. func[ADJ | (RIGHT << 2) | (UP << 3) | (NUNIT << 4)] = (internal::triangular_solve_matrix<Scalar,DenseIndex,OnTheRight,Lower|0, Conj, RowMajor,ColMajor>::run);
  78. func[NOTR | (LEFT << 2) | (LO << 3) | (NUNIT << 4)] = (internal::triangular_solve_matrix<Scalar,DenseIndex,OnTheLeft, Lower|0, false,ColMajor,ColMajor>::run);
  79. func[TR | (LEFT << 2) | (LO << 3) | (NUNIT << 4)] = (internal::triangular_solve_matrix<Scalar,DenseIndex,OnTheLeft, Upper|0, false,RowMajor,ColMajor>::run);
  80. func[ADJ | (LEFT << 2) | (LO << 3) | (NUNIT << 4)] = (internal::triangular_solve_matrix<Scalar,DenseIndex,OnTheLeft, Upper|0, Conj, RowMajor,ColMajor>::run);
  81. func[NOTR | (RIGHT << 2) | (LO << 3) | (NUNIT << 4)] = (internal::triangular_solve_matrix<Scalar,DenseIndex,OnTheRight,Lower|0, false,ColMajor,ColMajor>::run);
  82. func[TR | (RIGHT << 2) | (LO << 3) | (NUNIT << 4)] = (internal::triangular_solve_matrix<Scalar,DenseIndex,OnTheRight,Upper|0, false,RowMajor,ColMajor>::run);
  83. func[ADJ | (RIGHT << 2) | (LO << 3) | (NUNIT << 4)] = (internal::triangular_solve_matrix<Scalar,DenseIndex,OnTheRight,Upper|0, Conj, RowMajor,ColMajor>::run);
  84. func[NOTR | (LEFT << 2) | (UP << 3) | (UNIT << 4)] = (internal::triangular_solve_matrix<Scalar,DenseIndex,OnTheLeft, Upper|UnitDiag,false,ColMajor,ColMajor>::run);
  85. func[TR | (LEFT << 2) | (UP << 3) | (UNIT << 4)] = (internal::triangular_solve_matrix<Scalar,DenseIndex,OnTheLeft, Lower|UnitDiag,false,RowMajor,ColMajor>::run);
  86. func[ADJ | (LEFT << 2) | (UP << 3) | (UNIT << 4)] = (internal::triangular_solve_matrix<Scalar,DenseIndex,OnTheLeft, Lower|UnitDiag,Conj, RowMajor,ColMajor>::run);
  87. func[NOTR | (RIGHT << 2) | (UP << 3) | (UNIT << 4)] = (internal::triangular_solve_matrix<Scalar,DenseIndex,OnTheRight,Upper|UnitDiag,false,ColMajor,ColMajor>::run);
  88. func[TR | (RIGHT << 2) | (UP << 3) | (UNIT << 4)] = (internal::triangular_solve_matrix<Scalar,DenseIndex,OnTheRight,Lower|UnitDiag,false,RowMajor,ColMajor>::run);
  89. func[ADJ | (RIGHT << 2) | (UP << 3) | (UNIT << 4)] = (internal::triangular_solve_matrix<Scalar,DenseIndex,OnTheRight,Lower|UnitDiag,Conj, RowMajor,ColMajor>::run);
  90. func[NOTR | (LEFT << 2) | (LO << 3) | (UNIT << 4)] = (internal::triangular_solve_matrix<Scalar,DenseIndex,OnTheLeft, Lower|UnitDiag,false,ColMajor,ColMajor>::run);
  91. func[TR | (LEFT << 2) | (LO << 3) | (UNIT << 4)] = (internal::triangular_solve_matrix<Scalar,DenseIndex,OnTheLeft, Upper|UnitDiag,false,RowMajor,ColMajor>::run);
  92. func[ADJ | (LEFT << 2) | (LO << 3) | (UNIT << 4)] = (internal::triangular_solve_matrix<Scalar,DenseIndex,OnTheLeft, Upper|UnitDiag,Conj, RowMajor,ColMajor>::run);
  93. func[NOTR | (RIGHT << 2) | (LO << 3) | (UNIT << 4)] = (internal::triangular_solve_matrix<Scalar,DenseIndex,OnTheRight,Lower|UnitDiag,false,ColMajor,ColMajor>::run);
  94. func[TR | (RIGHT << 2) | (LO << 3) | (UNIT << 4)] = (internal::triangular_solve_matrix<Scalar,DenseIndex,OnTheRight,Upper|UnitDiag,false,RowMajor,ColMajor>::run);
  95. func[ADJ | (RIGHT << 2) | (LO << 3) | (UNIT << 4)] = (internal::triangular_solve_matrix<Scalar,DenseIndex,OnTheRight,Upper|UnitDiag,Conj, RowMajor,ColMajor>::run);
  96. init = true;
  97. }
  98. Scalar* a = reinterpret_cast<Scalar*>(pa);
  99. Scalar* b = reinterpret_cast<Scalar*>(pb);
  100. Scalar alpha = *reinterpret_cast<Scalar*>(palpha);
  101. int info = 0;
  102. if(SIDE(*side)==INVALID) info = 1;
  103. else if(UPLO(*uplo)==INVALID) info = 2;
  104. else if(OP(*opa)==INVALID) info = 3;
  105. else if(DIAG(*diag)==INVALID) info = 4;
  106. else if(*m<0) info = 5;
  107. else if(*n<0) info = 6;
  108. else if(*lda<std::max(1,(SIDE(*side)==LEFT)?*m:*n)) info = 9;
  109. else if(*ldb<std::max(1,*m)) info = 11;
  110. if(info)
  111. return xerbla_(SCALAR_SUFFIX_UP"TRSM ",&info,6);
  112. if(*m==0 || *n==0)
  113. return 0;
  114. int code = OP(*opa) | (SIDE(*side) << 2) | (UPLO(*uplo) << 3) | (DIAG(*diag) << 4);
  115. if(SIDE(*side)==LEFT)
  116. {
  117. internal::gemm_blocking_space<ColMajor,Scalar,Scalar,Dynamic,Dynamic,Dynamic,4> blocking(*m,*n,*m,1,false);
  118. func[code](*m, *n, a, *lda, b, *ldb, blocking);
  119. }
  120. else
  121. {
  122. internal::gemm_blocking_space<ColMajor,Scalar,Scalar,Dynamic,Dynamic,Dynamic,4> blocking(*m,*n,*n,1,false);
  123. func[code](*n, *m, a, *lda, b, *ldb, blocking);
  124. }
  125. if(alpha!=Scalar(1))
  126. matrix(b,*m,*n,*ldb) *= alpha;
  127. return 0;
  128. }
  129. // b = alpha*op(a)*b for side = 'L'or'l'
  130. // b = alpha*b*op(a) for side = 'R'or'r'
  131. int EIGEN_BLAS_FUNC(trmm)(char *side, char *uplo, char *opa, char *diag, int *m, int *n, RealScalar *palpha, RealScalar *pa, int *lda, RealScalar *pb, int *ldb)
  132. {
  133. // std::cerr << "in trmm " << *side << " " << *uplo << " " << *opa << " " << *diag << " " << *m << " " << *n << " " << *lda << " " << *ldb << " " << *palpha << "\n";
  134. typedef void (*functype)(DenseIndex, DenseIndex, DenseIndex, const Scalar *, DenseIndex, const Scalar *, DenseIndex, Scalar *, DenseIndex, const Scalar&, internal::level3_blocking<Scalar,Scalar>&);
  135. static functype func[32];
  136. static bool init = false;
  137. if(!init)
  138. {
  139. for(int k=0; k<32; ++k)
  140. func[k] = 0;
  141. func[NOTR | (LEFT << 2) | (UP << 3) | (NUNIT << 4)] = (internal::product_triangular_matrix_matrix<Scalar,DenseIndex,Upper|0, true, ColMajor,false,ColMajor,false,ColMajor>::run);
  142. func[TR | (LEFT << 2) | (UP << 3) | (NUNIT << 4)] = (internal::product_triangular_matrix_matrix<Scalar,DenseIndex,Lower|0, true, RowMajor,false,ColMajor,false,ColMajor>::run);
  143. func[ADJ | (LEFT << 2) | (UP << 3) | (NUNIT << 4)] = (internal::product_triangular_matrix_matrix<Scalar,DenseIndex,Lower|0, true, RowMajor,Conj, ColMajor,false,ColMajor>::run);
  144. func[NOTR | (RIGHT << 2) | (UP << 3) | (NUNIT << 4)] = (internal::product_triangular_matrix_matrix<Scalar,DenseIndex,Upper|0, false,ColMajor,false,ColMajor,false,ColMajor>::run);
  145. func[TR | (RIGHT << 2) | (UP << 3) | (NUNIT << 4)] = (internal::product_triangular_matrix_matrix<Scalar,DenseIndex,Lower|0, false,ColMajor,false,RowMajor,false,ColMajor>::run);
  146. func[ADJ | (RIGHT << 2) | (UP << 3) | (NUNIT << 4)] = (internal::product_triangular_matrix_matrix<Scalar,DenseIndex,Lower|0, false,ColMajor,false,RowMajor,Conj, ColMajor>::run);
  147. func[NOTR | (LEFT << 2) | (LO << 3) | (NUNIT << 4)] = (internal::product_triangular_matrix_matrix<Scalar,DenseIndex,Lower|0, true, ColMajor,false,ColMajor,false,ColMajor>::run);
  148. func[TR | (LEFT << 2) | (LO << 3) | (NUNIT << 4)] = (internal::product_triangular_matrix_matrix<Scalar,DenseIndex,Upper|0, true, RowMajor,false,ColMajor,false,ColMajor>::run);
  149. func[ADJ | (LEFT << 2) | (LO << 3) | (NUNIT << 4)] = (internal::product_triangular_matrix_matrix<Scalar,DenseIndex,Upper|0, true, RowMajor,Conj, ColMajor,false,ColMajor>::run);
  150. func[NOTR | (RIGHT << 2) | (LO << 3) | (NUNIT << 4)] = (internal::product_triangular_matrix_matrix<Scalar,DenseIndex,Lower|0, false,ColMajor,false,ColMajor,false,ColMajor>::run);
  151. func[TR | (RIGHT << 2) | (LO << 3) | (NUNIT << 4)] = (internal::product_triangular_matrix_matrix<Scalar,DenseIndex,Upper|0, false,ColMajor,false,RowMajor,false,ColMajor>::run);
  152. func[ADJ | (RIGHT << 2) | (LO << 3) | (NUNIT << 4)] = (internal::product_triangular_matrix_matrix<Scalar,DenseIndex,Upper|0, false,ColMajor,false,RowMajor,Conj, ColMajor>::run);
  153. func[NOTR | (LEFT << 2) | (UP << 3) | (UNIT << 4)] = (internal::product_triangular_matrix_matrix<Scalar,DenseIndex,Upper|UnitDiag,true, ColMajor,false,ColMajor,false,ColMajor>::run);
  154. func[TR | (LEFT << 2) | (UP << 3) | (UNIT << 4)] = (internal::product_triangular_matrix_matrix<Scalar,DenseIndex,Lower|UnitDiag,true, RowMajor,false,ColMajor,false,ColMajor>::run);
  155. func[ADJ | (LEFT << 2) | (UP << 3) | (UNIT << 4)] = (internal::product_triangular_matrix_matrix<Scalar,DenseIndex,Lower|UnitDiag,true, RowMajor,Conj, ColMajor,false,ColMajor>::run);
  156. func[NOTR | (RIGHT << 2) | (UP << 3) | (UNIT << 4)] = (internal::product_triangular_matrix_matrix<Scalar,DenseIndex,Upper|UnitDiag,false,ColMajor,false,ColMajor,false,ColMajor>::run);
  157. func[TR | (RIGHT << 2) | (UP << 3) | (UNIT << 4)] = (internal::product_triangular_matrix_matrix<Scalar,DenseIndex,Lower|UnitDiag,false,ColMajor,false,RowMajor,false,ColMajor>::run);
  158. func[ADJ | (RIGHT << 2) | (UP << 3) | (UNIT << 4)] = (internal::product_triangular_matrix_matrix<Scalar,DenseIndex,Lower|UnitDiag,false,ColMajor,false,RowMajor,Conj, ColMajor>::run);
  159. func[NOTR | (LEFT << 2) | (LO << 3) | (UNIT << 4)] = (internal::product_triangular_matrix_matrix<Scalar,DenseIndex,Lower|UnitDiag,true, ColMajor,false,ColMajor,false,ColMajor>::run);
  160. func[TR | (LEFT << 2) | (LO << 3) | (UNIT << 4)] = (internal::product_triangular_matrix_matrix<Scalar,DenseIndex,Upper|UnitDiag,true, RowMajor,false,ColMajor,false,ColMajor>::run);
  161. func[ADJ | (LEFT << 2) | (LO << 3) | (UNIT << 4)] = (internal::product_triangular_matrix_matrix<Scalar,DenseIndex,Upper|UnitDiag,true, RowMajor,Conj, ColMajor,false,ColMajor>::run);
  162. func[NOTR | (RIGHT << 2) | (LO << 3) | (UNIT << 4)] = (internal::product_triangular_matrix_matrix<Scalar,DenseIndex,Lower|UnitDiag,false,ColMajor,false,ColMajor,false,ColMajor>::run);
  163. func[TR | (RIGHT << 2) | (LO << 3) | (UNIT << 4)] = (internal::product_triangular_matrix_matrix<Scalar,DenseIndex,Upper|UnitDiag,false,ColMajor,false,RowMajor,false,ColMajor>::run);
  164. func[ADJ | (RIGHT << 2) | (LO << 3) | (UNIT << 4)] = (internal::product_triangular_matrix_matrix<Scalar,DenseIndex,Upper|UnitDiag,false,ColMajor,false,RowMajor,Conj, ColMajor>::run);
  165. init = true;
  166. }
  167. Scalar* a = reinterpret_cast<Scalar*>(pa);
  168. Scalar* b = reinterpret_cast<Scalar*>(pb);
  169. Scalar alpha = *reinterpret_cast<Scalar*>(palpha);
  170. int info = 0;
  171. if(SIDE(*side)==INVALID) info = 1;
  172. else if(UPLO(*uplo)==INVALID) info = 2;
  173. else if(OP(*opa)==INVALID) info = 3;
  174. else if(DIAG(*diag)==INVALID) info = 4;
  175. else if(*m<0) info = 5;
  176. else if(*n<0) info = 6;
  177. else if(*lda<std::max(1,(SIDE(*side)==LEFT)?*m:*n)) info = 9;
  178. else if(*ldb<std::max(1,*m)) info = 11;
  179. if(info)
  180. return xerbla_(SCALAR_SUFFIX_UP"TRMM ",&info,6);
  181. int code = OP(*opa) | (SIDE(*side) << 2) | (UPLO(*uplo) << 3) | (DIAG(*diag) << 4);
  182. if(*m==0 || *n==0)
  183. return 1;
  184. // FIXME find a way to avoid this copy
  185. Matrix<Scalar,Dynamic,Dynamic,ColMajor> tmp = matrix(b,*m,*n,*ldb);
  186. matrix(b,*m,*n,*ldb).setZero();
  187. if(SIDE(*side)==LEFT)
  188. {
  189. internal::gemm_blocking_space<ColMajor,Scalar,Scalar,Dynamic,Dynamic,Dynamic,4> blocking(*m,*n,*m,1,false);
  190. func[code](*m, *n, *m, a, *lda, tmp.data(), tmp.outerStride(), b, *ldb, alpha, blocking);
  191. }
  192. else
  193. {
  194. internal::gemm_blocking_space<ColMajor,Scalar,Scalar,Dynamic,Dynamic,Dynamic,4> blocking(*m,*n,*n,1,false);
  195. func[code](*m, *n, *n, tmp.data(), tmp.outerStride(), a, *lda, b, *ldb, alpha, blocking);
  196. }
  197. return 1;
  198. }
  199. // c = alpha*a*b + beta*c for side = 'L'or'l'
  200. // c = alpha*b*a + beta*c for side = 'R'or'r
  201. int EIGEN_BLAS_FUNC(symm)(char *side, char *uplo, int *m, int *n, RealScalar *palpha, RealScalar *pa, int *lda, RealScalar *pb, int *ldb, RealScalar *pbeta, RealScalar *pc, int *ldc)
  202. {
  203. // std::cerr << "in symm " << *side << " " << *uplo << " " << *m << "x" << *n << " lda:" << *lda << " ldb:" << *ldb << " ldc:" << *ldc << " alpha:" << *palpha << " beta:" << *pbeta << "\n";
  204. Scalar* a = reinterpret_cast<Scalar*>(pa);
  205. Scalar* b = reinterpret_cast<Scalar*>(pb);
  206. Scalar* c = reinterpret_cast<Scalar*>(pc);
  207. Scalar alpha = *reinterpret_cast<Scalar*>(palpha);
  208. Scalar beta = *reinterpret_cast<Scalar*>(pbeta);
  209. int info = 0;
  210. if(SIDE(*side)==INVALID) info = 1;
  211. else if(UPLO(*uplo)==INVALID) info = 2;
  212. else if(*m<0) info = 3;
  213. else if(*n<0) info = 4;
  214. else if(*lda<std::max(1,(SIDE(*side)==LEFT)?*m:*n)) info = 7;
  215. else if(*ldb<std::max(1,*m)) info = 9;
  216. else if(*ldc<std::max(1,*m)) info = 12;
  217. if(info)
  218. return xerbla_(SCALAR_SUFFIX_UP"SYMM ",&info,6);
  219. if(beta!=Scalar(1))
  220. {
  221. if(beta==Scalar(0)) matrix(c, *m, *n, *ldc).setZero();
  222. else matrix(c, *m, *n, *ldc) *= beta;
  223. }
  224. if(*m==0 || *n==0)
  225. {
  226. return 1;
  227. }
  228. #if ISCOMPLEX
  229. // FIXME add support for symmetric complex matrix
  230. int size = (SIDE(*side)==LEFT) ? (*m) : (*n);
  231. Matrix<Scalar,Dynamic,Dynamic,ColMajor> matA(size,size);
  232. if(UPLO(*uplo)==UP)
  233. {
  234. matA.triangularView<Upper>() = matrix(a,size,size,*lda);
  235. matA.triangularView<Lower>() = matrix(a,size,size,*lda).transpose();
  236. }
  237. else if(UPLO(*uplo)==LO)
  238. {
  239. matA.triangularView<Lower>() = matrix(a,size,size,*lda);
  240. matA.triangularView<Upper>() = matrix(a,size,size,*lda).transpose();
  241. }
  242. if(SIDE(*side)==LEFT)
  243. matrix(c, *m, *n, *ldc) += alpha * matA * matrix(b, *m, *n, *ldb);
  244. else if(SIDE(*side)==RIGHT)
  245. matrix(c, *m, *n, *ldc) += alpha * matrix(b, *m, *n, *ldb) * matA;
  246. #else
  247. if(SIDE(*side)==LEFT)
  248. if(UPLO(*uplo)==UP) internal::product_selfadjoint_matrix<Scalar, DenseIndex, RowMajor,true,false, ColMajor,false,false, ColMajor>::run(*m, *n, a, *lda, b, *ldb, c, *ldc, alpha);
  249. else if(UPLO(*uplo)==LO) internal::product_selfadjoint_matrix<Scalar, DenseIndex, ColMajor,true,false, ColMajor,false,false, ColMajor>::run(*m, *n, a, *lda, b, *ldb, c, *ldc, alpha);
  250. else return 0;
  251. else if(SIDE(*side)==RIGHT)
  252. if(UPLO(*uplo)==UP) internal::product_selfadjoint_matrix<Scalar, DenseIndex, ColMajor,false,false, RowMajor,true,false, ColMajor>::run(*m, *n, b, *ldb, a, *lda, c, *ldc, alpha);
  253. else if(UPLO(*uplo)==LO) internal::product_selfadjoint_matrix<Scalar, DenseIndex, ColMajor,false,false, ColMajor,true,false, ColMajor>::run(*m, *n, b, *ldb, a, *lda, c, *ldc, alpha);
  254. else return 0;
  255. else
  256. return 0;
  257. #endif
  258. return 0;
  259. }
  260. // c = alpha*a*a' + beta*c for op = 'N'or'n'
  261. // c = alpha*a'*a + beta*c for op = 'T'or't','C'or'c'
  262. int EIGEN_BLAS_FUNC(syrk)(char *uplo, char *op, int *n, int *k, RealScalar *palpha, RealScalar *pa, int *lda, RealScalar *pbeta, RealScalar *pc, int *ldc)
  263. {
  264. // std::cerr << "in syrk " << *uplo << " " << *op << " " << *n << " " << *k << " " << *palpha << " " << *lda << " " << *pbeta << " " << *ldc << "\n";
  265. #if !ISCOMPLEX
  266. typedef void (*functype)(DenseIndex, DenseIndex, const Scalar *, DenseIndex, const Scalar *, DenseIndex, Scalar *, DenseIndex, const Scalar&);
  267. static functype func[8];
  268. static bool init = false;
  269. if(!init)
  270. {
  271. for(int i=0; i<8; ++i)
  272. func[i] = 0;
  273. func[NOTR | (UP << 2)] = (internal::general_matrix_matrix_triangular_product<DenseIndex,Scalar,ColMajor,false,Scalar,RowMajor,ColMajor,Conj, Upper>::run);
  274. func[TR | (UP << 2)] = (internal::general_matrix_matrix_triangular_product<DenseIndex,Scalar,RowMajor,false,Scalar,ColMajor,ColMajor,Conj, Upper>::run);
  275. func[ADJ | (UP << 2)] = (internal::general_matrix_matrix_triangular_product<DenseIndex,Scalar,RowMajor,Conj, Scalar,ColMajor,ColMajor,false,Upper>::run);
  276. func[NOTR | (LO << 2)] = (internal::general_matrix_matrix_triangular_product<DenseIndex,Scalar,ColMajor,false,Scalar,RowMajor,ColMajor,Conj, Lower>::run);
  277. func[TR | (LO << 2)] = (internal::general_matrix_matrix_triangular_product<DenseIndex,Scalar,RowMajor,false,Scalar,ColMajor,ColMajor,Conj, Lower>::run);
  278. func[ADJ | (LO << 2)] = (internal::general_matrix_matrix_triangular_product<DenseIndex,Scalar,RowMajor,Conj, Scalar,ColMajor,ColMajor,false,Lower>::run);
  279. init = true;
  280. }
  281. #endif
  282. Scalar* a = reinterpret_cast<Scalar*>(pa);
  283. Scalar* c = reinterpret_cast<Scalar*>(pc);
  284. Scalar alpha = *reinterpret_cast<Scalar*>(palpha);
  285. Scalar beta = *reinterpret_cast<Scalar*>(pbeta);
  286. int info = 0;
  287. if(UPLO(*uplo)==INVALID) info = 1;
  288. else if(OP(*op)==INVALID) info = 2;
  289. else if(*n<0) info = 3;
  290. else if(*k<0) info = 4;
  291. else if(*lda<std::max(1,(OP(*op)==NOTR)?*n:*k)) info = 7;
  292. else if(*ldc<std::max(1,*n)) info = 10;
  293. if(info)
  294. return xerbla_(SCALAR_SUFFIX_UP"SYRK ",&info,6);
  295. if(beta!=Scalar(1))
  296. {
  297. if(UPLO(*uplo)==UP)
  298. if(beta==Scalar(0)) matrix(c, *n, *n, *ldc).triangularView<Upper>().setZero();
  299. else matrix(c, *n, *n, *ldc).triangularView<Upper>() *= beta;
  300. else
  301. if(beta==Scalar(0)) matrix(c, *n, *n, *ldc).triangularView<Lower>().setZero();
  302. else matrix(c, *n, *n, *ldc).triangularView<Lower>() *= beta;
  303. }
  304. if(*n==0 || *k==0)
  305. return 0;
  306. #if ISCOMPLEX
  307. // FIXME add support for symmetric complex matrix
  308. if(UPLO(*uplo)==UP)
  309. {
  310. if(OP(*op)==NOTR)
  311. matrix(c, *n, *n, *ldc).triangularView<Upper>() += alpha * matrix(a,*n,*k,*lda) * matrix(a,*n,*k,*lda).transpose();
  312. else
  313. matrix(c, *n, *n, *ldc).triangularView<Upper>() += alpha * matrix(a,*k,*n,*lda).transpose() * matrix(a,*k,*n,*lda);
  314. }
  315. else
  316. {
  317. if(OP(*op)==NOTR)
  318. matrix(c, *n, *n, *ldc).triangularView<Lower>() += alpha * matrix(a,*n,*k,*lda) * matrix(a,*n,*k,*lda).transpose();
  319. else
  320. matrix(c, *n, *n, *ldc).triangularView<Lower>() += alpha * matrix(a,*k,*n,*lda).transpose() * matrix(a,*k,*n,*lda);
  321. }
  322. #else
  323. int code = OP(*op) | (UPLO(*uplo) << 2);
  324. func[code](*n, *k, a, *lda, a, *lda, c, *ldc, alpha);
  325. #endif
  326. return 0;
  327. }
  328. // c = alpha*a*b' + alpha*b*a' + beta*c for op = 'N'or'n'
  329. // c = alpha*a'*b + alpha*b'*a + beta*c for op = 'T'or't'
  330. int EIGEN_BLAS_FUNC(syr2k)(char *uplo, char *op, int *n, int *k, RealScalar *palpha, RealScalar *pa, int *lda, RealScalar *pb, int *ldb, RealScalar *pbeta, RealScalar *pc, int *ldc)
  331. {
  332. Scalar* a = reinterpret_cast<Scalar*>(pa);
  333. Scalar* b = reinterpret_cast<Scalar*>(pb);
  334. Scalar* c = reinterpret_cast<Scalar*>(pc);
  335. Scalar alpha = *reinterpret_cast<Scalar*>(palpha);
  336. Scalar beta = *reinterpret_cast<Scalar*>(pbeta);
  337. // std::cerr << "in syr2k " << *uplo << " " << *op << " " << *n << " " << *k << " " << alpha << " " << *lda << " " << *ldb << " " << beta << " " << *ldc << "\n";
  338. int info = 0;
  339. if(UPLO(*uplo)==INVALID) info = 1;
  340. else if(OP(*op)==INVALID) info = 2;
  341. else if(*n<0) info = 3;
  342. else if(*k<0) info = 4;
  343. else if(*lda<std::max(1,(OP(*op)==NOTR)?*n:*k)) info = 7;
  344. else if(*ldb<std::max(1,(OP(*op)==NOTR)?*n:*k)) info = 9;
  345. else if(*ldc<std::max(1,*n)) info = 12;
  346. if(info)
  347. return xerbla_(SCALAR_SUFFIX_UP"SYR2K",&info,6);
  348. if(beta!=Scalar(1))
  349. {
  350. if(UPLO(*uplo)==UP)
  351. if(beta==Scalar(0)) matrix(c, *n, *n, *ldc).triangularView<Upper>().setZero();
  352. else matrix(c, *n, *n, *ldc).triangularView<Upper>() *= beta;
  353. else
  354. if(beta==Scalar(0)) matrix(c, *n, *n, *ldc).triangularView<Lower>().setZero();
  355. else matrix(c, *n, *n, *ldc).triangularView<Lower>() *= beta;
  356. }
  357. if(*k==0)
  358. return 1;
  359. if(OP(*op)==NOTR)
  360. {
  361. if(UPLO(*uplo)==UP)
  362. {
  363. matrix(c, *n, *n, *ldc).triangularView<Upper>()
  364. += alpha *matrix(a, *n, *k, *lda)*matrix(b, *n, *k, *ldb).transpose()
  365. + alpha*matrix(b, *n, *k, *ldb)*matrix(a, *n, *k, *lda).transpose();
  366. }
  367. else if(UPLO(*uplo)==LO)
  368. matrix(c, *n, *n, *ldc).triangularView<Lower>()
  369. += alpha*matrix(a, *n, *k, *lda)*matrix(b, *n, *k, *ldb).transpose()
  370. + alpha*matrix(b, *n, *k, *ldb)*matrix(a, *n, *k, *lda).transpose();
  371. }
  372. else if(OP(*op)==TR || OP(*op)==ADJ)
  373. {
  374. if(UPLO(*uplo)==UP)
  375. matrix(c, *n, *n, *ldc).triangularView<Upper>()
  376. += alpha*matrix(a, *k, *n, *lda).transpose()*matrix(b, *k, *n, *ldb)
  377. + alpha*matrix(b, *k, *n, *ldb).transpose()*matrix(a, *k, *n, *lda);
  378. else if(UPLO(*uplo)==LO)
  379. matrix(c, *n, *n, *ldc).triangularView<Lower>()
  380. += alpha*matrix(a, *k, *n, *lda).transpose()*matrix(b, *k, *n, *ldb)
  381. + alpha*matrix(b, *k, *n, *ldb).transpose()*matrix(a, *k, *n, *lda);
  382. }
  383. return 0;
  384. }
  385. #if ISCOMPLEX
  386. // c = alpha*a*b + beta*c for side = 'L'or'l'
  387. // c = alpha*b*a + beta*c for side = 'R'or'r
  388. int EIGEN_BLAS_FUNC(hemm)(char *side, char *uplo, int *m, int *n, RealScalar *palpha, RealScalar *pa, int *lda, RealScalar *pb, int *ldb, RealScalar *pbeta, RealScalar *pc, int *ldc)
  389. {
  390. Scalar* a = reinterpret_cast<Scalar*>(pa);
  391. Scalar* b = reinterpret_cast<Scalar*>(pb);
  392. Scalar* c = reinterpret_cast<Scalar*>(pc);
  393. Scalar alpha = *reinterpret_cast<Scalar*>(palpha);
  394. Scalar beta = *reinterpret_cast<Scalar*>(pbeta);
  395. // std::cerr << "in hemm " << *side << " " << *uplo << " " << *m << " " << *n << " " << alpha << " " << *lda << " " << beta << " " << *ldc << "\n";
  396. int info = 0;
  397. if(SIDE(*side)==INVALID) info = 1;
  398. else if(UPLO(*uplo)==INVALID) info = 2;
  399. else if(*m<0) info = 3;
  400. else if(*n<0) info = 4;
  401. else if(*lda<std::max(1,(SIDE(*side)==LEFT)?*m:*n)) info = 7;
  402. else if(*ldb<std::max(1,*m)) info = 9;
  403. else if(*ldc<std::max(1,*m)) info = 12;
  404. if(info)
  405. return xerbla_(SCALAR_SUFFIX_UP"HEMM ",&info,6);
  406. if(beta==Scalar(0)) matrix(c, *m, *n, *ldc).setZero();
  407. else if(beta!=Scalar(1)) matrix(c, *m, *n, *ldc) *= beta;
  408. if(*m==0 || *n==0)
  409. {
  410. return 1;
  411. }
  412. if(SIDE(*side)==LEFT)
  413. {
  414. if(UPLO(*uplo)==UP) internal::product_selfadjoint_matrix<Scalar,DenseIndex,RowMajor,true,Conj, ColMajor,false,false, ColMajor>
  415. ::run(*m, *n, a, *lda, b, *ldb, c, *ldc, alpha);
  416. else if(UPLO(*uplo)==LO) internal::product_selfadjoint_matrix<Scalar,DenseIndex,ColMajor,true,false, ColMajor,false,false, ColMajor>
  417. ::run(*m, *n, a, *lda, b, *ldb, c, *ldc, alpha);
  418. else return 0;
  419. }
  420. else if(SIDE(*side)==RIGHT)
  421. {
  422. if(UPLO(*uplo)==UP) matrix(c,*m,*n,*ldc) += alpha * matrix(b,*m,*n,*ldb) * matrix(a,*n,*n,*lda).selfadjointView<Upper>();/*internal::product_selfadjoint_matrix<Scalar,DenseIndex,ColMajor,false,false, RowMajor,true,Conj, ColMajor>
  423. ::run(*m, *n, b, *ldb, a, *lda, c, *ldc, alpha);*/
  424. else if(UPLO(*uplo)==LO) internal::product_selfadjoint_matrix<Scalar,DenseIndex,ColMajor,false,false, ColMajor,true,false, ColMajor>
  425. ::run(*m, *n, b, *ldb, a, *lda, c, *ldc, alpha);
  426. else return 0;
  427. }
  428. else
  429. {
  430. return 0;
  431. }
  432. return 0;
  433. }
  434. // c = alpha*a*conj(a') + beta*c for op = 'N'or'n'
  435. // c = alpha*conj(a')*a + beta*c for op = 'C'or'c'
  436. int EIGEN_BLAS_FUNC(herk)(char *uplo, char *op, int *n, int *k, RealScalar *palpha, RealScalar *pa, int *lda, RealScalar *pbeta, RealScalar *pc, int *ldc)
  437. {
  438. // std::cerr << "in herk " << *uplo << " " << *op << " " << *n << " " << *k << " " << *palpha << " " << *lda << " " << *pbeta << " " << *ldc << "\n";
  439. typedef void (*functype)(DenseIndex, DenseIndex, const Scalar *, DenseIndex, const Scalar *, DenseIndex, Scalar *, DenseIndex, const Scalar&);
  440. static functype func[8];
  441. static bool init = false;
  442. if(!init)
  443. {
  444. for(int i=0; i<8; ++i)
  445. func[i] = 0;
  446. func[NOTR | (UP << 2)] = (internal::general_matrix_matrix_triangular_product<DenseIndex,Scalar,ColMajor,false,Scalar,RowMajor,Conj, ColMajor,Upper>::run);
  447. func[ADJ | (UP << 2)] = (internal::general_matrix_matrix_triangular_product<DenseIndex,Scalar,RowMajor,Conj, Scalar,ColMajor,false,ColMajor,Upper>::run);
  448. func[NOTR | (LO << 2)] = (internal::general_matrix_matrix_triangular_product<DenseIndex,Scalar,ColMajor,false,Scalar,RowMajor,Conj, ColMajor,Lower>::run);
  449. func[ADJ | (LO << 2)] = (internal::general_matrix_matrix_triangular_product<DenseIndex,Scalar,RowMajor,Conj, Scalar,ColMajor,false,ColMajor,Lower>::run);
  450. init = true;
  451. }
  452. Scalar* a = reinterpret_cast<Scalar*>(pa);
  453. Scalar* c = reinterpret_cast<Scalar*>(pc);
  454. RealScalar alpha = *palpha;
  455. RealScalar beta = *pbeta;
  456. // std::cerr << "in herk " << *uplo << " " << *op << " " << *n << " " << *k << " " << alpha << " " << *lda << " " << beta << " " << *ldc << "\n";
  457. int info = 0;
  458. if(UPLO(*uplo)==INVALID) info = 1;
  459. else if((OP(*op)==INVALID) || (OP(*op)==TR)) info = 2;
  460. else if(*n<0) info = 3;
  461. else if(*k<0) info = 4;
  462. else if(*lda<std::max(1,(OP(*op)==NOTR)?*n:*k)) info = 7;
  463. else if(*ldc<std::max(1,*n)) info = 10;
  464. if(info)
  465. return xerbla_(SCALAR_SUFFIX_UP"HERK ",&info,6);
  466. int code = OP(*op) | (UPLO(*uplo) << 2);
  467. if(beta!=RealScalar(1))
  468. {
  469. if(UPLO(*uplo)==UP)
  470. if(beta==Scalar(0)) matrix(c, *n, *n, *ldc).triangularView<Upper>().setZero();
  471. else matrix(c, *n, *n, *ldc).triangularView<StrictlyUpper>() *= beta;
  472. else
  473. if(beta==Scalar(0)) matrix(c, *n, *n, *ldc).triangularView<Lower>().setZero();
  474. else matrix(c, *n, *n, *ldc).triangularView<StrictlyLower>() *= beta;
  475. if(beta!=Scalar(0))
  476. {
  477. matrix(c, *n, *n, *ldc).diagonal().real() *= beta;
  478. matrix(c, *n, *n, *ldc).diagonal().imag().setZero();
  479. }
  480. }
  481. if(*k>0 && alpha!=RealScalar(0))
  482. {
  483. func[code](*n, *k, a, *lda, a, *lda, c, *ldc, alpha);
  484. matrix(c, *n, *n, *ldc).diagonal().imag().setZero();
  485. }
  486. return 0;
  487. }
  488. // c = alpha*a*conj(b') + conj(alpha)*b*conj(a') + beta*c, for op = 'N'or'n'
  489. // c = alpha*conj(a')*b + conj(alpha)*conj(b')*a + beta*c, for op = 'C'or'c'
  490. int EIGEN_BLAS_FUNC(her2k)(char *uplo, char *op, int *n, int *k, RealScalar *palpha, RealScalar *pa, int *lda, RealScalar *pb, int *ldb, RealScalar *pbeta, RealScalar *pc, int *ldc)
  491. {
  492. Scalar* a = reinterpret_cast<Scalar*>(pa);
  493. Scalar* b = reinterpret_cast<Scalar*>(pb);
  494. Scalar* c = reinterpret_cast<Scalar*>(pc);
  495. Scalar alpha = *reinterpret_cast<Scalar*>(palpha);
  496. RealScalar beta = *pbeta;
  497. // std::cerr << "in her2k " << *uplo << " " << *op << " " << *n << " " << *k << " " << alpha << " " << *lda << " " << *ldb << " " << beta << " " << *ldc << "\n";
  498. int info = 0;
  499. if(UPLO(*uplo)==INVALID) info = 1;
  500. else if((OP(*op)==INVALID) || (OP(*op)==TR)) info = 2;
  501. else if(*n<0) info = 3;
  502. else if(*k<0) info = 4;
  503. else if(*lda<std::max(1,(OP(*op)==NOTR)?*n:*k)) info = 7;
  504. else if(*ldb<std::max(1,(OP(*op)==NOTR)?*n:*k)) info = 9;
  505. else if(*ldc<std::max(1,*n)) info = 12;
  506. if(info)
  507. return xerbla_(SCALAR_SUFFIX_UP"HER2K",&info,6);
  508. if(beta!=RealScalar(1))
  509. {
  510. if(UPLO(*uplo)==UP)
  511. if(beta==Scalar(0)) matrix(c, *n, *n, *ldc).triangularView<Upper>().setZero();
  512. else matrix(c, *n, *n, *ldc).triangularView<StrictlyUpper>() *= beta;
  513. else
  514. if(beta==Scalar(0)) matrix(c, *n, *n, *ldc).triangularView<Lower>().setZero();
  515. else matrix(c, *n, *n, *ldc).triangularView<StrictlyLower>() *= beta;
  516. if(beta!=Scalar(0))
  517. {
  518. matrix(c, *n, *n, *ldc).diagonal().real() *= beta;
  519. matrix(c, *n, *n, *ldc).diagonal().imag().setZero();
  520. }
  521. }
  522. else if(*k>0 && alpha!=Scalar(0))
  523. matrix(c, *n, *n, *ldc).diagonal().imag().setZero();
  524. if(*k==0)
  525. return 1;
  526. if(OP(*op)==NOTR)
  527. {
  528. if(UPLO(*uplo)==UP)
  529. {
  530. matrix(c, *n, *n, *ldc).triangularView<Upper>()
  531. += alpha *matrix(a, *n, *k, *lda)*matrix(b, *n, *k, *ldb).adjoint()
  532. + numext::conj(alpha)*matrix(b, *n, *k, *ldb)*matrix(a, *n, *k, *lda).adjoint();
  533. }
  534. else if(UPLO(*uplo)==LO)
  535. matrix(c, *n, *n, *ldc).triangularView<Lower>()
  536. += alpha*matrix(a, *n, *k, *lda)*matrix(b, *n, *k, *ldb).adjoint()
  537. + numext::conj(alpha)*matrix(b, *n, *k, *ldb)*matrix(a, *n, *k, *lda).adjoint();
  538. }
  539. else if(OP(*op)==ADJ)
  540. {
  541. if(UPLO(*uplo)==UP)
  542. matrix(c, *n, *n, *ldc).triangularView<Upper>()
  543. += alpha*matrix(a, *k, *n, *lda).adjoint()*matrix(b, *k, *n, *ldb)
  544. + numext::conj(alpha)*matrix(b, *k, *n, *ldb).adjoint()*matrix(a, *k, *n, *lda);
  545. else if(UPLO(*uplo)==LO)
  546. matrix(c, *n, *n, *ldc).triangularView<Lower>()
  547. += alpha*matrix(a, *k, *n, *lda).adjoint()*matrix(b, *k, *n, *ldb)
  548. + numext::conj(alpha)*matrix(b, *k, *n, *ldb).adjoint()*matrix(a, *k, *n, *lda);
  549. }
  550. return 1;
  551. }
  552. #endif // ISCOMPLEX