TriangularMatrixVector.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 Gael Guennebaud <gael.guennebaud@inria.fr>
5 //
6 // Eigen is free software; you can redistribute it and/or
7 // modify it under the terms of the GNU Lesser General Public
8 // License as published by the Free Software Foundation; either
9 // version 3 of the License, or (at your option) any later version.
10 //
11 // Alternatively, you can redistribute it and/or
12 // modify it under the terms of the GNU General Public License as
13 // published by the Free Software Foundation; either version 2 of
14 // the License, or (at your option) any later version.
15 //
16 // Eigen is distributed in the hope that it will be useful, but WITHOUT ANY
17 // WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
18 // FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the
19 // GNU General Public License for more details.
20 //
21 // You should have received a copy of the GNU Lesser General Public
22 // License and a copy of the GNU General Public License along with
23 // Eigen. If not, see <http://www.gnu.org/licenses/>.
24 
25 #ifndef EIGEN_TRIANGULARMATRIXVECTOR_H
26 #define EIGEN_TRIANGULARMATRIXVECTOR_H
27 
28 namespace Eigen {
29 
30 namespace internal {
31 
32 template<typename Index, int Mode, typename LhsScalar, bool ConjLhs, typename RhsScalar, bool ConjRhs, int StorageOrder, int Version=Specialized>
34 
35 template<typename Index, int Mode, typename LhsScalar, bool ConjLhs, typename RhsScalar, bool ConjRhs, int Version>
36 struct triangular_matrix_vector_product<Index,Mode,LhsScalar,ConjLhs,RhsScalar,ConjRhs,ColMajor,Version>
37 {
38  typedef typename scalar_product_traits<LhsScalar, RhsScalar>::ReturnType ResScalar;
39  enum {
40  IsLower = ((Mode&Lower)==Lower),
41  HasUnitDiag = (Mode & UnitDiag)==UnitDiag,
42  HasZeroDiag = (Mode & ZeroDiag)==ZeroDiag
43  };
44  static EIGEN_DONT_INLINE void run(Index _rows, Index _cols, const LhsScalar* _lhs, Index lhsStride,
45  const RhsScalar* _rhs, Index rhsIncr, ResScalar* _res, Index resIncr, ResScalar alpha)
46  {
47  static const Index PanelWidth = EIGEN_TUNE_TRIANGULAR_PANEL_WIDTH;
48  Index size = (std::min)(_rows,_cols);
49  Index rows = IsLower ? _rows : (std::min)(_rows,_cols);
50  Index cols = IsLower ? (std::min)(_rows,_cols) : _cols;
51 
52  typedef Map<const Matrix<LhsScalar,Dynamic,Dynamic,ColMajor>, 0, OuterStride<> > LhsMap;
53  const LhsMap lhs(_lhs,rows,cols,OuterStride<>(lhsStride));
54  typename conj_expr_if<ConjLhs,LhsMap>::type cjLhs(lhs);
55 
56  typedef Map<const Matrix<RhsScalar,Dynamic,1>, 0, InnerStride<> > RhsMap;
57  const RhsMap rhs(_rhs,cols,InnerStride<>(rhsIncr));
58  typename conj_expr_if<ConjRhs,RhsMap>::type cjRhs(rhs);
59 
60  typedef Map<Matrix<ResScalar,Dynamic,1> > ResMap;
61  ResMap res(_res,rows);
62 
63  for (Index pi=0; pi<size; pi+=PanelWidth)
64  {
65  Index actualPanelWidth = (std::min)(PanelWidth, size-pi);
66  for (Index k=0; k<actualPanelWidth; ++k)
67  {
68  Index i = pi + k;
69  Index s = IsLower ? ((HasUnitDiag||HasZeroDiag) ? i+1 : i ) : pi;
70  Index r = IsLower ? actualPanelWidth-k : k+1;
71  if ((!(HasUnitDiag||HasZeroDiag)) || (--r)>0)
72  res.segment(s,r) += (alpha * cjRhs.coeff(i)) * cjLhs.col(i).segment(s,r);
73  if (HasUnitDiag)
74  res.coeffRef(i) += alpha * cjRhs.coeff(i);
75  }
76  Index r = IsLower ? rows - pi - actualPanelWidth : pi;
77  if (r>0)
78  {
79  Index s = IsLower ? pi+actualPanelWidth : 0;
81  r, actualPanelWidth,
82  &lhs.coeffRef(s,pi), lhsStride,
83  &rhs.coeffRef(pi), rhsIncr,
84  &res.coeffRef(s), resIncr, alpha);
85  }
86  }
87  if((!IsLower) && cols>size)
88  {
90  rows, cols-size,
91  &lhs.coeffRef(0,size), lhsStride,
92  &rhs.coeffRef(size), rhsIncr,
93  _res, resIncr, alpha);
94  }
95  }
96 };
97 
98 template<typename Index, int Mode, typename LhsScalar, bool ConjLhs, typename RhsScalar, bool ConjRhs,int Version>
99 struct triangular_matrix_vector_product<Index,Mode,LhsScalar,ConjLhs,RhsScalar,ConjRhs,RowMajor,Version>
100 {
101  typedef typename scalar_product_traits<LhsScalar, RhsScalar>::ReturnType ResScalar;
102  enum {
103  IsLower = ((Mode&Lower)==Lower),
104  HasUnitDiag = (Mode & UnitDiag)==UnitDiag,
105  HasZeroDiag = (Mode & ZeroDiag)==ZeroDiag
106  };
107  static void run(Index _rows, Index _cols, const LhsScalar* _lhs, Index lhsStride,
108  const RhsScalar* _rhs, Index rhsIncr, ResScalar* _res, Index resIncr, ResScalar alpha)
109  {
110  static const Index PanelWidth = EIGEN_TUNE_TRIANGULAR_PANEL_WIDTH;
111  Index diagSize = (std::min)(_rows,_cols);
112  Index rows = IsLower ? _rows : diagSize;
113  Index cols = IsLower ? diagSize : _cols;
114 
115  typedef Map<const Matrix<LhsScalar,Dynamic,Dynamic,RowMajor>, 0, OuterStride<> > LhsMap;
116  const LhsMap lhs(_lhs,rows,cols,OuterStride<>(lhsStride));
117  typename conj_expr_if<ConjLhs,LhsMap>::type cjLhs(lhs);
118 
119  typedef Map<const Matrix<RhsScalar,Dynamic,1> > RhsMap;
120  const RhsMap rhs(_rhs,cols);
121  typename conj_expr_if<ConjRhs,RhsMap>::type cjRhs(rhs);
122 
123  typedef Map<Matrix<ResScalar,Dynamic,1>, 0, InnerStride<> > ResMap;
124  ResMap res(_res,rows,InnerStride<>(resIncr));
125 
126  for (Index pi=0; pi<diagSize; pi+=PanelWidth)
127  {
128  Index actualPanelWidth = (std::min)(PanelWidth, diagSize-pi);
129  for (Index k=0; k<actualPanelWidth; ++k)
130  {
131  Index i = pi + k;
132  Index s = IsLower ? pi : ((HasUnitDiag||HasZeroDiag) ? i+1 : i);
133  Index r = IsLower ? k+1 : actualPanelWidth-k;
134  if ((!(HasUnitDiag||HasZeroDiag)) || (--r)>0)
135  res.coeffRef(i) += alpha * (cjLhs.row(i).segment(s,r).cwiseProduct(cjRhs.segment(s,r).transpose())).sum();
136  if (HasUnitDiag)
137  res.coeffRef(i) += alpha * cjRhs.coeff(i);
138  }
139  Index r = IsLower ? pi : cols - pi - actualPanelWidth;
140  if (r>0)
141  {
142  Index s = IsLower ? 0 : pi + actualPanelWidth;
144  actualPanelWidth, r,
145  &lhs.coeffRef(pi,s), lhsStride,
146  &rhs.coeffRef(s), rhsIncr,
147  &res.coeffRef(pi), resIncr, alpha);
148  }
149  }
150  if(IsLower && rows>diagSize)
151  {
153  rows-diagSize, cols,
154  &lhs.coeffRef(diagSize,0), lhsStride,
155  &rhs.coeffRef(0), rhsIncr,
156  &res.coeffRef(diagSize), resIncr, alpha);
157  }
158  }
159 };
160 
161 /***************************************************************************
162 * Wrapper to product_triangular_vector
163 ***************************************************************************/
164 
165 template<int Mode, bool LhsIsTriangular, typename Lhs, typename Rhs>
166 struct traits<TriangularProduct<Mode,LhsIsTriangular,Lhs,false,Rhs,true> >
167  : traits<ProductBase<TriangularProduct<Mode,LhsIsTriangular,Lhs,false,Rhs,true>, Lhs, Rhs> >
168 {};
169 
170 template<int Mode, bool LhsIsTriangular, typename Lhs, typename Rhs>
171 struct traits<TriangularProduct<Mode,LhsIsTriangular,Lhs,true,Rhs,false> >
172  : traits<ProductBase<TriangularProduct<Mode,LhsIsTriangular,Lhs,true,Rhs,false>, Lhs, Rhs> >
173 {};
174 
175 
176 template<int StorageOrder>
177 struct trmv_selector;
178 
179 } // end namespace internal
180 
181 template<int Mode, typename Lhs, typename Rhs>
182 struct TriangularProduct<Mode,true,Lhs,false,Rhs,true>
183  : public ProductBase<TriangularProduct<Mode,true,Lhs,false,Rhs,true>, Lhs, Rhs >
184 {
185  EIGEN_PRODUCT_PUBLIC_INTERFACE(TriangularProduct)
186 
187  TriangularProduct(const Lhs& lhs, const Rhs& rhs) : Base(lhs,rhs) {}
188 
189  template<typename Dest> void scaleAndAddTo(Dest& dst, Scalar alpha) const
190  {
191  eigen_assert(dst.rows()==m_lhs.rows() && dst.cols()==m_rhs.cols());
192 
193  internal::trmv_selector<(int(internal::traits<Lhs>::Flags)&RowMajorBit) ? RowMajor : ColMajor>::run(*this, dst, alpha);
194  }
195 };
196 
197 template<int Mode, typename Lhs, typename Rhs>
198 struct TriangularProduct<Mode,false,Lhs,true,Rhs,false>
199  : public ProductBase<TriangularProduct<Mode,false,Lhs,true,Rhs,false>, Lhs, Rhs >
200 {
201  EIGEN_PRODUCT_PUBLIC_INTERFACE(TriangularProduct)
202 
203  TriangularProduct(const Lhs& lhs, const Rhs& rhs) : Base(lhs,rhs) {}
204 
205  template<typename Dest> void scaleAndAddTo(Dest& dst, Scalar alpha) const
206  {
207  eigen_assert(dst.rows()==m_lhs.rows() && dst.cols()==m_rhs.cols());
208 
209  typedef TriangularProduct<(Mode & (UnitDiag|ZeroDiag)) | ((Mode & Lower) ? Upper : Lower),true,Transpose<const Rhs>,false,Transpose<const Lhs>,true> TriangularProductTranspose;
210  Transpose<Dest> dstT(dst);
211  internal::trmv_selector<(int(internal::traits<Rhs>::Flags)&RowMajorBit) ? ColMajor : RowMajor>::run(
212  TriangularProductTranspose(m_rhs.transpose(),m_lhs.transpose()), dstT, alpha);
213  }
214 };
215 
216 namespace internal {
217 
218 // TODO: find a way to factorize this piece of code with gemv_selector since the logic is exactly the same.
219 
220 template<> struct trmv_selector<ColMajor>
221 {
222  template<int Mode, typename Lhs, typename Rhs, typename Dest>
223  static void run(const TriangularProduct<Mode,true,Lhs,false,Rhs,true>& prod, Dest& dest, typename TriangularProduct<Mode,true,Lhs,false,Rhs,true>::Scalar alpha)
224  {
225  typedef TriangularProduct<Mode,true,Lhs,false,Rhs,true> ProductType;
226  typedef typename ProductType::Index Index;
227  typedef typename ProductType::LhsScalar LhsScalar;
228  typedef typename ProductType::RhsScalar RhsScalar;
229  typedef typename ProductType::Scalar ResScalar;
230  typedef typename ProductType::RealScalar RealScalar;
231  typedef typename ProductType::ActualLhsType ActualLhsType;
232  typedef typename ProductType::ActualRhsType ActualRhsType;
233  typedef typename ProductType::LhsBlasTraits LhsBlasTraits;
234  typedef typename ProductType::RhsBlasTraits RhsBlasTraits;
235  typedef Map<Matrix<ResScalar,Dynamic,1>, Aligned> MappedDest;
236 
237  typename internal::add_const_on_value_type<ActualLhsType>::type actualLhs = LhsBlasTraits::extract(prod.lhs());
238  typename internal::add_const_on_value_type<ActualRhsType>::type actualRhs = RhsBlasTraits::extract(prod.rhs());
239 
240  ResScalar actualAlpha = alpha * LhsBlasTraits::extractScalarFactor(prod.lhs())
241  * RhsBlasTraits::extractScalarFactor(prod.rhs());
242 
243  enum {
244  // FIXME find a way to allow an inner stride on the result if packet_traits<Scalar>::size==1
245  // on, the other hand it is good for the cache to pack the vector anyways...
246  EvalToDestAtCompileTime = Dest::InnerStrideAtCompileTime==1,
248  MightCannotUseDest = (Dest::InnerStrideAtCompileTime!=1) || ComplexByReal
249  };
250 
251  gemv_static_vector_if<ResScalar,Dest::SizeAtCompileTime,Dest::MaxSizeAtCompileTime,MightCannotUseDest> static_dest;
252 
253  bool alphaIsCompatible = (!ComplexByReal) || (imag(actualAlpha)==RealScalar(0));
254  bool evalToDest = EvalToDestAtCompileTime && alphaIsCompatible;
255 
256  RhsScalar compatibleAlpha = get_factor<ResScalar,RhsScalar>::run(actualAlpha);
257 
258  ei_declare_aligned_stack_constructed_variable(ResScalar,actualDestPtr,dest.size(),
259  evalToDest ? dest.data() : static_dest.data());
260 
261  if(!evalToDest)
262  {
263  #ifdef EIGEN_DENSE_STORAGE_CTOR_PLUGIN
264  int size = dest.size();
265  EIGEN_DENSE_STORAGE_CTOR_PLUGIN
266  #endif
267  if(!alphaIsCompatible)
268  {
269  MappedDest(actualDestPtr, dest.size()).setZero();
270  compatibleAlpha = RhsScalar(1);
271  }
272  else
273  MappedDest(actualDestPtr, dest.size()) = dest;
274  }
275 
276  internal::triangular_matrix_vector_product
277  <Index,Mode,
278  LhsScalar, LhsBlasTraits::NeedToConjugate,
279  RhsScalar, RhsBlasTraits::NeedToConjugate,
280  ColMajor>
281  ::run(actualLhs.rows(),actualLhs.cols(),
282  actualLhs.data(),actualLhs.outerStride(),
283  actualRhs.data(),actualRhs.innerStride(),
284  actualDestPtr,1,compatibleAlpha);
285 
286  if (!evalToDest)
287  {
288  if(!alphaIsCompatible)
289  dest += actualAlpha * MappedDest(actualDestPtr, dest.size());
290  else
291  dest = MappedDest(actualDestPtr, dest.size());
292  }
293  }
294 };
295 
296 template<> struct trmv_selector<RowMajor>
297 {
298  template<int Mode, typename Lhs, typename Rhs, typename Dest>
299  static void run(const TriangularProduct<Mode,true,Lhs,false,Rhs,true>& prod, Dest& dest, typename TriangularProduct<Mode,true,Lhs,false,Rhs,true>::Scalar alpha)
300  {
301  typedef TriangularProduct<Mode,true,Lhs,false,Rhs,true> ProductType;
302  typedef typename ProductType::LhsScalar LhsScalar;
303  typedef typename ProductType::RhsScalar RhsScalar;
304  typedef typename ProductType::Scalar ResScalar;
305  typedef typename ProductType::Index Index;
306  typedef typename ProductType::ActualLhsType ActualLhsType;
307  typedef typename ProductType::ActualRhsType ActualRhsType;
308  typedef typename ProductType::_ActualRhsType _ActualRhsType;
309  typedef typename ProductType::LhsBlasTraits LhsBlasTraits;
310  typedef typename ProductType::RhsBlasTraits RhsBlasTraits;
311 
312  typename add_const<ActualLhsType>::type actualLhs = LhsBlasTraits::extract(prod.lhs());
313  typename add_const<ActualRhsType>::type actualRhs = RhsBlasTraits::extract(prod.rhs());
314 
315  ResScalar actualAlpha = alpha * LhsBlasTraits::extractScalarFactor(prod.lhs())
316  * RhsBlasTraits::extractScalarFactor(prod.rhs());
317 
318  enum {
319  DirectlyUseRhs = _ActualRhsType::InnerStrideAtCompileTime==1
320  };
321 
322  gemv_static_vector_if<RhsScalar,_ActualRhsType::SizeAtCompileTime,_ActualRhsType::MaxSizeAtCompileTime,!DirectlyUseRhs> static_rhs;
323 
324  ei_declare_aligned_stack_constructed_variable(RhsScalar,actualRhsPtr,actualRhs.size(),
325  DirectlyUseRhs ? const_cast<RhsScalar*>(actualRhs.data()) : static_rhs.data());
326 
327  if(!DirectlyUseRhs)
328  {
329  #ifdef EIGEN_DENSE_STORAGE_CTOR_PLUGIN
330  int size = actualRhs.size();
331  EIGEN_DENSE_STORAGE_CTOR_PLUGIN
332  #endif
333  Map<typename _ActualRhsType::PlainObject>(actualRhsPtr, actualRhs.size()) = actualRhs;
334  }
335 
336  internal::triangular_matrix_vector_product
337  <Index,Mode,
338  LhsScalar, LhsBlasTraits::NeedToConjugate,
339  RhsScalar, RhsBlasTraits::NeedToConjugate,
340  RowMajor>
341  ::run(actualLhs.rows(),actualLhs.cols(),
342  actualLhs.data(),actualLhs.outerStride(),
343  actualRhsPtr,1,
344  dest.data(),dest.innerStride(),
345  actualAlpha);
346  }
347 };
348 
349 } // end namespace internal
350 
351 } // end namespace Eigen
352 
353 #endif // EIGEN_TRIANGULARMATRIXVECTOR_H