SolveTriangular.h
Go to the documentation of this file.
00001 // This file is part of Eigen, a lightweight C++ template library
00002 // for linear algebra.
00003 //
00004 // Copyright (C) 2008-2009 Gael Guennebaud <gael.guennebaud@inria.fr>
00005 //
00006 // Eigen is free software; you can redistribute it and/or
00007 // modify it under the terms of the GNU Lesser General Public
00008 // License as published by the Free Software Foundation; either
00009 // version 3 of the License, or (at your option) any later version.
00010 //
00011 // Alternatively, you can redistribute it and/or
00012 // modify it under the terms of the GNU General Public License as
00013 // published by the Free Software Foundation; either version 2 of
00014 // the License, or (at your option) any later version.
00015 //
00016 // Eigen is distributed in the hope that it will be useful, but WITHOUT ANY
00017 // WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
00018 // FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the
00019 // GNU General Public License for more details.
00020 //
00021 // You should have received a copy of the GNU Lesser General Public
00022 // License and a copy of the GNU General Public License along with
00023 // Eigen. If not, see <http://www.gnu.org/licenses/>.
00024 
00025 #ifndef EIGEN_SOLVETRIANGULAR_H
00026 #define EIGEN_SOLVETRIANGULAR_H
00027 
00028 namespace internal {
00029 
00030 // Forward declarations:
00031 // The following two routines are implemented in the products/TriangularSolver*.h files
00032 template<typename LhsScalar, typename RhsScalar, typename Index, int Side, int Mode, bool Conjugate, int StorageOrder>
00033 struct triangular_solve_vector;
00034 
00035 template <typename Scalar, typename Index, int Side, int Mode, bool Conjugate, int TriStorageOrder, int OtherStorageOrder>
00036 struct triangular_solve_matrix;
00037 
00038 // small helper struct extracting some traits on the underlying solver operation
00039 template<typename Lhs, typename Rhs, int Side>
00040 class trsolve_traits
00041 {
00042   private:
00043     enum {
00044       RhsIsVectorAtCompileTime = (Side==OnTheLeft ? Rhs::ColsAtCompileTime : Rhs::RowsAtCompileTime)==1
00045     };
00046   public:
00047     enum {
00048       Unrolling   = (RhsIsVectorAtCompileTime && Rhs::SizeAtCompileTime != Dynamic && Rhs::SizeAtCompileTime <= 8)
00049                   ? CompleteUnrolling : NoUnrolling,
00050       RhsVectors  = RhsIsVectorAtCompileTime ? 1 : Dynamic
00051     };
00052 };
00053 
00054 template<typename Lhs, typename Rhs,
00055   int Side, // can be OnTheLeft/OnTheRight
00056   int Mode, // can be Upper/Lower | UnitDiag
00057   int Unrolling = trsolve_traits<Lhs,Rhs,Side>::Unrolling,
00058   int RhsVectors = trsolve_traits<Lhs,Rhs,Side>::RhsVectors
00059   >
00060 struct triangular_solver_selector;
00061 
00062 template<typename Lhs, typename Rhs, int Side, int Mode>
00063 struct triangular_solver_selector<Lhs,Rhs,Side,Mode,NoUnrolling,1>
00064 {
00065   typedef typename Lhs::Scalar LhsScalar;
00066   typedef typename Rhs::Scalar RhsScalar;
00067   typedef blas_traits<Lhs> LhsProductTraits;
00068   typedef typename LhsProductTraits::ExtractType ActualLhsType;
00069   typedef Map<Matrix<RhsScalar,Dynamic,1>, Aligned> MappedRhs;
00070   static void run(const Lhs& lhs, Rhs& rhs)
00071   {
00072     ActualLhsType actualLhs = LhsProductTraits::extract(lhs);
00073 
00074     // FIXME find a way to allow an inner stride if packet_traits<Scalar>::size==1
00075 
00076     bool useRhsDirectly = Rhs::InnerStrideAtCompileTime==1 || rhs.innerStride()==1;
00077 
00078     ei_declare_aligned_stack_constructed_variable(RhsScalar,actualRhs,rhs.size(),
00079                                                   (useRhsDirectly ? rhs.data() : 0));
00080                                                   
00081     if(!useRhsDirectly)
00082       MappedRhs(actualRhs,rhs.size()) = rhs;
00083 
00084     triangular_solve_vector<LhsScalar, RhsScalar, typename Lhs::Index, Side, Mode, LhsProductTraits::NeedToConjugate,
00085                             (int(Lhs::Flags) & RowMajorBit) ? RowMajor : ColMajor>
00086       ::run(actualLhs.cols(), actualLhs.data(), actualLhs.outerStride(), actualRhs);
00087 
00088     if(!useRhsDirectly)
00089       rhs = MappedRhs(actualRhs, rhs.size());
00090   }
00091 };
00092 
00093 // the rhs is a matrix
00094 template<typename Lhs, typename Rhs, int Side, int Mode>
00095 struct triangular_solver_selector<Lhs,Rhs,Side,Mode,NoUnrolling,Dynamic>
00096 {
00097   typedef typename Rhs::Scalar Scalar;
00098   typedef typename Rhs::Index Index;
00099   typedef blas_traits<Lhs> LhsProductTraits;
00100   typedef typename LhsProductTraits::DirectLinearAccessType ActualLhsType;
00101   static void run(const Lhs& lhs, Rhs& rhs)
00102   {
00103     const ActualLhsType actualLhs = LhsProductTraits::extract(lhs);
00104     triangular_solve_matrix<Scalar,Index,Side,Mode,LhsProductTraits::NeedToConjugate,(int(Lhs::Flags) & RowMajorBit) ? RowMajor : ColMajor,
00105                                (Rhs::Flags&RowMajorBit) ? RowMajor : ColMajor>
00106       ::run(lhs.rows(), Side==OnTheLeft? rhs.cols() : rhs.rows(), &actualLhs.coeffRef(0,0), actualLhs.outerStride(), &rhs.coeffRef(0,0), rhs.outerStride());
00107   }
00108 };
00109 
00110 /***************************************************************************
00111 * meta-unrolling implementation
00112 ***************************************************************************/
00113 
00114 template<typename Lhs, typename Rhs, int Mode, int Index, int Size,
00115          bool Stop = Index==Size>
00116 struct triangular_solver_unroller;
00117 
00118 template<typename Lhs, typename Rhs, int Mode, int Index, int Size>
00119 struct triangular_solver_unroller<Lhs,Rhs,Mode,Index,Size,false> {
00120   enum {
00121     IsLower = ((Mode&Lower)==Lower),
00122     I = IsLower ? Index : Size - Index - 1,
00123     S = IsLower ? 0     : I+1
00124   };
00125   static void run(const Lhs& lhs, Rhs& rhs)
00126   {
00127     if (Index>0)
00128       rhs.coeffRef(I) -= lhs.row(I).template segment<Index>(S).transpose()
00129                          .cwiseProduct(rhs.template segment<Index>(S)).sum();
00130 
00131     if(!(Mode & UnitDiag))
00132       rhs.coeffRef(I) /= lhs.coeff(I,I);
00133 
00134     triangular_solver_unroller<Lhs,Rhs,Mode,Index+1,Size>::run(lhs,rhs);
00135   }
00136 };
00137 
00138 template<typename Lhs, typename Rhs, int Mode, int Index, int Size>
00139 struct triangular_solver_unroller<Lhs,Rhs,Mode,Index,Size,true> {
00140   static void run(const Lhs&, Rhs&) {}
00141 };
00142 
00143 template<typename Lhs, typename Rhs, int Mode>
00144 struct triangular_solver_selector<Lhs,Rhs,OnTheLeft,Mode,CompleteUnrolling,1> {
00145   static void run(const Lhs& lhs, Rhs& rhs)
00146   { triangular_solver_unroller<Lhs,Rhs,Mode,0,Rhs::SizeAtCompileTime>::run(lhs,rhs); }
00147 };
00148 
00149 template<typename Lhs, typename Rhs, int Mode>
00150 struct triangular_solver_selector<Lhs,Rhs,OnTheRight,Mode,CompleteUnrolling,1> {
00151   static void run(const Lhs& lhs, Rhs& rhs)
00152   {
00153     Transpose<const Lhs> trLhs(lhs);
00154     Transpose<Rhs> trRhs(rhs);
00155     
00156     triangular_solver_unroller<Transpose<const Lhs>,Transpose<Rhs>,
00157                               ((Mode&Upper)==Upper ? Lower : Upper) | (Mode&UnitDiag),
00158                               0,Rhs::SizeAtCompileTime>::run(trLhs,trRhs);
00159   }
00160 };
00161 
00162 } // end namespace internal
00163 
00164 /***************************************************************************
00165 * TriangularView methods
00166 ***************************************************************************/
00167 
00175 template<typename MatrixType, unsigned int Mode>
00176 template<int Side, typename OtherDerived>
00177 void TriangularView<MatrixType,Mode>::solveInPlace(const MatrixBase<OtherDerived>& _other) const
00178 {
00179   OtherDerived& other = _other.const_cast_derived();
00180   eigen_assert(cols() == rows());
00181   eigen_assert( (Side==OnTheLeft && cols() == other.rows()) || (Side==OnTheRight && cols() == other.cols()) );
00182   eigen_assert(!(Mode & ZeroDiag));
00183   eigen_assert(Mode & (Upper|Lower));
00184 
00185   enum { copy = internal::traits<OtherDerived>::Flags & RowMajorBit  && OtherDerived::IsVectorAtCompileTime };
00186   typedef typename internal::conditional<copy,
00187     typename internal::plain_matrix_type_column_major<OtherDerived>::type, OtherDerived&>::type OtherCopy;
00188   OtherCopy otherCopy(other);
00189 
00190   internal::triangular_solver_selector<MatrixType, typename internal::remove_reference<OtherCopy>::type,
00191     Side, Mode>::run(nestedExpression(), otherCopy);
00192 
00193   if (copy)
00194     other = otherCopy;
00195 }
00196 
00218 template<typename Derived, unsigned int Mode>
00219 template<int Side, typename Other>
00220 const internal::triangular_solve_retval<Side,TriangularView<Derived,Mode>,Other>
00221 TriangularView<Derived,Mode>::solve(const MatrixBase<Other>& other) const
00222 {
00223   return internal::triangular_solve_retval<Side,TriangularView,Other>(*this, other.derived());
00224 }
00225 
00226 namespace internal {
00227 
00228 
00229 template<int Side, typename TriangularType, typename Rhs>
00230 struct traits<triangular_solve_retval<Side, TriangularType, Rhs> >
00231 {
00232   typedef typename internal::plain_matrix_type_column_major<Rhs>::type ReturnType;
00233 };
00234 
00235 template<int Side, typename TriangularType, typename Rhs> struct triangular_solve_retval
00236  : public ReturnByValue<triangular_solve_retval<Side, TriangularType, Rhs> >
00237 {
00238   typedef typename remove_all<typename Rhs::Nested>::type RhsNestedCleaned;
00239   typedef ReturnByValue<triangular_solve_retval> Base;
00240   typedef typename Base::Index Index;
00241 
00242   triangular_solve_retval(const TriangularType& tri, const Rhs& rhs)
00243     : m_triangularMatrix(tri), m_rhs(rhs)
00244   {}
00245 
00246   inline Index rows() const { return m_rhs.rows(); }
00247   inline Index cols() const { return m_rhs.cols(); }
00248 
00249   template<typename Dest> inline void evalTo(Dest& dst) const
00250   {
00251     if(!(is_same<RhsNestedCleaned,Dest>::value && extract_data(dst) == extract_data(m_rhs)))
00252       dst = m_rhs;
00253     m_triangularMatrix.template solveInPlace<Side>(dst);
00254   }
00255 
00256   protected:
00257     const TriangularType& m_triangularMatrix;
00258     const typename Rhs::Nested m_rhs;
00259 };
00260 
00261 } // namespace internal
00262 
00263 #endif // EIGEN_SOLVETRIANGULAR_H


re_vision
Author(s): Dorian Galvez-Lopez
autogenerated on Sun Jan 5 2014 11:32:45