level1_real_impl.h 4.2 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122
  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. // computes the sum of magnitudes of all vector elements or, for a complex vector x, the sum
  11. // res = |Rex1| + |Imx1| + |Rex2| + |Imx2| + ... + |Rexn| + |Imxn|, where x is a vector of order n
  12. RealScalar EIGEN_BLAS_FUNC(asum)(int *n, RealScalar *px, int *incx)
  13. {
  14. // std::cerr << "_asum " << *n << " " << *incx << "\n";
  15. Scalar* x = reinterpret_cast<Scalar*>(px);
  16. if(*n<=0) return 0;
  17. if(*incx==1) return make_vector(x,*n).cwiseAbs().sum();
  18. else return make_vector(x,*n,std::abs(*incx)).cwiseAbs().sum();
  19. }
  20. int EIGEN_CAT(i, EIGEN_BLAS_FUNC(amax))(int *n, RealScalar *px, int *incx)
  21. {
  22. if(*n<=0) return 0;
  23. Scalar* x = reinterpret_cast<Scalar*>(px);
  24. DenseIndex ret;
  25. if(*incx==1) make_vector(x,*n).cwiseAbs().maxCoeff(&ret);
  26. else make_vector(x,*n,std::abs(*incx)).cwiseAbs().maxCoeff(&ret);
  27. return int(ret)+1;
  28. }
  29. int EIGEN_CAT(i, EIGEN_BLAS_FUNC(amin))(int *n, RealScalar *px, int *incx)
  30. {
  31. if(*n<=0) return 0;
  32. Scalar* x = reinterpret_cast<Scalar*>(px);
  33. DenseIndex ret;
  34. if(*incx==1) make_vector(x,*n).cwiseAbs().minCoeff(&ret);
  35. else make_vector(x,*n,std::abs(*incx)).cwiseAbs().minCoeff(&ret);
  36. return int(ret)+1;
  37. }
  38. // computes a vector-vector dot product.
  39. Scalar EIGEN_BLAS_FUNC(dot)(int *n, RealScalar *px, int *incx, RealScalar *py, int *incy)
  40. {
  41. // std::cerr << "_dot " << *n << " " << *incx << " " << *incy << "\n";
  42. if(*n<=0) return 0;
  43. Scalar* x = reinterpret_cast<Scalar*>(px);
  44. Scalar* y = reinterpret_cast<Scalar*>(py);
  45. if(*incx==1 && *incy==1) return (make_vector(x,*n).cwiseProduct(make_vector(y,*n))).sum();
  46. else if(*incx>0 && *incy>0) return (make_vector(x,*n,*incx).cwiseProduct(make_vector(y,*n,*incy))).sum();
  47. else if(*incx<0 && *incy>0) return (make_vector(x,*n,-*incx).reverse().cwiseProduct(make_vector(y,*n,*incy))).sum();
  48. else if(*incx>0 && *incy<0) return (make_vector(x,*n,*incx).cwiseProduct(make_vector(y,*n,-*incy).reverse())).sum();
  49. else if(*incx<0 && *incy<0) return (make_vector(x,*n,-*incx).reverse().cwiseProduct(make_vector(y,*n,-*incy).reverse())).sum();
  50. else return 0;
  51. }
  52. // computes the Euclidean norm of a vector.
  53. // FIXME
  54. Scalar EIGEN_BLAS_FUNC(nrm2)(int *n, RealScalar *px, int *incx)
  55. {
  56. // std::cerr << "_nrm2 " << *n << " " << *incx << "\n";
  57. if(*n<=0) return 0;
  58. Scalar* x = reinterpret_cast<Scalar*>(px);
  59. if(*incx==1) return make_vector(x,*n).stableNorm();
  60. else return make_vector(x,*n,std::abs(*incx)).stableNorm();
  61. }
  62. int EIGEN_BLAS_FUNC(rot)(int *n, RealScalar *px, int *incx, RealScalar *py, int *incy, RealScalar *pc, RealScalar *ps)
  63. {
  64. // std::cerr << "_rot " << *n << " " << *incx << " " << *incy << "\n";
  65. if(*n<=0) return 0;
  66. Scalar* x = reinterpret_cast<Scalar*>(px);
  67. Scalar* y = reinterpret_cast<Scalar*>(py);
  68. Scalar c = *reinterpret_cast<Scalar*>(pc);
  69. Scalar s = *reinterpret_cast<Scalar*>(ps);
  70. StridedVectorType vx(make_vector(x,*n,std::abs(*incx)));
  71. StridedVectorType vy(make_vector(y,*n,std::abs(*incy)));
  72. Reverse<StridedVectorType> rvx(vx);
  73. Reverse<StridedVectorType> rvy(vy);
  74. if(*incx<0 && *incy>0) internal::apply_rotation_in_the_plane(rvx, vy, JacobiRotation<Scalar>(c,s));
  75. else if(*incx>0 && *incy<0) internal::apply_rotation_in_the_plane(vx, rvy, JacobiRotation<Scalar>(c,s));
  76. else internal::apply_rotation_in_the_plane(vx, vy, JacobiRotation<Scalar>(c,s));
  77. return 0;
  78. }
  79. /*
  80. // performs rotation of points in the modified plane.
  81. int EIGEN_BLAS_FUNC(rotm)(int *n, RealScalar *px, int *incx, RealScalar *py, int *incy, RealScalar *param)
  82. {
  83. Scalar* x = reinterpret_cast<Scalar*>(px);
  84. Scalar* y = reinterpret_cast<Scalar*>(py);
  85. // TODO
  86. return 0;
  87. }
  88. // computes the modified parameters for a Givens rotation.
  89. int EIGEN_BLAS_FUNC(rotmg)(RealScalar *d1, RealScalar *d2, RealScalar *x1, RealScalar *x2, RealScalar *param)
  90. {
  91. // TODO
  92. return 0;
  93. }
  94. */