Another set of UI messages fixes and tweaks! No functional changes.
[blender.git] / extern / Eigen3 / Eigen / src / Core / SolveTriangular.h
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 internal {
29
30 // Forward declarations:
31 // The following two routines are implemented in the products/TriangularSolver*.h files
32 template<typename LhsScalar, typename RhsScalar, typename Index, int Side, int Mode, bool Conjugate, int StorageOrder>
33 struct triangular_solve_vector;
34
35 template <typename Scalar, typename Index, int Side, int Mode, bool Conjugate, int TriStorageOrder, int OtherStorageOrder>
36 struct triangular_solve_matrix;
37
38 // small helper struct extracting some traits on the underlying solver operation
39 template<typename Lhs, typename Rhs, int Side>
40 class trsolve_traits
41 {
42   private:
43     enum {
44       RhsIsVectorAtCompileTime = (Side==OnTheLeft ? Rhs::ColsAtCompileTime : Rhs::RowsAtCompileTime)==1
45     };
46   public:
47     enum {
48       Unrolling   = (RhsIsVectorAtCompileTime && Rhs::SizeAtCompileTime != Dynamic && Rhs::SizeAtCompileTime <= 8)
49                   ? CompleteUnrolling : NoUnrolling,
50       RhsVectors  = RhsIsVectorAtCompileTime ? 1 : Dynamic
51     };
52 };
53
54 template<typename Lhs, typename Rhs,
55   int Side, // can be OnTheLeft/OnTheRight
56   int Mode, // can be Upper/Lower | UnitDiag
57   int Unrolling = trsolve_traits<Lhs,Rhs,Side>::Unrolling,
58   int RhsVectors = trsolve_traits<Lhs,Rhs,Side>::RhsVectors
59   >
60 struct triangular_solver_selector;
61
62 template<typename Lhs, typename Rhs, int Side, int Mode>
63 struct triangular_solver_selector<Lhs,Rhs,Side,Mode,NoUnrolling,1>
64 {
65   typedef typename Lhs::Scalar LhsScalar;
66   typedef typename Rhs::Scalar RhsScalar;
67   typedef blas_traits<Lhs> LhsProductTraits;
68   typedef typename LhsProductTraits::ExtractType ActualLhsType;
69   typedef Map<Matrix<RhsScalar,Dynamic,1>, Aligned> MappedRhs;
70   static void run(const Lhs& lhs, Rhs& rhs)
71   {
72     ActualLhsType actualLhs = LhsProductTraits::extract(lhs);
73
74     // FIXME find a way to allow an inner stride if packet_traits<Scalar>::size==1
75
76     bool useRhsDirectly = Rhs::InnerStrideAtCompileTime==1 || rhs.innerStride()==1;
77
78     ei_declare_aligned_stack_constructed_variable(RhsScalar,actualRhs,rhs.size(),
79                                                   (useRhsDirectly ? rhs.data() : 0));
80                                                   
81     if(!useRhsDirectly)
82       MappedRhs(actualRhs,rhs.size()) = rhs;
83
84     triangular_solve_vector<LhsScalar, RhsScalar, typename Lhs::Index, Side, Mode, LhsProductTraits::NeedToConjugate,
85                             (int(Lhs::Flags) & RowMajorBit) ? RowMajor : ColMajor>
86       ::run(actualLhs.cols(), actualLhs.data(), actualLhs.outerStride(), actualRhs);
87
88     if(!useRhsDirectly)
89       rhs = MappedRhs(actualRhs, rhs.size());
90   }
91 };
92
93 // the rhs is a matrix
94 template<typename Lhs, typename Rhs, int Side, int Mode>
95 struct triangular_solver_selector<Lhs,Rhs,Side,Mode,NoUnrolling,Dynamic>
96 {
97   typedef typename Rhs::Scalar Scalar;
98   typedef typename Rhs::Index Index;
99   typedef blas_traits<Lhs> LhsProductTraits;
100   typedef typename LhsProductTraits::DirectLinearAccessType ActualLhsType;
101   static void run(const Lhs& lhs, Rhs& rhs)
102   {
103     const ActualLhsType actualLhs = LhsProductTraits::extract(lhs);
104     triangular_solve_matrix<Scalar,Index,Side,Mode,LhsProductTraits::NeedToConjugate,(int(Lhs::Flags) & RowMajorBit) ? RowMajor : ColMajor,
105                                (Rhs::Flags&RowMajorBit) ? RowMajor : ColMajor>
106       ::run(lhs.rows(), Side==OnTheLeft? rhs.cols() : rhs.rows(), &actualLhs.coeffRef(0,0), actualLhs.outerStride(), &rhs.coeffRef(0,0), rhs.outerStride());
107   }
108 };
109
110 /***************************************************************************
111 * meta-unrolling implementation
112 ***************************************************************************/
113
114 template<typename Lhs, typename Rhs, int Mode, int Index, int Size,
115          bool Stop = Index==Size>
116 struct triangular_solver_unroller;
117
118 template<typename Lhs, typename Rhs, int Mode, int Index, int Size>
119 struct triangular_solver_unroller<Lhs,Rhs,Mode,Index,Size,false> {
120   enum {
121     IsLower = ((Mode&Lower)==Lower),
122     I = IsLower ? Index : Size - Index - 1,
123     S = IsLower ? 0     : I+1
124   };
125   static void run(const Lhs& lhs, Rhs& rhs)
126   {
127     if (Index>0)
128       rhs.coeffRef(I) -= lhs.row(I).template segment<Index>(S).transpose()
129                          .cwiseProduct(rhs.template segment<Index>(S)).sum();
130
131     if(!(Mode & UnitDiag))
132       rhs.coeffRef(I) /= lhs.coeff(I,I);
133
134     triangular_solver_unroller<Lhs,Rhs,Mode,Index+1,Size>::run(lhs,rhs);
135   }
136 };
137
138 template<typename Lhs, typename Rhs, int Mode, int Index, int Size>
139 struct triangular_solver_unroller<Lhs,Rhs,Mode,Index,Size,true> {
140   static void run(const Lhs&, Rhs&) {}
141 };
142
143 template<typename Lhs, typename Rhs, int Mode>
144 struct triangular_solver_selector<Lhs,Rhs,OnTheLeft,Mode,CompleteUnrolling,1> {
145   static void run(const Lhs& lhs, Rhs& rhs)
146   { triangular_solver_unroller<Lhs,Rhs,Mode,0,Rhs::SizeAtCompileTime>::run(lhs,rhs); }
147 };
148
149 template<typename Lhs, typename Rhs, int Mode>
150 struct triangular_solver_selector<Lhs,Rhs,OnTheRight,Mode,CompleteUnrolling,1> {
151   static void run(const Lhs& lhs, Rhs& rhs)
152   {
153     Transpose<const Lhs> trLhs(lhs);
154     Transpose<Rhs> trRhs(rhs);
155     
156     triangular_solver_unroller<Transpose<const Lhs>,Transpose<Rhs>,
157                               ((Mode&Upper)==Upper ? Lower : Upper) | (Mode&UnitDiag),
158                               0,Rhs::SizeAtCompileTime>::run(trLhs,trRhs);
159   }
160 };
161
162 } // end namespace internal
163
164 /***************************************************************************
165 * TriangularView methods
166 ***************************************************************************/
167
168 /** "in-place" version of TriangularView::solve() where the result is written in \a other
169   *
170   * \warning The parameter is only marked 'const' to make the C++ compiler accept a temporary expression here.
171   * This function will const_cast it, so constness isn't honored here.
172   *
173   * See TriangularView:solve() for the details.
174   */
175 template<typename MatrixType, unsigned int Mode>
176 template<int Side, typename OtherDerived>
177 void TriangularView<MatrixType,Mode>::solveInPlace(const MatrixBase<OtherDerived>& _other) const
178 {
179   OtherDerived& other = _other.const_cast_derived();
180   eigen_assert(cols() == rows());
181   eigen_assert( (Side==OnTheLeft && cols() == other.rows()) || (Side==OnTheRight && cols() == other.cols()) );
182   eigen_assert(!(Mode & ZeroDiag));
183   eigen_assert(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
197 /** \returns the product of the inverse of \c *this with \a other, \a *this being triangular.
198   *
199   * This function computes the inverse-matrix matrix product inverse(\c *this) * \a other if
200   * \a Side==OnTheLeft (the default), or the right-inverse-multiply  \a other * inverse(\c *this) if
201   * \a Side==OnTheRight.
202   *
203   * The matrix \c *this must be triangular and invertible (i.e., all the coefficients of the
204   * diagonal must be non zero). It works as a forward (resp. backward) substitution if \c *this
205   * is an upper (resp. lower) triangular matrix.
206   *
207   * Example: \include MatrixBase_marked.cpp
208   * Output: \verbinclude MatrixBase_marked.out
209   *
210   * This function returns an expression of the inverse-multiply and can works in-place if it is assigned
211   * to the same matrix or vector \a other.
212   *
213   * For users coming from BLAS, this function (and more specifically solveInPlace()) offer
214   * all the operations supported by the \c *TRSV and \c *TRSM BLAS routines.
215   *
216   * \sa TriangularView::solveInPlace()
217   */
218 template<typename Derived, unsigned int Mode>
219 template<int Side, typename Other>
220 const internal::triangular_solve_retval<Side,TriangularView<Derived,Mode>,Other>
221 TriangularView<Derived,Mode>::solve(const MatrixBase<Other>& other) const
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     const typename Rhs::Nested m_rhs;
259 };
260
261 } // namespace internal
262
263 #endif // EIGEN_SOLVETRIANGULAR_H