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.

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