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.

368 lines
11 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. // y = alpha*A*x + beta*y
  11. int STORMEIGEN_BLAS_FUNC(symv) (char *uplo, int *n, RealScalar *palpha, RealScalar *pa, int *lda, RealScalar *px, int *incx, RealScalar *pbeta, RealScalar *py, int *incy)
  12. {
  13. typedef void (*functype)(int, const Scalar*, int, const Scalar*, Scalar*, Scalar);
  14. static functype func[2];
  15. static bool init = false;
  16. if(!init)
  17. {
  18. for(int k=0; k<2; ++k)
  19. func[k] = 0;
  20. func[UP] = (internal::selfadjoint_matrix_vector_product<Scalar,int,ColMajor,Upper,false,false>::run);
  21. func[LO] = (internal::selfadjoint_matrix_vector_product<Scalar,int,ColMajor,Lower,false,false>::run);
  22. init = true;
  23. }
  24. Scalar* a = reinterpret_cast<Scalar*>(pa);
  25. Scalar* x = reinterpret_cast<Scalar*>(px);
  26. Scalar* y = reinterpret_cast<Scalar*>(py);
  27. Scalar alpha = *reinterpret_cast<Scalar*>(palpha);
  28. Scalar beta = *reinterpret_cast<Scalar*>(pbeta);
  29. // check arguments
  30. int info = 0;
  31. if(UPLO(*uplo)==INVALID) info = 1;
  32. else if(*n<0) info = 2;
  33. else if(*lda<std::max(1,*n)) info = 5;
  34. else if(*incx==0) info = 7;
  35. else if(*incy==0) info = 10;
  36. if(info)
  37. return xerbla_(SCALAR_SUFFIX_UP"SYMV ",&info,6);
  38. if(*n==0)
  39. return 0;
  40. Scalar* actual_x = get_compact_vector(x,*n,*incx);
  41. Scalar* actual_y = get_compact_vector(y,*n,*incy);
  42. if(beta!=Scalar(1))
  43. {
  44. if(beta==Scalar(0)) make_vector(actual_y, *n).setZero();
  45. else make_vector(actual_y, *n) *= beta;
  46. }
  47. int code = UPLO(*uplo);
  48. if(code>=2 || func[code]==0)
  49. return 0;
  50. func[code](*n, a, *lda, actual_x, actual_y, alpha);
  51. if(actual_x!=x) delete[] actual_x;
  52. if(actual_y!=y) delete[] copy_back(actual_y,y,*n,*incy);
  53. return 1;
  54. }
  55. // C := alpha*x*x' + C
  56. int STORMEIGEN_BLAS_FUNC(syr)(char *uplo, int *n, RealScalar *palpha, RealScalar *px, int *incx, RealScalar *pc, int *ldc)
  57. {
  58. // typedef void (*functype)(int, const Scalar *, int, Scalar *, int, Scalar);
  59. // static functype func[2];
  60. // static bool init = false;
  61. // if(!init)
  62. // {
  63. // for(int k=0; k<2; ++k)
  64. // func[k] = 0;
  65. //
  66. // func[UP] = (internal::selfadjoint_product<Scalar,ColMajor,ColMajor,false,UpperTriangular>::run);
  67. // func[LO] = (internal::selfadjoint_product<Scalar,ColMajor,ColMajor,false,LowerTriangular>::run);
  68. // init = true;
  69. // }
  70. typedef void (*functype)(int, Scalar*, int, const Scalar*, const Scalar*, const Scalar&);
  71. static functype func[2];
  72. static bool init = false;
  73. if(!init)
  74. {
  75. for(int k=0; k<2; ++k)
  76. func[k] = 0;
  77. func[UP] = (selfadjoint_rank1_update<Scalar,int,ColMajor,Upper,false,Conj>::run);
  78. func[LO] = (selfadjoint_rank1_update<Scalar,int,ColMajor,Lower,false,Conj>::run);
  79. init = true;
  80. }
  81. Scalar* x = reinterpret_cast<Scalar*>(px);
  82. Scalar* c = reinterpret_cast<Scalar*>(pc);
  83. Scalar alpha = *reinterpret_cast<Scalar*>(palpha);
  84. int info = 0;
  85. if(UPLO(*uplo)==INVALID) info = 1;
  86. else if(*n<0) info = 2;
  87. else if(*incx==0) info = 5;
  88. else if(*ldc<std::max(1,*n)) info = 7;
  89. if(info)
  90. return xerbla_(SCALAR_SUFFIX_UP"SYR ",&info,6);
  91. if(*n==0 || alpha==Scalar(0)) return 1;
  92. // if the increment is not 1, let's copy it to a temporary vector to enable vectorization
  93. Scalar* x_cpy = get_compact_vector(x,*n,*incx);
  94. int code = UPLO(*uplo);
  95. if(code>=2 || func[code]==0)
  96. return 0;
  97. func[code](*n, c, *ldc, x_cpy, x_cpy, alpha);
  98. if(x_cpy!=x) delete[] x_cpy;
  99. return 1;
  100. }
  101. // C := alpha*x*y' + alpha*y*x' + C
  102. int STORMEIGEN_BLAS_FUNC(syr2)(char *uplo, int *n, RealScalar *palpha, RealScalar *px, int *incx, RealScalar *py, int *incy, RealScalar *pc, int *ldc)
  103. {
  104. // typedef void (*functype)(int, const Scalar *, int, const Scalar *, int, Scalar *, int, Scalar);
  105. // static functype func[2];
  106. //
  107. // static bool init = false;
  108. // if(!init)
  109. // {
  110. // for(int k=0; k<2; ++k)
  111. // func[k] = 0;
  112. //
  113. // func[UP] = (internal::selfadjoint_product<Scalar,ColMajor,ColMajor,false,UpperTriangular>::run);
  114. // func[LO] = (internal::selfadjoint_product<Scalar,ColMajor,ColMajor,false,LowerTriangular>::run);
  115. //
  116. // init = true;
  117. // }
  118. typedef void (*functype)(int, Scalar*, int, const Scalar*, const Scalar*, Scalar);
  119. static functype func[2];
  120. static bool init = false;
  121. if(!init)
  122. {
  123. for(int k=0; k<2; ++k)
  124. func[k] = 0;
  125. func[UP] = (internal::rank2_update_selector<Scalar,int,Upper>::run);
  126. func[LO] = (internal::rank2_update_selector<Scalar,int,Lower>::run);
  127. init = true;
  128. }
  129. Scalar* x = reinterpret_cast<Scalar*>(px);
  130. Scalar* y = reinterpret_cast<Scalar*>(py);
  131. Scalar* c = reinterpret_cast<Scalar*>(pc);
  132. Scalar alpha = *reinterpret_cast<Scalar*>(palpha);
  133. int info = 0;
  134. if(UPLO(*uplo)==INVALID) info = 1;
  135. else if(*n<0) info = 2;
  136. else if(*incx==0) info = 5;
  137. else if(*incy==0) info = 7;
  138. else if(*ldc<std::max(1,*n)) info = 9;
  139. if(info)
  140. return xerbla_(SCALAR_SUFFIX_UP"SYR2 ",&info,6);
  141. if(alpha==Scalar(0))
  142. return 1;
  143. Scalar* x_cpy = get_compact_vector(x,*n,*incx);
  144. Scalar* y_cpy = get_compact_vector(y,*n,*incy);
  145. int code = UPLO(*uplo);
  146. if(code>=2 || func[code]==0)
  147. return 0;
  148. func[code](*n, c, *ldc, x_cpy, y_cpy, alpha);
  149. if(x_cpy!=x) delete[] x_cpy;
  150. if(y_cpy!=y) delete[] y_cpy;
  151. // int code = UPLO(*uplo);
  152. // if(code>=2 || func[code]==0)
  153. // return 0;
  154. // func[code](*n, a, *inca, b, *incb, c, *ldc, alpha);
  155. return 1;
  156. }
  157. /** DSBMV performs the matrix-vector operation
  158. *
  159. * y := alpha*A*x + beta*y,
  160. *
  161. * where alpha and beta are scalars, x and y are n element vectors and
  162. * A is an n by n symmetric band matrix, with k super-diagonals.
  163. */
  164. // int STORMEIGEN_BLAS_FUNC(sbmv)( char *uplo, int *n, int *k, RealScalar *alpha, RealScalar *a, int *lda,
  165. // RealScalar *x, int *incx, RealScalar *beta, RealScalar *y, int *incy)
  166. // {
  167. // return 1;
  168. // }
  169. /** DSPMV performs the matrix-vector operation
  170. *
  171. * y := alpha*A*x + beta*y,
  172. *
  173. * where alpha and beta are scalars, x and y are n element vectors and
  174. * A is an n by n symmetric matrix, supplied in packed form.
  175. *
  176. */
  177. // int STORMEIGEN_BLAS_FUNC(spmv)(char *uplo, int *n, RealScalar *alpha, RealScalar *ap, RealScalar *x, int *incx, RealScalar *beta, RealScalar *y, int *incy)
  178. // {
  179. // return 1;
  180. // }
  181. /** DSPR performs the symmetric rank 1 operation
  182. *
  183. * A := alpha*x*x' + A,
  184. *
  185. * where alpha is a real scalar, x is an n element vector and A is an
  186. * n by n symmetric matrix, supplied in packed form.
  187. */
  188. int STORMEIGEN_BLAS_FUNC(spr)(char *uplo, int *n, Scalar *palpha, Scalar *px, int *incx, Scalar *pap)
  189. {
  190. typedef void (*functype)(int, Scalar*, const Scalar*, Scalar);
  191. static functype func[2];
  192. static bool init = false;
  193. if(!init)
  194. {
  195. for(int k=0; k<2; ++k)
  196. func[k] = 0;
  197. func[UP] = (internal::selfadjoint_packed_rank1_update<Scalar,int,ColMajor,Upper,false,false>::run);
  198. func[LO] = (internal::selfadjoint_packed_rank1_update<Scalar,int,ColMajor,Lower,false,false>::run);
  199. init = true;
  200. }
  201. Scalar* x = reinterpret_cast<Scalar*>(px);
  202. Scalar* ap = reinterpret_cast<Scalar*>(pap);
  203. Scalar alpha = *reinterpret_cast<Scalar*>(palpha);
  204. int info = 0;
  205. if(UPLO(*uplo)==INVALID) info = 1;
  206. else if(*n<0) info = 2;
  207. else if(*incx==0) info = 5;
  208. if(info)
  209. return xerbla_(SCALAR_SUFFIX_UP"SPR ",&info,6);
  210. if(alpha==Scalar(0))
  211. return 1;
  212. Scalar* x_cpy = get_compact_vector(x, *n, *incx);
  213. int code = UPLO(*uplo);
  214. if(code>=2 || func[code]==0)
  215. return 0;
  216. func[code](*n, ap, x_cpy, alpha);
  217. if(x_cpy!=x) delete[] x_cpy;
  218. return 1;
  219. }
  220. /** DSPR2 performs the symmetric rank 2 operation
  221. *
  222. * A := alpha*x*y' + alpha*y*x' + A,
  223. *
  224. * where alpha is a scalar, x and y are n element vectors and A is an
  225. * n by n symmetric matrix, supplied in packed form.
  226. */
  227. int STORMEIGEN_BLAS_FUNC(spr2)(char *uplo, int *n, RealScalar *palpha, RealScalar *px, int *incx, RealScalar *py, int *incy, RealScalar *pap)
  228. {
  229. typedef void (*functype)(int, Scalar*, const Scalar*, const Scalar*, Scalar);
  230. static functype func[2];
  231. static bool init = false;
  232. if(!init)
  233. {
  234. for(int k=0; k<2; ++k)
  235. func[k] = 0;
  236. func[UP] = (internal::packed_rank2_update_selector<Scalar,int,Upper>::run);
  237. func[LO] = (internal::packed_rank2_update_selector<Scalar,int,Lower>::run);
  238. init = true;
  239. }
  240. Scalar* x = reinterpret_cast<Scalar*>(px);
  241. Scalar* y = reinterpret_cast<Scalar*>(py);
  242. Scalar* ap = reinterpret_cast<Scalar*>(pap);
  243. Scalar alpha = *reinterpret_cast<Scalar*>(palpha);
  244. int info = 0;
  245. if(UPLO(*uplo)==INVALID) info = 1;
  246. else if(*n<0) info = 2;
  247. else if(*incx==0) info = 5;
  248. else if(*incy==0) info = 7;
  249. if(info)
  250. return xerbla_(SCALAR_SUFFIX_UP"SPR2 ",&info,6);
  251. if(alpha==Scalar(0))
  252. return 1;
  253. Scalar* x_cpy = get_compact_vector(x, *n, *incx);
  254. Scalar* y_cpy = get_compact_vector(y, *n, *incy);
  255. int code = UPLO(*uplo);
  256. if(code>=2 || func[code]==0)
  257. return 0;
  258. func[code](*n, ap, x_cpy, y_cpy, alpha);
  259. if(x_cpy!=x) delete[] x_cpy;
  260. if(y_cpy!=y) delete[] y_cpy;
  261. return 1;
  262. }
  263. /** DGER performs the rank 1 operation
  264. *
  265. * A := alpha*x*y' + A,
  266. *
  267. * where alpha is a scalar, x is an m element vector, y is an n element
  268. * vector and A is an m by n matrix.
  269. */
  270. int STORMEIGEN_BLAS_FUNC(ger)(int *m, int *n, Scalar *palpha, Scalar *px, int *incx, Scalar *py, int *incy, Scalar *pa, int *lda)
  271. {
  272. Scalar* x = reinterpret_cast<Scalar*>(px);
  273. Scalar* y = reinterpret_cast<Scalar*>(py);
  274. Scalar* a = reinterpret_cast<Scalar*>(pa);
  275. Scalar alpha = *reinterpret_cast<Scalar*>(palpha);
  276. int info = 0;
  277. if(*m<0) info = 1;
  278. else if(*n<0) info = 2;
  279. else if(*incx==0) info = 5;
  280. else if(*incy==0) info = 7;
  281. else if(*lda<std::max(1,*m)) info = 9;
  282. if(info)
  283. return xerbla_(SCALAR_SUFFIX_UP"GER ",&info,6);
  284. if(alpha==Scalar(0))
  285. return 1;
  286. Scalar* x_cpy = get_compact_vector(x,*m,*incx);
  287. Scalar* y_cpy = get_compact_vector(y,*n,*incy);
  288. internal::general_rank1_update<Scalar,int,ColMajor,false,false>::run(*m, *n, a, *lda, x_cpy, y_cpy, alpha);
  289. if(x_cpy!=x) delete[] x_cpy;
  290. if(y_cpy!=y) delete[] y_cpy;
  291. return 1;
  292. }