• Home
  • Line#
  • Scopes#
  • Navigate#
  • Raw
  • Download
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 
21 namespace Eigen {
22 
23 namespace internal {
24 
25 struct constructor_without_unaligned_array_assert {};
26 
27 /** \internal
28   * Static array. If the MatrixOrArrayOptions require auto-alignment, the array will be automatically aligned:
29   * to 16 bytes boundary if the total size is a multiple of 16 bytes.
30   */
31 template <typename T, int Size, int MatrixOrArrayOptions,
32           int Alignment = (MatrixOrArrayOptions&DontAlign) ? 0
33                         : (((Size*sizeof(T))%16)==0) ? 16
34                         : 0 >
35 struct plain_array
36 {
37   T array[Size];
plain_arrayplain_array38   plain_array() {}
plain_arrayplain_array39   plain_array(constructor_without_unaligned_array_assert) {}
40 };
41 
42 #ifdef EIGEN_DISABLE_UNALIGNED_ARRAY_ASSERT
43   #define EIGEN_MAKE_UNALIGNED_ARRAY_ASSERT(sizemask)
44 #else
45   #define EIGEN_MAKE_UNALIGNED_ARRAY_ASSERT(sizemask) \
46     eigen_assert((reinterpret_cast<size_t>(array) & sizemask) == 0 \
47               && "this assertion is explained here: " \
48               "http://eigen.tuxfamily.org/dox-devel/TopicUnalignedArrayAssert.html" \
49               " **** READ THIS WEB PAGE !!! ****");
50 #endif
51 
52 template <typename T, int Size, int MatrixOrArrayOptions>
53 struct plain_array<T, Size, MatrixOrArrayOptions, 16>
54 {
55   EIGEN_USER_ALIGN16 T array[Size];
56   plain_array() { EIGEN_MAKE_UNALIGNED_ARRAY_ASSERT(0xf) }
57   plain_array(constructor_without_unaligned_array_assert) {}
58 };
59 
60 template <typename T, int MatrixOrArrayOptions, int Alignment>
61 struct plain_array<T, 0, MatrixOrArrayOptions, Alignment>
62 {
63   EIGEN_USER_ALIGN16 T array[1];
64   plain_array() {}
65   plain_array(constructor_without_unaligned_array_assert) {}
66 };
67 
68 } // end namespace internal
69 
70 /** \internal
71   *
72   * \class DenseStorage
73   * \ingroup Core_Module
74   *
75   * \brief Stores the data of a matrix
76   *
77   * This class stores the data of fixed-size, dynamic-size or mixed matrices
78   * in a way as compact as possible.
79   *
80   * \sa Matrix
81   */
82 template<typename T, int Size, int _Rows, int _Cols, int _Options> class DenseStorage;
83 
84 // purely fixed-size matrix
85 template<typename T, int Size, int _Rows, int _Cols, int _Options> class DenseStorage
86 {
87     internal::plain_array<T,Size,_Options> m_data;
88   public:
89     inline explicit DenseStorage() {}
90     inline DenseStorage(internal::constructor_without_unaligned_array_assert)
91       : m_data(internal::constructor_without_unaligned_array_assert()) {}
92     inline DenseStorage(DenseIndex,DenseIndex,DenseIndex) {}
93     inline void swap(DenseStorage& other) { std::swap(m_data,other.m_data); }
94     static inline DenseIndex rows(void) {return _Rows;}
95     static inline DenseIndex cols(void) {return _Cols;}
96     inline void conservativeResize(DenseIndex,DenseIndex,DenseIndex) {}
97     inline void resize(DenseIndex,DenseIndex,DenseIndex) {}
98     inline const T *data() const { return m_data.array; }
99     inline T *data() { return m_data.array; }
100 };
101 
102 // null matrix
103 template<typename T, int _Rows, int _Cols, int _Options> class DenseStorage<T, 0, _Rows, _Cols, _Options>
104 {
105   public:
106     inline explicit DenseStorage() {}
107     inline DenseStorage(internal::constructor_without_unaligned_array_assert) {}
108     inline DenseStorage(DenseIndex,DenseIndex,DenseIndex) {}
109     inline void swap(DenseStorage& ) {}
110     static inline DenseIndex rows(void) {return _Rows;}
111     static inline DenseIndex cols(void) {return _Cols;}
112     inline void conservativeResize(DenseIndex,DenseIndex,DenseIndex) {}
113     inline void resize(DenseIndex,DenseIndex,DenseIndex) {}
114     inline const T *data() const { return 0; }
115     inline T *data() { return 0; }
116 };
117 
118 // more specializations for null matrices; these are necessary to resolve ambiguities
119 template<typename T, int _Options> class DenseStorage<T, 0, Dynamic, Dynamic, _Options>
120 : public DenseStorage<T, 0, 0, 0, _Options> { };
121 
122 template<typename T, int _Rows, int _Options> class DenseStorage<T, 0, _Rows, Dynamic, _Options>
123 : public DenseStorage<T, 0, 0, 0, _Options> { };
124 
125 template<typename T, int _Cols, int _Options> class DenseStorage<T, 0, Dynamic, _Cols, _Options>
126 : public DenseStorage<T, 0, 0, 0, _Options> { };
127 
128 // dynamic-size matrix with fixed-size storage
129 template<typename T, int Size, int _Options> class DenseStorage<T, Size, Dynamic, Dynamic, _Options>
130 {
131     internal::plain_array<T,Size,_Options> m_data;
132     DenseIndex m_rows;
133     DenseIndex m_cols;
134   public:
135     inline explicit DenseStorage() : m_rows(0), m_cols(0) {}
136     inline DenseStorage(internal::constructor_without_unaligned_array_assert)
137       : m_data(internal::constructor_without_unaligned_array_assert()), m_rows(0), m_cols(0) {}
138     inline DenseStorage(DenseIndex, DenseIndex rows, DenseIndex cols) : m_rows(rows), m_cols(cols) {}
139     inline void swap(DenseStorage& other)
140     { std::swap(m_data,other.m_data); std::swap(m_rows,other.m_rows); std::swap(m_cols,other.m_cols); }
141     inline DenseIndex rows(void) const {return m_rows;}
142     inline DenseIndex cols(void) const {return m_cols;}
143     inline void conservativeResize(DenseIndex, DenseIndex rows, DenseIndex cols) { m_rows = rows; m_cols = cols; }
144     inline void resize(DenseIndex, DenseIndex rows, DenseIndex cols) { m_rows = rows; m_cols = cols; }
145     inline const T *data() const { return m_data.array; }
146     inline T *data() { return m_data.array; }
147 };
148 
149 // dynamic-size matrix with fixed-size storage and fixed width
150 template<typename T, int Size, int _Cols, int _Options> class DenseStorage<T, Size, Dynamic, _Cols, _Options>
151 {
152     internal::plain_array<T,Size,_Options> m_data;
153     DenseIndex m_rows;
154   public:
155     inline explicit DenseStorage() : m_rows(0) {}
156     inline DenseStorage(internal::constructor_without_unaligned_array_assert)
157       : m_data(internal::constructor_without_unaligned_array_assert()), m_rows(0) {}
158     inline DenseStorage(DenseIndex, DenseIndex rows, DenseIndex) : m_rows(rows) {}
159     inline void swap(DenseStorage& other) { std::swap(m_data,other.m_data); std::swap(m_rows,other.m_rows); }
160     inline DenseIndex rows(void) const {return m_rows;}
161     inline DenseIndex cols(void) const {return _Cols;}
162     inline void conservativeResize(DenseIndex, DenseIndex rows, DenseIndex) { m_rows = rows; }
163     inline void resize(DenseIndex, DenseIndex rows, DenseIndex) { m_rows = rows; }
164     inline const T *data() const { return m_data.array; }
165     inline T *data() { return m_data.array; }
166 };
167 
168 // dynamic-size matrix with fixed-size storage and fixed height
169 template<typename T, int Size, int _Rows, int _Options> class DenseStorage<T, Size, _Rows, Dynamic, _Options>
170 {
171     internal::plain_array<T,Size,_Options> m_data;
172     DenseIndex m_cols;
173   public:
174     inline explicit DenseStorage() : m_cols(0) {}
175     inline DenseStorage(internal::constructor_without_unaligned_array_assert)
176       : m_data(internal::constructor_without_unaligned_array_assert()), m_cols(0) {}
177     inline DenseStorage(DenseIndex, DenseIndex, DenseIndex cols) : m_cols(cols) {}
178     inline void swap(DenseStorage& other) { std::swap(m_data,other.m_data); std::swap(m_cols,other.m_cols); }
179     inline DenseIndex rows(void) const {return _Rows;}
180     inline DenseIndex cols(void) const {return m_cols;}
181     inline void conservativeResize(DenseIndex, DenseIndex, DenseIndex cols) { m_cols = cols; }
182     inline void resize(DenseIndex, DenseIndex, DenseIndex cols) { m_cols = cols; }
183     inline const T *data() const { return m_data.array; }
184     inline T *data() { return m_data.array; }
185 };
186 
187 // purely dynamic matrix.
188 template<typename T, int _Options> class DenseStorage<T, Dynamic, Dynamic, Dynamic, _Options>
189 {
190     T *m_data;
191     DenseIndex m_rows;
192     DenseIndex m_cols;
193   public:
194     inline explicit DenseStorage() : m_data(0), m_rows(0), m_cols(0) {}
195     inline DenseStorage(internal::constructor_without_unaligned_array_assert)
196        : m_data(0), m_rows(0), m_cols(0) {}
197     inline DenseStorage(DenseIndex size, DenseIndex rows, DenseIndex cols)
198       : m_data(internal::conditional_aligned_new_auto<T,(_Options&DontAlign)==0>(size)), m_rows(rows), m_cols(cols)
199     { EIGEN_INTERNAL_DENSE_STORAGE_CTOR_PLUGIN }
200     inline ~DenseStorage() { internal::conditional_aligned_delete_auto<T,(_Options&DontAlign)==0>(m_data, m_rows*m_cols); }
201     inline void swap(DenseStorage& other)
202     { std::swap(m_data,other.m_data); std::swap(m_rows,other.m_rows); std::swap(m_cols,other.m_cols); }
203     inline DenseIndex rows(void) const {return m_rows;}
204     inline DenseIndex cols(void) const {return m_cols;}
205     inline void conservativeResize(DenseIndex size, DenseIndex rows, DenseIndex cols)
206     {
207       m_data = internal::conditional_aligned_realloc_new_auto<T,(_Options&DontAlign)==0>(m_data, size, m_rows*m_cols);
208       m_rows = rows;
209       m_cols = cols;
210     }
211     void resize(DenseIndex size, DenseIndex rows, DenseIndex cols)
212     {
213       if(size != m_rows*m_cols)
214       {
215         internal::conditional_aligned_delete_auto<T,(_Options&DontAlign)==0>(m_data, m_rows*m_cols);
216         if (size)
217           m_data = internal::conditional_aligned_new_auto<T,(_Options&DontAlign)==0>(size);
218         else
219           m_data = 0;
220         EIGEN_INTERNAL_DENSE_STORAGE_CTOR_PLUGIN
221       }
222       m_rows = rows;
223       m_cols = cols;
224     }
225     inline const T *data() const { return m_data; }
226     inline T *data() { return m_data; }
227 };
228 
229 // matrix with dynamic width and fixed height (so that matrix has dynamic size).
230 template<typename T, int _Rows, int _Options> class DenseStorage<T, Dynamic, _Rows, Dynamic, _Options>
231 {
232     T *m_data;
233     DenseIndex m_cols;
234   public:
235     inline explicit DenseStorage() : m_data(0), m_cols(0) {}
236     inline DenseStorage(internal::constructor_without_unaligned_array_assert) : m_data(0), m_cols(0) {}
237     inline DenseStorage(DenseIndex size, DenseIndex, DenseIndex cols) : m_data(internal::conditional_aligned_new_auto<T,(_Options&DontAlign)==0>(size)), m_cols(cols)
238     { EIGEN_INTERNAL_DENSE_STORAGE_CTOR_PLUGIN }
239     inline ~DenseStorage() { internal::conditional_aligned_delete_auto<T,(_Options&DontAlign)==0>(m_data, _Rows*m_cols); }
240     inline void swap(DenseStorage& other) { std::swap(m_data,other.m_data); std::swap(m_cols,other.m_cols); }
241     static inline DenseIndex rows(void) {return _Rows;}
242     inline DenseIndex cols(void) const {return m_cols;}
243     inline void conservativeResize(DenseIndex size, DenseIndex, DenseIndex cols)
244     {
245       m_data = internal::conditional_aligned_realloc_new_auto<T,(_Options&DontAlign)==0>(m_data, size, _Rows*m_cols);
246       m_cols = cols;
247     }
248     EIGEN_STRONG_INLINE void resize(DenseIndex size, DenseIndex, DenseIndex cols)
249     {
250       if(size != _Rows*m_cols)
251       {
252         internal::conditional_aligned_delete_auto<T,(_Options&DontAlign)==0>(m_data, _Rows*m_cols);
253         if (size)
254           m_data = internal::conditional_aligned_new_auto<T,(_Options&DontAlign)==0>(size);
255         else
256           m_data = 0;
257         EIGEN_INTERNAL_DENSE_STORAGE_CTOR_PLUGIN
258       }
259       m_cols = cols;
260     }
261     inline const T *data() const { return m_data; }
262     inline T *data() { return m_data; }
263 };
264 
265 // matrix with dynamic height and fixed width (so that matrix has dynamic size).
266 template<typename T, int _Cols, int _Options> class DenseStorage<T, Dynamic, Dynamic, _Cols, _Options>
267 {
268     T *m_data;
269     DenseIndex m_rows;
270   public:
271     inline explicit DenseStorage() : m_data(0), m_rows(0) {}
272     inline DenseStorage(internal::constructor_without_unaligned_array_assert) : m_data(0), m_rows(0) {}
273     inline DenseStorage(DenseIndex size, DenseIndex rows, DenseIndex) : m_data(internal::conditional_aligned_new_auto<T,(_Options&DontAlign)==0>(size)), m_rows(rows)
274     { EIGEN_INTERNAL_DENSE_STORAGE_CTOR_PLUGIN }
275     inline ~DenseStorage() { internal::conditional_aligned_delete_auto<T,(_Options&DontAlign)==0>(m_data, _Cols*m_rows); }
276     inline void swap(DenseStorage& other) { std::swap(m_data,other.m_data); std::swap(m_rows,other.m_rows); }
277     inline DenseIndex rows(void) const {return m_rows;}
278     static inline DenseIndex cols(void) {return _Cols;}
279     inline void conservativeResize(DenseIndex size, DenseIndex rows, DenseIndex)
280     {
281       m_data = internal::conditional_aligned_realloc_new_auto<T,(_Options&DontAlign)==0>(m_data, size, m_rows*_Cols);
282       m_rows = rows;
283     }
284     EIGEN_STRONG_INLINE void resize(DenseIndex size, DenseIndex rows, DenseIndex)
285     {
286       if(size != m_rows*_Cols)
287       {
288         internal::conditional_aligned_delete_auto<T,(_Options&DontAlign)==0>(m_data, _Cols*m_rows);
289         if (size)
290           m_data = internal::conditional_aligned_new_auto<T,(_Options&DontAlign)==0>(size);
291         else
292           m_data = 0;
293         EIGEN_INTERNAL_DENSE_STORAGE_CTOR_PLUGIN
294       }
295       m_rows = rows;
296     }
297     inline const T *data() const { return m_data; }
298     inline T *data() { return m_data; }
299 };
300 
301 } // end namespace Eigen
302 
303 #endif // EIGEN_MATRIX_H
304