diff options
author | Gael Guennebaud <g.gael@free.fr> | 2018-09-21 09:36:05 +0000 |
---|---|---|
committer | Gael Guennebaud <g.gael@free.fr> | 2018-09-21 09:36:05 +0000 |
commit | 3ec29859146f33e323ebd6d41e46acda63dffcef (patch) | |
tree | 4cbd7d59080ac24c0268aa50a50ded50dbeb5a13 | |
parent | 651e5d4866fb82f30e548029c22834b18897c116 (diff) | |
parent | f0ef3467de0b9a726024e7d5efb3aff007ad6f93 (diff) |
Merged indexing cleanup (pull request PR-506)
-rw-r--r-- | Eigen/src/Core/ArithmeticSequence.h | 41 | ||||
-rw-r--r-- | Eigen/src/Core/util/IndexedViewHelper.h | 41 | ||||
-rw-r--r-- | Eigen/src/Core/util/SymbolicIndex.h | 2 | ||||
-rw-r--r-- | Eigen/src/SVD/BDCSVD.h | 10 | ||||
-rw-r--r-- | Eigen/src/SparseCore/SparseMatrix.h | 4 | ||||
-rw-r--r-- | doc/examples/nullary_indexing.cpp | 8 | ||||
-rw-r--r-- | test/indexed_view.cpp | 10 | ||||
-rw-r--r-- | test/symbolic_index.cpp | 20 |
8 files changed, 80 insertions, 56 deletions
diff --git a/Eigen/src/Core/ArithmeticSequence.h b/Eigen/src/Core/ArithmeticSequence.h index 941028893..db6da0001 100644 --- a/Eigen/src/Core/ArithmeticSequence.h +++ b/Eigen/src/Core/ArithmeticSequence.h @@ -329,9 +329,9 @@ seq(const symbolic::BaseExpr<FirstTypeDerived> &f, const symbolic::BaseExpr<Last * \sa lastN(SizeType), seqN(FirstType,SizeType), seq(FirstType,LastType,IncrType) */ template<typename SizeType,typename IncrType> auto lastN(SizeType size, IncrType incr) --> decltype(seqN(Eigen::placeholders::last-(size-fix<1>())*incr, size, incr)) +-> decltype(seqN(Eigen::last-(size-fix<1>())*incr, size, incr)) { - return seqN(Eigen::placeholders::last-(size-fix<1>())*incr, size, incr); + return seqN(Eigen::last-(size-fix<1>())*incr, size, incr); } /** \cpp11 @@ -342,9 +342,9 @@ auto lastN(SizeType size, IncrType incr) * \sa lastN(SizeType,IncrType, seqN(FirstType,SizeType), seq(FirstType,LastType) */ template<typename SizeType> auto lastN(SizeType size) --> decltype(seqN(Eigen::placeholders::last+fix<1>()-size, size)) +-> decltype(seqN(Eigen::last+fix<1>()-size, size)) { - return seqN(Eigen::placeholders::last+fix<1>()-size, size); + return seqN(Eigen::last+fix<1>()-size, size); } #endif @@ -375,6 +375,39 @@ struct get_compile_time_incr<ArithmeticSequence<FirstType,SizeType,IncrType> > { } // end namespace internal +/** \namespace Eigen::indexing + * \ingroup Core_Module + * + * The sole purpose of this namespace is to be able to import all functions + * and symbols that are expected to be used within operator() for indexing + * and slicing. If you already imported the whole Eigen namespace: + * \code using namespace Eigen; \endcode + * then you are already all set. Otherwise, if you don't want/cannot import + * the whole Eigen namespace, the following line: + * \code using namespace Eigen::indexing; \endcode + * is equivalent to: + * \code + using Eigen::all; + using Eigen::seq; + using Eigen::seqN; + using Eigen::lastN; // c++11 only + using Eigen::last; + using Eigen::lastp1; + using Eigen::fix; + \endcode + */ +namespace indexing { + using Eigen::all; + using Eigen::seq; + using Eigen::seqN; + #if EIGEN_HAS_CXX11 + using Eigen::lastN; + #endif + using Eigen::last; + using Eigen::lastp1; + using Eigen::fix; +} + } // end namespace Eigen #endif // EIGEN_ARITHMETIC_SEQUENCE_H diff --git a/Eigen/src/Core/util/IndexedViewHelper.h b/Eigen/src/Core/util/IndexedViewHelper.h index 5d7cf8e91..40e16fdb4 100644 --- a/Eigen/src/Core/util/IndexedViewHelper.h +++ b/Eigen/src/Core/util/IndexedViewHelper.h @@ -13,13 +13,6 @@ namespace Eigen { -/** \namespace Eigen::placeholders - * \ingroup Core_Module - * - * Namespace containing symbolic placeholder and identifiers - */ -namespace placeholders { - namespace internal { struct symbolic_last_tag {}; } @@ -35,36 +28,35 @@ struct symbolic_last_tag {}; * A typical usage example would be: * \code * using namespace Eigen; - * using Eigen::placeholders::last; + * using Eigen::last; * VectorXd v(n); * v(seq(2,last-2)).setOnes(); * \endcode * * \sa end */ -static const symbolic::SymbolExpr<internal::symbolic_last_tag> last; +static const symbolic::SymbolExpr<internal::symbolic_last_tag> last; // PLEASE use Eigen::last instead of Eigen::placeholders::last -/** \var end +/** \var lastp1 * \ingroup Core_Module * - * Can be used as a parameter to Eigen::seq and Eigen::seqN functions to symbolically reference the last+1 element/row/columns - * of the underlying vector or matrix once passed to DenseBase::operator()(const RowIndices&, const ColIndices&). + * Can be used as a parameter to Eigen::seq and Eigen::seqN functions to symbolically + * reference the last+1 element/row/columns of the underlying vector or matrix once + * passed to DenseBase::operator()(const RowIndices&, const ColIndices&). * * This symbolic placeholder support standard arithmetic operation. - * It is essentially an alias to last+1 + * It is essentially an alias to last+fix<1>. * * \sa last */ #ifdef EIGEN_PARSED_BY_DOXYGEN -static const auto end = last+1; +static const auto lastp1 = last+fix<1>; #else // Using a FixedExpr<1> expression is important here to make sure the compiler // can fully optimize the computation starting indices with zero overhead. -static const symbolic::AddExpr<symbolic::SymbolExpr<internal::symbolic_last_tag>,symbolic::ValueExpr<Eigen::internal::FixedInt<1> > > end(last+fix<1>()); +static const symbolic::AddExpr<symbolic::SymbolExpr<internal::symbolic_last_tag>,symbolic::ValueExpr<Eigen::internal::FixedInt<1> > > lastp1(last+fix<1>()); #endif -} // end namespace placeholders - namespace internal { // Replace symbolic last/end "keywords" by their true runtime value @@ -76,7 +68,7 @@ FixedInt<N> eval_expr_given_size(FixedInt<N> x, Index /*size*/) { return x; } template<typename Derived> Index eval_expr_given_size(const symbolic::BaseExpr<Derived> &x, Index size) { - return x.derived().eval(placeholders::last=size-1); + return x.derived().eval(last=size-1); } // Extract increment/step at compile time @@ -172,14 +164,21 @@ template<int Size> struct get_compile_time_incr<AllRange<Size> > { } // end namespace internal -namespace placeholders { - /** \var all * \ingroup Core_Module * Can be used as a parameter to DenseBase::operator()(const RowIndices&, const ColIndices&) to index all rows or columns */ -static const Eigen::internal::all_t all; +static const Eigen::internal::all_t all; // PLEASE use Eigen::all instead of Eigen::placeholders::all + + +namespace placeholders { + typedef symbolic::SymbolExpr<internal::symbolic_last_tag> last_t; + typedef symbolic::AddExpr<symbolic::SymbolExpr<internal::symbolic_last_tag>,symbolic::ValueExpr<Eigen::internal::FixedInt<1> > > end_t; + typedef Eigen::internal::all_t all_t; + EIGEN_DEPRECATED static const all_t all = Eigen::all; // PLEASE use Eigen::all instead of Eigen::placeholders::all + EIGEN_DEPRECATED static const last_t last = Eigen::last; // PLEASE use Eigen::last instead of Eigen::placeholders::last + EIGEN_DEPRECATED static const end_t end = Eigen::lastp1; // PLEASE use Eigen::lastp1 instead of Eigen::placeholders::end } } // end namespace Eigen diff --git a/Eigen/src/Core/util/SymbolicIndex.h b/Eigen/src/Core/util/SymbolicIndex.h index 41d477056..17cf46f05 100644 --- a/Eigen/src/Core/util/SymbolicIndex.h +++ b/Eigen/src/Core/util/SymbolicIndex.h @@ -35,7 +35,7 @@ namespace Eigen { * std::cout << expr98.eval(x=6) << "\n"; * \endcode * - * It is currently only used internally to define and manipulate the placeholders::last and placeholders::end symbols in Eigen::seq and Eigen::seqN. + * It is currently only used internally to define and manipulate the Eigen::last and Eigen::lastp1 symbols in Eigen::seq and Eigen::seqN. * */ namespace symbolic { diff --git a/Eigen/src/SVD/BDCSVD.h b/Eigen/src/SVD/BDCSVD.h index 11df14918..4daa9dd21 100644 --- a/Eigen/src/SVD/BDCSVD.h +++ b/Eigen/src/SVD/BDCSVD.h @@ -945,7 +945,7 @@ void BDCSVD<MatrixType>::perturbCol0 zhat.setZero(); return; } - Index last = perm(m-1); + Index lastIdx = perm(m-1); // The offset permits to skip deflated entries while computing zhat for (Index k = 0; k < n; ++k) { @@ -955,12 +955,12 @@ void BDCSVD<MatrixType>::perturbCol0 { // see equation (3.6) RealScalar dk = diag(k); - RealScalar prod = (singVals(last) + dk) * (mus(last) + (shifts(last) - dk)); + RealScalar prod = (singVals(lastIdx) + dk) * (mus(lastIdx) + (shifts(lastIdx) - dk)); #ifdef EIGEN_BDCSVD_SANITY_CHECKS if(prod<0) { std::cout << "k = " << k << " ; z(k)=" << col0(k) << ", diag(k)=" << dk << "\n"; - std::cout << "prod = " << "(" << singVals(last) << " + " << dk << ") * (" << mus(last) << " + (" << shifts(last) << " - " << dk << "))" << "\n"; - std::cout << " = " << singVals(last) + dk << " * " << mus(last) + (shifts(last) - dk) << "\n"; + std::cout << "prod = " << "(" << singVals(lastIdx) << " + " << dk << ") * (" << mus(lastIdx) << " + (" << shifts(lastIdx) << " - " << dk << "))" << "\n"; + std::cout << " = " << singVals(lastIdx) + dk << " * " << mus(lastIdx) + (shifts(lastIdx) - dk) << "\n"; } assert(prod>=0); #endif @@ -1000,7 +1000,7 @@ void BDCSVD<MatrixType>::perturbCol0 } } #ifdef EIGEN_BDCSVD_DEBUG_VERBOSE - std::cout << "zhat(" << k << ") = sqrt( " << prod << ") ; " << (singVals(last) + dk) << " * " << mus(last) + shifts(last) << " - " << dk << "\n"; + std::cout << "zhat(" << k << ") = sqrt( " << prod << ") ; " << (singVals(lastIdx) + dk) << " * " << mus(lastIdx) + shifts(lastIdx) << " - " << dk << "\n"; #endif RealScalar tmp = sqrt(prod); #ifdef EIGEN_BDCSVD_SANITY_CHECKS diff --git a/Eigen/src/SparseCore/SparseMatrix.h b/Eigen/src/SparseCore/SparseMatrix.h index 8f77194b6..8bfa5f6b8 100644 --- a/Eigen/src/SparseCore/SparseMatrix.h +++ b/Eigen/src/SparseCore/SparseMatrix.h @@ -604,9 +604,9 @@ class SparseMatrix m_outerIndex = newOuterIndex; if (outerChange > 0) { - StorageIndex last = m_outerSize == 0 ? 0 : m_outerIndex[m_outerSize]; + StorageIndex lastIdx = m_outerSize == 0 ? 0 : m_outerIndex[m_outerSize]; for(Index i=m_outerSize; i<m_outerSize+outerChange+1; i++) - m_outerIndex[i] = last; + m_outerIndex[i] = lastIdx; } m_outerSize += outerChange; } diff --git a/doc/examples/nullary_indexing.cpp b/doc/examples/nullary_indexing.cpp index e27c3585a..ca1745628 100644 --- a/doc/examples/nullary_indexing.cpp +++ b/doc/examples/nullary_indexing.cpp @@ -30,7 +30,7 @@ public: // [function] template <class ArgType, class RowIndexType, class ColIndexType> CwiseNullaryOp<indexing_functor<ArgType,RowIndexType,ColIndexType>, typename indexing_functor<ArgType,RowIndexType,ColIndexType>::MatrixType> -indexing(const Eigen::MatrixBase<ArgType>& arg, const RowIndexType& row_indices, const ColIndexType& col_indices) +mat_indexing(const Eigen::MatrixBase<ArgType>& arg, const RowIndexType& row_indices, const ColIndexType& col_indices) { typedef indexing_functor<ArgType,RowIndexType,ColIndexType> Func; typedef typename Func::MatrixType MatrixType; @@ -45,7 +45,7 @@ int main() Eigen::MatrixXi A = Eigen::MatrixXi::Random(4,4); Array3i ri(1,2,1); ArrayXi ci(6); ci << 3,2,1,0,0,2; - Eigen::MatrixXi B = indexing(A, ri, ci); + Eigen::MatrixXi B = mat_indexing(A, ri, ci); std::cout << "A =" << std::endl; std::cout << A << std::endl << std::endl; std::cout << "A([" << ri.transpose() << "], [" << ci.transpose() << "]) =" << std::endl; @@ -53,11 +53,11 @@ int main() std::cout << "[main1]\n"; std::cout << "[main2]\n"; - B = indexing(A, ri+1, ci); + B = mat_indexing(A, ri+1, ci); std::cout << "A(ri+1,ci) =" << std::endl; std::cout << B << std::endl << std::endl; #if __cplusplus >= 201103L - B = indexing(A, ArrayXi::LinSpaced(13,0,12).unaryExpr([](int x){return x%4;}), ArrayXi::LinSpaced(4,0,3)); + B = mat_indexing(A, ArrayXi::LinSpaced(13,0,12).unaryExpr([](int x){return x%4;}), ArrayXi::LinSpaced(4,0,3)); std::cout << "A(ArrayXi::LinSpaced(13,0,12).unaryExpr([](int x){return x%4;}), ArrayXi::LinSpaced(4,0,3)) =" << std::endl; std::cout << B << std::endl << std::endl; #endif diff --git a/test/indexed_view.cpp b/test/indexed_view.cpp index e6ffe5a74..9a284cf0a 100644 --- a/test/indexed_view.cpp +++ b/test/indexed_view.cpp @@ -82,10 +82,6 @@ enum DummyEnum { XX=0, YY=1 }; void check_indexed_view() { - using Eigen::placeholders::all; - using Eigen::placeholders::last; - using Eigen::placeholders::end; - Index n = 10; ArrayXd a = ArrayXd::LinSpaced(n,0,n-1); @@ -239,7 +235,7 @@ void check_indexed_view() VERIFY_IS_APPROX( A(seq(n-1,2,-2), seqN(n-1-6,4)), A(seq(last,2,-2), seqN(last-6,4)) ); VERIFY_IS_APPROX( A(seq(n-1-6,n-1-2), seqN(n-1-6,4)), A(seq(last-6,last-2), seqN(6+last-6-6,4)) ); VERIFY_IS_APPROX( A(seq((n-1)/2,(n)/2+3), seqN(2,4)), A(seq(last/2,(last+1)/2+3), seqN(last+2-last,4)) ); - VERIFY_IS_APPROX( A(seq(n-2,2,-2), seqN(n-8,4)), A(seq(end-2,2,-2), seqN(end-8,4)) ); + VERIFY_IS_APPROX( A(seq(n-2,2,-2), seqN(n-8,4)), A(seq(lastp1-2,2,-2), seqN(lastp1-8,4)) ); // Check all combinations of seq: VERIFY_IS_APPROX( A(seq(1,n-1-2,2), seq(1,n-1-2,2)), A(seq(1,last-2,2), seq(1,last-2,fix<2>)) ); @@ -249,7 +245,7 @@ void check_indexed_view() VERIFY_IS_APPROX( A(seq(n-1-5,n-1-2), seq(n-1-5,n-1-2)), A(seq(last-5,last-2), seq(last-5,last-2)) ); VERIFY_IS_APPROX( A.col(A.cols()-1), A(all,last) ); - VERIFY_IS_APPROX( A(A.rows()-2, A.cols()/2), A(last-1, end/2) ); + VERIFY_IS_APPROX( A(A.rows()-2, A.cols()/2), A(last-1, lastp1/2) ); VERIFY_IS_APPROX( a(a.size()-2), a(last-1) ); VERIFY_IS_APPROX( a(a.size()/2), a((last+1)/2) ); @@ -272,7 +268,7 @@ void check_indexed_view() VERIFY( is_same_eq(a.head(4), a(seq(0,3))) ); VERIFY( is_same_eq(a.tail(4), a(seqN(last-3,4))) ); - VERIFY( is_same_eq(a.tail(4), a(seq(end-4,last))) ); + VERIFY( is_same_eq(a.tail(4), a(seq(lastp1-4,last))) ); VERIFY( is_same_eq(a.segment<4>(3), a(seqN(3,fix<4>))) ); } diff --git a/test/symbolic_index.cpp b/test/symbolic_index.cpp index 02db2eb51..bccf1b884 100644 --- a/test/symbolic_index.cpp +++ b/test/symbolic_index.cpp @@ -54,7 +54,6 @@ is_same_type(const T1&, const T2&) template<typename T1,typename T2> bool is_same_symb(const T1& a, const T2& b, Index size) { - using Eigen::placeholders::last; return a.eval(last=size-1) == b.eval(last=size-1); } @@ -72,14 +71,11 @@ void check_isnot_symbolic(const T&) { void check_symbolic_index() { - using Eigen::placeholders::last; - using Eigen::placeholders::end; - check_is_symbolic(last); - check_is_symbolic(end); + check_is_symbolic(lastp1); check_is_symbolic(last+1); - check_is_symbolic(last-end); - check_is_symbolic(2*last-end/2); + check_is_symbolic(last-lastp1); + check_is_symbolic(2*last-lastp1/2); check_isnot_symbolic(fix<3>()); Index size=100; @@ -93,14 +89,14 @@ void check_symbolic_index() VERIFY( is_same_type( fix<9>()|fix<2>(), fix<9|2>() ) ); VERIFY( is_same_type( fix<9>()/2, int(9/2) ) ); - VERIFY( is_same_symb( end-1, last, size) ); - VERIFY( is_same_symb( end-fix<1>, last, size) ); + VERIFY( is_same_symb( lastp1-1, last, size) ); + VERIFY( is_same_symb( lastp1-fix<1>, last, size) ); VERIFY_IS_EQUAL( ( (last*5-2)/3 ).eval(last=size-1), ((size-1)*5-2)/3 ); VERIFY_IS_EQUAL( ( (last*fix<5>-fix<2>)/fix<3> ).eval(last=size-1), ((size-1)*5-2)/3 ); - VERIFY_IS_EQUAL( ( -last*end ).eval(last=size-1), -(size-1)*size ); - VERIFY_IS_EQUAL( ( end-3*last ).eval(last=size-1), size- 3*(size-1) ); - VERIFY_IS_EQUAL( ( (end-3*last)/end ).eval(last=size-1), (size- 3*(size-1))/size ); + VERIFY_IS_EQUAL( ( -last*lastp1 ).eval(last=size-1), -(size-1)*size ); + VERIFY_IS_EQUAL( ( lastp1-3*last ).eval(last=size-1), size- 3*(size-1) ); + VERIFY_IS_EQUAL( ( (lastp1-3*last)/lastp1 ).eval(last=size-1), (size- 3*(size-1))/size ); #if EIGEN_HAS_CXX14 { |