1
// This file is part of Eigen, a lightweight C++ template library
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>
8
// Eigen is free software; you can redistribute it and/or
9
// modify it under the terms of the GNU Lesser General Public
10
// License as published by the Free Software Foundation; either
11
// version 3 of the License, or (at your option) any later version.
13
// Alternatively, you can redistribute it and/or
14
// modify it under the terms of the GNU General Public License as
15
// published by the Free Software Foundation; either version 2 of
16
// the License, or (at your option) any later version.
18
// Eigen is distributed in the hope that it will be useful, but WITHOUT ANY
19
// WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
20
// FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the
21
// GNU General Public License for more details.
23
// You should have received a copy of the GNU Lesser General Public
24
// License and a copy of the GNU General Public License along with
25
// Eigen. If not, see <http://www.gnu.org/licenses/>.
27
#ifndef EIGEN_MATRIXSTORAGE_H
28
#define EIGEN_MATRIXSTORAGE_H
30
#ifdef EIGEN_DENSE_STORAGE_CTOR_PLUGIN
31
#define EIGEN_INTERNAL_DENSE_STORAGE_CTOR_PLUGIN EIGEN_DENSE_STORAGE_CTOR_PLUGIN;
33
#define EIGEN_INTERNAL_DENSE_STORAGE_CTOR_PLUGIN
38
struct constructor_without_unaligned_array_assert {};
41
* Static array. If the MatrixOrArrayOptions require auto-alignment, the array will be automatically aligned:
42
* to 16 bytes boundary if the total size is a multiple of 16 bytes.
44
template <typename T, int Size, int MatrixOrArrayOptions,
45
int Alignment = (MatrixOrArrayOptions&DontAlign) ? 0
46
: (((Size*sizeof(T))%16)==0) ? 16
52
plain_array(constructor_without_unaligned_array_assert) {}
55
#ifdef EIGEN_DISABLE_UNALIGNED_ARRAY_ASSERT
56
#define EIGEN_MAKE_UNALIGNED_ARRAY_ASSERT(sizemask)
58
#define EIGEN_MAKE_UNALIGNED_ARRAY_ASSERT(sizemask) \
59
eigen_assert((reinterpret_cast<size_t>(array) & sizemask) == 0 \
60
&& "this assertion is explained here: " \
61
"http://eigen.tuxfamily.org/dox-devel/TopicUnalignedArrayAssert.html" \
62
" **** READ THIS WEB PAGE !!! ****");
65
template <typename T, int Size, int MatrixOrArrayOptions>
66
struct plain_array<T, Size, MatrixOrArrayOptions, 16>
68
EIGEN_USER_ALIGN16 T array[Size];
69
plain_array() { EIGEN_MAKE_UNALIGNED_ARRAY_ASSERT(0xf) }
70
plain_array(constructor_without_unaligned_array_assert) {}
73
template <typename T, int MatrixOrArrayOptions, int Alignment>
74
struct plain_array<T, 0, MatrixOrArrayOptions, Alignment>
76
EIGEN_USER_ALIGN16 T array[1];
78
plain_array(constructor_without_unaligned_array_assert) {}
81
} // end namespace internal
86
* \ingroup Core_Module
88
* \brief Stores the data of a matrix
90
* This class stores the data of fixed-size, dynamic-size or mixed matrices
91
* in a way as compact as possible.
95
template<typename T, int Size, int _Rows, int _Cols, int _Options> class DenseStorage;
97
// purely fixed-size matrix
98
template<typename T, int Size, int _Rows, int _Cols, int _Options> class DenseStorage
100
internal::plain_array<T,Size,_Options> m_data;
102
inline explicit DenseStorage() {}
103
inline DenseStorage(internal::constructor_without_unaligned_array_assert)
104
: m_data(internal::constructor_without_unaligned_array_assert()) {}
105
inline DenseStorage(DenseIndex,DenseIndex,DenseIndex) {}
106
inline void swap(DenseStorage& other) { std::swap(m_data,other.m_data); }
107
inline static DenseIndex rows(void) {return _Rows;}
108
inline static DenseIndex cols(void) {return _Cols;}
109
inline void conservativeResize(DenseIndex,DenseIndex,DenseIndex) {}
110
inline void resize(DenseIndex,DenseIndex,DenseIndex) {}
111
inline const T *data() const { return m_data.array; }
112
inline T *data() { return m_data.array; }
116
template<typename T, int _Rows, int _Cols, int _Options> class DenseStorage<T, 0, _Rows, _Cols, _Options>
119
inline explicit DenseStorage() {}
120
inline DenseStorage(internal::constructor_without_unaligned_array_assert) {}
121
inline DenseStorage(DenseIndex,DenseIndex,DenseIndex) {}
122
inline void swap(DenseStorage& ) {}
123
inline static DenseIndex rows(void) {return _Rows;}
124
inline static DenseIndex cols(void) {return _Cols;}
125
inline void conservativeResize(DenseIndex,DenseIndex,DenseIndex) {}
126
inline void resize(DenseIndex,DenseIndex,DenseIndex) {}
127
inline const T *data() const { return 0; }
128
inline T *data() { return 0; }
131
// dynamic-size matrix with fixed-size storage
132
template<typename T, int Size, int _Options> class DenseStorage<T, Size, Dynamic, Dynamic, _Options>
134
internal::plain_array<T,Size,_Options> m_data;
138
inline explicit DenseStorage() : m_rows(0), m_cols(0) {}
139
inline DenseStorage(internal::constructor_without_unaligned_array_assert)
140
: m_data(internal::constructor_without_unaligned_array_assert()), m_rows(0), m_cols(0) {}
141
inline DenseStorage(DenseIndex, DenseIndex rows, DenseIndex cols) : m_rows(rows), m_cols(cols) {}
142
inline void swap(DenseStorage& other)
143
{ std::swap(m_data,other.m_data); std::swap(m_rows,other.m_rows); std::swap(m_cols,other.m_cols); }
144
inline DenseIndex rows(void) const {return m_rows;}
145
inline DenseIndex cols(void) const {return m_cols;}
146
inline void conservativeResize(DenseIndex, DenseIndex rows, DenseIndex cols) { m_rows = rows; m_cols = cols; }
147
inline void resize(DenseIndex, DenseIndex rows, DenseIndex cols) { m_rows = rows; m_cols = cols; }
148
inline const T *data() const { return m_data.array; }
149
inline T *data() { return m_data.array; }
152
// dynamic-size matrix with fixed-size storage and fixed width
153
template<typename T, int Size, int _Cols, int _Options> class DenseStorage<T, Size, Dynamic, _Cols, _Options>
155
internal::plain_array<T,Size,_Options> m_data;
158
inline explicit DenseStorage() : m_rows(0) {}
159
inline DenseStorage(internal::constructor_without_unaligned_array_assert)
160
: m_data(internal::constructor_without_unaligned_array_assert()), m_rows(0) {}
161
inline DenseStorage(DenseIndex, DenseIndex rows, DenseIndex) : m_rows(rows) {}
162
inline void swap(DenseStorage& other) { std::swap(m_data,other.m_data); std::swap(m_rows,other.m_rows); }
163
inline DenseIndex rows(void) const {return m_rows;}
164
inline DenseIndex cols(void) const {return _Cols;}
165
inline void conservativeResize(DenseIndex, DenseIndex rows, DenseIndex) { m_rows = rows; }
166
inline void resize(DenseIndex, DenseIndex rows, DenseIndex) { m_rows = rows; }
167
inline const T *data() const { return m_data.array; }
168
inline T *data() { return m_data.array; }
171
// dynamic-size matrix with fixed-size storage and fixed height
172
template<typename T, int Size, int _Rows, int _Options> class DenseStorage<T, Size, _Rows, Dynamic, _Options>
174
internal::plain_array<T,Size,_Options> m_data;
177
inline explicit DenseStorage() : m_cols(0) {}
178
inline DenseStorage(internal::constructor_without_unaligned_array_assert)
179
: m_data(internal::constructor_without_unaligned_array_assert()), m_cols(0) {}
180
inline DenseStorage(DenseIndex, DenseIndex, DenseIndex cols) : m_cols(cols) {}
181
inline void swap(DenseStorage& other) { std::swap(m_data,other.m_data); std::swap(m_cols,other.m_cols); }
182
inline DenseIndex rows(void) const {return _Rows;}
183
inline DenseIndex cols(void) const {return m_cols;}
184
inline void conservativeResize(DenseIndex, DenseIndex, DenseIndex cols) { m_cols = cols; }
185
inline void resize(DenseIndex, DenseIndex, DenseIndex cols) { m_cols = cols; }
186
inline const T *data() const { return m_data.array; }
187
inline T *data() { return m_data.array; }
190
// purely dynamic matrix.
191
template<typename T, int _Options> class DenseStorage<T, Dynamic, Dynamic, Dynamic, _Options>
197
inline explicit DenseStorage() : m_data(0), m_rows(0), m_cols(0) {}
198
inline DenseStorage(internal::constructor_without_unaligned_array_assert)
199
: m_data(0), m_rows(0), m_cols(0) {}
200
inline DenseStorage(DenseIndex size, DenseIndex rows, DenseIndex cols)
201
: m_data(internal::conditional_aligned_new_auto<T,(_Options&DontAlign)==0>(size)), m_rows(rows), m_cols(cols)
202
{ EIGEN_INTERNAL_DENSE_STORAGE_CTOR_PLUGIN }
203
inline ~DenseStorage() { internal::conditional_aligned_delete_auto<T,(_Options&DontAlign)==0>(m_data, m_rows*m_cols); }
204
inline void swap(DenseStorage& other)
205
{ std::swap(m_data,other.m_data); std::swap(m_rows,other.m_rows); std::swap(m_cols,other.m_cols); }
206
inline DenseIndex rows(void) const {return m_rows;}
207
inline DenseIndex cols(void) const {return m_cols;}
208
inline void conservativeResize(DenseIndex size, DenseIndex rows, DenseIndex cols)
210
m_data = internal::conditional_aligned_realloc_new_auto<T,(_Options&DontAlign)==0>(m_data, size, m_rows*m_cols);
214
void resize(DenseIndex size, DenseIndex rows, DenseIndex cols)
216
if(size != m_rows*m_cols)
218
internal::conditional_aligned_delete_auto<T,(_Options&DontAlign)==0>(m_data, m_rows*m_cols);
220
m_data = internal::conditional_aligned_new_auto<T,(_Options&DontAlign)==0>(size);
223
EIGEN_INTERNAL_DENSE_STORAGE_CTOR_PLUGIN
228
inline const T *data() const { return m_data; }
229
inline T *data() { return m_data; }
232
// matrix with dynamic width and fixed height (so that matrix has dynamic size).
233
template<typename T, int _Rows, int _Options> class DenseStorage<T, Dynamic, _Rows, Dynamic, _Options>
238
inline explicit DenseStorage() : m_data(0), m_cols(0) {}
239
inline DenseStorage(internal::constructor_without_unaligned_array_assert) : m_data(0), m_cols(0) {}
240
inline DenseStorage(DenseIndex size, DenseIndex, DenseIndex cols) : m_data(internal::conditional_aligned_new_auto<T,(_Options&DontAlign)==0>(size)), m_cols(cols)
241
{ EIGEN_INTERNAL_DENSE_STORAGE_CTOR_PLUGIN }
242
inline ~DenseStorage() { internal::conditional_aligned_delete_auto<T,(_Options&DontAlign)==0>(m_data, _Rows*m_cols); }
243
inline void swap(DenseStorage& other) { std::swap(m_data,other.m_data); std::swap(m_cols,other.m_cols); }
244
inline static DenseIndex rows(void) {return _Rows;}
245
inline DenseIndex cols(void) const {return m_cols;}
246
inline void conservativeResize(DenseIndex size, DenseIndex, DenseIndex cols)
248
m_data = internal::conditional_aligned_realloc_new_auto<T,(_Options&DontAlign)==0>(m_data, size, _Rows*m_cols);
251
EIGEN_STRONG_INLINE void resize(DenseIndex size, DenseIndex, DenseIndex cols)
253
if(size != _Rows*m_cols)
255
internal::conditional_aligned_delete_auto<T,(_Options&DontAlign)==0>(m_data, _Rows*m_cols);
257
m_data = internal::conditional_aligned_new_auto<T,(_Options&DontAlign)==0>(size);
260
EIGEN_INTERNAL_DENSE_STORAGE_CTOR_PLUGIN
264
inline const T *data() const { return m_data; }
265
inline T *data() { return m_data; }
268
// matrix with dynamic height and fixed width (so that matrix has dynamic size).
269
template<typename T, int _Cols, int _Options> class DenseStorage<T, Dynamic, Dynamic, _Cols, _Options>
274
inline explicit DenseStorage() : m_data(0), m_rows(0) {}
275
inline DenseStorage(internal::constructor_without_unaligned_array_assert) : m_data(0), m_rows(0) {}
276
inline DenseStorage(DenseIndex size, DenseIndex rows, DenseIndex) : m_data(internal::conditional_aligned_new_auto<T,(_Options&DontAlign)==0>(size)), m_rows(rows)
277
{ EIGEN_INTERNAL_DENSE_STORAGE_CTOR_PLUGIN }
278
inline ~DenseStorage() { internal::conditional_aligned_delete_auto<T,(_Options&DontAlign)==0>(m_data, _Cols*m_rows); }
279
inline void swap(DenseStorage& other) { std::swap(m_data,other.m_data); std::swap(m_rows,other.m_rows); }
280
inline DenseIndex rows(void) const {return m_rows;}
281
inline static DenseIndex cols(void) {return _Cols;}
282
inline void conservativeResize(DenseIndex size, DenseIndex rows, DenseIndex)
284
m_data = internal::conditional_aligned_realloc_new_auto<T,(_Options&DontAlign)==0>(m_data, size, m_rows*_Cols);
287
EIGEN_STRONG_INLINE void resize(DenseIndex size, DenseIndex rows, DenseIndex)
289
if(size != m_rows*_Cols)
291
internal::conditional_aligned_delete_auto<T,(_Options&DontAlign)==0>(m_data, _Cols*m_rows);
293
m_data = internal::conditional_aligned_new_auto<T,(_Options&DontAlign)==0>(size);
296
EIGEN_INTERNAL_DENSE_STORAGE_CTOR_PLUGIN
300
inline const T *data() const { return m_data; }
301
inline T *data() { return m_data; }
304
#endif // EIGEN_MATRIX_H