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
15 {
16 // std::cerr << "_asum " << *n << " " << *incx << "\n";
17 
18  Scalar* x = reinterpret_cast<Scalar*>(px);
19 
20  if(*n<=0) return 0;
21 
22  if(*incx==1) return make_vector(x,*n).cwiseAbs().sum();
23  else return make_vector(x,*n,std::abs(*incx)).cwiseAbs().sum();
24 }
25 
26 // computes a vector-vector dot product.
28 {
29 // std::cerr << "_dot " << *n << " " << *incx << " " << *incy << "\n";
30 
31  if(*n<=0) return 0;
32 
33  Scalar* x = reinterpret_cast<Scalar*>(px);
34  Scalar* y = reinterpret_cast<Scalar*>(py);
35 
36  if(*incx==1 && *incy==1) return (make_vector(x,*n).cwiseProduct(make_vector(y,*n))).sum();
37  else if(*incx>0 && *incy>0) return (make_vector(x,*n,*incx).cwiseProduct(make_vector(y,*n,*incy))).sum();
38  else if(*incx<0 && *incy>0) return (make_vector(x,*n,-*incx).reverse().cwiseProduct(make_vector(y,*n,*incy))).sum();
39  else if(*incx>0 && *incy<0) return (make_vector(x,*n,*incx).cwiseProduct(make_vector(y,*n,-*incy).reverse())).sum();
40  else if(*incx<0 && *incy<0) return (make_vector(x,*n,-*incx).reverse().cwiseProduct(make_vector(y,*n,-*incy).reverse())).sum();
41  else return 0;
42 }
43 
44 // computes the Euclidean norm of a vector.
45 // FIXME
47 {
48 // std::cerr << "_nrm2 " << *n << " " << *incx << "\n";
49  if(*n<=0) return 0;
50 
51  Scalar* x = reinterpret_cast<Scalar*>(px);
52 
53  if(*incx==1) return make_vector(x,*n).stableNorm();
54  else return make_vector(x,*n,std::abs(*incx)).stableNorm();
55 }
56 
58 {
59 // std::cerr << "_rot " << *n << " " << *incx << " " << *incy << "\n";
60  if(*n<=0) return 0;
61 
62  Scalar* x = reinterpret_cast<Scalar*>(px);
63  Scalar* y = reinterpret_cast<Scalar*>(py);
64  Scalar c = *reinterpret_cast<Scalar*>(pc);
65  Scalar s = *reinterpret_cast<Scalar*>(ps);
66 
69 
72 
73  if(*incx<0 && *incy>0) internal::apply_rotation_in_the_plane(rvx, vy, JacobiRotation<Scalar>(c,s));
76 
77 
78  return 0;
79 }
80 
81 /*
82 // performs rotation of points in the modified plane.
83 int EIGEN_BLAS_FUNC(rotm)(int *n, RealScalar *px, int *incx, RealScalar *py, int *incy, RealScalar *param)
84 {
85  Scalar* x = reinterpret_cast<Scalar*>(px);
86  Scalar* y = reinterpret_cast<Scalar*>(py);
87 
88  // TODO
89 
90  return 0;
91 }
92 
93 // computes the modified parameters for a Givens rotation.
94 int EIGEN_BLAS_FUNC(rotmg)(RealScalar *d1, RealScalar *d2, RealScalar *x1, RealScalar *x2, RealScalar *param)
95 {
96  // TODO
97 
98  return 0;
99 }
100 */
int EIGEN_BLAS_FUNC() rot(int *n, RealScalar *px, int *incx, RealScalar *py, int *incy, RealScalar *pc, RealScalar *ps)
SCALAR Scalar
Definition: bench_gemm.cpp:33
Reverse< StridedVectorType > rvy(vy)
Scalar * y
RealScalar EIGEN_BLAS_FUNC() asum(int *n, RealScalar *px, int *incx)
A matrix or vector expression mapping an existing array of data.
Definition: Map.h:94
int RealScalar int RealScalar int RealScalar * pc
RealScalar RealScalar int * incx
int n
Scalar Scalar * c
Definition: benchVecAdd.cpp:17
Rotation given by a cosine-sine pair.
#define EIGEN_BLAS_FUNC(X)
int RealScalar int RealScalar int RealScalar RealScalar * ps
StridedVectorType vy(make_vector(y,*n, std::abs(*incy)))
Scalar EIGEN_BLAS_FUNC() dot(int *n, RealScalar *px, int *incx, RealScalar *py, int *incy)
int RealScalar int RealScalar * py
RealScalar RealScalar * px
RealScalar s
const mpreal sum(const mpreal tab[], const unsigned long int n, int &status, mp_rnd_t mode=mpreal::get_default_rnd())
Definition: mpreal.h:2381
Scalar EIGEN_BLAS_FUNC() nrm2(int *n, RealScalar *px, int *incx)
NumTraits< Scalar >::Real RealScalar
Definition: bench_gemm.cpp:34
Reverse< StridedVectorType > rvx(vx)
void reverse(const MatrixType &m)
Map< Matrix< T, Dynamic, 1 >, 0, InnerStride< Dynamic > > make_vector(T *data, int size, int incr)
void apply_rotation_in_the_plane(DenseBase< VectorX > &xpr_x, DenseBase< VectorY > &xpr_y, const JacobiRotation< OtherScalar > &j)
Definition: Jacobi.h:432
StridedVectorType vx(make_vector(x,*n, std::abs(*incx)))
Expression of the reverse of a vector or matrix.
Definition: Reverse.h:63
set noclip points set clip one set noclip two set bar set border lt lw set xdata set ydata set zdata set x2data set y2data set boxwidth set dummy x
#define abs(x)
Definition: datatypes.h:17
int RealScalar int RealScalar int * incy


gtsam
Author(s):
autogenerated on Sat May 8 2021 02:42:30