level1_cplx_impl.h
Go to the documentation of this file.
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 
10 #include "common.h"
11 
14  inline RealScalar operator()(const Scalar &a) const { return Eigen::numext::norm1(a); }
15 };
16 namespace Eigen {
17 namespace internal {
18 template <>
21 };
22 } // namespace internal
23 } // namespace Eigen
24 
25 // computes the sum of magnitudes of all vector elements or, for a complex vector x, the sum
26 // res = |Rex1| + |Imx1| + |Rex2| + |Imx2| + ... + |Rexn| + |Imxn|, where x is a vector of order n
28  // std::cerr << "__asum " << *n << " " << *incx << "\n";
29  Complex *x = reinterpret_cast<Complex *>(px);
30 
31  if (*n <= 0) return 0;
32 
33  if (*incx == 1)
34  return make_vector(x, *n).unaryExpr<scalar_norm1_op>().sum();
35  else
36  return make_vector(x, *n, std::abs(*incx)).unaryExpr<scalar_norm1_op>().sum();
37 }
38 
39 extern "C" int EIGEN_CAT(i, EIGEN_BLAS_FUNC_NAME(amax))(int *n, RealScalar *px, int *incx) {
40  if (*n <= 0) return 0;
41  Scalar *x = reinterpret_cast<Scalar *>(px);
42 
44  if (*incx == 1)
45  make_vector(x, *n).unaryExpr<scalar_norm1_op>().maxCoeff(&ret);
46  else
47  make_vector(x, *n, std::abs(*incx)).unaryExpr<scalar_norm1_op>().maxCoeff(&ret);
48  return int(ret) + 1;
49 }
50 
51 extern "C" int EIGEN_CAT(i, EIGEN_BLAS_FUNC_NAME(amin))(int *n, RealScalar *px, int *incx) {
52  if (*n <= 0) return 0;
53  Scalar *x = reinterpret_cast<Scalar *>(px);
54 
56  if (*incx == 1)
57  make_vector(x, *n).unaryExpr<scalar_norm1_op>().minCoeff(&ret);
58  else
59  make_vector(x, *n, std::abs(*incx)).unaryExpr<scalar_norm1_op>().minCoeff(&ret);
60  return int(ret) + 1;
61 }
62 
63 // computes a dot product of a conjugated vector with another vector.
65  // std::cerr << "_dotc " << *n << " " << *incx << " " << *incy << "\n";
66  Scalar *res = reinterpret_cast<Scalar *>(pres);
67 
68  if (*n <= 0) {
69  *res = Scalar(0);
70  return;
71  }
72 
73  Scalar *x = reinterpret_cast<Scalar *>(px);
74  Scalar *y = reinterpret_cast<Scalar *>(py);
75 
76  if (*incx == 1 && *incy == 1)
77  *res = (make_vector(x, *n).dot(make_vector(y, *n)));
78  else if (*incx > 0 && *incy > 0)
79  *res = (make_vector(x, *n, *incx).dot(make_vector(y, *n, *incy)));
80  else if (*incx < 0 && *incy > 0)
81  *res = (make_vector(x, *n, -*incx).reverse().dot(make_vector(y, *n, *incy)));
82  else if (*incx > 0 && *incy < 0)
83  *res = (make_vector(x, *n, *incx).dot(make_vector(y, *n, -*incy).reverse()));
84  else if (*incx < 0 && *incy < 0)
85  *res = (make_vector(x, *n, -*incx).reverse().dot(make_vector(y, *n, -*incy).reverse()));
86 }
87 
88 // computes a vector-vector dot product without complex conjugation.
90  Scalar *res = reinterpret_cast<Scalar *>(pres);
91 
92  if (*n <= 0) {
93  *res = Scalar(0);
94  return;
95  }
96 
97  Scalar *x = reinterpret_cast<Scalar *>(px);
98  Scalar *y = reinterpret_cast<Scalar *>(py);
99 
100  if (*incx == 1 && *incy == 1)
101  *res = (make_vector(x, *n).cwiseProduct(make_vector(y, *n))).sum();
102  else if (*incx > 0 && *incy > 0)
103  *res = (make_vector(x, *n, *incx).cwiseProduct(make_vector(y, *n, *incy))).sum();
104  else if (*incx < 0 && *incy > 0)
105  *res = (make_vector(x, *n, -*incx).reverse().cwiseProduct(make_vector(y, *n, *incy))).sum();
106  else if (*incx > 0 && *incy < 0)
107  *res = (make_vector(x, *n, *incx).cwiseProduct(make_vector(y, *n, -*incy).reverse())).sum();
108  else if (*incx < 0 && *incy < 0)
109  *res = (make_vector(x, *n, -*incx).reverse().cwiseProduct(make_vector(y, *n, -*incy).reverse())).sum();
110 }
111 
113  // std::cerr << "__nrm2 " << *n << " " << *incx << "\n";
114  if (*n <= 0) return 0;
115 
116  Scalar *x = reinterpret_cast<Scalar *>(px);
117 
118  if (*incx == 1) return make_vector(x, *n).stableNorm();
119 
120  return make_vector(x, *n, *incx).stableNorm();
121 }
122 
124 (int *n, RealScalar *px, int *incx, RealScalar *py, int *incy, RealScalar *pc, RealScalar *ps) {
125  if (*n <= 0) return;
126 
127  Scalar *x = reinterpret_cast<Scalar *>(px);
128  Scalar *y = reinterpret_cast<Scalar *>(py);
131 
134 
137 
138  // TODO implement mixed real-scalar rotations
139  if (*incx < 0 && *incy > 0)
141  else if (*incx > 0 && *incy < 0)
143  else
145 }
146 
148  if (*n <= 0) return;
149 
150  Scalar *x = reinterpret_cast<Scalar *>(px);
152 
153  // std::cerr << "__scal " << *n << " " << alpha << " " << *incx << "\n";
154 
155  if (*incx == 1)
156  make_vector(x, *n) *= alpha;
157  else
158  make_vector(x, *n, std::abs(*incx)) *= alpha;
159 }
AnnoyingScalar abs(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:135
int i
Definition: BiCGSTAB_step_by_step.cpp:9
cout<< "Here is the matrix m:"<< endl<< m<< endl;Matrix< ptrdiff_t, 3, 1 > res
Definition: PartialRedux_count.cpp:3
void reverse(const MatrixType &m)
Definition: array_reverse.cpp:17
SCALAR Scalar
Definition: bench_gemm.cpp:45
NumTraits< Scalar >::Real RealScalar
Definition: bench_gemm.cpp:46
#define REAL_SCALAR_SUFFIX
Definition: blas/complex_double.cpp:13
Rotation given by a cosine-sine pair.
Definition: Jacobi.h:38
A matrix or vector expression mapping an existing array of data.
Definition: Map.h:96
Expression of the reverse of a vector or matrix.
Definition: Reverse.h:65
std::complex< RealScalar > Complex
Definition: common.h:71
#define EIGEN_BLAS_FUNC_NAME(X)
Definition: common.h:150
RealScalar RealScalar * px
Definition: level1_cplx_impl.h:27
Scalar * y
Definition: level1_cplx_impl.h:128
if incx return make_vector(x, *n).unaryExpr< scalar_norm1_op >().sum()
RealScalar * palpha
Definition: level1_cplx_impl.h:147
RealScalar s
Definition: level1_cplx_impl.h:130
return int(ret)+1
int RealScalar int RealScalar int * incy
Definition: level1_cplx_impl.h:124
StridedVectorType vy(make_vector(y, *n, std::abs(*incy)))
Scalar * x
Definition: level1_cplx_impl.h:41
RealScalar c
Definition: level1_cplx_impl.h:129
int RealScalar int RealScalar int RealScalar RealScalar * ps
Definition: level1_cplx_impl.h:124
Eigen::DenseIndex ret
Definition: level1_cplx_impl.h:43
StridedVectorType vx(make_vector(x, *n, std::abs(*incx)))
Eigen::Reverse< StridedVectorType > rvy(vy)
int * n
Definition: level1_cplx_impl.h:124
Eigen::Reverse< StridedVectorType > rvx(vx)
int RealScalar int RealScalar * py
Definition: level1_cplx_impl.h:124
RealScalar alpha
Definition: level1_cplx_impl.h:151
EIGEN_BLAS_FUNC() dotuw(int *n, RealScalar *px, int *incx, RealScalar *py, int *incy, RealScalar *pres)
Definition: level1_cplx_impl.h:89
EIGEN_BLAS_FUNC(EIGEN_CAT(REAL_SCALAR_SUFFIX, scal))(int *n
EIGEN_BLAS_FUNC() dotcw(int *n, RealScalar *px, int *incx, RealScalar *py, int *incy, RealScalar *pres)
Definition: level1_cplx_impl.h:64
int RealScalar int RealScalar int RealScalar * pc
Definition: level1_cplx_impl.h:124
RealScalar RealScalar int * incx
Definition: level1_cplx_impl.h:27
RealScalar EIGEN_CAT(REAL_SCALAR_SUFFIX, EIGEN_BLAS_FUNC_NAME(asum))(int *n
EIGEN_BLAS_FUNC() scal(int *n, RealScalar *palpha, RealScalar *px, int *incx)
Definition: level1_impl.h:105
RealScalar EIGEN_BLAS_FUNC_NAME() asum(int *n, Scalar *px, int *incx)
Definition: level1_real_impl.h:14
EIGEN_BLAS_FUNC() rot(int *n, Scalar *px, int *incx, Scalar *py, int *incy, Scalar *pc, Scalar *ps)
Definition: level1_real_impl.h:88
Scalar EIGEN_BLAS_FUNC_NAME() nrm2(int *n, Scalar *px, int *incx)
Definition: level1_real_impl.h:76
const Scalar * a
Definition: level2_cplx_impl.h:32
EIGEN_DEVICE_FUNC void apply_rotation_in_the_plane(DenseBase< VectorX > &xpr_x, DenseBase< VectorY > &xpr_y, const JacobiRotation< OtherScalar > &j)
Definition: Jacobi.h:400
Namespace containing all symbols from the Eigen library.
Definition: bench_norm.cpp:70
EIGEN_DEFAULT_DENSE_INDEX_TYPE DenseIndex
Definition: Meta.h:75
Definition: Eigen_Colamd.h:49
#define amin(a, b)
Definition: oomph_metis_from_parmetis_3.1.1/macros.h:22
#define amax(a, b)
Definition: oomph_metis_from_parmetis_3.1.1/macros.h:21
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition: NumTraits.h:217
Definition: XprHelper.h:205
@ PacketAccess
Definition: XprHelper.h:206
@ Cost
Definition: XprHelper.h:206
Definition: level1_cplx_impl.h:12
RealScalar result_type
Definition: level1_cplx_impl.h:13
RealScalar operator()(const Scalar &a) const
Definition: level1_cplx_impl.h:14