DenseStorage.h
1// This file is part of Eigen, a lightweight C++ template library
2// for linear algebra.
3//
4// Copyright (C) 2008 Gael Guennebaud <gael.guennebaud@inria.fr>
5// Copyright (C) 2006-2009 Benoit Jacob <jacob.benoit.1@gmail.com>
6// Copyright (C) 2010 Hauke Heibel <hauke.heibel@gmail.com>
7//
8// This Source Code Form is subject to the terms of the Mozilla
9// Public License v. 2.0. If a copy of the MPL was not distributed
10// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
11
12#ifndef EIGEN_MATRIXSTORAGE_H
13#define EIGEN_MATRIXSTORAGE_H
14
15#ifdef EIGEN_DENSE_STORAGE_CTOR_PLUGIN
16 #define EIGEN_INTERNAL_DENSE_STORAGE_CTOR_PLUGIN EIGEN_DENSE_STORAGE_CTOR_PLUGIN;
17#else
18 #define EIGEN_INTERNAL_DENSE_STORAGE_CTOR_PLUGIN
19#endif
20
21namespace Eigen {
22
23namespace internal {
24
25struct constructor_without_unaligned_array_assert {};
26
31template <typename T, int Size, int MatrixOrArrayOptions,
32 int Alignment = (MatrixOrArrayOptions&DontAlign) ? 0
33 : (((Size*sizeof(T))%16)==0) ? 16
34 : 0 >
35struct plain_array
36{
37 T array[Size];
38 plain_array() {}
39 plain_array(constructor_without_unaligned_array_assert) {}
40};
41
42#if defined(EIGEN_DISABLE_UNALIGNED_ARRAY_ASSERT)
43 #define EIGEN_MAKE_UNALIGNED_ARRAY_ASSERT(sizemask)
44#elif EIGEN_GNUC_AT_LEAST(4,7)
45 // GCC 4.7 is too aggressive in its optimizations and remove the alignement test based on the fact the array is declared to be aligned.
46 // See this bug report: http://gcc.gnu.org/bugzilla/show_bug.cgi?id=53900
47 // Hiding the origin of the array pointer behind a function argument seems to do the trick even if the function is inlined:
48 template<typename PtrType>
49 EIGEN_ALWAYS_INLINE PtrType eigen_unaligned_array_assert_workaround_gcc47(PtrType array) { return array; }
50 #define EIGEN_MAKE_UNALIGNED_ARRAY_ASSERT(sizemask) \
51 eigen_assert((reinterpret_cast<size_t>(eigen_unaligned_array_assert_workaround_gcc47(array)) & sizemask) == 0 \
52 && "this assertion is explained here: " \
53 "http://eigen.tuxfamily.org/dox-devel/group__TopicUnalignedArrayAssert.html" \
54 " **** READ THIS WEB PAGE !!! ****");
55#else
56 #define EIGEN_MAKE_UNALIGNED_ARRAY_ASSERT(sizemask) \
57 eigen_assert((reinterpret_cast<size_t>(array) & sizemask) == 0 \
58 && "this assertion is explained here: " \
59 "http://eigen.tuxfamily.org/dox-devel/group__TopicUnalignedArrayAssert.html" \
60 " **** READ THIS WEB PAGE !!! ****");
61#endif
62
63template <typename T, int Size, int MatrixOrArrayOptions>
64struct plain_array<T, Size, MatrixOrArrayOptions, 16>
65{
66 EIGEN_USER_ALIGN16 T array[Size];
67 plain_array() { EIGEN_MAKE_UNALIGNED_ARRAY_ASSERT(0xf) }
68 plain_array(constructor_without_unaligned_array_assert) {}
69};
70
71template <typename T, int MatrixOrArrayOptions, int Alignment>
72struct plain_array<T, 0, MatrixOrArrayOptions, Alignment>
73{
74 EIGEN_USER_ALIGN16 T array[1];
75 plain_array() {}
76 plain_array(constructor_without_unaligned_array_assert) {}
77};
78
79} // end namespace internal
80
93template<typename T, int Size, int _Rows, int _Cols, int _Options> class DenseStorage;
94
95// purely fixed-size matrix
96template<typename T, int Size, int _Rows, int _Cols, int _Options> class DenseStorage
97{
98 internal::plain_array<T,Size,_Options> m_data;
99 public:
100 inline explicit DenseStorage() {}
101 inline DenseStorage(internal::constructor_without_unaligned_array_assert)
102 : m_data(internal::constructor_without_unaligned_array_assert()) {}
103 inline DenseStorage(DenseIndex,DenseIndex,DenseIndex) {}
104 inline void swap(DenseStorage& other) { std::swap(m_data,other.m_data); }
105 static inline DenseIndex rows(void) {return _Rows;}
106 static inline DenseIndex cols(void) {return _Cols;}
107 inline void conservativeResize(DenseIndex,DenseIndex,DenseIndex) {}
108 inline void resize(DenseIndex,DenseIndex,DenseIndex) {}
109 inline const T *data() const { return m_data.array; }
110 inline T *data() { return m_data.array; }
111};
112
113// null matrix
114template<typename T, int _Rows, int _Cols, int _Options> class DenseStorage<T, 0, _Rows, _Cols, _Options>
115{
116 public:
117 inline explicit DenseStorage() {}
118 inline DenseStorage(internal::constructor_without_unaligned_array_assert) {}
119 inline DenseStorage(DenseIndex,DenseIndex,DenseIndex) {}
120 inline void swap(DenseStorage& ) {}
121 static inline DenseIndex rows(void) {return _Rows;}
122 static inline DenseIndex cols(void) {return _Cols;}
123 inline void conservativeResize(DenseIndex,DenseIndex,DenseIndex) {}
124 inline void resize(DenseIndex,DenseIndex,DenseIndex) {}
125 inline const T *data() const { return 0; }
126 inline T *data() { return 0; }
127};
128
129// more specializations for null matrices; these are necessary to resolve ambiguities
130template<typename T, int _Options> class DenseStorage<T, 0, Dynamic, Dynamic, _Options>
131: public DenseStorage<T, 0, 0, 0, _Options> { };
132
133template<typename T, int _Rows, int _Options> class DenseStorage<T, 0, _Rows, Dynamic, _Options>
134: public DenseStorage<T, 0, 0, 0, _Options> { };
135
136template<typename T, int _Cols, int _Options> class DenseStorage<T, 0, Dynamic, _Cols, _Options>
137: public DenseStorage<T, 0, 0, 0, _Options> { };
138
139// dynamic-size matrix with fixed-size storage
140template<typename T, int Size, int _Options> class DenseStorage<T, Size, Dynamic, Dynamic, _Options>
141{
142 internal::plain_array<T,Size,_Options> m_data;
143 DenseIndex m_rows;
144 DenseIndex m_cols;
145 public:
146 inline explicit DenseStorage() : m_rows(0), m_cols(0) {}
147 inline DenseStorage(internal::constructor_without_unaligned_array_assert)
148 : m_data(internal::constructor_without_unaligned_array_assert()), m_rows(0), m_cols(0) {}
149 inline DenseStorage(DenseIndex, DenseIndex rows, DenseIndex cols) : m_rows(rows), m_cols(cols) {}
150 inline void swap(DenseStorage& other)
151 { std::swap(m_data,other.m_data); std::swap(m_rows,other.m_rows); std::swap(m_cols,other.m_cols); }
152 inline DenseIndex rows(void) const {return m_rows;}
153 inline DenseIndex cols(void) const {return m_cols;}
154 inline void conservativeResize(DenseIndex, DenseIndex rows, DenseIndex cols) { m_rows = rows; m_cols = cols; }
155 inline void resize(DenseIndex, DenseIndex rows, DenseIndex cols) { m_rows = rows; m_cols = cols; }
156 inline const T *data() const { return m_data.array; }
157 inline T *data() { return m_data.array; }
158};
159
160// dynamic-size matrix with fixed-size storage and fixed width
161template<typename T, int Size, int _Cols, int _Options> class DenseStorage<T, Size, Dynamic, _Cols, _Options>
162{
163 internal::plain_array<T,Size,_Options> m_data;
164 DenseIndex m_rows;
165 public:
166 inline explicit DenseStorage() : m_rows(0) {}
167 inline DenseStorage(internal::constructor_without_unaligned_array_assert)
168 : m_data(internal::constructor_without_unaligned_array_assert()), m_rows(0) {}
169 inline DenseStorage(DenseIndex, DenseIndex rows, DenseIndex) : m_rows(rows) {}
170 inline void swap(DenseStorage& other) { std::swap(m_data,other.m_data); std::swap(m_rows,other.m_rows); }
171 inline DenseIndex rows(void) const {return m_rows;}
172 inline DenseIndex cols(void) const {return _Cols;}
173 inline void conservativeResize(DenseIndex, DenseIndex rows, DenseIndex) { m_rows = rows; }
174 inline void resize(DenseIndex, DenseIndex rows, DenseIndex) { m_rows = rows; }
175 inline const T *data() const { return m_data.array; }
176 inline T *data() { return m_data.array; }
177};
178
179// dynamic-size matrix with fixed-size storage and fixed height
180template<typename T, int Size, int _Rows, int _Options> class DenseStorage<T, Size, _Rows, Dynamic, _Options>
181{
182 internal::plain_array<T,Size,_Options> m_data;
183 DenseIndex m_cols;
184 public:
185 inline explicit DenseStorage() : m_cols(0) {}
186 inline DenseStorage(internal::constructor_without_unaligned_array_assert)
187 : m_data(internal::constructor_without_unaligned_array_assert()), m_cols(0) {}
188 inline DenseStorage(DenseIndex, DenseIndex, DenseIndex cols) : m_cols(cols) {}
189 inline void swap(DenseStorage& other) { std::swap(m_data,other.m_data); std::swap(m_cols,other.m_cols); }
190 inline DenseIndex rows(void) const {return _Rows;}
191 inline DenseIndex cols(void) const {return m_cols;}
192 inline void conservativeResize(DenseIndex, DenseIndex, DenseIndex cols) { m_cols = cols; }
193 inline void resize(DenseIndex, DenseIndex, DenseIndex cols) { m_cols = cols; }
194 inline const T *data() const { return m_data.array; }
195 inline T *data() { return m_data.array; }
196};
197
198// purely dynamic matrix.
199template<typename T, int _Options> class DenseStorage<T, Dynamic, Dynamic, Dynamic, _Options>
200{
201 T *m_data;
202 DenseIndex m_rows;
203 DenseIndex m_cols;
204 public:
205 inline explicit DenseStorage() : m_data(0), m_rows(0), m_cols(0) {}
206 inline DenseStorage(internal::constructor_without_unaligned_array_assert)
207 : m_data(0), m_rows(0), m_cols(0) {}
208 inline DenseStorage(DenseIndex size, DenseIndex rows, DenseIndex cols)
209 : m_data(internal::conditional_aligned_new_auto<T,(_Options&DontAlign)==0>(size)), m_rows(rows), m_cols(cols)
210 { EIGEN_INTERNAL_DENSE_STORAGE_CTOR_PLUGIN }
211 inline ~DenseStorage() { internal::conditional_aligned_delete_auto<T,(_Options&DontAlign)==0>(m_data, m_rows*m_cols); }
212 inline void swap(DenseStorage& other)
213 { std::swap(m_data,other.m_data); std::swap(m_rows,other.m_rows); std::swap(m_cols,other.m_cols); }
214 inline DenseIndex rows(void) const {return m_rows;}
215 inline DenseIndex cols(void) const {return m_cols;}
216 inline void conservativeResize(DenseIndex size, DenseIndex rows, DenseIndex cols)
217 {
218 m_data = internal::conditional_aligned_realloc_new_auto<T,(_Options&DontAlign)==0>(m_data, size, m_rows*m_cols);
219 m_rows = rows;
220 m_cols = cols;
221 }
222 void resize(DenseIndex size, DenseIndex rows, DenseIndex cols)
223 {
224 if(size != m_rows*m_cols)
225 {
226 internal::conditional_aligned_delete_auto<T,(_Options&DontAlign)==0>(m_data, m_rows*m_cols);
227 if (size)
228 m_data = internal::conditional_aligned_new_auto<T,(_Options&DontAlign)==0>(size);
229 else
230 m_data = 0;
231 EIGEN_INTERNAL_DENSE_STORAGE_CTOR_PLUGIN
232 }
233 m_rows = rows;
234 m_cols = cols;
235 }
236 inline const T *data() const { return m_data; }
237 inline T *data() { return m_data; }
238};
239
240// matrix with dynamic width and fixed height (so that matrix has dynamic size).
241template<typename T, int _Rows, int _Options> class DenseStorage<T, Dynamic, _Rows, Dynamic, _Options>
242{
243 T *m_data;
244 DenseIndex m_cols;
245 public:
246 inline explicit DenseStorage() : m_data(0), m_cols(0) {}
247 inline DenseStorage(internal::constructor_without_unaligned_array_assert) : m_data(0), m_cols(0) {}
248 inline DenseStorage(DenseIndex size, DenseIndex, DenseIndex cols) : m_data(internal::conditional_aligned_new_auto<T,(_Options&DontAlign)==0>(size)), m_cols(cols)
249 { EIGEN_INTERNAL_DENSE_STORAGE_CTOR_PLUGIN }
250 inline ~DenseStorage() { internal::conditional_aligned_delete_auto<T,(_Options&DontAlign)==0>(m_data, _Rows*m_cols); }
251 inline void swap(DenseStorage& other) { std::swap(m_data,other.m_data); std::swap(m_cols,other.m_cols); }
252 static inline DenseIndex rows(void) {return _Rows;}
253 inline DenseIndex cols(void) const {return m_cols;}
254 inline void conservativeResize(DenseIndex size, DenseIndex, DenseIndex cols)
255 {
256 m_data = internal::conditional_aligned_realloc_new_auto<T,(_Options&DontAlign)==0>(m_data, size, _Rows*m_cols);
257 m_cols = cols;
258 }
259 EIGEN_STRONG_INLINE void resize(DenseIndex size, DenseIndex, DenseIndex cols)
260 {
261 if(size != _Rows*m_cols)
262 {
263 internal::conditional_aligned_delete_auto<T,(_Options&DontAlign)==0>(m_data, _Rows*m_cols);
264 if (size)
265 m_data = internal::conditional_aligned_new_auto<T,(_Options&DontAlign)==0>(size);
266 else
267 m_data = 0;
268 EIGEN_INTERNAL_DENSE_STORAGE_CTOR_PLUGIN
269 }
270 m_cols = cols;
271 }
272 inline const T *data() const { return m_data; }
273 inline T *data() { return m_data; }
274};
275
276// matrix with dynamic height and fixed width (so that matrix has dynamic size).
277template<typename T, int _Cols, int _Options> class DenseStorage<T, Dynamic, Dynamic, _Cols, _Options>
278{
279 T *m_data;
280 DenseIndex m_rows;
281 public:
282 inline explicit DenseStorage() : m_data(0), m_rows(0) {}
283 inline DenseStorage(internal::constructor_without_unaligned_array_assert) : m_data(0), m_rows(0) {}
284 inline DenseStorage(DenseIndex size, DenseIndex rows, DenseIndex) : m_data(internal::conditional_aligned_new_auto<T,(_Options&DontAlign)==0>(size)), m_rows(rows)
285 { EIGEN_INTERNAL_DENSE_STORAGE_CTOR_PLUGIN }
286 inline ~DenseStorage() { internal::conditional_aligned_delete_auto<T,(_Options&DontAlign)==0>(m_data, _Cols*m_rows); }
287 inline void swap(DenseStorage& other) { std::swap(m_data,other.m_data); std::swap(m_rows,other.m_rows); }
288 inline DenseIndex rows(void) const {return m_rows;}
289 static inline DenseIndex cols(void) {return _Cols;}
290 inline void conservativeResize(DenseIndex size, DenseIndex rows, DenseIndex)
291 {
292 m_data = internal::conditional_aligned_realloc_new_auto<T,(_Options&DontAlign)==0>(m_data, size, m_rows*_Cols);
293 m_rows = rows;
294 }
295 EIGEN_STRONG_INLINE void resize(DenseIndex size, DenseIndex rows, DenseIndex)
296 {
297 if(size != m_rows*_Cols)
298 {
299 internal::conditional_aligned_delete_auto<T,(_Options&DontAlign)==0>(m_data, _Cols*m_rows);
300 if (size)
301 m_data = internal::conditional_aligned_new_auto<T,(_Options&DontAlign)==0>(size);
302 else
303 m_data = 0;
304 EIGEN_INTERNAL_DENSE_STORAGE_CTOR_PLUGIN
305 }
306 m_rows = rows;
307 }
308 inline const T *data() const { return m_data; }
309 inline T *data() { return m_data; }
310};
311
312} // end namespace Eigen
313
314#endif // EIGEN_MATRIX_H
Definition LDLT.h:18