stdvector.cpp File Reference
#include "main.h"
#include <Eigen/StdVector>
#include <Eigen/Geometry>

Functions

template<typename MatrixType >
void check_stdvector_matrix (const MatrixType &m)
 
template<typename TransformType >
void check_stdvector_transform (const TransformType &)
 
template<typename QuaternionType >
void check_stdvector_quaternion (const QuaternionType &)
 
void std_vector_gcc_warning ()
 
 EIGEN_DECLARE_TEST (stdvector)
 

Function Documentation

◆ check_stdvector_matrix()

template<typename MatrixType >
void check_stdvector_matrix ( const MatrixType m)
15  {
16  Index rows = m.rows();
17  Index cols = m.cols();
18  MatrixType x = MatrixType::Random(rows, cols), y = MatrixType::Random(rows, cols);
19  std::vector<MatrixType, Eigen::aligned_allocator<MatrixType> > v(10, MatrixType::Zero(rows, cols)), w(20, y);
20  v[5] = x;
21  w[6] = v[5];
22  VERIFY_IS_APPROX(w[6], v[5]);
23  v = w;
24  for (int i = 0; i < 20; i++) {
25  VERIFY_IS_APPROX(w[i], v[i]);
26  }
27 
28  v.resize(21);
29  v[20] = x;
30  VERIFY_IS_APPROX(v[20], x);
31  v.resize(22, y);
32  VERIFY_IS_APPROX(v[21], y);
33  v.push_back(x);
34  VERIFY_IS_APPROX(v[22], x);
35  VERIFY((std::uintptr_t) & (v[22]) == (std::uintptr_t) & (v[21]) + sizeof(MatrixType));
36 
37  // do a lot of push_back such that the vector gets internally resized
38  // (with memory reallocation)
39  MatrixType* ref = &w[0];
40  for (int i = 0; i < 30 || ((ref == &w[0]) && i < 300); ++i) v.push_back(w[i % w.size()]);
41  for (unsigned int i = 23; i < v.size(); ++i) {
42  VERIFY(v[i] == w[(i - 23) % w.size()]);
43  }
44 }
Array< int, Dynamic, 1 > v
Definition: Array_initializer_list_vector_cxx11.cpp:1
int i
Definition: BiCGSTAB_step_by_step.cpp:9
RowVector3d w
Definition: Matrix_resize_int.cpp:3
int rows
Definition: Tutorial_commainit_02.cpp:1
int cols
Definition: Tutorial_commainit_02.cpp:1
MatrixXf MatrixType
Definition: benchmark-blocking-sizes.cpp:52
#define VERIFY_IS_APPROX(a, b)
Definition: integer_types.cpp:13
Scalar * y
Definition: level1_cplx_impl.h:128
int * m
Definition: level2_cplx_impl.h:294
#define VERIFY(a)
Definition: main.h:362
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:83
double Zero
Definition: pseudosolid_node_update_elements.cc:35
list x
Definition: plotDoE.py:28

References cols, i, m, rows, v, VERIFY, VERIFY_IS_APPROX, w, plotDoE::x, y, and oomph::PseudoSolidHelper::Zero.

Referenced by EIGEN_DECLARE_TEST().

◆ check_stdvector_quaternion()

template<typename QuaternionType >
void check_stdvector_quaternion ( const QuaternionType &  )
78  {
79  typedef typename QuaternionType::Coefficients Coefficients;
80  QuaternionType x(Coefficients::Random()), y(Coefficients::Random()), qi = QuaternionType::Identity();
81  std::vector<QuaternionType, Eigen::aligned_allocator<QuaternionType> > v(10, qi), w(20, y);
82  v[5] = x;
83  w[6] = v[5];
84  VERIFY_IS_APPROX(w[6], v[5]);
85  v = w;
86  for (int i = 0; i < 20; i++) {
87  VERIFY_IS_APPROX(w[i], v[i]);
88  }
89 
90  v.resize(21);
91  v[20] = x;
92  VERIFY_IS_APPROX(v[20], x);
93  v.resize(22, y);
94  VERIFY_IS_APPROX(v[21], y);
95  v.push_back(x);
96  VERIFY_IS_APPROX(v[22], x);
97  VERIFY((std::uintptr_t) & (v[22]) == (std::uintptr_t) & (v[21]) + sizeof(QuaternionType));
98 
99  // do a lot of push_back such that the vector gets internally resized
100  // (with memory reallocation)
101  QuaternionType* ref = &w[0];
102  for (int i = 0; i < 30 || ((ref == &w[0]) && i < 300); ++i) v.push_back(w[i % w.size()]);
103  for (unsigned int i = 23; i < v.size(); ++i) {
104  VERIFY(v[i].coeffs() == w[(i - 23) % w.size()].coeffs());
105  }
106 }

References i, v, VERIFY, VERIFY_IS_APPROX, w, plotDoE::x, and y.

Referenced by EIGEN_DECLARE_TEST().

◆ check_stdvector_transform()

template<typename TransformType >
void check_stdvector_transform ( const TransformType &  )
47  {
48  typedef typename TransformType::MatrixType MatrixType;
49  TransformType x(MatrixType::Random()), y(MatrixType::Random());
50  std::vector<TransformType, Eigen::aligned_allocator<TransformType> > v(10), w(20, y);
51  v[5] = x;
52  w[6] = v[5];
53  VERIFY_IS_APPROX(w[6], v[5]);
54  v = w;
55  for (int i = 0; i < 20; i++) {
56  VERIFY_IS_APPROX(w[i], v[i]);
57  }
58 
59  v.resize(21);
60  v[20] = x;
61  VERIFY_IS_APPROX(v[20], x);
62  v.resize(22, y);
63  VERIFY_IS_APPROX(v[21], y);
64  v.push_back(x);
65  VERIFY_IS_APPROX(v[22], x);
66  VERIFY((std::uintptr_t) & (v[22]) == (std::uintptr_t) & (v[21]) + sizeof(TransformType));
67 
68  // do a lot of push_back such that the vector gets internally resized
69  // (with memory reallocation)
70  TransformType* ref = &w[0];
71  for (int i = 0; i < 30 || ((ref == &w[0]) && i < 300); ++i) v.push_back(w[i % w.size()]);
72  for (unsigned int i = 23; i < v.size(); ++i) {
73  VERIFY(v[i].matrix() == w[(i - 23) % w.size()].matrix());
74  }
75 }
Eigen::Map< Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic, Eigen::ColMajor >, 0, Eigen::OuterStride<> > matrix(T *data, int rows, int cols, int stride)
Definition: common.h:85

References i, matrix(), v, VERIFY, VERIFY_IS_APPROX, w, plotDoE::x, and y.

Referenced by EIGEN_DECLARE_TEST().

◆ EIGEN_DECLARE_TEST()

EIGEN_DECLARE_TEST ( stdvector  )
117  {
118  // some non vectorizable fixed sizes
122 
123  // some vectorizable fixed sizes
128 
129  // some dynamic sizes
130  CALL_SUBTEST_3(check_stdvector_matrix(MatrixXd(1, 1)));
131  CALL_SUBTEST_3(check_stdvector_matrix(VectorXd(20)));
132  CALL_SUBTEST_3(check_stdvector_matrix(RowVectorXf(20)));
133  CALL_SUBTEST_3(check_stdvector_matrix(MatrixXcf(10, 10)));
134 
135  // some Transform
139  // CALL_SUBTEST(heck_stdvector_transform(Projective4d()));
140 
141  // some Quaternion
144 }
Transform< float, 2, Projective > Projective2f
Definition: Transform.h:698
Transform< double, 3, Projective > Projective3d
Definition: Transform.h:704
Quaternion< double > Quaterniond
Definition: Eigen/Eigen/src/Geometry/Quaternion.h:387
Quaternion< float > Quaternionf
Definition: Eigen/Eigen/src/Geometry/Quaternion.h:384
Transform< float, 3, Projective > Projective3f
Definition: Transform.h:700
#define CALL_SUBTEST_3(FUNC)
Definition: split_test_helper.h:16
#define CALL_SUBTEST_1(FUNC)
Definition: split_test_helper.h:4
#define CALL_SUBTEST_5(FUNC)
Definition: split_test_helper.h:28
#define CALL_SUBTEST_2(FUNC)
Definition: split_test_helper.h:10
#define CALL_SUBTEST_4(FUNC)
Definition: split_test_helper.h:22
void check_stdvector_matrix(const MatrixType &m)
Definition: stdvector.cpp:15
void check_stdvector_transform(const TransformType &)
Definition: stdvector.cpp:47
void check_stdvector_quaternion(const QuaternionType &)
Definition: stdvector.cpp:78

References CALL_SUBTEST_1, CALL_SUBTEST_2, CALL_SUBTEST_3, CALL_SUBTEST_4, CALL_SUBTEST_5, check_stdvector_matrix(), check_stdvector_quaternion(), and check_stdvector_transform().

◆ std_vector_gcc_warning()

void std_vector_gcc_warning ( )
111  {
112  typedef Eigen::Vector3f T;
113  std::vector<T, Eigen::aligned_allocator<T> > v;
114  v.push_back(T(1.0f, 2.0f, 3.0f));
115 }
Eigen::Triplet< double > T
Definition: EigenUnitTest.cpp:11

References v.