aboutsummaryrefslogtreecommitdiffhomepage
diff options
context:
space:
mode:
-rw-r--r--Eigen/src/StlSupport/StdDeque.h18
-rw-r--r--Eigen/src/StlSupport/StdList.h18
-rw-r--r--test/CMakeLists.txt2
-rw-r--r--test/stddeque_overload.cpp159
-rw-r--r--test/stdlist_overload.cpp192
5 files changed, 363 insertions, 26 deletions
diff --git a/Eigen/src/StlSupport/StdDeque.h b/Eigen/src/StlSupport/StdDeque.h
index 25930cb85..cf1fedf92 100644
--- a/Eigen/src/StlSupport/StdDeque.h
+++ b/Eigen/src/StlSupport/StdDeque.h
@@ -13,32 +13,24 @@
#include "details.h"
-// Define the explicit instantiation (e.g. necessary for the Intel compiler)
-#if EIGEN_COMP_GNUC || EIGEN_COMP_ICC
- #define EIGEN_EXPLICIT_STL_DEQUE_INSTANTIATION(...) template class std::deque<__VA_ARGS__, EIGEN_ALIGNED_ALLOCATOR<__VA_ARGS__> >;
-#else
- #define EIGEN_EXPLICIT_STL_DEQUE_INSTANTIATION(...)
-#endif
-
/**
* This section contains a convenience MACRO which allows an easy specialization of
* std::deque such that for data types with alignment issues the correct allocator
* is used automatically.
*/
#define EIGEN_DEFINE_STL_DEQUE_SPECIALIZATION(...) \
-EIGEN_EXPLICIT_STL_DEQUE_INSTANTIATION(__VA_ARGS__) \
namespace std \
{ \
- template<typename _Ay> \
- class deque<__VA_ARGS__, _Ay> \
+ template<> \
+ class deque<__VA_ARGS__, std::allocator<__VA_ARGS__> > \
: public deque<__VA_ARGS__, EIGEN_ALIGNED_ALLOCATOR<__VA_ARGS__> > \
{ \
typedef deque<__VA_ARGS__, EIGEN_ALIGNED_ALLOCATOR<__VA_ARGS__> > deque_base; \
public: \
typedef __VA_ARGS__ value_type; \
- typedef typename deque_base::allocator_type allocator_type; \
- typedef typename deque_base::size_type size_type; \
- typedef typename deque_base::iterator iterator; \
+ typedef deque_base::allocator_type allocator_type; \
+ typedef deque_base::size_type size_type; \
+ typedef deque_base::iterator iterator; \
explicit deque(const allocator_type& a = allocator_type()) : deque_base(a) {} \
template<typename InputIterator> \
deque(InputIterator first, InputIterator last, const allocator_type& a = allocator_type()) : deque_base(first, last, a) {} \
diff --git a/Eigen/src/StlSupport/StdList.h b/Eigen/src/StlSupport/StdList.h
index 7412b50aa..e1eba4985 100644
--- a/Eigen/src/StlSupport/StdList.h
+++ b/Eigen/src/StlSupport/StdList.h
@@ -12,32 +12,24 @@
#include "details.h"
-// Define the explicit instantiation (e.g. necessary for the Intel compiler)
-#if EIGEN_COMP_GNUC || EIGEN_COMP_ICC
- #define EIGEN_EXPLICIT_STL_LIST_INSTANTIATION(...) template class std::list<__VA_ARGS__, EIGEN_ALIGNED_ALLOCATOR<__VA_ARGS__> >;
-#else
- #define EIGEN_EXPLICIT_STL_LIST_INSTANTIATION(...)
-#endif
-
/**
* This section contains a convenience MACRO which allows an easy specialization of
* std::list such that for data types with alignment issues the correct allocator
* is used automatically.
*/
#define EIGEN_DEFINE_STL_LIST_SPECIALIZATION(...) \
-EIGEN_EXPLICIT_STL_LIST_INSTANTIATION(__VA_ARGS__) \
namespace std \
{ \
- template<typename _Ay> \
- class list<__VA_ARGS__, _Ay> \
+ template<> \
+ class list<__VA_ARGS__, std::allocator<__VA_ARGS__> > \
: public list<__VA_ARGS__, EIGEN_ALIGNED_ALLOCATOR<__VA_ARGS__> > \
{ \
typedef list<__VA_ARGS__, EIGEN_ALIGNED_ALLOCATOR<__VA_ARGS__> > list_base; \
public: \
typedef __VA_ARGS__ value_type; \
- typedef typename list_base::allocator_type allocator_type; \
- typedef typename list_base::size_type size_type; \
- typedef typename list_base::iterator iterator; \
+ typedef list_base::allocator_type allocator_type; \
+ typedef list_base::size_type size_type; \
+ typedef list_base::iterator iterator; \
explicit list(const allocator_type& a = allocator_type()) : list_base(a) {} \
template<typename InputIterator> \
list(InputIterator first, InputIterator last, const allocator_type& a = allocator_type()) : list_base(first, last, a) {} \
diff --git a/test/CMakeLists.txt b/test/CMakeLists.txt
index bbebf29cd..4420e0c51 100644
--- a/test/CMakeLists.txt
+++ b/test/CMakeLists.txt
@@ -226,7 +226,9 @@ ei_add_test(geo_homogeneous)
ei_add_test(stdvector)
ei_add_test(stdvector_overload)
ei_add_test(stdlist)
+ei_add_test(stdlist_overload)
ei_add_test(stddeque)
+ei_add_test(stddeque_overload)
ei_add_test(sparse_basic)
ei_add_test(sparse_block)
ei_add_test(sparse_vector)
diff --git a/test/stddeque_overload.cpp b/test/stddeque_overload.cpp
new file mode 100644
index 000000000..d887e35ba
--- /dev/null
+++ b/test/stddeque_overload.cpp
@@ -0,0 +1,159 @@
+// This file is part of Eigen, a lightweight C++ template library
+// for linear algebra.
+//
+// Copyright (C) 2008 Benoit Jacob <jacob.benoit.1@gmail.com>
+// Copyright (C) 2010 Hauke Heibel <hauke.heibel@gmail.com>
+//
+// 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/.
+
+#include "main.h"
+
+#include <Eigen/StdDeque>
+#include <Eigen/Geometry>
+
+EIGEN_DEFINE_STL_DEQUE_SPECIALIZATION(Vector4f)
+
+EIGEN_DEFINE_STL_DEQUE_SPECIALIZATION(Matrix2f)
+EIGEN_DEFINE_STL_DEQUE_SPECIALIZATION(Matrix4f)
+EIGEN_DEFINE_STL_DEQUE_SPECIALIZATION(Matrix4d)
+
+EIGEN_DEFINE_STL_DEQUE_SPECIALIZATION(Affine3f)
+EIGEN_DEFINE_STL_DEQUE_SPECIALIZATION(Affine3d)
+
+EIGEN_DEFINE_STL_DEQUE_SPECIALIZATION(Quaternionf)
+EIGEN_DEFINE_STL_DEQUE_SPECIALIZATION(Quaterniond)
+
+template<typename MatrixType>
+void check_stddeque_matrix(const MatrixType& m)
+{
+ typename MatrixType::Index rows = m.rows();
+ typename MatrixType::Index cols = m.cols();
+ MatrixType x = MatrixType::Random(rows,cols), y = MatrixType::Random(rows,cols);
+ std::deque<MatrixType> v(10, MatrixType(rows,cols)), w(20, y);
+ v[5] = x;
+ w[6] = v[5];
+ VERIFY_IS_APPROX(w[6], v[5]);
+ v = w;
+ for(int i = 0; i < 20; i++)
+ {
+ VERIFY_IS_APPROX(w[i], v[i]);
+ }
+
+ v.resize(21);
+ v[20] = x;
+ VERIFY_IS_APPROX(v[20], x);
+ v.resize(22,y);
+ VERIFY_IS_APPROX(v[21], y);
+ v.push_back(x);
+ VERIFY_IS_APPROX(v[22], x);
+ VERIFY((size_t)&(v[22]) == (size_t)&(v[21]) + sizeof(MatrixType));
+
+ // do a lot of push_back such that the deque gets internally resized
+ // (with memory reallocation)
+ MatrixType* ref = &w[0];
+ for(int i=0; i<30 || ((ref==&w[0]) && i<300); ++i)
+ v.push_back(w[i%w.size()]);
+ for(unsigned int i=23; i<v.size(); ++i)
+ {
+ VERIFY(v[i]==w[(i-23)%w.size()]);
+ }
+}
+
+template<typename TransformType>
+void check_stddeque_transform(const TransformType&)
+{
+ typedef typename TransformType::MatrixType MatrixType;
+ TransformType x(MatrixType::Random()), y(MatrixType::Random());
+ std::deque<TransformType> v(10), w(20, y);
+ v[5] = x;
+ w[6] = v[5];
+ VERIFY_IS_APPROX(w[6], v[5]);
+ v = w;
+ for(int i = 0; i < 20; i++)
+ {
+ VERIFY_IS_APPROX(w[i], v[i]);
+ }
+
+ v.resize(21);
+ v[20] = x;
+ VERIFY_IS_APPROX(v[20], x);
+ v.resize(22,y);
+ VERIFY_IS_APPROX(v[21], y);
+ v.push_back(x);
+ VERIFY_IS_APPROX(v[22], x);
+
+ // do a lot of push_back such that the deque gets internally resized
+ // (with memory reallocation)
+ TransformType* ref = &w[0];
+ for(int i=0; i<30 || ((ref==&w[0]) && i<300); ++i)
+ v.push_back(w[i%w.size()]);
+ for(unsigned int i=23; i<v.size(); ++i)
+ {
+ VERIFY(v[i].matrix()==w[(i-23)%w.size()].matrix());
+ }
+}
+
+template<typename QuaternionType>
+void check_stddeque_quaternion(const QuaternionType&)
+{
+ typedef typename QuaternionType::Coefficients Coefficients;
+ QuaternionType x(Coefficients::Random()), y(Coefficients::Random());
+ std::deque<QuaternionType> v(10), w(20, y);
+ v[5] = x;
+ w[6] = v[5];
+ VERIFY_IS_APPROX(w[6], v[5]);
+ v = w;
+ for(int i = 0; i < 20; i++)
+ {
+ VERIFY_IS_APPROX(w[i], v[i]);
+ }
+
+ v.resize(21);
+ v[20] = x;
+ VERIFY_IS_APPROX(v[20], x);
+ v.resize(22,y);
+ VERIFY_IS_APPROX(v[21], y);
+ v.push_back(x);
+ VERIFY_IS_APPROX(v[22], x);
+
+ // do a lot of push_back such that the deque gets internally resized
+ // (with memory reallocation)
+ QuaternionType* ref = &w[0];
+ for(int i=0; i<30 || ((ref==&w[0]) && i<300); ++i)
+ v.push_back(w[i%w.size()]);
+ for(unsigned int i=23; i<v.size(); ++i)
+ {
+ VERIFY(v[i].coeffs()==w[(i-23)%w.size()].coeffs());
+ }
+}
+
+void test_stddeque_overload()
+{
+ // some non vectorizable fixed sizes
+ CALL_SUBTEST_1(check_stddeque_matrix(Vector2f()));
+ CALL_SUBTEST_1(check_stddeque_matrix(Matrix3f()));
+ CALL_SUBTEST_2(check_stddeque_matrix(Matrix3d()));
+
+ // some vectorizable fixed sizes
+ CALL_SUBTEST_1(check_stddeque_matrix(Matrix2f()));
+ CALL_SUBTEST_1(check_stddeque_matrix(Vector4f()));
+ CALL_SUBTEST_1(check_stddeque_matrix(Matrix4f()));
+ CALL_SUBTEST_2(check_stddeque_matrix(Matrix4d()));
+
+ // some dynamic sizes
+ CALL_SUBTEST_3(check_stddeque_matrix(MatrixXd(1,1)));
+ CALL_SUBTEST_3(check_stddeque_matrix(VectorXd(20)));
+ CALL_SUBTEST_3(check_stddeque_matrix(RowVectorXf(20)));
+ CALL_SUBTEST_3(check_stddeque_matrix(MatrixXcf(10,10)));
+
+ // some Transform
+ CALL_SUBTEST_4(check_stddeque_transform(Affine2f())); // does not need the specialization (2+1)^2 = 9
+ CALL_SUBTEST_4(check_stddeque_transform(Affine3f()));
+ CALL_SUBTEST_4(check_stddeque_transform(Affine3d()));
+
+ // some Quaternion
+ CALL_SUBTEST_5(check_stddeque_quaternion(Quaternionf()));
+ CALL_SUBTEST_5(check_stddeque_quaternion(Quaterniond()));
+}
diff --git a/test/stdlist_overload.cpp b/test/stdlist_overload.cpp
new file mode 100644
index 000000000..bb910bd43
--- /dev/null
+++ b/test/stdlist_overload.cpp
@@ -0,0 +1,192 @@
+// This file is part of Eigen, a lightweight C++ template library
+// for linear algebra.
+//
+// Copyright (C) 2008 Benoit Jacob <jacob.benoit.1@gmail.com>
+// Copyright (C) 2010 Hauke Heibel <hauke.heibel@gmail.com>
+//
+// 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/.
+
+#include "main.h"
+
+#include <Eigen/StdList>
+#include <Eigen/Geometry>
+
+EIGEN_DEFINE_STL_LIST_SPECIALIZATION(Vector4f)
+
+EIGEN_DEFINE_STL_LIST_SPECIALIZATION(Matrix2f)
+EIGEN_DEFINE_STL_LIST_SPECIALIZATION(Matrix4f)
+EIGEN_DEFINE_STL_LIST_SPECIALIZATION(Matrix4d)
+
+EIGEN_DEFINE_STL_LIST_SPECIALIZATION(Affine3f)
+EIGEN_DEFINE_STL_LIST_SPECIALIZATION(Affine3d)
+
+EIGEN_DEFINE_STL_LIST_SPECIALIZATION(Quaternionf)
+EIGEN_DEFINE_STL_LIST_SPECIALIZATION(Quaterniond)
+
+template <class Container, class Position>
+typename Container::iterator get(Container & c, Position position)
+{
+ typename Container::iterator it = c.begin();
+ std::advance(it, position);
+ return it;
+}
+
+template <class Container, class Position, class Value>
+void set(Container & c, Position position, const Value & value)
+{
+ typename Container::iterator it = c.begin();
+ std::advance(it, position);
+ *it = value;
+}
+
+template<typename MatrixType>
+void check_stdlist_matrix(const MatrixType& m)
+{
+ typename MatrixType::Index rows = m.rows();
+ typename MatrixType::Index cols = m.cols();
+ MatrixType x = MatrixType::Random(rows,cols), y = MatrixType::Random(rows,cols);
+ std::list<MatrixType> v(10, MatrixType(rows,cols)), w(20, y);
+ typename std::list<MatrixType>::iterator itv = get(v, 5);
+ typename std::list<MatrixType>::iterator itw = get(w, 6);
+ *itv = x;
+ *itw = *itv;
+ VERIFY_IS_APPROX(*itw, *itv);
+ v = w;
+ itv = v.begin();
+ itw = w.begin();
+ for(int i = 0; i < 20; i++)
+ {
+ VERIFY_IS_APPROX(*itw, *itv);
+ ++itv;
+ ++itw;
+ }
+
+ v.resize(21);
+ set(v, 20, x);
+ VERIFY_IS_APPROX(*get(v, 20), x);
+ v.resize(22,y);
+ VERIFY_IS_APPROX(*get(v, 21), y);
+ v.push_back(x);
+ VERIFY_IS_APPROX(*get(v, 22), x);
+
+ // do a lot of push_back such that the list gets internally resized
+ // (with memory reallocation)
+ MatrixType* ref = &(*get(w, 0));
+ for(int i=0; i<30 || ((ref==&(*get(w, 0))) && i<300); ++i)
+ v.push_back(*get(w, i%w.size()));
+ for(unsigned int i=23; i<v.size(); ++i)
+ {
+ VERIFY((*get(v, i))==(*get(w, (i-23)%w.size())));
+ }
+}
+
+template<typename TransformType>
+void check_stdlist_transform(const TransformType&)
+{
+ typedef typename TransformType::MatrixType MatrixType;
+ TransformType x(MatrixType::Random()), y(MatrixType::Random());
+ std::list<TransformType> v(10), w(20, y);
+ typename std::list<TransformType>::iterator itv = get(v, 5);
+ typename std::list<TransformType>::iterator itw = get(w, 6);
+ *itv = x;
+ *itw = *itv;
+ VERIFY_IS_APPROX(*itw, *itv);
+ v = w;
+ itv = v.begin();
+ itw = w.begin();
+ for(int i = 0; i < 20; i++)
+ {
+ VERIFY_IS_APPROX(*itw, *itv);
+ ++itv;
+ ++itw;
+ }
+
+ v.resize(21);
+ set(v, 20, x);
+ VERIFY_IS_APPROX(*get(v, 20), x);
+ v.resize(22,y);
+ VERIFY_IS_APPROX(*get(v, 21), y);
+ v.push_back(x);
+ VERIFY_IS_APPROX(*get(v, 22), x);
+
+ // do a lot of push_back such that the list gets internally resized
+ // (with memory reallocation)
+ TransformType* ref = &(*get(w, 0));
+ for(int i=0; i<30 || ((ref==&(*get(w, 0))) && i<300); ++i)
+ v.push_back(*get(w, i%w.size()));
+ for(unsigned int i=23; i<v.size(); ++i)
+ {
+ VERIFY(get(v, i)->matrix()==get(w, (i-23)%w.size())->matrix());
+ }
+}
+
+template<typename QuaternionType>
+void check_stdlist_quaternion(const QuaternionType&)
+{
+ typedef typename QuaternionType::Coefficients Coefficients;
+ QuaternionType x(Coefficients::Random()), y(Coefficients::Random());
+ std::list<QuaternionType> v(10), w(20, y);
+ typename std::list<QuaternionType>::iterator itv = get(v, 5);
+ typename std::list<QuaternionType>::iterator itw = get(w, 6);
+ *itv = x;
+ *itw = *itv;
+ VERIFY_IS_APPROX(*itw, *itv);
+ v = w;
+ itv = v.begin();
+ itw = w.begin();
+ for(int i = 0; i < 20; i++)
+ {
+ VERIFY_IS_APPROX(*itw, *itv);
+ ++itv;
+ ++itw;
+ }
+
+ v.resize(21);
+ set(v, 20, x);
+ VERIFY_IS_APPROX(*get(v, 20), x);
+ v.resize(22,y);
+ VERIFY_IS_APPROX(*get(v, 21), y);
+ v.push_back(x);
+ VERIFY_IS_APPROX(*get(v, 22), x);
+
+ // do a lot of push_back such that the list gets internally resized
+ // (with memory reallocation)
+ QuaternionType* ref = &(*get(w, 0));
+ for(int i=0; i<30 || ((ref==&(*get(w, 0))) && i<300); ++i)
+ v.push_back(*get(w, i%w.size()));
+ for(unsigned int i=23; i<v.size(); ++i)
+ {
+ VERIFY(get(v, i)->coeffs()==get(w, (i-23)%w.size())->coeffs());
+ }
+}
+
+void test_stdlist_overload()
+{
+ // some non vectorizable fixed sizes
+ CALL_SUBTEST_1(check_stdlist_matrix(Vector2f()));
+ CALL_SUBTEST_1(check_stdlist_matrix(Matrix3f()));
+ CALL_SUBTEST_2(check_stdlist_matrix(Matrix3d()));
+
+ // some vectorizable fixed sizes
+ CALL_SUBTEST_1(check_stdlist_matrix(Matrix2f()));
+ CALL_SUBTEST_1(check_stdlist_matrix(Vector4f()));
+ CALL_SUBTEST_1(check_stdlist_matrix(Matrix4f()));
+ CALL_SUBTEST_2(check_stdlist_matrix(Matrix4d()));
+
+ // some dynamic sizes
+ CALL_SUBTEST_3(check_stdlist_matrix(MatrixXd(1,1)));
+ CALL_SUBTEST_3(check_stdlist_matrix(VectorXd(20)));
+ CALL_SUBTEST_3(check_stdlist_matrix(RowVectorXf(20)));
+ CALL_SUBTEST_3(check_stdlist_matrix(MatrixXcf(10,10)));
+
+ // some Transform
+ CALL_SUBTEST_4(check_stdlist_transform(Affine2f())); // does not need the specialization (2+1)^2 = 9
+ CALL_SUBTEST_4(check_stdlist_transform(Affine3f()));
+ CALL_SUBTEST_4(check_stdlist_transform(Affine3d()));
+
+ // some Quaternion
+ CALL_SUBTEST_5(check_stdlist_quaternion(Quaternionf()));
+ CALL_SUBTEST_5(check_stdlist_quaternion(Quaterniond()));
+}