blob: e2e73680e02042cd07b6beb00a959b03bf558f3c [file] [log] [blame]
Benoit Jacob335d3bc2009-01-09 23:26:45 +00001// This file is part of Eigen, a lightweight C++ template library
2// for linear algebra. Eigen itself is part of the KDE project.
3//
Gael Guennebaudf645d1f2009-01-20 16:50:47 +00004// Copyright (C) 2008 Benoit Jacob <jacob.benoit.1@gmail.com>
Benoit Jacob335d3bc2009-01-09 23:26:45 +00005//
6// Eigen is free software; you can redistribute it and/or
7// modify it under the terms of the GNU Lesser General Public
8// License as published by the Free Software Foundation; either
9// version 3 of the License, or (at your option) any later version.
10//
11// Alternatively, you can redistribute it and/or
12// modify it under the terms of the GNU General Public License as
13// published by the Free Software Foundation; either version 2 of
14// the License, or (at your option) any later version.
15//
16// Eigen is distributed in the hope that it will be useful, but WITHOUT ANY
17// WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
18// FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the
19// GNU General Public License for more details.
20//
21// You should have received a copy of the GNU Lesser General Public
22// License and a copy of the GNU General Public License along with
23// Eigen. If not, see <http://www.gnu.org/licenses/>.
24
Gael Guennebaud3009d792009-02-07 11:16:15 +000025#include "main.h"
Benoit Jacobeac79b62009-05-09 03:41:17 +000026#include <Eigen/StdVector>
Benoit Jacob4d44ca22009-01-12 16:06:04 +000027#include <Eigen/Geometry>
Benoit Jacob335d3bc2009-01-09 23:26:45 +000028
29template<typename MatrixType>
Benoit Jacob4d44ca22009-01-12 16:06:04 +000030void check_stdvector_matrix(const MatrixType& m)
Benoit Jacob335d3bc2009-01-09 23:26:45 +000031{
Benoit Jacob0c1ef2f2009-01-10 13:10:23 +000032 int rows = m.rows();
33 int cols = m.cols();
34 MatrixType x = MatrixType::Random(rows,cols), y = MatrixType::Random(rows,cols);
Gael Guennebaud469b8aa2009-05-03 15:39:37 +000035 std::vector<MatrixType,Eigen::aligned_allocator<MatrixType> > v(10, MatrixType(rows,cols)), w(20, y);
Benoit Jacob335d3bc2009-01-09 23:26:45 +000036 v[5] = x;
37 w[6] = v[5];
38 VERIFY_IS_APPROX(w[6], v[5]);
39 v = w;
40 for(int i = 0; i < 20; i++)
41 {
42 VERIFY_IS_APPROX(w[i], v[i]);
43 }
Gael Guennebaud9e8f4372009-01-11 21:59:04 +000044
Benoit Jacob335d3bc2009-01-09 23:26:45 +000045 v.resize(21);
Benoit Jacob5f43a422009-01-21 17:10:23 +000046 v[20] = x;
Benoit Jacob335d3bc2009-01-09 23:26:45 +000047 VERIFY_IS_APPROX(v[20], x);
48 v.resize(22,y);
Gael Guennebaud9e8f4372009-01-11 21:59:04 +000049 VERIFY_IS_APPROX(v[21], y);
Benoit Jacob335d3bc2009-01-09 23:26:45 +000050 v.push_back(x);
51 VERIFY_IS_APPROX(v[22], x);
Benoit Jacob0c1ef2f2009-01-10 13:10:23 +000052 VERIFY((size_t)&(v[22]) == (size_t)&(v[21]) + sizeof(MatrixType));
Gael Guennebaud9e8f4372009-01-11 21:59:04 +000053
54 // do a lot of push_back such that the vector gets internally resized
55 // (with memory reallocation)
56 MatrixType* ref = &w[0];
Benoit Jacob336ad582009-01-12 13:41:40 +000057 for(int i=0; i<30 || ((ref==&w[0]) && i<300); ++i)
Gael Guennebaud9e8f4372009-01-11 21:59:04 +000058 v.push_back(w[i%w.size()]);
Benoit Jacob336ad582009-01-12 13:41:40 +000059 for(unsigned int i=23; i<v.size(); ++i)
Gael Guennebaud9e8f4372009-01-11 21:59:04 +000060 {
61 VERIFY(v[i]==w[(i-23)%w.size()]);
62 }
Benoit Jacob335d3bc2009-01-09 23:26:45 +000063}
64
Benoit Jacob4d44ca22009-01-12 16:06:04 +000065template<typename TransformType>
66void check_stdvector_transform(const TransformType&)
67{
68 typedef typename TransformType::MatrixType MatrixType;
69 TransformType x(MatrixType::Random()), y(MatrixType::Random());
Gael Guennebaud469b8aa2009-05-03 15:39:37 +000070 std::vector<TransformType,Eigen::aligned_allocator<TransformType> > v(10), w(20, y);
Benoit Jacob4d44ca22009-01-12 16:06:04 +000071 v[5] = x;
72 w[6] = v[5];
73 VERIFY_IS_APPROX(w[6], v[5]);
74 v = w;
75 for(int i = 0; i < 20; i++)
76 {
77 VERIFY_IS_APPROX(w[i], v[i]);
78 }
79
80 v.resize(21);
81 v[20] = x;
82 VERIFY_IS_APPROX(v[20], x);
83 v.resize(22,y);
84 VERIFY_IS_APPROX(v[21], y);
85 v.push_back(x);
86 VERIFY_IS_APPROX(v[22], x);
87 VERIFY((size_t)&(v[22]) == (size_t)&(v[21]) + sizeof(TransformType));
88
89 // do a lot of push_back such that the vector gets internally resized
90 // (with memory reallocation)
91 TransformType* ref = &w[0];
92 for(int i=0; i<30 || ((ref==&w[0]) && i<300); ++i)
93 v.push_back(w[i%w.size()]);
94 for(unsigned int i=23; i<v.size(); ++i)
95 {
96 VERIFY(v[i].matrix()==w[(i-23)%w.size()].matrix());
97 }
98}
99
100template<typename QuaternionType>
101void check_stdvector_quaternion(const QuaternionType&)
102{
103 typedef typename QuaternionType::Coefficients Coefficients;
104 QuaternionType x(Coefficients::Random()), y(Coefficients::Random());
Gael Guennebaud469b8aa2009-05-03 15:39:37 +0000105 std::vector<QuaternionType,Eigen::aligned_allocator<QuaternionType> > v(10), w(20, y);
Benoit Jacob4d44ca22009-01-12 16:06:04 +0000106 v[5] = x;
107 w[6] = v[5];
108 VERIFY_IS_APPROX(w[6], v[5]);
109 v = w;
110 for(int i = 0; i < 20; i++)
111 {
112 VERIFY_IS_APPROX(w[i], v[i]);
113 }
114
115 v.resize(21);
116 v[20] = x;
117 VERIFY_IS_APPROX(v[20], x);
118 v.resize(22,y);
119 VERIFY_IS_APPROX(v[21], y);
120 v.push_back(x);
121 VERIFY_IS_APPROX(v[22], x);
122 VERIFY((size_t)&(v[22]) == (size_t)&(v[21]) + sizeof(QuaternionType));
123
124 // do a lot of push_back such that the vector gets internally resized
125 // (with memory reallocation)
126 QuaternionType* ref = &w[0];
127 for(int i=0; i<30 || ((ref==&w[0]) && i<300); ++i)
128 v.push_back(w[i%w.size()]);
129 for(unsigned int i=23; i<v.size(); ++i)
130 {
131 VERIFY(v[i].coeffs()==w[(i-23)%w.size()].coeffs());
132 }
133}
134
Benoit Jacob335d3bc2009-01-09 23:26:45 +0000135void test_stdvector()
136{
137 // some non vectorizable fixed sizes
Benoit Jacob4d44ca22009-01-12 16:06:04 +0000138 CALL_SUBTEST(check_stdvector_matrix(Vector2f()));
139 CALL_SUBTEST(check_stdvector_matrix(Matrix3f()));
140 CALL_SUBTEST(check_stdvector_matrix(Matrix3d()));
Benoit Jacob335d3bc2009-01-09 23:26:45 +0000141
142 // some vectorizable fixed sizes
Benoit Jacob4d44ca22009-01-12 16:06:04 +0000143 CALL_SUBTEST(check_stdvector_matrix(Matrix2f()));
144 CALL_SUBTEST(check_stdvector_matrix(Vector4f()));
145 CALL_SUBTEST(check_stdvector_matrix(Matrix4f()));
146 CALL_SUBTEST(check_stdvector_matrix(Matrix4d()));
Benoit Jacob0c1ef2f2009-01-10 13:10:23 +0000147
148 // some dynamic sizes
Benoit Jacob4d44ca22009-01-12 16:06:04 +0000149 CALL_SUBTEST(check_stdvector_matrix(MatrixXd(1,1)));
150 CALL_SUBTEST(check_stdvector_matrix(VectorXd(20)));
151 CALL_SUBTEST(check_stdvector_matrix(RowVectorXf(20)));
152 CALL_SUBTEST(check_stdvector_matrix(MatrixXcf(10,10)));
153
154 // some Transform
155 CALL_SUBTEST(check_stdvector_transform(Transform2f()));
156 CALL_SUBTEST(check_stdvector_transform(Transform3f()));
157 CALL_SUBTEST(check_stdvector_transform(Transform3d()));
158 //CALL_SUBTEST(check_stdvector_transform(Transform4d()));
159
160 // some Quaternion
161 CALL_SUBTEST(check_stdvector_quaternion(Quaternionf()));
162 CALL_SUBTEST(check_stdvector_quaternion(Quaternionf()));
Benoit Jacob335d3bc2009-01-09 23:26:45 +0000163}