1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2008 Benoit Jacob <jacob.benoit.1@gmail.com>
5 // Copyright (C) 2010 Hauke Heibel <hauke.heibel@gmail.com>
6 //
7 // This Source Code Form is subject to the terms of the Mozilla
8 // Public License v. 2.0. If a copy of the MPL was not distributed
9 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
10
11 #include "main.h"
12
13 #include <Eigen/StdVector>
14 #include <Eigen/Geometry>
15
16 EIGEN_DEFINE_STL_VECTOR_SPECIALIZATION(Vector4f)
17
EIGEN_DEFINE_STL_VECTOR_SPECIALIZATION(Matrix2f)18 EIGEN_DEFINE_STL_VECTOR_SPECIALIZATION(Matrix2f)
19 EIGEN_DEFINE_STL_VECTOR_SPECIALIZATION(Matrix4f)
20 EIGEN_DEFINE_STL_VECTOR_SPECIALIZATION(Matrix4d)
21
22 EIGEN_DEFINE_STL_VECTOR_SPECIALIZATION(Affine3f)
23 EIGEN_DEFINE_STL_VECTOR_SPECIALIZATION(Affine3d)
24
25 EIGEN_DEFINE_STL_VECTOR_SPECIALIZATION(Quaternionf)
26 EIGEN_DEFINE_STL_VECTOR_SPECIALIZATION(Quaterniond)
27
28 template<typename MatrixType>
29 void check_stdvector_matrix(const MatrixType& m)
30 {
31 typename MatrixType::Index rows = m.rows();
32 typename MatrixType::Index cols = m.cols();
33 MatrixType x = MatrixType::Random(rows,cols), y = MatrixType::Random(rows,cols);
34 std::vector<MatrixType> v(10, MatrixType(rows,cols)), w(20, y);
35 v[5] = x;
36 w[6] = v[5];
37 VERIFY_IS_APPROX(w[6], v[5]);
38 v = w;
39 for(int i = 0; i < 20; i++)
40 {
41 VERIFY_IS_APPROX(w[i], v[i]);
42 }
43
44 v.resize(21);
45 v[20] = x;
46 VERIFY_IS_APPROX(v[20], x);
47 v.resize(22,y);
48 VERIFY_IS_APPROX(v[21], y);
49 v.push_back(x);
50 VERIFY_IS_APPROX(v[22], x);
51 VERIFY((size_t)&(v[22]) == (size_t)&(v[21]) + sizeof(MatrixType));
52
53 // do a lot of push_back such that the vector gets internally resized
54 // (with memory reallocation)
55 MatrixType* ref = &w[0];
56 for(int i=0; i<30 || ((ref==&w[0]) && i<300); ++i)
57 v.push_back(w[i%w.size()]);
58 for(unsigned int i=23; i<v.size(); ++i)
59 {
60 VERIFY(v[i]==w[(i-23)%w.size()]);
61 }
62 }
63
64 template<typename TransformType>
check_stdvector_transform(const TransformType &)65 void check_stdvector_transform(const TransformType&)
66 {
67 typedef typename TransformType::MatrixType MatrixType;
68 TransformType x(MatrixType::Random()), y(MatrixType::Random());
69 std::vector<TransformType> v(10), w(20, y);
70 v[5] = x;
71 w[6] = v[5];
72 VERIFY_IS_APPROX(w[6], v[5]);
73 v = w;
74 for(int i = 0; i < 20; i++)
75 {
76 VERIFY_IS_APPROX(w[i], v[i]);
77 }
78
79 v.resize(21);
80 v[20] = x;
81 VERIFY_IS_APPROX(v[20], x);
82 v.resize(22,y);
83 VERIFY_IS_APPROX(v[21], y);
84 v.push_back(x);
85 VERIFY_IS_APPROX(v[22], x);
86 VERIFY((size_t)&(v[22]) == (size_t)&(v[21]) + sizeof(TransformType));
87
88 // do a lot of push_back such that the vector gets internally resized
89 // (with memory reallocation)
90 TransformType* ref = &w[0];
91 for(int i=0; i<30 || ((ref==&w[0]) && i<300); ++i)
92 v.push_back(w[i%w.size()]);
93 for(unsigned int i=23; i<v.size(); ++i)
94 {
95 VERIFY(v[i].matrix()==w[(i-23)%w.size()].matrix());
96 }
97 }
98
99 template<typename QuaternionType>
check_stdvector_quaternion(const QuaternionType &)100 void check_stdvector_quaternion(const QuaternionType&)
101 {
102 typedef typename QuaternionType::Coefficients Coefficients;
103 QuaternionType x(Coefficients::Random()), y(Coefficients::Random());
104 std::vector<QuaternionType> v(10), w(20, y);
105 v[5] = x;
106 w[6] = v[5];
107 VERIFY_IS_APPROX(w[6], v[5]);
108 v = w;
109 for(int i = 0; i < 20; i++)
110 {
111 VERIFY_IS_APPROX(w[i], v[i]);
112 }
113
114 v.resize(21);
115 v[20] = x;
116 VERIFY_IS_APPROX(v[20], x);
117 v.resize(22,y);
118 VERIFY_IS_APPROX(v[21], y);
119 v.push_back(x);
120 VERIFY_IS_APPROX(v[22], x);
121 VERIFY((size_t)&(v[22]) == (size_t)&(v[21]) + sizeof(QuaternionType));
122
123 // do a lot of push_back such that the vector gets internally resized
124 // (with memory reallocation)
125 QuaternionType* ref = &w[0];
126 for(int i=0; i<30 || ((ref==&w[0]) && i<300); ++i)
127 v.push_back(w[i%w.size()]);
128 for(unsigned int i=23; i<v.size(); ++i)
129 {
130 VERIFY(v[i].coeffs()==w[(i-23)%w.size()].coeffs());
131 }
132 }
133
test_stdvector_overload()134 void test_stdvector_overload()
135 {
136 // some non vectorizable fixed sizes
137 CALL_SUBTEST_1(check_stdvector_matrix(Vector2f()));
138 CALL_SUBTEST_1(check_stdvector_matrix(Matrix3f()));
139 CALL_SUBTEST_2(check_stdvector_matrix(Matrix3d()));
140
141 // some vectorizable fixed sizes
142 CALL_SUBTEST_1(check_stdvector_matrix(Matrix2f()));
143 CALL_SUBTEST_1(check_stdvector_matrix(Vector4f()));
144 CALL_SUBTEST_1(check_stdvector_matrix(Matrix4f()));
145 CALL_SUBTEST_2(check_stdvector_matrix(Matrix4d()));
146
147 // some dynamic sizes
148 CALL_SUBTEST_3(check_stdvector_matrix(MatrixXd(1,1)));
149 CALL_SUBTEST_3(check_stdvector_matrix(VectorXd(20)));
150 CALL_SUBTEST_3(check_stdvector_matrix(RowVectorXf(20)));
151 CALL_SUBTEST_3(check_stdvector_matrix(MatrixXcf(10,10)));
152
153 // some Transform
154 CALL_SUBTEST_4(check_stdvector_transform(Affine2f())); // does not need the specialization (2+1)^2 = 9
155 CALL_SUBTEST_4(check_stdvector_transform(Affine3f()));
156 CALL_SUBTEST_4(check_stdvector_transform(Affine3d()));
157
158 // some Quaternion
159 CALL_SUBTEST_5(check_stdvector_quaternion(Quaternionf()));
160 CALL_SUBTEST_5(check_stdvector_quaternion(Quaterniond()));
161 }
162