// This file is part of Eigen, a lightweight C++ template library // for linear algebra. // // Copyright (C) 2009 Gael Guennebaud // // This Source Code Form is subject to the terms of the Mozilla // Public License v. 2.0. If a copy of the MPL was not distributed // with this file, You can obtain one at http://mozilla.org/MPL/2.0/. #ifndef EIGEN_TRIANGULARMATRIXVECTOR_H #define EIGEN_TRIANGULARMATRIXVECTOR_H namespace Eigen { namespace internal { template struct triangular_matrix_vector_product; template struct triangular_matrix_vector_product { typedef typename scalar_product_traits::ReturnType ResScalar; enum { IsLower = ((Mode&Lower)==Lower), HasUnitDiag = (Mode & UnitDiag)==UnitDiag, HasZeroDiag = (Mode & ZeroDiag)==ZeroDiag }; static EIGEN_DONT_INLINE void run(Index _rows, Index _cols, const LhsScalar* _lhs, Index lhsStride, const RhsScalar* _rhs, Index rhsIncr, ResScalar* _res, Index resIncr, const ResScalar& alpha); }; template EIGEN_DONT_INLINE void triangular_matrix_vector_product ::run(Index _rows, Index _cols, const LhsScalar* _lhs, Index lhsStride, const RhsScalar* _rhs, Index rhsIncr, ResScalar* _res, Index resIncr, const ResScalar& alpha) { static const Index PanelWidth = EIGEN_TUNE_TRIANGULAR_PANEL_WIDTH; Index size = (std::min)(_rows,_cols); Index rows = IsLower ? _rows : (std::min)(_rows,_cols); Index cols = IsLower ? (std::min)(_rows,_cols) : _cols; typedef Map, 0, OuterStride<> > LhsMap; const LhsMap lhs(_lhs,rows,cols,OuterStride<>(lhsStride)); typename conj_expr_if::type cjLhs(lhs); typedef Map, 0, InnerStride<> > RhsMap; const RhsMap rhs(_rhs,cols,InnerStride<>(rhsIncr)); typename conj_expr_if::type cjRhs(rhs); typedef Map > ResMap; ResMap res(_res,rows); for (Index pi=0; pi0) res.segment(s,r) += (alpha * cjRhs.coeff(i)) * cjLhs.col(i).segment(s,r); if (HasUnitDiag) res.coeffRef(i) += alpha * cjRhs.coeff(i); } Index r = IsLower ? rows - pi - actualPanelWidth : pi; if (r>0) { Index s = IsLower ? pi+actualPanelWidth : 0; general_matrix_vector_product::run( r, actualPanelWidth, &lhs.coeffRef(s,pi), lhsStride, &rhs.coeffRef(pi), rhsIncr, &res.coeffRef(s), resIncr, alpha); } } if((!IsLower) && cols>size) { general_matrix_vector_product::run( rows, cols-size, &lhs.coeffRef(0,size), lhsStride, &rhs.coeffRef(size), rhsIncr, _res, resIncr, alpha); } } template struct triangular_matrix_vector_product { typedef typename scalar_product_traits::ReturnType ResScalar; enum { IsLower = ((Mode&Lower)==Lower), HasUnitDiag = (Mode & UnitDiag)==UnitDiag, HasZeroDiag = (Mode & ZeroDiag)==ZeroDiag }; static EIGEN_DONT_INLINE void run(Index _rows, Index _cols, const LhsScalar* _lhs, Index lhsStride, const RhsScalar* _rhs, Index rhsIncr, ResScalar* _res, Index resIncr, const ResScalar& alpha); }; template EIGEN_DONT_INLINE void triangular_matrix_vector_product ::run(Index _rows, Index _cols, const LhsScalar* _lhs, Index lhsStride, const RhsScalar* _rhs, Index rhsIncr, ResScalar* _res, Index resIncr, const ResScalar& alpha) { static const Index PanelWidth = EIGEN_TUNE_TRIANGULAR_PANEL_WIDTH; Index diagSize = (std::min)(_rows,_cols); Index rows = IsLower ? _rows : diagSize; Index cols = IsLower ? diagSize : _cols; typedef Map, 0, OuterStride<> > LhsMap; const LhsMap lhs(_lhs,rows,cols,OuterStride<>(lhsStride)); typename conj_expr_if::type cjLhs(lhs); typedef Map > RhsMap; const RhsMap rhs(_rhs,cols); typename conj_expr_if::type cjRhs(rhs); typedef Map, 0, InnerStride<> > ResMap; ResMap res(_res,rows,InnerStride<>(resIncr)); for (Index pi=0; pi0) res.coeffRef(i) += alpha * (cjLhs.row(i).segment(s,r).cwiseProduct(cjRhs.segment(s,r).transpose())).sum(); if (HasUnitDiag) res.coeffRef(i) += alpha * cjRhs.coeff(i); } Index r = IsLower ? pi : cols - pi - actualPanelWidth; if (r>0) { Index s = IsLower ? 0 : pi + actualPanelWidth; general_matrix_vector_product::run( actualPanelWidth, r, &lhs.coeffRef(pi,s), lhsStride, &rhs.coeffRef(s), rhsIncr, &res.coeffRef(pi), resIncr, alpha); } } if(IsLower && rows>diagSize) { general_matrix_vector_product::run( rows-diagSize, cols, &lhs.coeffRef(diagSize,0), lhsStride, &rhs.coeffRef(0), rhsIncr, &res.coeffRef(diagSize), resIncr, alpha); } } /*************************************************************************** * Wrapper to product_triangular_vector ***************************************************************************/ template struct traits > : traits, Lhs, Rhs> > {}; template struct traits > : traits, Lhs, Rhs> > {}; template struct trmv_selector; } // end namespace internal template struct TriangularProduct : public ProductBase, Lhs, Rhs > { EIGEN_PRODUCT_PUBLIC_INTERFACE(TriangularProduct) TriangularProduct(const Lhs& lhs, const Rhs& rhs) : Base(lhs,rhs) {} template void scaleAndAddTo(Dest& dst, const Scalar& alpha) const { eigen_assert(dst.rows()==m_lhs.rows() && dst.cols()==m_rhs.cols()); internal::trmv_selector<(int(internal::traits::Flags)&RowMajorBit) ? RowMajor : ColMajor>::run(*this, dst, alpha); } }; template struct TriangularProduct : public ProductBase, Lhs, Rhs > { EIGEN_PRODUCT_PUBLIC_INTERFACE(TriangularProduct) TriangularProduct(const Lhs& lhs, const Rhs& rhs) : Base(lhs,rhs) {} template void scaleAndAddTo(Dest& dst, const Scalar& alpha) const { eigen_assert(dst.rows()==m_lhs.rows() && dst.cols()==m_rhs.cols()); typedef TriangularProduct<(Mode & (UnitDiag|ZeroDiag)) | ((Mode & Lower) ? Upper : Lower),true,Transpose,false,Transpose,true> TriangularProductTranspose; Transpose dstT(dst); internal::trmv_selector<(int(internal::traits::Flags)&RowMajorBit) ? ColMajor : RowMajor>::run( TriangularProductTranspose(m_rhs.transpose(),m_lhs.transpose()), dstT, alpha); } }; namespace internal { // TODO: find a way to factorize this piece of code with gemv_selector since the logic is exactly the same. template<> struct trmv_selector { template static void run(const TriangularProduct& prod, Dest& dest, const typename TriangularProduct::Scalar& alpha) { typedef TriangularProduct ProductType; typedef typename ProductType::Index Index; typedef typename ProductType::LhsScalar LhsScalar; typedef typename ProductType::RhsScalar RhsScalar; typedef typename ProductType::Scalar ResScalar; typedef typename ProductType::RealScalar RealScalar; typedef typename ProductType::ActualLhsType ActualLhsType; typedef typename ProductType::ActualRhsType ActualRhsType; typedef typename ProductType::LhsBlasTraits LhsBlasTraits; typedef typename ProductType::RhsBlasTraits RhsBlasTraits; typedef Map, Aligned> MappedDest; typename internal::add_const_on_value_type::type actualLhs = LhsBlasTraits::extract(prod.lhs()); typename internal::add_const_on_value_type::type actualRhs = RhsBlasTraits::extract(prod.rhs()); ResScalar actualAlpha = alpha * LhsBlasTraits::extractScalarFactor(prod.lhs()) * RhsBlasTraits::extractScalarFactor(prod.rhs()); enum { // FIXME find a way to allow an inner stride on the result if packet_traits::size==1 // on, the other hand it is good for the cache to pack the vector anyways... EvalToDestAtCompileTime = Dest::InnerStrideAtCompileTime==1, ComplexByReal = (NumTraits::IsComplex) && (!NumTraits::IsComplex), MightCannotUseDest = (Dest::InnerStrideAtCompileTime!=1) || ComplexByReal }; gemv_static_vector_if static_dest; bool alphaIsCompatible = (!ComplexByReal) || (numext::imag(actualAlpha)==RealScalar(0)); bool evalToDest = EvalToDestAtCompileTime && alphaIsCompatible; RhsScalar compatibleAlpha = get_factor::run(actualAlpha); ei_declare_aligned_stack_constructed_variable(ResScalar,actualDestPtr,dest.size(), evalToDest ? dest.data() : static_dest.data()); if(!evalToDest) { #ifdef EIGEN_DENSE_STORAGE_CTOR_PLUGIN Index size = dest.size(); EIGEN_DENSE_STORAGE_CTOR_PLUGIN #endif if(!alphaIsCompatible) { MappedDest(actualDestPtr, dest.size()).setZero(); compatibleAlpha = RhsScalar(1); } else MappedDest(actualDestPtr, dest.size()) = dest; } internal::triangular_matrix_vector_product ::run(actualLhs.rows(),actualLhs.cols(), actualLhs.data(),actualLhs.outerStride(), actualRhs.data(),actualRhs.innerStride(), actualDestPtr,1,compatibleAlpha); if (!evalToDest) { if(!alphaIsCompatible) dest += actualAlpha * MappedDest(actualDestPtr, dest.size()); else dest = MappedDest(actualDestPtr, dest.size()); } } }; template<> struct trmv_selector { template static void run(const TriangularProduct& prod, Dest& dest, const typename TriangularProduct::Scalar& alpha) { typedef TriangularProduct ProductType; typedef typename ProductType::LhsScalar LhsScalar; typedef typename ProductType::RhsScalar RhsScalar; typedef typename ProductType::Scalar ResScalar; typedef typename ProductType::Index Index; typedef typename ProductType::ActualLhsType ActualLhsType; typedef typename ProductType::ActualRhsType ActualRhsType; typedef typename ProductType::_ActualRhsType _ActualRhsType; typedef typename ProductType::LhsBlasTraits LhsBlasTraits; typedef typename ProductType::RhsBlasTraits RhsBlasTraits; typename add_const::type actualLhs = LhsBlasTraits::extract(prod.lhs()); typename add_const::type actualRhs = RhsBlasTraits::extract(prod.rhs()); ResScalar actualAlpha = alpha * LhsBlasTraits::extractScalarFactor(prod.lhs()) * RhsBlasTraits::extractScalarFactor(prod.rhs()); enum { DirectlyUseRhs = _ActualRhsType::InnerStrideAtCompileTime==1 }; gemv_static_vector_if static_rhs; ei_declare_aligned_stack_constructed_variable(RhsScalar,actualRhsPtr,actualRhs.size(), DirectlyUseRhs ? const_cast(actualRhs.data()) : static_rhs.data()); if(!DirectlyUseRhs) { #ifdef EIGEN_DENSE_STORAGE_CTOR_PLUGIN int size = actualRhs.size(); EIGEN_DENSE_STORAGE_CTOR_PLUGIN #endif Map(actualRhsPtr, actualRhs.size()) = actualRhs; } internal::triangular_matrix_vector_product ::run(actualLhs.rows(),actualLhs.cols(), actualLhs.data(),actualLhs.outerStride(), actualRhsPtr,1, dest.data(),dest.innerStride(), actualAlpha); } }; } // end namespace internal } // end namespace Eigen #endif // EIGEN_TRIANGULARMATRIXVECTOR_H