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.

133 lines
4.8 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. struct scalar_norm1_op {
  11. typedef RealScalar result_type;
  12. EIGEN_EMPTY_STRUCT_CTOR(scalar_norm1_op)
  13. inline RealScalar operator() (const Scalar& a) const { return numext::norm1(a); }
  14. };
  15. namespace StormEigen {
  16. namespace internal {
  17. template<> struct functor_traits<scalar_norm1_op >
  18. {
  19. enum { Cost = 3 * NumTraits<Scalar>::AddCost, PacketAccess = 0 };
  20. };
  21. }
  22. }
  23. // computes the sum of magnitudes of all vector elements or, for a complex vector x, the sum
  24. // res = |Rex1| + |Imx1| + |Rex2| + |Imx2| + ... + |Rexn| + |Imxn|, where x is a vector of order n
  25. RealScalar EIGEN_CAT(EIGEN_CAT(REAL_SCALAR_SUFFIX,SCALAR_SUFFIX),asum_)(int *n, RealScalar *px, int *incx)
  26. {
  27. // std::cerr << "__asum " << *n << " " << *incx << "\n";
  28. Complex* x = reinterpret_cast<Complex*>(px);
  29. if(*n<=0) return 0;
  30. if(*incx==1) return make_vector(x,*n).unaryExpr<scalar_norm1_op>().sum();
  31. else return make_vector(x,*n,std::abs(*incx)).unaryExpr<scalar_norm1_op>().sum();
  32. }
  33. // computes a dot product of a conjugated vector with another vector.
  34. int EIGEN_BLAS_FUNC(dotcw)(int *n, RealScalar *px, int *incx, RealScalar *py, int *incy, RealScalar* pres)
  35. {
  36. // std::cerr << "_dotc " << *n << " " << *incx << " " << *incy << "\n";
  37. Scalar* res = reinterpret_cast<Scalar*>(pres);
  38. if(*n<=0)
  39. {
  40. *res = Scalar(0);
  41. return 0;
  42. }
  43. Scalar* x = reinterpret_cast<Scalar*>(px);
  44. Scalar* y = reinterpret_cast<Scalar*>(py);
  45. if(*incx==1 && *incy==1) *res = (make_vector(x,*n).dot(make_vector(y,*n)));
  46. else if(*incx>0 && *incy>0) *res = (make_vector(x,*n,*incx).dot(make_vector(y,*n,*incy)));
  47. else if(*incx<0 && *incy>0) *res = (make_vector(x,*n,-*incx).reverse().dot(make_vector(y,*n,*incy)));
  48. else if(*incx>0 && *incy<0) *res = (make_vector(x,*n,*incx).dot(make_vector(y,*n,-*incy).reverse()));
  49. else if(*incx<0 && *incy<0) *res = (make_vector(x,*n,-*incx).reverse().dot(make_vector(y,*n,-*incy).reverse()));
  50. return 0;
  51. }
  52. // computes a vector-vector dot product without complex conjugation.
  53. int EIGEN_BLAS_FUNC(dotuw)(int *n, RealScalar *px, int *incx, RealScalar *py, int *incy, RealScalar* pres)
  54. {
  55. Scalar* res = reinterpret_cast<Scalar*>(pres);
  56. if(*n<=0)
  57. {
  58. *res = Scalar(0);
  59. return 0;
  60. }
  61. Scalar* x = reinterpret_cast<Scalar*>(px);
  62. Scalar* y = reinterpret_cast<Scalar*>(py);
  63. if(*incx==1 && *incy==1) *res = (make_vector(x,*n).cwiseProduct(make_vector(y,*n))).sum();
  64. else if(*incx>0 && *incy>0) *res = (make_vector(x,*n,*incx).cwiseProduct(make_vector(y,*n,*incy))).sum();
  65. else if(*incx<0 && *incy>0) *res = (make_vector(x,*n,-*incx).reverse().cwiseProduct(make_vector(y,*n,*incy))).sum();
  66. else if(*incx>0 && *incy<0) *res = (make_vector(x,*n,*incx).cwiseProduct(make_vector(y,*n,-*incy).reverse())).sum();
  67. else if(*incx<0 && *incy<0) *res = (make_vector(x,*n,-*incx).reverse().cwiseProduct(make_vector(y,*n,-*incy).reverse())).sum();
  68. return 0;
  69. }
  70. RealScalar EIGEN_CAT(EIGEN_CAT(REAL_SCALAR_SUFFIX,SCALAR_SUFFIX),nrm2_)(int *n, RealScalar *px, int *incx)
  71. {
  72. // std::cerr << "__nrm2 " << *n << " " << *incx << "\n";
  73. if(*n<=0) return 0;
  74. Scalar* x = reinterpret_cast<Scalar*>(px);
  75. if(*incx==1)
  76. return make_vector(x,*n).stableNorm();
  77. return make_vector(x,*n,*incx).stableNorm();
  78. }
  79. int EIGEN_CAT(EIGEN_CAT(SCALAR_SUFFIX,REAL_SCALAR_SUFFIX),rot_)(int *n, RealScalar *px, int *incx, RealScalar *py, int *incy, RealScalar *pc, RealScalar *ps)
  80. {
  81. if(*n<=0) return 0;
  82. Scalar* x = reinterpret_cast<Scalar*>(px);
  83. Scalar* y = reinterpret_cast<Scalar*>(py);
  84. RealScalar c = *pc;
  85. RealScalar s = *ps;
  86. StridedVectorType vx(make_vector(x,*n,std::abs(*incx)));
  87. StridedVectorType vy(make_vector(y,*n,std::abs(*incy)));
  88. Reverse<StridedVectorType> rvx(vx);
  89. Reverse<StridedVectorType> rvy(vy);
  90. // TODO implement mixed real-scalar rotations
  91. if(*incx<0 && *incy>0) internal::apply_rotation_in_the_plane(rvx, vy, JacobiRotation<Scalar>(c,s));
  92. else if(*incx>0 && *incy<0) internal::apply_rotation_in_the_plane(vx, rvy, JacobiRotation<Scalar>(c,s));
  93. else internal::apply_rotation_in_the_plane(vx, vy, JacobiRotation<Scalar>(c,s));
  94. return 0;
  95. }
  96. int EIGEN_CAT(EIGEN_CAT(SCALAR_SUFFIX,REAL_SCALAR_SUFFIX),scal_)(int *n, RealScalar *palpha, RealScalar *px, int *incx)
  97. {
  98. if(*n<=0) return 0;
  99. Scalar* x = reinterpret_cast<Scalar*>(px);
  100. RealScalar alpha = *palpha;
  101. // std::cerr << "__scal " << *n << " " << alpha << " " << *incx << "\n";
  102. if(*incx==1) make_vector(x,*n) *= alpha;
  103. else make_vector(x,*n,std::abs(*incx)) *= alpha;
  104. return 0;
  105. }