Eigen  3.2.7
 All Classes Namespaces Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
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 <[email protected]>
5 // Copyright (C) 2006-2009 Benoit Jacob <[email protected]>
6 // Copyright (C) 2010 Hauke Heibel <[email protected]>
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 
21 namespace Eigen {
22 
23 namespace internal {
24 
25 struct constructor_without_unaligned_array_assert {};
26 
27 template<typename T, int Size> void check_static_allocation_size()
28 {
29  // if EIGEN_STACK_ALLOCATION_LIMIT is defined to 0, then no limit
30  #if EIGEN_STACK_ALLOCATION_LIMIT
31  EIGEN_STATIC_ASSERT(Size * sizeof(T) <= EIGEN_STACK_ALLOCATION_LIMIT, OBJECT_ALLOCATED_ON_STACK_IS_TOO_BIG);
32  #endif
33 }
34 
39 template <typename T, int Size, int MatrixOrArrayOptions,
40  int Alignment = (MatrixOrArrayOptions&DontAlign) ? 0
41  : (((Size*sizeof(T))%16)==0) ? 16
42  : 0 >
43 struct plain_array
44 {
45  T array[Size];
46 
47  plain_array()
48  {
49  check_static_allocation_size<T,Size>();
50  }
51 
52  plain_array(constructor_without_unaligned_array_assert)
53  {
54  check_static_allocation_size<T,Size>();
55  }
56 };
57 
58 #if defined(EIGEN_DISABLE_UNALIGNED_ARRAY_ASSERT)
59  #define EIGEN_MAKE_UNALIGNED_ARRAY_ASSERT(sizemask)
60 #elif EIGEN_GNUC_AT_LEAST(4,7)
61  // 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.
62  // See this bug report: http://gcc.gnu.org/bugzilla/show_bug.cgi?id=53900
63  // Hiding the origin of the array pointer behind a function argument seems to do the trick even if the function is inlined:
64  template<typename PtrType>
65  EIGEN_ALWAYS_INLINE PtrType eigen_unaligned_array_assert_workaround_gcc47(PtrType array) { return array; }
66  #define EIGEN_MAKE_UNALIGNED_ARRAY_ASSERT(sizemask) \
67  eigen_assert((reinterpret_cast<size_t>(eigen_unaligned_array_assert_workaround_gcc47(array)) & sizemask) == 0 \
68  && "this assertion is explained here: " \
69  "http://eigen.tuxfamily.org/dox-devel/group__TopicUnalignedArrayAssert.html" \
70  " **** READ THIS WEB PAGE !!! ****");
71 #else
72  #define EIGEN_MAKE_UNALIGNED_ARRAY_ASSERT(sizemask) \
73  eigen_assert((reinterpret_cast<size_t>(array) & sizemask) == 0 \
74  && "this assertion is explained here: " \
75  "http://eigen.tuxfamily.org/dox-devel/group__TopicUnalignedArrayAssert.html" \
76  " **** READ THIS WEB PAGE !!! ****");
77 #endif
78 
79 template <typename T, int Size, int MatrixOrArrayOptions>
80 struct plain_array<T, Size, MatrixOrArrayOptions, 16>
81 {
82  EIGEN_USER_ALIGN16 T array[Size];
83 
84  plain_array()
85  {
86  EIGEN_MAKE_UNALIGNED_ARRAY_ASSERT(0xf);
87  check_static_allocation_size<T,Size>();
88  }
89 
90  plain_array(constructor_without_unaligned_array_assert)
91  {
92  check_static_allocation_size<T,Size>();
93  }
94 };
95 
96 template <typename T, int MatrixOrArrayOptions, int Alignment>
97 struct plain_array<T, 0, MatrixOrArrayOptions, Alignment>
98 {
99  EIGEN_USER_ALIGN16 T array[1];
100  plain_array() {}
101  plain_array(constructor_without_unaligned_array_assert) {}
102 };
103 
104 } // end namespace internal
105 
118 template<typename T, int Size, int _Rows, int _Cols, int _Options> class DenseStorage;
119 
120 // purely fixed-size matrix
121 template<typename T, int Size, int _Rows, int _Cols, int _Options> class DenseStorage
122 {
123  internal::plain_array<T,Size,_Options> m_data;
124  public:
125  DenseStorage() {}
126  DenseStorage(internal::constructor_without_unaligned_array_assert)
127  : m_data(internal::constructor_without_unaligned_array_assert()) {}
128  DenseStorage(const DenseStorage& other) : m_data(other.m_data) {}
129  DenseStorage& operator=(const DenseStorage& other)
130  {
131  if (this != &other) m_data = other.m_data;
132  return *this;
133  }
134  DenseStorage(DenseIndex,DenseIndex,DenseIndex) {}
135  void swap(DenseStorage& other) { std::swap(m_data,other.m_data); }
136  static DenseIndex rows(void) {return _Rows;}
137  static DenseIndex cols(void) {return _Cols;}
138  void conservativeResize(DenseIndex,DenseIndex,DenseIndex) {}
139  void resize(DenseIndex,DenseIndex,DenseIndex) {}
140  const T *data() const { return m_data.array; }
141  T *data() { return m_data.array; }
142 };
143 
144 // null matrix
145 template<typename T, int _Rows, int _Cols, int _Options> class DenseStorage<T, 0, _Rows, _Cols, _Options>
146 {
147  public:
148  DenseStorage() {}
149  DenseStorage(internal::constructor_without_unaligned_array_assert) {}
150  DenseStorage(const DenseStorage&) {}
151  DenseStorage& operator=(const DenseStorage&) { return *this; }
152  DenseStorage(DenseIndex,DenseIndex,DenseIndex) {}
153  void swap(DenseStorage& ) {}
154  static DenseIndex rows(void) {return _Rows;}
155  static DenseIndex cols(void) {return _Cols;}
156  void conservativeResize(DenseIndex,DenseIndex,DenseIndex) {}
157  void resize(DenseIndex,DenseIndex,DenseIndex) {}
158  const T *data() const { return 0; }
159  T *data() { return 0; }
160 };
161 
162 // more specializations for null matrices; these are necessary to resolve ambiguities
163 template<typename T, int _Options> class DenseStorage<T, 0, Dynamic, Dynamic, _Options>
164 : public DenseStorage<T, 0, 0, 0, _Options> { };
165 
166 template<typename T, int _Rows, int _Options> class DenseStorage<T, 0, _Rows, Dynamic, _Options>
167 : public DenseStorage<T, 0, 0, 0, _Options> { };
168 
169 template<typename T, int _Cols, int _Options> class DenseStorage<T, 0, Dynamic, _Cols, _Options>
170 : public DenseStorage<T, 0, 0, 0, _Options> { };
171 
172 // dynamic-size matrix with fixed-size storage
173 template<typename T, int Size, int _Options> class DenseStorage<T, Size, Dynamic, Dynamic, _Options>
174 {
175  internal::plain_array<T,Size,_Options> m_data;
176  DenseIndex m_rows;
177  DenseIndex m_cols;
178  public:
179  DenseStorage() : m_rows(0), m_cols(0) {}
180  DenseStorage(internal::constructor_without_unaligned_array_assert)
181  : m_data(internal::constructor_without_unaligned_array_assert()), m_rows(0), m_cols(0) {}
182  DenseStorage(const DenseStorage& other) : m_data(other.m_data), m_rows(other.m_rows), m_cols(other.m_cols) {}
183  DenseStorage& operator=(const DenseStorage& other)
184  {
185  if (this != &other)
186  {
187  m_data = other.m_data;
188  m_rows = other.m_rows;
189  m_cols = other.m_cols;
190  }
191  return *this;
192  }
193  DenseStorage(DenseIndex, DenseIndex nbRows, DenseIndex nbCols) : m_rows(nbRows), m_cols(nbCols) {}
194  void swap(DenseStorage& other)
195  { std::swap(m_data,other.m_data); std::swap(m_rows,other.m_rows); std::swap(m_cols,other.m_cols); }
196  DenseIndex rows() const {return m_rows;}
197  DenseIndex cols() const {return m_cols;}
198  void conservativeResize(DenseIndex, DenseIndex nbRows, DenseIndex nbCols) { m_rows = nbRows; m_cols = nbCols; }
199  void resize(DenseIndex, DenseIndex nbRows, DenseIndex nbCols) { m_rows = nbRows; m_cols = nbCols; }
200  const T *data() const { return m_data.array; }
201  T *data() { return m_data.array; }
202 };
203 
204 // dynamic-size matrix with fixed-size storage and fixed width
205 template<typename T, int Size, int _Cols, int _Options> class DenseStorage<T, Size, Dynamic, _Cols, _Options>
206 {
207  internal::plain_array<T,Size,_Options> m_data;
208  DenseIndex m_rows;
209  public:
210  DenseStorage() : m_rows(0) {}
211  DenseStorage(internal::constructor_without_unaligned_array_assert)
212  : m_data(internal::constructor_without_unaligned_array_assert()), m_rows(0) {}
213  DenseStorage(const DenseStorage& other) : m_data(other.m_data), m_rows(other.m_rows) {}
214  DenseStorage& operator=(const DenseStorage& other)
215  {
216  if (this != &other)
217  {
218  m_data = other.m_data;
219  m_rows = other.m_rows;
220  }
221  return *this;
222  }
223  DenseStorage(DenseIndex, DenseIndex nbRows, DenseIndex) : m_rows(nbRows) {}
224  void swap(DenseStorage& other) { std::swap(m_data,other.m_data); std::swap(m_rows,other.m_rows); }
225  DenseIndex rows(void) const {return m_rows;}
226  DenseIndex cols(void) const {return _Cols;}
227  void conservativeResize(DenseIndex, DenseIndex nbRows, DenseIndex) { m_rows = nbRows; }
228  void resize(DenseIndex, DenseIndex nbRows, DenseIndex) { m_rows = nbRows; }
229  const T *data() const { return m_data.array; }
230  T *data() { return m_data.array; }
231 };
232 
233 // dynamic-size matrix with fixed-size storage and fixed height
234 template<typename T, int Size, int _Rows, int _Options> class DenseStorage<T, Size, _Rows, Dynamic, _Options>
235 {
236  internal::plain_array<T,Size,_Options> m_data;
237  DenseIndex m_cols;
238  public:
239  DenseStorage() : m_cols(0) {}
240  DenseStorage(internal::constructor_without_unaligned_array_assert)
241  : m_data(internal::constructor_without_unaligned_array_assert()), m_cols(0) {}
242  DenseStorage(const DenseStorage& other) : m_data(other.m_data), m_cols(other.m_cols) {}
243  DenseStorage& operator=(const DenseStorage& other)
244  {
245  if (this != &other)
246  {
247  m_data = other.m_data;
248  m_cols = other.m_cols;
249  }
250  return *this;
251  }
252  DenseStorage(DenseIndex, DenseIndex, DenseIndex nbCols) : m_cols(nbCols) {}
253  void swap(DenseStorage& other) { std::swap(m_data,other.m_data); std::swap(m_cols,other.m_cols); }
254  DenseIndex rows(void) const {return _Rows;}
255  DenseIndex cols(void) const {return m_cols;}
256  void conservativeResize(DenseIndex, DenseIndex, DenseIndex nbCols) { m_cols = nbCols; }
257  void resize(DenseIndex, DenseIndex, DenseIndex nbCols) { m_cols = nbCols; }
258  const T *data() const { return m_data.array; }
259  T *data() { return m_data.array; }
260 };
261 
262 // purely dynamic matrix.
263 template<typename T, int _Options> class DenseStorage<T, Dynamic, Dynamic, Dynamic, _Options>
264 {
265  T *m_data;
266  DenseIndex m_rows;
267  DenseIndex m_cols;
268  public:
269  DenseStorage() : m_data(0), m_rows(0), m_cols(0) {}
270  DenseStorage(internal::constructor_without_unaligned_array_assert)
271  : m_data(0), m_rows(0), m_cols(0) {}
272  DenseStorage(DenseIndex size, DenseIndex nbRows, DenseIndex nbCols)
273  : m_data(internal::conditional_aligned_new_auto<T,(_Options&DontAlign)==0>(size)), m_rows(nbRows), m_cols(nbCols)
274  { EIGEN_INTERNAL_DENSE_STORAGE_CTOR_PLUGIN }
275 #ifdef EIGEN_HAVE_RVALUE_REFERENCES
276  DenseStorage(DenseStorage&& other)
277  : m_data(std::move(other.m_data))
278  , m_rows(std::move(other.m_rows))
279  , m_cols(std::move(other.m_cols))
280  {
281  other.m_data = nullptr;
282  }
283  DenseStorage& operator=(DenseStorage&& other)
284  {
285  using std::swap;
286  swap(m_data, other.m_data);
287  swap(m_rows, other.m_rows);
288  swap(m_cols, other.m_cols);
289  return *this;
290  }
291 #endif
292  ~DenseStorage() { internal::conditional_aligned_delete_auto<T,(_Options&DontAlign)==0>(m_data, m_rows*m_cols); }
293  void swap(DenseStorage& other)
294  { std::swap(m_data,other.m_data); std::swap(m_rows,other.m_rows); std::swap(m_cols,other.m_cols); }
295  DenseIndex rows(void) const {return m_rows;}
296  DenseIndex cols(void) const {return m_cols;}
297  void conservativeResize(DenseIndex size, DenseIndex nbRows, DenseIndex nbCols)
298  {
299  m_data = internal::conditional_aligned_realloc_new_auto<T,(_Options&DontAlign)==0>(m_data, size, m_rows*m_cols);
300  m_rows = nbRows;
301  m_cols = nbCols;
302  }
303  void resize(DenseIndex size, DenseIndex nbRows, DenseIndex nbCols)
304  {
305  if(size != m_rows*m_cols)
306  {
307  internal::conditional_aligned_delete_auto<T,(_Options&DontAlign)==0>(m_data, m_rows*m_cols);
308  if (size)
309  m_data = internal::conditional_aligned_new_auto<T,(_Options&DontAlign)==0>(size);
310  else
311  m_data = 0;
312  EIGEN_INTERNAL_DENSE_STORAGE_CTOR_PLUGIN
313  }
314  m_rows = nbRows;
315  m_cols = nbCols;
316  }
317  const T *data() const { return m_data; }
318  T *data() { return m_data; }
319  private:
320  DenseStorage(const DenseStorage&);
321  DenseStorage& operator=(const DenseStorage&);
322 };
323 
324 // matrix with dynamic width and fixed height (so that matrix has dynamic size).
325 template<typename T, int _Rows, int _Options> class DenseStorage<T, Dynamic, _Rows, Dynamic, _Options>
326 {
327  T *m_data;
328  DenseIndex m_cols;
329  public:
330  DenseStorage() : m_data(0), m_cols(0) {}
331  DenseStorage(internal::constructor_without_unaligned_array_assert) : m_data(0), m_cols(0) {}
332  DenseStorage(DenseIndex size, DenseIndex, DenseIndex nbCols) : m_data(internal::conditional_aligned_new_auto<T,(_Options&DontAlign)==0>(size)), m_cols(nbCols)
333  { EIGEN_INTERNAL_DENSE_STORAGE_CTOR_PLUGIN }
334 #ifdef EIGEN_HAVE_RVALUE_REFERENCES
335  DenseStorage(DenseStorage&& other)
336  : m_data(std::move(other.m_data))
337  , m_cols(std::move(other.m_cols))
338  {
339  other.m_data = nullptr;
340  }
341  DenseStorage& operator=(DenseStorage&& other)
342  {
343  using std::swap;
344  swap(m_data, other.m_data);
345  swap(m_cols, other.m_cols);
346  return *this;
347  }
348 #endif
349  ~DenseStorage() { internal::conditional_aligned_delete_auto<T,(_Options&DontAlign)==0>(m_data, _Rows*m_cols); }
350  void swap(DenseStorage& other) { std::swap(m_data,other.m_data); std::swap(m_cols,other.m_cols); }
351  static DenseIndex rows(void) {return _Rows;}
352  DenseIndex cols(void) const {return m_cols;}
353  void conservativeResize(DenseIndex size, DenseIndex, DenseIndex nbCols)
354  {
355  m_data = internal::conditional_aligned_realloc_new_auto<T,(_Options&DontAlign)==0>(m_data, size, _Rows*m_cols);
356  m_cols = nbCols;
357  }
358  EIGEN_STRONG_INLINE void resize(DenseIndex size, DenseIndex, DenseIndex nbCols)
359  {
360  if(size != _Rows*m_cols)
361  {
362  internal::conditional_aligned_delete_auto<T,(_Options&DontAlign)==0>(m_data, _Rows*m_cols);
363  if (size)
364  m_data = internal::conditional_aligned_new_auto<T,(_Options&DontAlign)==0>(size);
365  else
366  m_data = 0;
367  EIGEN_INTERNAL_DENSE_STORAGE_CTOR_PLUGIN
368  }
369  m_cols = nbCols;
370  }
371  const T *data() const { return m_data; }
372  T *data() { return m_data; }
373  private:
374  DenseStorage(const DenseStorage&);
375  DenseStorage& operator=(const DenseStorage&);
376 };
377 
378 // matrix with dynamic height and fixed width (so that matrix has dynamic size).
379 template<typename T, int _Cols, int _Options> class DenseStorage<T, Dynamic, Dynamic, _Cols, _Options>
380 {
381  T *m_data;
382  DenseIndex m_rows;
383  public:
384  DenseStorage() : m_data(0), m_rows(0) {}
385  DenseStorage(internal::constructor_without_unaligned_array_assert) : m_data(0), m_rows(0) {}
386  DenseStorage(DenseIndex size, DenseIndex nbRows, DenseIndex) : m_data(internal::conditional_aligned_new_auto<T,(_Options&DontAlign)==0>(size)), m_rows(nbRows)
387  { EIGEN_INTERNAL_DENSE_STORAGE_CTOR_PLUGIN }
388 #ifdef EIGEN_HAVE_RVALUE_REFERENCES
389  DenseStorage(DenseStorage&& other)
390  : m_data(std::move(other.m_data))
391  , m_rows(std::move(other.m_rows))
392  {
393  other.m_data = nullptr;
394  }
395  DenseStorage& operator=(DenseStorage&& other)
396  {
397  using std::swap;
398  swap(m_data, other.m_data);
399  swap(m_rows, other.m_rows);
400  return *this;
401  }
402 #endif
403  ~DenseStorage() { internal::conditional_aligned_delete_auto<T,(_Options&DontAlign)==0>(m_data, _Cols*m_rows); }
404  void swap(DenseStorage& other) { std::swap(m_data,other.m_data); std::swap(m_rows,other.m_rows); }
405  DenseIndex rows(void) const {return m_rows;}
406  static DenseIndex cols(void) {return _Cols;}
407  void conservativeResize(DenseIndex size, DenseIndex nbRows, DenseIndex)
408  {
409  m_data = internal::conditional_aligned_realloc_new_auto<T,(_Options&DontAlign)==0>(m_data, size, m_rows*_Cols);
410  m_rows = nbRows;
411  }
412  EIGEN_STRONG_INLINE void resize(DenseIndex size, DenseIndex nbRows, DenseIndex)
413  {
414  if(size != m_rows*_Cols)
415  {
416  internal::conditional_aligned_delete_auto<T,(_Options&DontAlign)==0>(m_data, _Cols*m_rows);
417  if (size)
418  m_data = internal::conditional_aligned_new_auto<T,(_Options&DontAlign)==0>(size);
419  else
420  m_data = 0;
421  EIGEN_INTERNAL_DENSE_STORAGE_CTOR_PLUGIN
422  }
423  m_rows = nbRows;
424  }
425  const T *data() const { return m_data; }
426  T *data() { return m_data; }
427  private:
428  DenseStorage(const DenseStorage&);
429  DenseStorage& operator=(const DenseStorage&);
430 };
431 
432 } // end namespace Eigen
433 
434 #endif // EIGEN_MATRIX_H
Definition: Constants.h:270