level1_real_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 
12 // computes the sum of magnitudes of all vector elements or, for a complex vector x, the sum
13 // res = |Rex1| + |Imx1| + |Rex2| + |Imx2| + ... + |Rexn| + |Imxn|, where x is a vector of order n
14 extern "C" RealScalar EIGEN_BLAS_FUNC_NAME(asum)(int *n, Scalar *px, int *incx) {
15  // std::cerr << "_asum " << *n << " " << *incx << "\n";
16 
17  Scalar *x = reinterpret_cast<Scalar *>(px);
18 
19  if (*n <= 0) return 0;
20 
21  if (*incx == 1)
22  return make_vector(x, *n).cwiseAbs().sum();
23  else
24  return make_vector(x, *n, std::abs(*incx)).cwiseAbs().sum();
25 }
26 
27 extern "C" int EIGEN_CAT(i, EIGEN_BLAS_FUNC_NAME(amax))(int *n, Scalar *px, int *incx) {
28  if (*n <= 0) return 0;
29  Scalar *x = reinterpret_cast<Scalar *>(px);
30 
32  if (*incx == 1)
33  make_vector(x, *n).cwiseAbs().maxCoeff(&ret);
34  else
35  make_vector(x, *n, std::abs(*incx)).cwiseAbs().maxCoeff(&ret);
36  return int(ret) + 1;
37 }
38 
39 extern "C" int EIGEN_CAT(i, EIGEN_BLAS_FUNC_NAME(amin))(int *n, Scalar *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).cwiseAbs().minCoeff(&ret);
46  else
47  make_vector(x, *n, std::abs(*incx)).cwiseAbs().minCoeff(&ret);
48  return int(ret) + 1;
49 }
50 
51 // computes a vector-vector dot product.
52 extern "C" Scalar EIGEN_BLAS_FUNC_NAME(dot)(int *n, Scalar *px, int *incx, Scalar *py, int *incy) {
53  // std::cerr << "_dot " << *n << " " << *incx << " " << *incy << "\n";
54 
55  if (*n <= 0) return 0;
56 
57  Scalar *x = reinterpret_cast<Scalar *>(px);
58  Scalar *y = reinterpret_cast<Scalar *>(py);
59 
60  if (*incx == 1 && *incy == 1)
61  return (make_vector(x, *n).cwiseProduct(make_vector(y, *n))).sum();
62  else if (*incx > 0 && *incy > 0)
63  return (make_vector(x, *n, *incx).cwiseProduct(make_vector(y, *n, *incy))).sum();
64  else if (*incx < 0 && *incy > 0)
65  return (make_vector(x, *n, -*incx).reverse().cwiseProduct(make_vector(y, *n, *incy))).sum();
66  else if (*incx > 0 && *incy < 0)
67  return (make_vector(x, *n, *incx).cwiseProduct(make_vector(y, *n, -*incy).reverse())).sum();
68  else if (*incx < 0 && *incy < 0)
69  return (make_vector(x, *n, -*incx).reverse().cwiseProduct(make_vector(y, *n, -*incy).reverse())).sum();
70  else
71  return 0;
72 }
73 
74 // computes the Euclidean norm of a vector.
75 // FIXME
76 extern "C" Scalar EIGEN_BLAS_FUNC_NAME(nrm2)(int *n, Scalar *px, int *incx) {
77  // std::cerr << "_nrm2 " << *n << " " << *incx << "\n";
78  if (*n <= 0) return 0;
79 
80  Scalar *x = reinterpret_cast<Scalar *>(px);
81 
82  if (*incx == 1)
83  return make_vector(x, *n).stableNorm();
84  else
85  return make_vector(x, *n, std::abs(*incx)).stableNorm();
86 }
87 
88 EIGEN_BLAS_FUNC(rot)(int *n, Scalar *px, int *incx, Scalar *py, int *incy, Scalar *pc, Scalar *ps) {
89  // std::cerr << "_rot " << *n << " " << *incx << " " << *incy << "\n";
90  if (*n <= 0) return;
91 
92  Scalar *x = reinterpret_cast<Scalar *>(px);
93  Scalar *y = reinterpret_cast<Scalar *>(py);
94  Scalar c = *reinterpret_cast<Scalar *>(pc);
95  Scalar s = *reinterpret_cast<Scalar *>(ps);
96 
99 
102 
103  if (*incx < 0 && *incy > 0)
105  else if (*incx > 0 && *incy < 0)
107  else
109 }
110 
111 /*
112 // performs rotation of points in the modified plane.
113 EIGEN_BLAS_FUNC(rotm)(int *n, Scalar *px, int *incx, Scalar *py, int *incy, Scalar *param)
114 {
115  Scalar* x = reinterpret_cast<Scalar*>(px);
116  Scalar* y = reinterpret_cast<Scalar*>(py);
117 
118  // TODO
119 
120  return 0;
121 }
122 
123 // computes the modified parameters for a Givens rotation.
124 EIGEN_BLAS_FUNC(rotmg)(Scalar *d1, Scalar *d2, Scalar *x1, Scalar *x2, Scalar *param)
125 {
126  // TODO
127 
128  return 0;
129 }
130 */
AnnoyingScalar abs(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:135
int i
Definition: BiCGSTAB_step_by_step.cpp:9
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
SCALAR Scalar
Definition: bench_gemm.cpp:45
NumTraits< Scalar >::Real RealScalar
Definition: bench_gemm.cpp:46
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
#define EIGEN_BLAS_FUNC_NAME(X)
Definition: common.h:150
Scalar * y
Definition: level1_cplx_impl.h:128
RealScalar s
Definition: level1_cplx_impl.h:130
int RealScalar int RealScalar int * incy
Definition: level1_cplx_impl.h:124
StridedVectorType vy(make_vector(y, *n, std::abs(*incy)))
int RealScalar int RealScalar int RealScalar RealScalar * ps
Definition: level1_cplx_impl.h:124
StridedVectorType vx(make_vector(x, *n, std::abs(*incx)))
Eigen::Reverse< StridedVectorType > rvy(vy)
Eigen::Reverse< StridedVectorType > rvx(vx)
int RealScalar int RealScalar * py
Definition: level1_cplx_impl.h:124
EIGEN_BLAS_FUNC(EIGEN_CAT(REAL_SCALAR_SUFFIX, scal))(int *n
int RealScalar int RealScalar int RealScalar * pc
Definition: level1_cplx_impl.h:124
int EIGEN_CAT(i, EIGEN_BLAS_FUNC_NAME(amax))(int *n
RealScalar EIGEN_BLAS_FUNC_NAME() asum(int *n, Scalar *px, int *incx)
Definition: level1_real_impl.h:14
Scalar * x
Definition: level1_real_impl.h:29
Scalar EIGEN_BLAS_FUNC_NAME() dot(int *n, Scalar *px, int *incx, Scalar *py, int *incy)
Definition: level1_real_impl.h:52
return int(ret)+1
Eigen::DenseIndex ret
Definition: level1_real_impl.h:31
EIGEN_BLAS_FUNC() rot(int *n, Scalar *px, int *incx, Scalar *py, int *incy, Scalar *pc, Scalar *ps)
Definition: level1_real_impl.h:88
int Scalar int * incx
Definition: level1_real_impl.h:27
Scalar EIGEN_BLAS_FUNC_NAME() nrm2(int *n, Scalar *px, int *incx)
Definition: level1_real_impl.h:76
if incx make_vector(x, *n).cwiseAbs().maxCoeff(&ret)
int Scalar * px
Definition: level1_real_impl.h:27
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
EIGEN_DEFAULT_DENSE_INDEX_TYPE DenseIndex
Definition: Meta.h:75
int c
Definition: calibrate.py:100
#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