File size: 4,956 Bytes
90f0b29 |
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 |
// This file is part of Eigen, a lightweight C++ template library
// for linear algebra.
//
// Copyright (C) 2009-2010 Gael Guennebaud <gael.guennebaud@inria.fr>
//
// This Source Code Form is subject to the terms of the Mozilla
// Public License v. 2.0. If a copy of the MPL was not distributed
// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
#include "common.h"
struct scalar_norm1_op {
typedef RealScalar result_type;
EIGEN_EMPTY_STRUCT_CTOR(scalar_norm1_op)
inline RealScalar operator() (const Scalar& a) const { return numext::norm1(a); }
};
namespace Eigen {
namespace internal {
template<> struct functor_traits<scalar_norm1_op >
{
enum { Cost = 3 * NumTraits<Scalar>::AddCost, PacketAccess = 0 };
};
}
}
// computes the sum of magnitudes of all vector elements or, for a complex vector x, the sum
// res = |Rex1| + |Imx1| + |Rex2| + |Imx2| + ... + |Rexn| + |Imxn|, where x is a vector of order n
RealScalar EIGEN_CAT(EIGEN_CAT(REAL_SCALAR_SUFFIX,SCALAR_SUFFIX),asum_)(int *n, RealScalar *px, int *incx)
{
// std::cerr << "__asum " << *n << " " << *incx << "\n";
Complex* x = reinterpret_cast<Complex*>(px);
if(*n<=0) return 0;
if(*incx==1) return make_vector(x,*n).unaryExpr<scalar_norm1_op>().sum();
else return make_vector(x,*n,std::abs(*incx)).unaryExpr<scalar_norm1_op>().sum();
}
// computes a dot product of a conjugated vector with another vector.
int EIGEN_BLAS_FUNC(dotcw)(int *n, RealScalar *px, int *incx, RealScalar *py, int *incy, RealScalar* pres)
{
// std::cerr << "_dotc " << *n << " " << *incx << " " << *incy << "\n";
Scalar* res = reinterpret_cast<Scalar*>(pres);
if(*n<=0)
{
*res = Scalar(0);
return 0;
}
Scalar* x = reinterpret_cast<Scalar*>(px);
Scalar* y = reinterpret_cast<Scalar*>(py);
if(*incx==1 && *incy==1) *res = (make_vector(x,*n).dot(make_vector(y,*n)));
else if(*incx>0 && *incy>0) *res = (make_vector(x,*n,*incx).dot(make_vector(y,*n,*incy)));
else if(*incx<0 && *incy>0) *res = (make_vector(x,*n,-*incx).reverse().dot(make_vector(y,*n,*incy)));
else if(*incx>0 && *incy<0) *res = (make_vector(x,*n,*incx).dot(make_vector(y,*n,-*incy).reverse()));
else if(*incx<0 && *incy<0) *res = (make_vector(x,*n,-*incx).reverse().dot(make_vector(y,*n,-*incy).reverse()));
return 0;
}
// computes a vector-vector dot product without complex conjugation.
int EIGEN_BLAS_FUNC(dotuw)(int *n, RealScalar *px, int *incx, RealScalar *py, int *incy, RealScalar* pres)
{
Scalar* res = reinterpret_cast<Scalar*>(pres);
if(*n<=0)
{
*res = Scalar(0);
return 0;
}
Scalar* x = reinterpret_cast<Scalar*>(px);
Scalar* y = reinterpret_cast<Scalar*>(py);
if(*incx==1 && *incy==1) *res = (make_vector(x,*n).cwiseProduct(make_vector(y,*n))).sum();
else if(*incx>0 && *incy>0) *res = (make_vector(x,*n,*incx).cwiseProduct(make_vector(y,*n,*incy))).sum();
else if(*incx<0 && *incy>0) *res = (make_vector(x,*n,-*incx).reverse().cwiseProduct(make_vector(y,*n,*incy))).sum();
else if(*incx>0 && *incy<0) *res = (make_vector(x,*n,*incx).cwiseProduct(make_vector(y,*n,-*incy).reverse())).sum();
else if(*incx<0 && *incy<0) *res = (make_vector(x,*n,-*incx).reverse().cwiseProduct(make_vector(y,*n,-*incy).reverse())).sum();
return 0;
}
RealScalar EIGEN_CAT(EIGEN_CAT(REAL_SCALAR_SUFFIX,SCALAR_SUFFIX),nrm2_)(int *n, RealScalar *px, int *incx)
{
// std::cerr << "__nrm2 " << *n << " " << *incx << "\n";
if(*n<=0) return 0;
Scalar* x = reinterpret_cast<Scalar*>(px);
if(*incx==1)
return make_vector(x,*n).stableNorm();
return make_vector(x,*n,*incx).stableNorm();
}
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)
{
if(*n<=0) return 0;
Scalar* x = reinterpret_cast<Scalar*>(px);
Scalar* y = reinterpret_cast<Scalar*>(py);
RealScalar c = *pc;
RealScalar s = *ps;
StridedVectorType vx(make_vector(x,*n,std::abs(*incx)));
StridedVectorType vy(make_vector(y,*n,std::abs(*incy)));
Reverse<StridedVectorType> rvx(vx);
Reverse<StridedVectorType> rvy(vy);
// TODO implement mixed real-scalar rotations
if(*incx<0 && *incy>0) internal::apply_rotation_in_the_plane(rvx, vy, JacobiRotation<Scalar>(c,s));
else if(*incx>0 && *incy<0) internal::apply_rotation_in_the_plane(vx, rvy, JacobiRotation<Scalar>(c,s));
else internal::apply_rotation_in_the_plane(vx, vy, JacobiRotation<Scalar>(c,s));
return 0;
}
int EIGEN_CAT(EIGEN_CAT(SCALAR_SUFFIX,REAL_SCALAR_SUFFIX),scal_)(int *n, RealScalar *palpha, RealScalar *px, int *incx)
{
if(*n<=0) return 0;
Scalar* x = reinterpret_cast<Scalar*>(px);
RealScalar alpha = *palpha;
// std::cerr << "__scal " << *n << " " << alpha << " " << *incx << "\n";
if(*incx==1) make_vector(x,*n) *= alpha;
else make_vector(x,*n,std::abs(*incx)) *= alpha;
return 0;
}
|