SelfadjointMatrixVector.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) 2008-2009 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 #ifndef EIGEN_SELFADJOINT_MATRIX_VECTOR_H
11 #define EIGEN_SELFADJOINT_MATRIX_VECTOR_H
12 
13 #include "../InternalHeaderCheck.h"
14 
15 namespace Eigen {
16 
17 namespace internal {
18 
19 /* Optimized selfadjoint matrix * vector product:
20  * This algorithm processes 2 columns at once that allows to both reduce
21  * the number of load/stores of the result by a factor 2 and to reduce
22  * the instruction dependency.
23  */
24 
25 template<typename Scalar, typename Index, int StorageOrder, int UpLo, bool ConjugateLhs, bool ConjugateRhs, int Version=Specialized>
26 struct selfadjoint_matrix_vector_product;
27 
28 template<typename Scalar, typename Index, int StorageOrder, int UpLo, bool ConjugateLhs, bool ConjugateRhs, int Version>
29 struct selfadjoint_matrix_vector_product
30 
31 {
33 void run(
34  Index size,
35  const Scalar* lhs, Index lhsStride,
36  const Scalar* rhs,
37  Scalar* res,
38  Scalar alpha);
39 };
40 
41 template<typename Scalar, typename Index, int StorageOrder, int UpLo, bool ConjugateLhs, bool ConjugateRhs, int Version>
43 void selfadjoint_matrix_vector_product<Scalar,Index,StorageOrder,UpLo,ConjugateLhs,ConjugateRhs,Version>::run(
44  Index size,
45  const Scalar* lhs, Index lhsStride,
46  const Scalar* rhs,
47  Scalar* res,
48  Scalar alpha)
49 {
50  typedef typename packet_traits<Scalar>::type Packet;
51  typedef typename NumTraits<Scalar>::Real RealScalar;
52  const Index PacketSize = sizeof(Packet)/sizeof(Scalar);
53 
54  enum {
55  IsRowMajor = StorageOrder==RowMajor ? 1 : 0,
56  IsLower = UpLo == Lower ? 1 : 0,
57  FirstTriangular = IsRowMajor == IsLower
58  };
59 
60  conj_helper<Scalar,Scalar,NumTraits<Scalar>::IsComplex && logical_xor(ConjugateLhs, IsRowMajor), ConjugateRhs> cj0;
61  conj_helper<Scalar,Scalar,NumTraits<Scalar>::IsComplex && logical_xor(ConjugateLhs, !IsRowMajor), ConjugateRhs> cj1;
62  conj_helper<RealScalar,Scalar,false, ConjugateRhs> cjd;
63 
64  conj_helper<Packet,Packet,NumTraits<Scalar>::IsComplex && logical_xor(ConjugateLhs, IsRowMajor), ConjugateRhs> pcj0;
65  conj_helper<Packet,Packet,NumTraits<Scalar>::IsComplex && logical_xor(ConjugateLhs, !IsRowMajor), ConjugateRhs> pcj1;
66 
67  Scalar cjAlpha = ConjugateRhs ? numext::conj(alpha) : alpha;
68 
69  Index bound = numext::maxi(Index(0), size-8) & 0xfffffffe;
70  if (FirstTriangular)
71  bound = size - bound;
72 
73  for (Index j=FirstTriangular ? bound : 0;
74  j<(FirstTriangular ? size : bound);j+=2)
75  {
76  const Scalar* EIGEN_RESTRICT A0 = lhs + j*lhsStride;
77  const Scalar* EIGEN_RESTRICT A1 = lhs + (j+1)*lhsStride;
78 
79  Scalar t0 = cjAlpha * rhs[j];
80  Packet ptmp0 = pset1<Packet>(t0);
81  Scalar t1 = cjAlpha * rhs[j+1];
82  Packet ptmp1 = pset1<Packet>(t1);
83 
84  Scalar t2(0);
85  Packet ptmp2 = pset1<Packet>(t2);
86  Scalar t3(0);
87  Packet ptmp3 = pset1<Packet>(t3);
88 
89  Index starti = FirstTriangular ? 0 : j+2;
90  Index endi = FirstTriangular ? j : size;
91  Index alignedStart = (starti) + internal::first_default_aligned(&res[starti], endi-starti);
92  Index alignedEnd = alignedStart + ((endi-alignedStart)/(PacketSize))*(PacketSize);
93 
94  res[j] += cjd.pmul(numext::real(A0[j]), t0);
95  res[j+1] += cjd.pmul(numext::real(A1[j+1]), t1);
96  if(FirstTriangular)
97  {
98  res[j] += cj0.pmul(A1[j], t1);
99  t3 += cj1.pmul(A1[j], rhs[j]);
100  }
101  else
102  {
103  res[j+1] += cj0.pmul(A0[j+1],t0);
104  t2 += cj1.pmul(A0[j+1], rhs[j+1]);
105  }
106 
107  for (Index i=starti; i<alignedStart; ++i)
108  {
109  res[i] += cj0.pmul(A0[i], t0) + cj0.pmul(A1[i],t1);
110  t2 += cj1.pmul(A0[i], rhs[i]);
111  t3 += cj1.pmul(A1[i], rhs[i]);
112  }
113  // Yes this an optimization for gcc 4.3 and 4.4 (=> huge speed up)
114  // gcc 4.2 does this optimization automatically.
115  const Scalar* EIGEN_RESTRICT a0It = A0 + alignedStart;
116  const Scalar* EIGEN_RESTRICT a1It = A1 + alignedStart;
117  const Scalar* EIGEN_RESTRICT rhsIt = rhs + alignedStart;
118  Scalar* EIGEN_RESTRICT resIt = res + alignedStart;
119  for (Index i=alignedStart; i<alignedEnd; i+=PacketSize)
120  {
121  Packet A0i = ploadu<Packet>(a0It); a0It += PacketSize;
122  Packet A1i = ploadu<Packet>(a1It); a1It += PacketSize;
123  Packet Bi = ploadu<Packet>(rhsIt); rhsIt += PacketSize; // FIXME should be aligned in most cases
124  Packet Xi = pload <Packet>(resIt);
125 
126  Xi = pcj0.pmadd(A0i,ptmp0, pcj0.pmadd(A1i,ptmp1,Xi));
127  ptmp2 = pcj1.pmadd(A0i, Bi, ptmp2);
128  ptmp3 = pcj1.pmadd(A1i, Bi, ptmp3);
129  pstore(resIt,Xi); resIt += PacketSize;
130  }
131  for (Index i=alignedEnd; i<endi; i++)
132  {
133  res[i] += cj0.pmul(A0[i], t0) + cj0.pmul(A1[i],t1);
134  t2 += cj1.pmul(A0[i], rhs[i]);
135  t3 += cj1.pmul(A1[i], rhs[i]);
136  }
137 
138  res[j] += alpha * (t2 + predux(ptmp2));
139  res[j+1] += alpha * (t3 + predux(ptmp3));
140  }
141  for (Index j=FirstTriangular ? 0 : bound;j<(FirstTriangular ? bound : size);j++)
142  {
143  const Scalar* EIGEN_RESTRICT A0 = lhs + j*lhsStride;
144 
145  Scalar t1 = cjAlpha * rhs[j];
146  Scalar t2(0);
147  res[j] += cjd.pmul(numext::real(A0[j]), t1);
148  for (Index i=FirstTriangular ? 0 : j+1; i<(FirstTriangular ? j : size); i++)
149  {
150  res[i] += cj0.pmul(A0[i], t1);
151  t2 += cj1.pmul(A0[i], rhs[i]);
152  }
153  res[j] += alpha * t2;
154  }
155 }
156 
157 } // end namespace internal
158 
159 
163 namespace internal {
164 
165 template<typename Lhs, int LhsMode, typename Rhs>
166 struct selfadjoint_product_impl<Lhs,LhsMode,false,Rhs,0,true>
167 {
168  typedef typename Product<Lhs,Rhs>::Scalar Scalar;
169 
170  typedef internal::blas_traits<Lhs> LhsBlasTraits;
171  typedef typename LhsBlasTraits::DirectLinearAccessType ActualLhsType;
172  typedef internal::remove_all_t<ActualLhsType> ActualLhsTypeCleaned;
173 
174  typedef internal::blas_traits<Rhs> RhsBlasTraits;
175  typedef typename RhsBlasTraits::DirectLinearAccessType ActualRhsType;
176  typedef internal::remove_all_t<ActualRhsType> ActualRhsTypeCleaned;
177 
178  enum { LhsUpLo = LhsMode&(Upper|Lower) };
179 
180  template<typename Dest>
181  static EIGEN_DEVICE_FUNC
182  void run(Dest& dest, const Lhs &a_lhs, const Rhs &a_rhs, const Scalar& alpha)
183  {
184  typedef typename Dest::Scalar ResScalar;
185  typedef typename Rhs::Scalar RhsScalar;
186  typedef Map<Matrix<ResScalar,Dynamic,1>, plain_enum_min(AlignedMax,internal::packet_traits<ResScalar>::size)> MappedDest;
187 
188  eigen_assert(dest.rows()==a_lhs.rows() && dest.cols()==a_rhs.cols());
189 
190  add_const_on_value_type_t<ActualLhsType> lhs = LhsBlasTraits::extract(a_lhs);
191  add_const_on_value_type_t<ActualRhsType> rhs = RhsBlasTraits::extract(a_rhs);
192 
193  Scalar actualAlpha = alpha * LhsBlasTraits::extractScalarFactor(a_lhs)
194  * RhsBlasTraits::extractScalarFactor(a_rhs);
195 
196  enum {
197  EvalToDest = (Dest::InnerStrideAtCompileTime==1),
198  UseRhs = (ActualRhsTypeCleaned::InnerStrideAtCompileTime==1)
199  };
200 
201  internal::gemv_static_vector_if<ResScalar,Dest::SizeAtCompileTime,Dest::MaxSizeAtCompileTime,!EvalToDest> static_dest;
202  internal::gemv_static_vector_if<RhsScalar,ActualRhsTypeCleaned::SizeAtCompileTime,ActualRhsTypeCleaned::MaxSizeAtCompileTime,!UseRhs> static_rhs;
203 
204  ei_declare_aligned_stack_constructed_variable(ResScalar,actualDestPtr,dest.size(),
205  EvalToDest ? dest.data() : static_dest.data());
206 
207  ei_declare_aligned_stack_constructed_variable(RhsScalar,actualRhsPtr,rhs.size(),
208  UseRhs ? const_cast<RhsScalar*>(rhs.data()) : static_rhs.data());
209 
210  if(!EvalToDest)
211  {
212  #ifdef EIGEN_DENSE_STORAGE_CTOR_PLUGIN
213  Index size = dest.size();
214  EIGEN_DENSE_STORAGE_CTOR_PLUGIN
215  #endif
216  MappedDest(actualDestPtr, dest.size()) = dest;
217  }
218 
219  if(!UseRhs)
220  {
221  #ifdef EIGEN_DENSE_STORAGE_CTOR_PLUGIN
222  Index size = rhs.size();
223  EIGEN_DENSE_STORAGE_CTOR_PLUGIN
224  #endif
225  Map<typename ActualRhsTypeCleaned::PlainObject>(actualRhsPtr, rhs.size()) = rhs;
226  }
227 
228 
229  internal::selfadjoint_matrix_vector_product<Scalar, Index, (internal::traits<ActualLhsTypeCleaned>::Flags&RowMajorBit) ? RowMajor : ColMajor,
230  int(LhsUpLo), bool(LhsBlasTraits::NeedToConjugate), bool(RhsBlasTraits::NeedToConjugate)>::run
231  (
232  lhs.rows(), // size
233  &lhs.coeffRef(0,0), lhs.outerStride(), // lhs info
234  actualRhsPtr, // rhs info
235  actualDestPtr, // result info
236  actualAlpha // scale factor
237  );
238 
239  if(!EvalToDest)
240  dest = MappedDest(actualDestPtr, dest.size());
241  }
242 };
243 
244 template<typename Lhs, typename Rhs, int RhsMode>
245 struct selfadjoint_product_impl<Lhs,0,true,Rhs,RhsMode,false>
246 {
247  typedef typename Product<Lhs,Rhs>::Scalar Scalar;
248  enum { RhsUpLo = RhsMode&(Upper|Lower) };
249 
250  template<typename Dest>
251  static void run(Dest& dest, const Lhs &a_lhs, const Rhs &a_rhs, const Scalar& alpha)
252  {
253  // let's simply transpose the product
254  Transpose<Dest> destT(dest);
255  selfadjoint_product_impl<Transpose<const Rhs>, int(RhsUpLo)==Upper ? Lower : Upper, false,
256  Transpose<const Lhs>, 0, true>::run(destT, a_rhs.transpose(), a_lhs.transpose(), alpha);
257  }
258 };
259 
260 } // end namespace internal
261 
262 } // end namespace Eigen
263 
264 #endif // EIGEN_SELFADJOINT_MATRIX_VECTOR_H
RealReturnType real() const
#define EIGEN_RESTRICT
Definition: Macros.h:1055
#define EIGEN_DEVICE_FUNC
Definition: Macros.h:883
#define EIGEN_DONT_INLINE
Definition: Macros.h:844
#define eigen_assert(x)
Definition: Macros.h:902
#define ei_declare_aligned_stack_constructed_variable(TYPE, NAME, SIZE, BUFFER)
Definition: Memory.h:847
cout<< "Here is the matrix m:"<< endl<< m<< endl;Matrix< ptrdiff_t, 3, 1 > res
@ Lower
Definition: Constants.h:211
@ Upper
Definition: Constants.h:213
@ AlignedMax
Definition: Constants.h:254
@ ColMajor
Definition: Constants.h:321
@ RowMajor
Definition: Constants.h:323
const unsigned int RowMajorBit
Definition: Constants.h:68
constexpr int plain_enum_min(A a, B b)
Definition: Meta.h:516
void pstore(Scalar *to, const Packet &from)
unpacket_traits< Packet >::type predux(const Packet &a)
constexpr bool logical_xor(bool a, bool b)
Definition: Meta.h:574
static Index first_default_aligned(const DenseBase< Derived > &m)
EIGEN_ALWAYS_INLINE T maxi(const T &x, const T &y)
: InteropHeaders
Definition: Core:139
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:82
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_conjugate_op< typename Derived::Scalar >, const Derived > conj(const Eigen::ArrayBase< Derived > &x)
std::ptrdiff_t j