blob: 233268d1d0573b24a4e899d7dcdf0cf52e7e0505 [file] [log] [blame]
Benoit Jacob3958e7f2008-12-31 00:05:22 +00001// This file is part of Eigen, a lightweight C++ template library
Benoit Jacob6347b1d2009-05-22 20:25:33 +02002// for linear algebra.
Benoit Jacob3958e7f2008-12-31 00:05:22 +00003//
4// Copyright (C) 2008 Benoit Jacob <jacob.benoit.1@gmail.com>
5//
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
25#include "main.h"
26
Benoit Jacobd4157782009-10-05 10:11:11 -040027struct TestNew1
Benoit Jacob3958e7f2008-12-31 00:05:22 +000028{
29 MatrixXd m; // good: m will allocate its own array, taking care of alignment.
Benoit Jacobd4157782009-10-05 10:11:11 -040030 TestNew1() : m(20,20) {}
Benoit Jacob3958e7f2008-12-31 00:05:22 +000031};
32
Benoit Jacobd4157782009-10-05 10:11:11 -040033struct TestNew2
Benoit Jacob3958e7f2008-12-31 00:05:22 +000034{
Benoit Jacobd4157782009-10-05 10:11:11 -040035 Matrix3d m; // good: m's size isn't a multiple of 16 bytes, so m doesn't have to be 16-byte aligned,
36 // 8-byte alignment is good enough here, which we'll get automatically
Benoit Jacob3958e7f2008-12-31 00:05:22 +000037};
38
Benoit Jacobd4157782009-10-05 10:11:11 -040039struct TestNew3
Benoit Jacob3958e7f2008-12-31 00:05:22 +000040{
Benoit Jacobd4157782009-10-05 10:11:11 -040041 Vector2f m; // good: m's size isn't a multiple of 16 bytes, so m doesn't have to be 16-byte aligned
Benoit Jacob3958e7f2008-12-31 00:05:22 +000042};
43
Benoit Jacobd4157782009-10-05 10:11:11 -040044struct TestNew4
Benoit Jacob3958e7f2008-12-31 00:05:22 +000045{
Benoit Jacobeb7dcbb2009-01-08 15:37:13 +000046 EIGEN_MAKE_ALIGNED_OPERATOR_NEW
Benoit Jacob3958e7f2008-12-31 00:05:22 +000047 Vector2d m;
48 float f; // make the struct have sizeof%16!=0 to make it a little more tricky when we allow an array of 2 such objects
49};
50
Benoit Jacobd4157782009-10-05 10:11:11 -040051struct TestNew5
Benoit Jacob3958e7f2008-12-31 00:05:22 +000052{
Benoit Jacobeb7dcbb2009-01-08 15:37:13 +000053 EIGEN_MAKE_ALIGNED_OPERATOR_NEW
Benoit Jacobd4157782009-10-05 10:11:11 -040054 float f; // try the f at first -- the EIGEN_ALIGN16 attribute of m should make that still work
Benoit Jacob3958e7f2008-12-31 00:05:22 +000055 Matrix4f m;
56};
57
Benoit Jacobd4157782009-10-05 10:11:11 -040058struct TestNew6
Benoit Jacob15ca6652009-01-04 15:26:32 +000059{
Benoit Jacob2db43422009-01-06 18:07:16 +000060 Matrix<float,2,2,DontAlign> m; // good: no alignment requested
Benoit Jacob15ca6652009-01-04 15:26:32 +000061 float f;
62};
63
Benoit Jacob1c29d702009-01-06 03:16:50 +000064template<bool Align> struct Depends
65{
Benoit Jacobeb7dcbb2009-01-08 15:37:13 +000066 EIGEN_MAKE_ALIGNED_OPERATOR_NEW_IF(Align)
Benoit Jacob1c29d702009-01-06 03:16:50 +000067 Vector2d m;
68 float f;
69};
70
Benoit Jacob3958e7f2008-12-31 00:05:22 +000071template<typename T>
72void check_unalignedassert_good()
73{
74 T *x, *y;
75 x = new T;
76 delete x;
77 y = new T[2];
78 delete[] y;
79}
80
Benoit Jacob95bda5e2009-05-03 13:50:56 +000081#if EIGEN_ALIGN
Benoit Jacob3958e7f2008-12-31 00:05:22 +000082template<typename T>
Benoit Jacobd4157782009-10-05 10:11:11 -040083void construct_at_boundary(int boundary)
Benoit Jacob3958e7f2008-12-31 00:05:22 +000084{
Benoit Jacobd4157782009-10-05 10:11:11 -040085 char buf[sizeof(T)+256];
86 size_t _buf = reinterpret_cast<size_t>(buf);
87 _buf += (16 - (_buf % 16)); // make 16-byte aligned
88 _buf += boundary; // make exact boundary-aligned
89 T *x = ::new(reinterpret_cast<void*>(_buf)) T;
Benoit Jacob3958e7f2008-12-31 00:05:22 +000090 x->~T();
91}
Benoit Jacobf81479d2009-02-04 16:55:38 +000092#endif
Benoit Jacob3958e7f2008-12-31 00:05:22 +000093
94void unalignedassert()
95{
Benoit Jacobbb1cc0d2009-10-05 10:55:42 -040096 construct_at_boundary<Vector2f>(4);
Benoit Jacobd4157782009-10-05 10:11:11 -040097 construct_at_boundary<Vector3f>(4);
98 construct_at_boundary<Vector4f>(16);
99 construct_at_boundary<Matrix2f>(16);
100 construct_at_boundary<Matrix3f>(4);
101 construct_at_boundary<Matrix4f>(16);
102
103 construct_at_boundary<Vector2d>(16);
Benoit Jacobbb1cc0d2009-10-05 10:55:42 -0400104 construct_at_boundary<Vector3d>(4);
Benoit Jacobd4157782009-10-05 10:11:11 -0400105 construct_at_boundary<Vector4d>(16);
106 construct_at_boundary<Matrix2d>(16);
Benoit Jacobbb1cc0d2009-10-05 10:55:42 -0400107 construct_at_boundary<Matrix3d>(4);
Benoit Jacobd4157782009-10-05 10:11:11 -0400108 construct_at_boundary<Matrix4d>(16);
109
Benoit Jacobbb1cc0d2009-10-05 10:55:42 -0400110 construct_at_boundary<Vector2cf>(16);
111 construct_at_boundary<Vector3cf>(4);
112 construct_at_boundary<Vector2cd>(16);
113 construct_at_boundary<Vector3cd>(16);
114
Benoit Jacobd4157782009-10-05 10:11:11 -0400115 check_unalignedassert_good<TestNew1>();
116 check_unalignedassert_good<TestNew2>();
117 check_unalignedassert_good<TestNew3>();
Benoit Jacobf81479d2009-02-04 16:55:38 +0000118
Benoit Jacobd4157782009-10-05 10:11:11 -0400119 check_unalignedassert_good<TestNew4>();
120 check_unalignedassert_good<TestNew5>();
121 check_unalignedassert_good<TestNew6>();
Benoit Jacob1c29d702009-01-06 03:16:50 +0000122 check_unalignedassert_good<Depends<true> >();
Benoit Jacobd4157782009-10-05 10:11:11 -0400123
Benoit Jacob95bda5e2009-05-03 13:50:56 +0000124#if EIGEN_ALIGN
Benoit Jacobd4157782009-10-05 10:11:11 -0400125 VERIFY_RAISES_ASSERT(construct_at_boundary<Vector4f>(8));
Benoit Jacobd4157782009-10-05 10:11:11 -0400126 VERIFY_RAISES_ASSERT(construct_at_boundary<Matrix4f>(8));
127 VERIFY_RAISES_ASSERT(construct_at_boundary<Vector2d>(8));
Benoit Jacobd4157782009-10-05 10:11:11 -0400128 VERIFY_RAISES_ASSERT(construct_at_boundary<Vector4d>(8));
129 VERIFY_RAISES_ASSERT(construct_at_boundary<Matrix2d>(8));
Benoit Jacobd4157782009-10-05 10:11:11 -0400130 VERIFY_RAISES_ASSERT(construct_at_boundary<Matrix4d>(8));
Benoit Jacobbb1cc0d2009-10-05 10:55:42 -0400131 VERIFY_RAISES_ASSERT(construct_at_boundary<Vector2cf>(8));
132 VERIFY_RAISES_ASSERT(construct_at_boundary<Vector2cd>(8));
Benoit Jacobf81479d2009-02-04 16:55:38 +0000133#endif
Benoit Jacob3958e7f2008-12-31 00:05:22 +0000134}
135
136void test_unalignedassert()
137{
138 CALL_SUBTEST(unalignedassert());
139}