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.

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