SolveTriangular.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 // 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_SOLVETRIANGULAR_H
26 #define EIGEN_SOLVETRIANGULAR_H
27 
28 namespace Eigen {
29 
30 namespace internal {
31 
32 // Forward declarations:
33 // The following two routines are implemented in the products/TriangularSolver*.h files
34 template<typename LhsScalar, typename RhsScalar, typename Index, int Side, int Mode, bool Conjugate, int StorageOrder>
35 struct triangular_solve_vector;
36 
37 template <typename Scalar, typename Index, int Side, int Mode, bool Conjugate, int TriStorageOrder, int OtherStorageOrder>
38 struct triangular_solve_matrix;
39 
40 // small helper struct extracting some traits on the underlying solver operation
41 template<typename Lhs, typename Rhs, int Side>
42 class trsolve_traits
43 {
44  private:
45  enum {
46  RhsIsVectorAtCompileTime = (Side==OnTheLeft ? Rhs::ColsAtCompileTime : Rhs::RowsAtCompileTime)==1
47  };
48  public:
49  enum {
50  Unrolling = (RhsIsVectorAtCompileTime && Rhs::SizeAtCompileTime != Dynamic && Rhs::SizeAtCompileTime <= 8)
52  RhsVectors = RhsIsVectorAtCompileTime ? 1 : Dynamic
53  };
54 };
55 
56 template<typename Lhs, typename Rhs,
57  int Side, // can be OnTheLeft/OnTheRight
58  int Mode, // can be Upper/Lower | UnitDiag
59  int Unrolling = trsolve_traits<Lhs,Rhs,Side>::Unrolling,
60  int RhsVectors = trsolve_traits<Lhs,Rhs,Side>::RhsVectors
61  >
62 struct triangular_solver_selector;
63 
64 template<typename Lhs, typename Rhs, int Side, int Mode>
65 struct triangular_solver_selector<Lhs,Rhs,Side,Mode,NoUnrolling,1>
66 {
67  typedef typename Lhs::Scalar LhsScalar;
68  typedef typename Rhs::Scalar RhsScalar;
69  typedef blas_traits<Lhs> LhsProductTraits;
70  typedef typename LhsProductTraits::ExtractType ActualLhsType;
71  typedef Map<Matrix<RhsScalar,Dynamic,1>, Aligned> MappedRhs;
72  static void run(const Lhs& lhs, Rhs& rhs)
73  {
74  ActualLhsType actualLhs = LhsProductTraits::extract(lhs);
75 
76  // FIXME find a way to allow an inner stride if packet_traits<Scalar>::size==1
77 
78  bool useRhsDirectly = Rhs::InnerStrideAtCompileTime==1 || rhs.innerStride()==1;
79 
80  ei_declare_aligned_stack_constructed_variable(RhsScalar,actualRhs,rhs.size(),
81  (useRhsDirectly ? rhs.data() : 0));
82 
83  if(!useRhsDirectly)
84  MappedRhs(actualRhs,rhs.size()) = rhs;
85 
86  triangular_solve_vector<LhsScalar, RhsScalar, typename Lhs::Index, Side, Mode, LhsProductTraits::NeedToConjugate,
87  (int(Lhs::Flags) & RowMajorBit) ? RowMajor : ColMajor>
88  ::run(actualLhs.cols(), actualLhs.data(), actualLhs.outerStride(), actualRhs);
89 
90  if(!useRhsDirectly)
91  rhs = MappedRhs(actualRhs, rhs.size());
92  }
93 };
94 
95 // the rhs is a matrix
96 template<typename Lhs, typename Rhs, int Side, int Mode>
97 struct triangular_solver_selector<Lhs,Rhs,Side,Mode,NoUnrolling,Dynamic>
98 {
99  typedef typename Rhs::Scalar Scalar;
100  typedef typename Rhs::Index Index;
101  typedef blas_traits<Lhs> LhsProductTraits;
102  typedef typename LhsProductTraits::DirectLinearAccessType ActualLhsType;
103  static void run(const Lhs& lhs, Rhs& rhs)
104  {
105  typename internal::add_const_on_value_type<ActualLhsType>::type actualLhs = LhsProductTraits::extract(lhs);
106  triangular_solve_matrix<Scalar,Index,Side,Mode,LhsProductTraits::NeedToConjugate,(int(Lhs::Flags) & RowMajorBit) ? RowMajor : ColMajor,
107  (Rhs::Flags&RowMajorBit) ? RowMajor : ColMajor>
108  ::run(lhs.rows(), Side==OnTheLeft? rhs.cols() : rhs.rows(), &actualLhs.coeffRef(0,0), actualLhs.outerStride(), &rhs.coeffRef(0,0), rhs.outerStride());
109  }
110 };
111 
112 /***************************************************************************
113 * meta-unrolling implementation
114 ***************************************************************************/
115 
116 template<typename Lhs, typename Rhs, int Mode, int Index, int Size,
117  bool Stop = Index==Size>
118 struct triangular_solver_unroller;
119 
120 template<typename Lhs, typename Rhs, int Mode, int Index, int Size>
121 struct triangular_solver_unroller<Lhs,Rhs,Mode,Index,Size,false> {
122  enum {
123  IsLower = ((Mode&Lower)==Lower),
124  I = IsLower ? Index : Size - Index - 1,
125  S = IsLower ? 0 : I+1
126  };
127  static void run(const Lhs& lhs, Rhs& rhs)
128  {
129  if (Index>0)
130  rhs.coeffRef(I) -= lhs.row(I).template segment<Index>(S).transpose()
131  .cwiseProduct(rhs.template segment<Index>(S)).sum();
132 
133  if(!(Mode & UnitDiag))
134  rhs.coeffRef(I) /= lhs.coeff(I,I);
135 
136  triangular_solver_unroller<Lhs,Rhs,Mode,Index+1,Size>::run(lhs,rhs);
137  }
138 };
139 
140 template<typename Lhs, typename Rhs, int Mode, int Index, int Size>
141 struct triangular_solver_unroller<Lhs,Rhs,Mode,Index,Size,true> {
142  static void run(const Lhs&, Rhs&) {}
143 };
144 
145 template<typename Lhs, typename Rhs, int Mode>
146 struct triangular_solver_selector<Lhs,Rhs,OnTheLeft,Mode,CompleteUnrolling,1> {
147  static void run(const Lhs& lhs, Rhs& rhs)
148  { triangular_solver_unroller<Lhs,Rhs,Mode,0,Rhs::SizeAtCompileTime>::run(lhs,rhs); }
149 };
150 
151 template<typename Lhs, typename Rhs, int Mode>
152 struct triangular_solver_selector<Lhs,Rhs,OnTheRight,Mode,CompleteUnrolling,1> {
153  static void run(const Lhs& lhs, Rhs& rhs)
154  {
155  Transpose<const Lhs> trLhs(lhs);
156  Transpose<Rhs> trRhs(rhs);
157 
158  triangular_solver_unroller<Transpose<const Lhs>,Transpose<Rhs>,
159  ((Mode&Upper)==Upper ? Lower : Upper) | (Mode&UnitDiag),
160  0,Rhs::SizeAtCompileTime>::run(trLhs,trRhs);
161  }
162 };
163 
164 } // end namespace internal
165 
166 /***************************************************************************
167 * TriangularView methods
168 ***************************************************************************/
169 
177 template<typename MatrixType, unsigned int Mode>
178 template<int Side, typename OtherDerived>
180 {
181  OtherDerived& other = _other.const_cast_derived();
182  eigen_assert( cols() == rows() && ((Side==OnTheLeft && cols() == other.rows()) || (Side==OnTheRight && cols() == other.cols())) );
183  eigen_assert((!(Mode & ZeroDiag)) && bool(Mode & (Upper|Lower)));
184 
185  enum { copy = internal::traits<OtherDerived>::Flags & RowMajorBit && OtherDerived::IsVectorAtCompileTime };
186  typedef typename internal::conditional<copy,
187  typename internal::plain_matrix_type_column_major<OtherDerived>::type, OtherDerived&>::type OtherCopy;
188  OtherCopy otherCopy(other);
189 
190  internal::triangular_solver_selector<MatrixType, typename internal::remove_reference<OtherCopy>::type,
191  Side, Mode>::run(nestedExpression(), otherCopy);
192 
193  if (copy)
194  other = otherCopy;
195 }
196 
218 template<typename Derived, unsigned int Mode>
219 template<int Side, typename Other>
220 const internal::triangular_solve_retval<Side,TriangularView<Derived,Mode>,Other>
222 {
223  return internal::triangular_solve_retval<Side,TriangularView,Other>(*this, other.derived());
224 }
225 
226 namespace internal {
227 
228 
229 template<int Side, typename TriangularType, typename Rhs>
230 struct traits<triangular_solve_retval<Side, TriangularType, Rhs> >
231 {
232  typedef typename internal::plain_matrix_type_column_major<Rhs>::type ReturnType;
233 };
234 
235 template<int Side, typename TriangularType, typename Rhs> struct triangular_solve_retval
236  : public ReturnByValue<triangular_solve_retval<Side, TriangularType, Rhs> >
237 {
238  typedef typename remove_all<typename Rhs::Nested>::type RhsNestedCleaned;
239  typedef ReturnByValue<triangular_solve_retval> Base;
240  typedef typename Base::Index Index;
241 
242  triangular_solve_retval(const TriangularType& tri, const Rhs& rhs)
243  : m_triangularMatrix(tri), m_rhs(rhs)
244  {}
245 
246  inline Index rows() const { return m_rhs.rows(); }
247  inline Index cols() const { return m_rhs.cols(); }
248 
249  template<typename Dest> inline void evalTo(Dest& dst) const
250  {
251  if(!(is_same<RhsNestedCleaned,Dest>::value && extract_data(dst) == extract_data(m_rhs)))
252  dst = m_rhs;
253  m_triangularMatrix.template solveInPlace<Side>(dst);
254  }
255 
256  protected:
257  const TriangularType& m_triangularMatrix;
258  typename Rhs::Nested m_rhs;
259 };
260 
261 } // namespace internal
262 
263 } // end namespace Eigen
264 
265 #endif // EIGEN_SOLVETRIANGULAR_H