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.

457 lines
19 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 "common.h"
  10. int EIGEN_BLAS_FUNC(gemv)(char *opa, int *m, int *n, RealScalar *palpha, RealScalar *pa, int *lda, RealScalar *pb, int *incb, RealScalar *pbeta, RealScalar *pc, int *incc)
  11. {
  12. typedef void (*functype)(int, int, const Scalar *, int, const Scalar *, int , Scalar *, int, Scalar);
  13. static functype func[4];
  14. static bool init = false;
  15. if(!init)
  16. {
  17. for(int k=0; k<4; ++k)
  18. func[k] = 0;
  19. func[NOTR] = (internal::general_matrix_vector_product<int,Scalar,ColMajor,false,Scalar,false>::run);
  20. func[TR ] = (internal::general_matrix_vector_product<int,Scalar,RowMajor,false,Scalar,false>::run);
  21. func[ADJ ] = (internal::general_matrix_vector_product<int,Scalar,RowMajor,Conj, Scalar,false>::run);
  22. init = true;
  23. }
  24. Scalar* a = reinterpret_cast<Scalar*>(pa);
  25. Scalar* b = reinterpret_cast<Scalar*>(pb);
  26. Scalar* c = reinterpret_cast<Scalar*>(pc);
  27. Scalar alpha = *reinterpret_cast<Scalar*>(palpha);
  28. Scalar beta = *reinterpret_cast<Scalar*>(pbeta);
  29. // check arguments
  30. int info = 0;
  31. if(OP(*opa)==INVALID) info = 1;
  32. else if(*m<0) info = 2;
  33. else if(*n<0) info = 3;
  34. else if(*lda<std::max(1,*m)) info = 6;
  35. else if(*incb==0) info = 8;
  36. else if(*incc==0) info = 11;
  37. if(info)
  38. return xerbla_(SCALAR_SUFFIX_UP"GEMV ",&info,6);
  39. if(*m==0 || *n==0 || (alpha==Scalar(0) && beta==Scalar(1)))
  40. return 0;
  41. int actual_m = *m;
  42. int actual_n = *n;
  43. if(OP(*opa)!=NOTR)
  44. std::swap(actual_m,actual_n);
  45. Scalar* actual_b = get_compact_vector(b,actual_n,*incb);
  46. Scalar* actual_c = get_compact_vector(c,actual_m,*incc);
  47. if(beta!=Scalar(1))
  48. {
  49. if(beta==Scalar(0)) vector(actual_c, actual_m).setZero();
  50. else vector(actual_c, actual_m) *= beta;
  51. }
  52. int code = OP(*opa);
  53. func[code](actual_m, actual_n, a, *lda, actual_b, 1, actual_c, 1, alpha);
  54. if(actual_b!=b) delete[] actual_b;
  55. if(actual_c!=c) delete[] copy_back(actual_c,c,actual_m,*incc);
  56. return 1;
  57. }
  58. int EIGEN_BLAS_FUNC(trsv)(char *uplo, char *opa, char *diag, int *n, RealScalar *pa, int *lda, RealScalar *pb, int *incb)
  59. {
  60. typedef void (*functype)(int, const Scalar *, int, Scalar *);
  61. static functype func[16];
  62. static bool init = false;
  63. if(!init)
  64. {
  65. for(int k=0; k<16; ++k)
  66. func[k] = 0;
  67. func[NOTR | (UP << 2) | (NUNIT << 3)] = (internal::triangular_solve_vector<Scalar,Scalar,int,OnTheLeft, Upper|0, false,ColMajor>::run);
  68. func[TR | (UP << 2) | (NUNIT << 3)] = (internal::triangular_solve_vector<Scalar,Scalar,int,OnTheLeft, Lower|0, false,RowMajor>::run);
  69. func[ADJ | (UP << 2) | (NUNIT << 3)] = (internal::triangular_solve_vector<Scalar,Scalar,int,OnTheLeft, Lower|0, Conj, RowMajor>::run);
  70. func[NOTR | (LO << 2) | (NUNIT << 3)] = (internal::triangular_solve_vector<Scalar,Scalar,int,OnTheLeft, Lower|0, false,ColMajor>::run);
  71. func[TR | (LO << 2) | (NUNIT << 3)] = (internal::triangular_solve_vector<Scalar,Scalar,int,OnTheLeft, Upper|0, false,RowMajor>::run);
  72. func[ADJ | (LO << 2) | (NUNIT << 3)] = (internal::triangular_solve_vector<Scalar,Scalar,int,OnTheLeft, Upper|0, Conj, RowMajor>::run);
  73. func[NOTR | (UP << 2) | (UNIT << 3)] = (internal::triangular_solve_vector<Scalar,Scalar,int,OnTheLeft, Upper|UnitDiag,false,ColMajor>::run);
  74. func[TR | (UP << 2) | (UNIT << 3)] = (internal::triangular_solve_vector<Scalar,Scalar,int,OnTheLeft, Lower|UnitDiag,false,RowMajor>::run);
  75. func[ADJ | (UP << 2) | (UNIT << 3)] = (internal::triangular_solve_vector<Scalar,Scalar,int,OnTheLeft, Lower|UnitDiag,Conj, RowMajor>::run);
  76. func[NOTR | (LO << 2) | (UNIT << 3)] = (internal::triangular_solve_vector<Scalar,Scalar,int,OnTheLeft, Lower|UnitDiag,false,ColMajor>::run);
  77. func[TR | (LO << 2) | (UNIT << 3)] = (internal::triangular_solve_vector<Scalar,Scalar,int,OnTheLeft, Upper|UnitDiag,false,RowMajor>::run);
  78. func[ADJ | (LO << 2) | (UNIT << 3)] = (internal::triangular_solve_vector<Scalar,Scalar,int,OnTheLeft, Upper|UnitDiag,Conj, RowMajor>::run);
  79. init = true;
  80. }
  81. Scalar* a = reinterpret_cast<Scalar*>(pa);
  82. Scalar* b = reinterpret_cast<Scalar*>(pb);
  83. int info = 0;
  84. if(UPLO(*uplo)==INVALID) info = 1;
  85. else if(OP(*opa)==INVALID) info = 2;
  86. else if(DIAG(*diag)==INVALID) info = 3;
  87. else if(*n<0) info = 4;
  88. else if(*lda<std::max(1,*n)) info = 6;
  89. else if(*incb==0) info = 8;
  90. if(info)
  91. return xerbla_(SCALAR_SUFFIX_UP"TRSV ",&info,6);
  92. Scalar* actual_b = get_compact_vector(b,*n,*incb);
  93. int code = OP(*opa) | (UPLO(*uplo) << 2) | (DIAG(*diag) << 3);
  94. func[code](*n, a, *lda, actual_b);
  95. if(actual_b!=b) delete[] copy_back(actual_b,b,*n,*incb);
  96. return 0;
  97. }
  98. int EIGEN_BLAS_FUNC(trmv)(char *uplo, char *opa, char *diag, int *n, RealScalar *pa, int *lda, RealScalar *pb, int *incb)
  99. {
  100. typedef void (*functype)(int, int, const Scalar *, int, const Scalar *, int, Scalar *, int, Scalar);
  101. static functype func[16];
  102. static bool init = false;
  103. if(!init)
  104. {
  105. for(int k=0; k<16; ++k)
  106. func[k] = 0;
  107. func[NOTR | (UP << 2) | (NUNIT << 3)] = (internal::triangular_matrix_vector_product<int,Upper|0, Scalar,false,Scalar,false,ColMajor>::run);
  108. func[TR | (UP << 2) | (NUNIT << 3)] = (internal::triangular_matrix_vector_product<int,Lower|0, Scalar,false,Scalar,false,RowMajor>::run);
  109. func[ADJ | (UP << 2) | (NUNIT << 3)] = (internal::triangular_matrix_vector_product<int,Lower|0, Scalar,Conj, Scalar,false,RowMajor>::run);
  110. func[NOTR | (LO << 2) | (NUNIT << 3)] = (internal::triangular_matrix_vector_product<int,Lower|0, Scalar,false,Scalar,false,ColMajor>::run);
  111. func[TR | (LO << 2) | (NUNIT << 3)] = (internal::triangular_matrix_vector_product<int,Upper|0, Scalar,false,Scalar,false,RowMajor>::run);
  112. func[ADJ | (LO << 2) | (NUNIT << 3)] = (internal::triangular_matrix_vector_product<int,Upper|0, Scalar,Conj, Scalar,false,RowMajor>::run);
  113. func[NOTR | (UP << 2) | (UNIT << 3)] = (internal::triangular_matrix_vector_product<int,Upper|UnitDiag,Scalar,false,Scalar,false,ColMajor>::run);
  114. func[TR | (UP << 2) | (UNIT << 3)] = (internal::triangular_matrix_vector_product<int,Lower|UnitDiag,Scalar,false,Scalar,false,RowMajor>::run);
  115. func[ADJ | (UP << 2) | (UNIT << 3)] = (internal::triangular_matrix_vector_product<int,Lower|UnitDiag,Scalar,Conj, Scalar,false,RowMajor>::run);
  116. func[NOTR | (LO << 2) | (UNIT << 3)] = (internal::triangular_matrix_vector_product<int,Lower|UnitDiag,Scalar,false,Scalar,false,ColMajor>::run);
  117. func[TR | (LO << 2) | (UNIT << 3)] = (internal::triangular_matrix_vector_product<int,Upper|UnitDiag,Scalar,false,Scalar,false,RowMajor>::run);
  118. func[ADJ | (LO << 2) | (UNIT << 3)] = (internal::triangular_matrix_vector_product<int,Upper|UnitDiag,Scalar,Conj, Scalar,false,RowMajor>::run);
  119. init = true;
  120. }
  121. Scalar* a = reinterpret_cast<Scalar*>(pa);
  122. Scalar* b = reinterpret_cast<Scalar*>(pb);
  123. int info = 0;
  124. if(UPLO(*uplo)==INVALID) info = 1;
  125. else if(OP(*opa)==INVALID) info = 2;
  126. else if(DIAG(*diag)==INVALID) info = 3;
  127. else if(*n<0) info = 4;
  128. else if(*lda<std::max(1,*n)) info = 6;
  129. else if(*incb==0) info = 8;
  130. if(info)
  131. return xerbla_(SCALAR_SUFFIX_UP"TRMV ",&info,6);
  132. if(*n==0)
  133. return 1;
  134. Scalar* actual_b = get_compact_vector(b,*n,*incb);
  135. Matrix<Scalar,Dynamic,1> res(*n);
  136. res.setZero();
  137. int code = OP(*opa) | (UPLO(*uplo) << 2) | (DIAG(*diag) << 3);
  138. if(code>=16 || func[code]==0)
  139. return 0;
  140. func[code](*n, *n, a, *lda, actual_b, 1, res.data(), 1, Scalar(1));
  141. copy_back(res.data(),b,*n,*incb);
  142. if(actual_b!=b) delete[] actual_b;
  143. return 0;
  144. }
  145. /** GBMV performs one of the matrix-vector operations
  146. *
  147. * y := alpha*A*x + beta*y, or y := alpha*A'*x + beta*y,
  148. *
  149. * where alpha and beta are scalars, x and y are vectors and A is an
  150. * m by n band matrix, with kl sub-diagonals and ku super-diagonals.
  151. */
  152. int EIGEN_BLAS_FUNC(gbmv)(char *trans, int *m, int *n, int *kl, int *ku, RealScalar *palpha, RealScalar *pa, int *lda,
  153. RealScalar *px, int *incx, RealScalar *pbeta, RealScalar *py, int *incy)
  154. {
  155. Scalar* a = reinterpret_cast<Scalar*>(pa);
  156. Scalar* x = reinterpret_cast<Scalar*>(px);
  157. Scalar* y = reinterpret_cast<Scalar*>(py);
  158. Scalar alpha = *reinterpret_cast<Scalar*>(palpha);
  159. Scalar beta = *reinterpret_cast<Scalar*>(pbeta);
  160. int coeff_rows = *kl+*ku+1;
  161. int info = 0;
  162. if(OP(*trans)==INVALID) info = 1;
  163. else if(*m<0) info = 2;
  164. else if(*n<0) info = 3;
  165. else if(*kl<0) info = 4;
  166. else if(*ku<0) info = 5;
  167. else if(*lda<coeff_rows) info = 8;
  168. else if(*incx==0) info = 10;
  169. else if(*incy==0) info = 13;
  170. if(info)
  171. return xerbla_(SCALAR_SUFFIX_UP"GBMV ",&info,6);
  172. if(*m==0 || *n==0 || (alpha==Scalar(0) && beta==Scalar(1)))
  173. return 0;
  174. int actual_m = *m;
  175. int actual_n = *n;
  176. if(OP(*trans)!=NOTR)
  177. std::swap(actual_m,actual_n);
  178. Scalar* actual_x = get_compact_vector(x,actual_n,*incx);
  179. Scalar* actual_y = get_compact_vector(y,actual_m,*incy);
  180. if(beta!=Scalar(1))
  181. {
  182. if(beta==Scalar(0)) vector(actual_y, actual_m).setZero();
  183. else vector(actual_y, actual_m) *= beta;
  184. }
  185. MatrixType mat_coeffs(a,coeff_rows,*n,*lda);
  186. int nb = std::min(*n,(*m)+(*ku));
  187. for(int j=0; j<nb; ++j)
  188. {
  189. int start = std::max(0,j - *ku);
  190. int end = std::min((*m)-1,j + *kl);
  191. int len = end - start + 1;
  192. int offset = (*ku) - j + start;
  193. if(OP(*trans)==NOTR)
  194. vector(actual_y+start,len) += (alpha*actual_x[j]) * mat_coeffs.col(j).segment(offset,len);
  195. else if(OP(*trans)==TR)
  196. actual_y[j] += alpha * ( mat_coeffs.col(j).segment(offset,len).transpose() * vector(actual_x+start,len) ).value();
  197. else
  198. actual_y[j] += alpha * ( mat_coeffs.col(j).segment(offset,len).adjoint() * vector(actual_x+start,len) ).value();
  199. }
  200. if(actual_x!=x) delete[] actual_x;
  201. if(actual_y!=y) delete[] copy_back(actual_y,y,actual_m,*incy);
  202. return 0;
  203. }
  204. #if 0
  205. /** TBMV performs one of the matrix-vector operations
  206. *
  207. * x := A*x, or x := A'*x,
  208. *
  209. * where x is an n element vector and A is an n by n unit, or non-unit,
  210. * upper or lower triangular band matrix, with ( k + 1 ) diagonals.
  211. */
  212. int EIGEN_BLAS_FUNC(tbmv)(char *uplo, char *opa, char *diag, int *n, int *k, RealScalar *pa, int *lda, RealScalar *px, int *incx)
  213. {
  214. Scalar* a = reinterpret_cast<Scalar*>(pa);
  215. Scalar* x = reinterpret_cast<Scalar*>(px);
  216. int coeff_rows = *k + 1;
  217. int info = 0;
  218. if(UPLO(*uplo)==INVALID) info = 1;
  219. else if(OP(*opa)==INVALID) info = 2;
  220. else if(DIAG(*diag)==INVALID) info = 3;
  221. else if(*n<0) info = 4;
  222. else if(*k<0) info = 5;
  223. else if(*lda<coeff_rows) info = 7;
  224. else if(*incx==0) info = 9;
  225. if(info)
  226. return xerbla_(SCALAR_SUFFIX_UP"TBMV ",&info,6);
  227. if(*n==0)
  228. return 0;
  229. int actual_n = *n;
  230. Scalar* actual_x = get_compact_vector(x,actual_n,*incx);
  231. MatrixType mat_coeffs(a,coeff_rows,*n,*lda);
  232. int ku = UPLO(*uplo)==UPPER ? *k : 0;
  233. int kl = UPLO(*uplo)==LOWER ? *k : 0;
  234. for(int j=0; j<*n; ++j)
  235. {
  236. int start = std::max(0,j - ku);
  237. int end = std::min((*m)-1,j + kl);
  238. int len = end - start + 1;
  239. int offset = (ku) - j + start;
  240. if(OP(*trans)==NOTR)
  241. vector(actual_y+start,len) += (alpha*actual_x[j]) * mat_coeffs.col(j).segment(offset,len);
  242. else if(OP(*trans)==TR)
  243. actual_y[j] += alpha * ( mat_coeffs.col(j).segment(offset,len).transpose() * vector(actual_x+start,len) ).value();
  244. else
  245. actual_y[j] += alpha * ( mat_coeffs.col(j).segment(offset,len).adjoint() * vector(actual_x+start,len) ).value();
  246. }
  247. if(actual_x!=x) delete[] actual_x;
  248. if(actual_y!=y) delete[] copy_back(actual_y,y,actual_m,*incy);
  249. return 0;
  250. }
  251. #endif
  252. /** DTBSV solves one of the systems of equations
  253. *
  254. * A*x = b, or A'*x = b,
  255. *
  256. * where b and x are n element vectors and A is an n by n unit, or
  257. * non-unit, upper or lower triangular band matrix, with ( k + 1 )
  258. * diagonals.
  259. *
  260. * No test for singularity or near-singularity is included in this
  261. * routine. Such tests must be performed before calling this routine.
  262. */
  263. int EIGEN_BLAS_FUNC(tbsv)(char *uplo, char *op, char *diag, int *n, int *k, RealScalar *pa, int *lda, RealScalar *px, int *incx)
  264. {
  265. typedef void (*functype)(int, int, const Scalar *, int, Scalar *);
  266. static functype func[16];
  267. static bool init = false;
  268. if(!init)
  269. {
  270. for(int k=0; k<16; ++k)
  271. func[k] = 0;
  272. func[NOTR | (UP << 2) | (NUNIT << 3)] = (internal::band_solve_triangular_selector<int,Upper|0, Scalar,false,Scalar,ColMajor>::run);
  273. func[TR | (UP << 2) | (NUNIT << 3)] = (internal::band_solve_triangular_selector<int,Lower|0, Scalar,false,Scalar,RowMajor>::run);
  274. func[ADJ | (UP << 2) | (NUNIT << 3)] = (internal::band_solve_triangular_selector<int,Lower|0, Scalar,Conj, Scalar,RowMajor>::run);
  275. func[NOTR | (LO << 2) | (NUNIT << 3)] = (internal::band_solve_triangular_selector<int,Lower|0, Scalar,false,Scalar,ColMajor>::run);
  276. func[TR | (LO << 2) | (NUNIT << 3)] = (internal::band_solve_triangular_selector<int,Upper|0, Scalar,false,Scalar,RowMajor>::run);
  277. func[ADJ | (LO << 2) | (NUNIT << 3)] = (internal::band_solve_triangular_selector<int,Upper|0, Scalar,Conj, Scalar,RowMajor>::run);
  278. func[NOTR | (UP << 2) | (UNIT << 3)] = (internal::band_solve_triangular_selector<int,Upper|UnitDiag,Scalar,false,Scalar,ColMajor>::run);
  279. func[TR | (UP << 2) | (UNIT << 3)] = (internal::band_solve_triangular_selector<int,Lower|UnitDiag,Scalar,false,Scalar,RowMajor>::run);
  280. func[ADJ | (UP << 2) | (UNIT << 3)] = (internal::band_solve_triangular_selector<int,Lower|UnitDiag,Scalar,Conj, Scalar,RowMajor>::run);
  281. func[NOTR | (LO << 2) | (UNIT << 3)] = (internal::band_solve_triangular_selector<int,Lower|UnitDiag,Scalar,false,Scalar,ColMajor>::run);
  282. func[TR | (LO << 2) | (UNIT << 3)] = (internal::band_solve_triangular_selector<int,Upper|UnitDiag,Scalar,false,Scalar,RowMajor>::run);
  283. func[ADJ | (LO << 2) | (UNIT << 3)] = (internal::band_solve_triangular_selector<int,Upper|UnitDiag,Scalar,Conj, Scalar,RowMajor>::run);
  284. init = true;
  285. }
  286. Scalar* a = reinterpret_cast<Scalar*>(pa);
  287. Scalar* x = reinterpret_cast<Scalar*>(px);
  288. int coeff_rows = *k+1;
  289. int info = 0;
  290. if(UPLO(*uplo)==INVALID) info = 1;
  291. else if(OP(*op)==INVALID) info = 2;
  292. else if(DIAG(*diag)==INVALID) info = 3;
  293. else if(*n<0) info = 4;
  294. else if(*k<0) info = 5;
  295. else if(*lda<coeff_rows) info = 7;
  296. else if(*incx==0) info = 9;
  297. if(info)
  298. return xerbla_(SCALAR_SUFFIX_UP"TBSV ",&info,6);
  299. if(*n==0 || (*k==0 && DIAG(*diag)==UNIT))
  300. return 0;
  301. int actual_n = *n;
  302. Scalar* actual_x = get_compact_vector(x,actual_n,*incx);
  303. int code = OP(*op) | (UPLO(*uplo) << 2) | (DIAG(*diag) << 3);
  304. if(code>=16 || func[code]==0)
  305. return 0;
  306. func[code](*n, *k, a, *lda, actual_x);
  307. if(actual_x!=x) delete[] copy_back(actual_x,x,actual_n,*incx);
  308. return 0;
  309. }
  310. /** DTPMV performs one of the matrix-vector operations
  311. *
  312. * x := A*x, or x := A'*x,
  313. *
  314. * where x is an n element vector and A is an n by n unit, or non-unit,
  315. * upper or lower triangular matrix, supplied in packed form.
  316. */
  317. // int EIGEN_BLAS_FUNC(tpmv)(char *uplo, char *trans, char *diag, int *n, RealScalar *ap, RealScalar *x, int *incx)
  318. // {
  319. // return 1;
  320. // }
  321. /** DTPSV solves one of the systems of equations
  322. *
  323. * A*x = b, or A'*x = b,
  324. *
  325. * where b and x are n element vectors and A is an n by n unit, or
  326. * non-unit, upper or lower triangular matrix, supplied in packed form.
  327. *
  328. * No test for singularity or near-singularity is included in this
  329. * routine. Such tests must be performed before calling this routine.
  330. */
  331. // int EIGEN_BLAS_FUNC(tpsv)(char *uplo, char *trans, char *diag, int *n, RealScalar *ap, RealScalar *x, int *incx)
  332. // {
  333. // return 1;
  334. // }
  335. /** DGER performs the rank 1 operation
  336. *
  337. * A := alpha*x*y' + A,
  338. *
  339. * where alpha is a scalar, x is an m element vector, y is an n element
  340. * vector and A is an m by n matrix.
  341. */
  342. int EIGEN_BLAS_FUNC(ger)(int *m, int *n, Scalar *palpha, Scalar *px, int *incx, Scalar *py, int *incy, Scalar *pa, int *lda)
  343. {
  344. Scalar* x = reinterpret_cast<Scalar*>(px);
  345. Scalar* y = reinterpret_cast<Scalar*>(py);
  346. Scalar* a = reinterpret_cast<Scalar*>(pa);
  347. Scalar alpha = *reinterpret_cast<Scalar*>(palpha);
  348. int info = 0;
  349. if(*m<0) info = 1;
  350. else if(*n<0) info = 2;
  351. else if(*incx==0) info = 5;
  352. else if(*incy==0) info = 7;
  353. else if(*lda<std::max(1,*m)) info = 9;
  354. if(info)
  355. return xerbla_(SCALAR_SUFFIX_UP"GER ",&info,6);
  356. if(alpha==Scalar(0))
  357. return 1;
  358. Scalar* x_cpy = get_compact_vector(x,*m,*incx);
  359. Scalar* y_cpy = get_compact_vector(y,*n,*incy);
  360. // TODO perform direct calls to underlying implementation
  361. matrix(a,*m,*n,*lda) += alpha * vector(x_cpy,*m) * vector(y_cpy,*n).adjoint();
  362. if(x_cpy!=x) delete[] x_cpy;
  363. if(y_cpy!=y) delete[] y_cpy;
  364. return 1;
  365. }