level1_impl.h 3.8 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144
  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(axpy)(const int *n, const RealScalar *palpha, const RealScalar *px, const int *incx, RealScalar *py, const int *incy)
  11. {
  12. const Scalar* x = reinterpret_cast<const Scalar*>(px);
  13. Scalar* y = reinterpret_cast<Scalar*>(py);
  14. Scalar alpha = *reinterpret_cast<const Scalar*>(palpha);
  15. if(*n<=0) return 0;
  16. if(*incx==1 && *incy==1) make_vector(y,*n) += alpha * make_vector(x,*n);
  17. else if(*incx>0 && *incy>0) make_vector(y,*n,*incy) += alpha * make_vector(x,*n,*incx);
  18. else if(*incx>0 && *incy<0) make_vector(y,*n,-*incy).reverse() += alpha * make_vector(x,*n,*incx);
  19. else if(*incx<0 && *incy>0) make_vector(y,*n,*incy) += alpha * make_vector(x,*n,-*incx).reverse();
  20. else if(*incx<0 && *incy<0) make_vector(y,*n,-*incy).reverse() += alpha * make_vector(x,*n,-*incx).reverse();
  21. return 0;
  22. }
  23. int EIGEN_BLAS_FUNC(copy)(int *n, RealScalar *px, int *incx, RealScalar *py, int *incy)
  24. {
  25. if(*n<=0) return 0;
  26. Scalar* x = reinterpret_cast<Scalar*>(px);
  27. Scalar* y = reinterpret_cast<Scalar*>(py);
  28. // be careful, *incx==0 is allowed !!
  29. if(*incx==1 && *incy==1)
  30. make_vector(y,*n) = make_vector(x,*n);
  31. else
  32. {
  33. if(*incx<0) x = x - (*n-1)*(*incx);
  34. if(*incy<0) y = y - (*n-1)*(*incy);
  35. for(int i=0;i<*n;++i)
  36. {
  37. *y = *x;
  38. x += *incx;
  39. y += *incy;
  40. }
  41. }
  42. return 0;
  43. }
  44. int EIGEN_BLAS_FUNC(rotg)(RealScalar *pa, RealScalar *pb, RealScalar *pc, RealScalar *ps)
  45. {
  46. using std::sqrt;
  47. using std::abs;
  48. Scalar& a = *reinterpret_cast<Scalar*>(pa);
  49. Scalar& b = *reinterpret_cast<Scalar*>(pb);
  50. RealScalar* c = pc;
  51. Scalar* s = reinterpret_cast<Scalar*>(ps);
  52. #if !ISCOMPLEX
  53. Scalar r,z;
  54. Scalar aa = abs(a);
  55. Scalar ab = abs(b);
  56. if((aa+ab)==Scalar(0))
  57. {
  58. *c = 1;
  59. *s = 0;
  60. r = 0;
  61. z = 0;
  62. }
  63. else
  64. {
  65. r = sqrt(a*a + b*b);
  66. Scalar amax = aa>ab ? a : b;
  67. r = amax>0 ? r : -r;
  68. *c = a/r;
  69. *s = b/r;
  70. z = 1;
  71. if (aa > ab) z = *s;
  72. if (ab > aa && *c!=RealScalar(0))
  73. z = Scalar(1)/ *c;
  74. }
  75. *pa = r;
  76. *pb = z;
  77. #else
  78. Scalar alpha;
  79. RealScalar norm,scale;
  80. if(abs(a)==RealScalar(0))
  81. {
  82. *c = RealScalar(0);
  83. *s = Scalar(1);
  84. a = b;
  85. }
  86. else
  87. {
  88. scale = abs(a) + abs(b);
  89. norm = scale*sqrt((numext::abs2(a/scale)) + (numext::abs2(b/scale)));
  90. alpha = a/abs(a);
  91. *c = abs(a)/norm;
  92. *s = alpha*numext::conj(b)/norm;
  93. a = alpha*norm;
  94. }
  95. #endif
  96. // JacobiRotation<Scalar> r;
  97. // r.makeGivens(a,b);
  98. // *c = r.c();
  99. // *s = r.s();
  100. return 0;
  101. }
  102. int EIGEN_BLAS_FUNC(scal)(int *n, RealScalar *palpha, RealScalar *px, int *incx)
  103. {
  104. if(*n<=0) return 0;
  105. Scalar* x = reinterpret_cast<Scalar*>(px);
  106. Scalar alpha = *reinterpret_cast<Scalar*>(palpha);
  107. if(*incx==1) make_vector(x,*n) *= alpha;
  108. else make_vector(x,*n,std::abs(*incx)) *= alpha;
  109. return 0;
  110. }
  111. int EIGEN_BLAS_FUNC(swap)(int *n, RealScalar *px, int *incx, RealScalar *py, int *incy)
  112. {
  113. if(*n<=0) return 0;
  114. Scalar* x = reinterpret_cast<Scalar*>(px);
  115. Scalar* y = reinterpret_cast<Scalar*>(py);
  116. if(*incx==1 && *incy==1) make_vector(y,*n).swap(make_vector(x,*n));
  117. else if(*incx>0 && *incy>0) make_vector(y,*n,*incy).swap(make_vector(x,*n,*incx));
  118. else if(*incx>0 && *incy<0) make_vector(y,*n,-*incy).reverse().swap(make_vector(x,*n,*incx));
  119. else if(*incx<0 && *incy>0) make_vector(y,*n,*incy).swap(make_vector(x,*n,-*incx).reverse());
  120. else if(*incx<0 && *incy<0) make_vector(y,*n,-*incy).reverse().swap(make_vector(x,*n,-*incx).reverse());
  121. return 1;
  122. }