Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-01-18 09:56:18

0001 // This file is part of Eigen, a lightweight C++ template library
0002 // for linear algebra.
0003 //
0004 // Copyright (C) 2008-2009 Gael Guennebaud <gael.guennebaud@inria.fr>
0005 //
0006 // This Source Code Form is subject to the terms of the Mozilla
0007 // Public License v. 2.0. If a copy of the MPL was not distributed
0008 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
0009 
0010 #ifndef EIGEN_SOLVETRIANGULAR_H
0011 #define EIGEN_SOLVETRIANGULAR_H
0012 
0013 namespace Eigen {
0014 
0015 namespace internal {
0016 
0017 // Forward declarations:
0018 // The following two routines are implemented in the products/TriangularSolver*.h files
0019 template<typename LhsScalar, typename RhsScalar, typename Index, int Side, int Mode, bool Conjugate, int StorageOrder>
0020 struct triangular_solve_vector;
0021 
0022 template <typename Scalar, typename Index, int Side, int Mode, bool Conjugate, int TriStorageOrder, int OtherStorageOrder, int OtherInnerStride>
0023 struct triangular_solve_matrix;
0024 
0025 // small helper struct extracting some traits on the underlying solver operation
0026 template<typename Lhs, typename Rhs, int Side>
0027 class trsolve_traits
0028 {
0029   private:
0030     enum {
0031       RhsIsVectorAtCompileTime = (Side==OnTheLeft ? Rhs::ColsAtCompileTime : Rhs::RowsAtCompileTime)==1
0032     };
0033   public:
0034     enum {
0035       Unrolling   = (RhsIsVectorAtCompileTime && Rhs::SizeAtCompileTime != Dynamic && Rhs::SizeAtCompileTime <= 8)
0036                   ? CompleteUnrolling : NoUnrolling,
0037       RhsVectors  = RhsIsVectorAtCompileTime ? 1 : Dynamic
0038     };
0039 };
0040 
0041 template<typename Lhs, typename Rhs,
0042   int Side, // can be OnTheLeft/OnTheRight
0043   int Mode, // can be Upper/Lower | UnitDiag
0044   int Unrolling = trsolve_traits<Lhs,Rhs,Side>::Unrolling,
0045   int RhsVectors = trsolve_traits<Lhs,Rhs,Side>::RhsVectors
0046   >
0047 struct triangular_solver_selector;
0048 
0049 template<typename Lhs, typename Rhs, int Side, int Mode>
0050 struct triangular_solver_selector<Lhs,Rhs,Side,Mode,NoUnrolling,1>
0051 {
0052   typedef typename Lhs::Scalar LhsScalar;
0053   typedef typename Rhs::Scalar RhsScalar;
0054   typedef blas_traits<Lhs> LhsProductTraits;
0055   typedef typename LhsProductTraits::ExtractType ActualLhsType;
0056   typedef Map<Matrix<RhsScalar,Dynamic,1>, Aligned> MappedRhs;
0057   static EIGEN_DEVICE_FUNC void run(const Lhs& lhs, Rhs& rhs)
0058   {
0059     ActualLhsType actualLhs = LhsProductTraits::extract(lhs);
0060 
0061     // FIXME find a way to allow an inner stride if packet_traits<Scalar>::size==1
0062 
0063     bool useRhsDirectly = Rhs::InnerStrideAtCompileTime==1 || rhs.innerStride()==1;
0064 
0065     ei_declare_aligned_stack_constructed_variable(RhsScalar,actualRhs,rhs.size(),
0066                                                   (useRhsDirectly ? rhs.data() : 0));
0067 
0068     if(!useRhsDirectly)
0069       MappedRhs(actualRhs,rhs.size()) = rhs;
0070 
0071     triangular_solve_vector<LhsScalar, RhsScalar, Index, Side, Mode, LhsProductTraits::NeedToConjugate,
0072                             (int(Lhs::Flags) & RowMajorBit) ? RowMajor : ColMajor>
0073       ::run(actualLhs.cols(), actualLhs.data(), actualLhs.outerStride(), actualRhs);
0074 
0075     if(!useRhsDirectly)
0076       rhs = MappedRhs(actualRhs, rhs.size());
0077   }
0078 };
0079 
0080 // the rhs is a matrix
0081 template<typename Lhs, typename Rhs, int Side, int Mode>
0082 struct triangular_solver_selector<Lhs,Rhs,Side,Mode,NoUnrolling,Dynamic>
0083 {
0084   typedef typename Rhs::Scalar Scalar;
0085   typedef blas_traits<Lhs> LhsProductTraits;
0086   typedef typename LhsProductTraits::DirectLinearAccessType ActualLhsType;
0087 
0088   static EIGEN_DEVICE_FUNC void run(const Lhs& lhs, Rhs& rhs)
0089   {
0090     typename internal::add_const_on_value_type<ActualLhsType>::type actualLhs = LhsProductTraits::extract(lhs);
0091 
0092     const Index size = lhs.rows();
0093     const Index othersize = Side==OnTheLeft? rhs.cols() : rhs.rows();
0094 
0095     typedef internal::gemm_blocking_space<(Rhs::Flags&RowMajorBit) ? RowMajor : ColMajor,Scalar,Scalar,
0096               Rhs::MaxRowsAtCompileTime, Rhs::MaxColsAtCompileTime, Lhs::MaxRowsAtCompileTime,4> BlockingType;
0097 
0098     BlockingType blocking(rhs.rows(), rhs.cols(), size, 1, false);
0099 
0100     triangular_solve_matrix<Scalar,Index,Side,Mode,LhsProductTraits::NeedToConjugate,(int(Lhs::Flags) & RowMajorBit) ? RowMajor : ColMajor,
0101                                (Rhs::Flags&RowMajorBit) ? RowMajor : ColMajor, Rhs::InnerStrideAtCompileTime>
0102       ::run(size, othersize, &actualLhs.coeffRef(0,0), actualLhs.outerStride(), &rhs.coeffRef(0,0), rhs.innerStride(), rhs.outerStride(), blocking);
0103   }
0104 };
0105 
0106 /***************************************************************************
0107 * meta-unrolling implementation
0108 ***************************************************************************/
0109 
0110 template<typename Lhs, typename Rhs, int Mode, int LoopIndex, int Size,
0111          bool Stop = LoopIndex==Size>
0112 struct triangular_solver_unroller;
0113 
0114 template<typename Lhs, typename Rhs, int Mode, int LoopIndex, int Size>
0115 struct triangular_solver_unroller<Lhs,Rhs,Mode,LoopIndex,Size,false> {
0116   enum {
0117     IsLower = ((Mode&Lower)==Lower),
0118     DiagIndex  = IsLower ? LoopIndex : Size - LoopIndex - 1,
0119     StartIndex = IsLower ? 0         : DiagIndex+1
0120   };
0121   static EIGEN_DEVICE_FUNC void run(const Lhs& lhs, Rhs& rhs)
0122   {
0123     if (LoopIndex>0)
0124       rhs.coeffRef(DiagIndex) -= lhs.row(DiagIndex).template segment<LoopIndex>(StartIndex).transpose()
0125                                 .cwiseProduct(rhs.template segment<LoopIndex>(StartIndex)).sum();
0126 
0127     if(!(Mode & UnitDiag))
0128       rhs.coeffRef(DiagIndex) /= lhs.coeff(DiagIndex,DiagIndex);
0129 
0130     triangular_solver_unroller<Lhs,Rhs,Mode,LoopIndex+1,Size>::run(lhs,rhs);
0131   }
0132 };
0133 
0134 template<typename Lhs, typename Rhs, int Mode, int LoopIndex, int Size>
0135 struct triangular_solver_unroller<Lhs,Rhs,Mode,LoopIndex,Size,true> {
0136   static EIGEN_DEVICE_FUNC void run(const Lhs&, Rhs&) {}
0137 };
0138 
0139 template<typename Lhs, typename Rhs, int Mode>
0140 struct triangular_solver_selector<Lhs,Rhs,OnTheLeft,Mode,CompleteUnrolling,1> {
0141   static EIGEN_DEVICE_FUNC void run(const Lhs& lhs, Rhs& rhs)
0142   { triangular_solver_unroller<Lhs,Rhs,Mode,0,Rhs::SizeAtCompileTime>::run(lhs,rhs); }
0143 };
0144 
0145 template<typename Lhs, typename Rhs, int Mode>
0146 struct triangular_solver_selector<Lhs,Rhs,OnTheRight,Mode,CompleteUnrolling,1> {
0147   static EIGEN_DEVICE_FUNC void run(const Lhs& lhs, Rhs& rhs)
0148   {
0149     Transpose<const Lhs> trLhs(lhs);
0150     Transpose<Rhs> trRhs(rhs);
0151 
0152     triangular_solver_unroller<Transpose<const Lhs>,Transpose<Rhs>,
0153                               ((Mode&Upper)==Upper ? Lower : Upper) | (Mode&UnitDiag),
0154                               0,Rhs::SizeAtCompileTime>::run(trLhs,trRhs);
0155   }
0156 };
0157 
0158 } // end namespace internal
0159 
0160 /***************************************************************************
0161 * TriangularView methods
0162 ***************************************************************************/
0163 
0164 #ifndef EIGEN_PARSED_BY_DOXYGEN
0165 template<typename MatrixType, unsigned int Mode>
0166 template<int Side, typename OtherDerived>
0167 EIGEN_DEVICE_FUNC void TriangularViewImpl<MatrixType,Mode,Dense>::solveInPlace(const MatrixBase<OtherDerived>& _other) const
0168 {
0169   OtherDerived& other = _other.const_cast_derived();
0170   eigen_assert( derived().cols() == derived().rows() && ((Side==OnTheLeft && derived().cols() == other.rows()) || (Side==OnTheRight && derived().cols() == other.cols())) );
0171   eigen_assert((!(int(Mode) & int(ZeroDiag))) && bool(int(Mode) & (int(Upper) | int(Lower))));
0172   // If solving for a 0x0 matrix, nothing to do, simply return.
0173   if (derived().cols() == 0)
0174     return;
0175 
0176   enum { copy = (internal::traits<OtherDerived>::Flags & RowMajorBit)  && OtherDerived::IsVectorAtCompileTime && OtherDerived::SizeAtCompileTime!=1};
0177   typedef typename internal::conditional<copy,
0178     typename internal::plain_matrix_type_column_major<OtherDerived>::type, OtherDerived&>::type OtherCopy;
0179   OtherCopy otherCopy(other);
0180 
0181   internal::triangular_solver_selector<MatrixType, typename internal::remove_reference<OtherCopy>::type,
0182     Side, Mode>::run(derived().nestedExpression(), otherCopy);
0183 
0184   if (copy)
0185     other = otherCopy;
0186 }
0187 
0188 template<typename Derived, unsigned int Mode>
0189 template<int Side, typename Other>
0190 const internal::triangular_solve_retval<Side,TriangularView<Derived,Mode>,Other>
0191 TriangularViewImpl<Derived,Mode,Dense>::solve(const MatrixBase<Other>& other) const
0192 {
0193   return internal::triangular_solve_retval<Side,TriangularViewType,Other>(derived(), other.derived());
0194 }
0195 #endif
0196 
0197 namespace internal {
0198 
0199 
0200 template<int Side, typename TriangularType, typename Rhs>
0201 struct traits<triangular_solve_retval<Side, TriangularType, Rhs> >
0202 {
0203   typedef typename internal::plain_matrix_type_column_major<Rhs>::type ReturnType;
0204 };
0205 
0206 template<int Side, typename TriangularType, typename Rhs> struct triangular_solve_retval
0207  : public ReturnByValue<triangular_solve_retval<Side, TriangularType, Rhs> >
0208 {
0209   typedef typename remove_all<typename Rhs::Nested>::type RhsNestedCleaned;
0210   typedef ReturnByValue<triangular_solve_retval> Base;
0211 
0212   triangular_solve_retval(const TriangularType& tri, const Rhs& rhs)
0213     : m_triangularMatrix(tri), m_rhs(rhs)
0214   {}
0215 
0216   inline EIGEN_CONSTEXPR Index rows() const EIGEN_NOEXCEPT { return m_rhs.rows(); }
0217   inline EIGEN_CONSTEXPR Index cols() const EIGEN_NOEXCEPT { return m_rhs.cols(); }
0218 
0219   template<typename Dest> inline void evalTo(Dest& dst) const
0220   {
0221     if(!is_same_dense(dst,m_rhs))
0222       dst = m_rhs;
0223     m_triangularMatrix.template solveInPlace<Side>(dst);
0224   }
0225 
0226   protected:
0227     const TriangularType& m_triangularMatrix;
0228     typename Rhs::Nested m_rhs;
0229 };
0230 
0231 } // namespace internal
0232 
0233 } // end namespace Eigen
0234 
0235 #endif // EIGEN_SOLVETRIANGULAR_H