level2_real_impl.h 10 KB

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