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.

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