Eigen  3.2.10
 
Loading...
Searching...
No Matches
Block.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-2010 Benoit Jacob <jacob.benoit.1@gmail.com>
6//
7// This Source Code Form is subject to the terms of the Mozilla
8// Public License v. 2.0. If a copy of the MPL was not distributed
9// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
10
11#ifndef EIGEN_BLOCK_H
12#define EIGEN_BLOCK_H
13
14namespace Eigen {
15
47
48namespace internal {
49template<typename XprType, int BlockRows, int BlockCols, bool InnerPanel>
50struct traits<Block<XprType, BlockRows, BlockCols, InnerPanel> > : traits<XprType>
51{
52 typedef typename traits<XprType>::Scalar Scalar;
53 typedef typename traits<XprType>::StorageKind StorageKind;
54 typedef typename traits<XprType>::XprKind XprKind;
55 typedef typename nested<XprType>::type XprTypeNested;
56 typedef typename remove_reference<XprTypeNested>::type _XprTypeNested;
57 enum{
58 MatrixRows = traits<XprType>::RowsAtCompileTime,
59 MatrixCols = traits<XprType>::ColsAtCompileTime,
60 RowsAtCompileTime = MatrixRows == 0 ? 0 : BlockRows,
61 ColsAtCompileTime = MatrixCols == 0 ? 0 : BlockCols,
62 MaxRowsAtCompileTime = BlockRows==0 ? 0
63 : RowsAtCompileTime != Dynamic ? int(RowsAtCompileTime)
64 : int(traits<XprType>::MaxRowsAtCompileTime),
65 MaxColsAtCompileTime = BlockCols==0 ? 0
66 : ColsAtCompileTime != Dynamic ? int(ColsAtCompileTime)
67 : int(traits<XprType>::MaxColsAtCompileTime),
68 XprTypeIsRowMajor = (int(traits<XprType>::Flags)&RowMajorBit) != 0,
69 IsRowMajor = (MaxRowsAtCompileTime==1&&MaxColsAtCompileTime!=1) ? 1
70 : (MaxColsAtCompileTime==1&&MaxRowsAtCompileTime!=1) ? 0
71 : XprTypeIsRowMajor,
72 HasSameStorageOrderAsXprType = (IsRowMajor == XprTypeIsRowMajor),
73 InnerSize = IsRowMajor ? int(ColsAtCompileTime) : int(RowsAtCompileTime),
74 InnerStrideAtCompileTime = HasSameStorageOrderAsXprType
75 ? int(inner_stride_at_compile_time<XprType>::ret)
76 : int(outer_stride_at_compile_time<XprType>::ret),
77 OuterStrideAtCompileTime = HasSameStorageOrderAsXprType
78 ? int(outer_stride_at_compile_time<XprType>::ret)
79 : int(inner_stride_at_compile_time<XprType>::ret),
80 MaskPacketAccessBit = (InnerSize == Dynamic || (InnerSize % packet_traits<Scalar>::size) == 0)
81 && (InnerStrideAtCompileTime == 1)
82 ? PacketAccessBit : 0,
83 MaskAlignedBit = (InnerPanel && (OuterStrideAtCompileTime!=Dynamic) && (((OuterStrideAtCompileTime * int(sizeof(Scalar))) % 16) == 0)) ? AlignedBit : 0,
84 FlagsLinearAccessBit = (RowsAtCompileTime == 1 || ColsAtCompileTime == 1 || (InnerPanel && (traits<XprType>::Flags&LinearAccessBit))) ? LinearAccessBit : 0,
85 FlagsLvalueBit = is_lvalue<XprType>::value ? LvalueBit : 0,
86 FlagsRowMajorBit = IsRowMajor ? RowMajorBit : 0,
87 Flags0 = traits<XprType>::Flags & ( (HereditaryBits & ~RowMajorBit) |
89 MaskPacketAccessBit |
90 MaskAlignedBit),
91 Flags = Flags0 | FlagsLinearAccessBit | FlagsLvalueBit | FlagsRowMajorBit
92 };
93};
94
95template<typename XprType, int BlockRows=Dynamic, int BlockCols=Dynamic, bool InnerPanel = false,
96 bool HasDirectAccess = internal::has_direct_access<XprType>::ret> class BlockImpl_dense;
97
98} // end namespace internal
99
100template<typename XprType, int BlockRows, int BlockCols, bool InnerPanel, typename StorageKind> class BlockImpl;
101
102template<typename XprType, int BlockRows, int BlockCols, bool InnerPanel> class Block
103 : public BlockImpl<XprType, BlockRows, BlockCols, InnerPanel, typename internal::traits<XprType>::StorageKind>
104{
105 typedef BlockImpl<XprType, BlockRows, BlockCols, InnerPanel, typename internal::traits<XprType>::StorageKind> Impl;
106 public:
107 //typedef typename Impl::Base Base;
108 typedef Impl Base;
109 EIGEN_GENERIC_PUBLIC_INTERFACE(Block)
110 EIGEN_INHERIT_ASSIGNMENT_OPERATORS(Block)
111
112
114 inline Block(XprType& xpr, Index i) : Impl(xpr,i)
115 {
116 eigen_assert( (i>=0) && (
117 ((BlockRows==1) && (BlockCols==XprType::ColsAtCompileTime) && i<xpr.rows())
118 ||((BlockRows==XprType::RowsAtCompileTime) && (BlockCols==1) && i<xpr.cols())));
119 }
120
123 inline Block(XprType& xpr, Index a_startRow, Index a_startCol)
124 : Impl(xpr, a_startRow, a_startCol)
125 {
126 EIGEN_STATIC_ASSERT(RowsAtCompileTime!=Dynamic && ColsAtCompileTime!=Dynamic,THIS_METHOD_IS_ONLY_FOR_FIXED_SIZE)
127 eigen_assert(a_startRow >= 0 && BlockRows >= 1 && a_startRow + BlockRows <= xpr.rows()
128 && a_startCol >= 0 && BlockCols >= 1 && a_startCol + BlockCols <= xpr.cols());
129 }
130
133 inline Block(XprType& xpr,
134 Index a_startRow, Index a_startCol,
135 Index blockRows, Index blockCols)
136 : Impl(xpr, a_startRow, a_startCol, blockRows, blockCols)
137 {
138 eigen_assert((RowsAtCompileTime==Dynamic || RowsAtCompileTime==blockRows)
139 && (ColsAtCompileTime==Dynamic || ColsAtCompileTime==blockCols));
140 eigen_assert(a_startRow >= 0 && blockRows >= 0 && a_startRow <= xpr.rows() - blockRows
141 && a_startCol >= 0 && blockCols >= 0 && a_startCol <= xpr.cols() - blockCols);
142 }
143};
144
145// The generic default implementation for dense block simplu forward to the internal::BlockImpl_dense
146// that must be specialized for direct and non-direct access...
147template<typename XprType, int BlockRows, int BlockCols, bool InnerPanel>
148class BlockImpl<XprType, BlockRows, BlockCols, InnerPanel, Dense>
149 : public internal::BlockImpl_dense<XprType, BlockRows, BlockCols, InnerPanel>
150{
151 typedef internal::BlockImpl_dense<XprType, BlockRows, BlockCols, InnerPanel> Impl;
152 typedef typename XprType::Index Index;
153 public:
154 typedef Impl Base;
155 EIGEN_INHERIT_ASSIGNMENT_OPERATORS(BlockImpl)
156 inline BlockImpl(XprType& xpr, Index i) : Impl(xpr,i) {}
157 inline BlockImpl(XprType& xpr, Index a_startRow, Index a_startCol) : Impl(xpr, a_startRow, a_startCol) {}
158 inline BlockImpl(XprType& xpr, Index a_startRow, Index a_startCol, Index blockRows, Index blockCols)
159 : Impl(xpr, a_startRow, a_startCol, blockRows, blockCols) {}
160};
161
162namespace internal {
163
165template<typename XprType, int BlockRows, int BlockCols, bool InnerPanel, bool HasDirectAccess> class BlockImpl_dense
166 : public internal::dense_xpr_base<Block<XprType, BlockRows, BlockCols, InnerPanel> >::type
167{
168 typedef Block<XprType, BlockRows, BlockCols, InnerPanel> BlockType;
169 public:
170
171 typedef typename internal::dense_xpr_base<BlockType>::type Base;
172 EIGEN_DENSE_PUBLIC_INTERFACE(BlockType)
173 EIGEN_INHERIT_ASSIGNMENT_OPERATORS(BlockImpl_dense)
174
175 class InnerIterator;
176
179 inline BlockImpl_dense(XprType& xpr, Index i)
180 : m_xpr(xpr),
181 // It is a row if and only if BlockRows==1 and BlockCols==XprType::ColsAtCompileTime,
182 // and it is a column if and only if BlockRows==XprType::RowsAtCompileTime and BlockCols==1,
183 // all other cases are invalid.
184 // The case a 1x1 matrix seems ambiguous, but the result is the same anyway.
185 m_startRow( (BlockRows==1) && (BlockCols==XprType::ColsAtCompileTime) ? i : 0),
186 m_startCol( (BlockRows==XprType::RowsAtCompileTime) && (BlockCols==1) ? i : 0),
187 m_blockRows(BlockRows==1 ? 1 : xpr.rows()),
188 m_blockCols(BlockCols==1 ? 1 : xpr.cols())
189 {}
190
193 inline BlockImpl_dense(XprType& xpr, Index a_startRow, Index a_startCol)
194 : m_xpr(xpr), m_startRow(a_startRow), m_startCol(a_startCol),
195 m_blockRows(BlockRows), m_blockCols(BlockCols)
196 {}
197
200 inline BlockImpl_dense(XprType& xpr,
201 Index a_startRow, Index a_startCol,
202 Index blockRows, Index blockCols)
203 : m_xpr(xpr), m_startRow(a_startRow), m_startCol(a_startCol),
204 m_blockRows(blockRows), m_blockCols(blockCols)
205 {}
206
207 inline Index rows() const { return m_blockRows.value(); }
208 inline Index cols() const { return m_blockCols.value(); }
209
210 inline Scalar& coeffRef(Index rowId, Index colId)
211 {
212 EIGEN_STATIC_ASSERT_LVALUE(XprType)
213 return m_xpr.const_cast_derived()
214 .coeffRef(rowId + m_startRow.value(), colId + m_startCol.value());
215 }
216
217 inline const Scalar& coeffRef(Index rowId, Index colId) const
218 {
219 return m_xpr.derived()
220 .coeffRef(rowId + m_startRow.value(), colId + m_startCol.value());
221 }
222
223 EIGEN_STRONG_INLINE const CoeffReturnType coeff(Index rowId, Index colId) const
224 {
225 return m_xpr.coeff(rowId + m_startRow.value(), colId + m_startCol.value());
226 }
227
228 inline Scalar& coeffRef(Index index)
229 {
230 EIGEN_STATIC_ASSERT_LVALUE(XprType)
231 return m_xpr.const_cast_derived()
232 .coeffRef(m_startRow.value() + (RowsAtCompileTime == 1 ? 0 : index),
233 m_startCol.value() + (RowsAtCompileTime == 1 ? index : 0));
234 }
235
236 inline const Scalar& coeffRef(Index index) const
237 {
238 return m_xpr.const_cast_derived()
239 .coeffRef(m_startRow.value() + (RowsAtCompileTime == 1 ? 0 : index),
240 m_startCol.value() + (RowsAtCompileTime == 1 ? index : 0));
241 }
242
243 inline const CoeffReturnType coeff(Index index) const
244 {
245 return m_xpr
246 .coeff(m_startRow.value() + (RowsAtCompileTime == 1 ? 0 : index),
247 m_startCol.value() + (RowsAtCompileTime == 1 ? index : 0));
248 }
249
250 template<int LoadMode>
251 inline PacketScalar packet(Index rowId, Index colId) const
252 {
253 return m_xpr.template packet<Unaligned>
254 (rowId + m_startRow.value(), colId + m_startCol.value());
255 }
256
257 template<int LoadMode>
258 inline void writePacket(Index rowId, Index colId, const PacketScalar& val)
259 {
260 m_xpr.const_cast_derived().template writePacket<Unaligned>
261 (rowId + m_startRow.value(), colId + m_startCol.value(), val);
262 }
263
264 template<int LoadMode>
265 inline PacketScalar packet(Index index) const
266 {
267 return m_xpr.template packet<Unaligned>
268 (m_startRow.value() + (RowsAtCompileTime == 1 ? 0 : index),
269 m_startCol.value() + (RowsAtCompileTime == 1 ? index : 0));
270 }
271
272 template<int LoadMode>
273 inline void writePacket(Index index, const PacketScalar& val)
274 {
275 m_xpr.const_cast_derived().template writePacket<Unaligned>
276 (m_startRow.value() + (RowsAtCompileTime == 1 ? 0 : index),
277 m_startCol.value() + (RowsAtCompileTime == 1 ? index : 0), val);
278 }
279
280 #ifdef EIGEN_PARSED_BY_DOXYGEN
282 inline const Scalar* data() const;
283 inline Index innerStride() const;
284 inline Index outerStride() const;
285 #endif
286
287 const typename internal::remove_all<typename XprType::Nested>::type& nestedExpression() const
288 {
289 return m_xpr;
290 }
291
292 Index startRow() const
293 {
294 return m_startRow.value();
295 }
296
297 Index startCol() const
298 {
299 return m_startCol.value();
300 }
301
302 protected:
303
304 const typename XprType::Nested m_xpr;
305 const internal::variable_if_dynamic<Index, XprType::RowsAtCompileTime == 1 ? 0 : Dynamic> m_startRow;
306 const internal::variable_if_dynamic<Index, XprType::ColsAtCompileTime == 1 ? 0 : Dynamic> m_startCol;
307 const internal::variable_if_dynamic<Index, RowsAtCompileTime> m_blockRows;
308 const internal::variable_if_dynamic<Index, ColsAtCompileTime> m_blockCols;
309};
310
312template<typename XprType, int BlockRows, int BlockCols, bool InnerPanel>
313class BlockImpl_dense<XprType,BlockRows,BlockCols, InnerPanel,true>
314 : public MapBase<Block<XprType, BlockRows, BlockCols, InnerPanel> >
315{
316 typedef Block<XprType, BlockRows, BlockCols, InnerPanel> BlockType;
317 public:
318
319 typedef MapBase<BlockType> Base;
320 EIGEN_DENSE_PUBLIC_INTERFACE(BlockType)
321 EIGEN_INHERIT_ASSIGNMENT_OPERATORS(BlockImpl_dense)
322
323
325 inline BlockImpl_dense(XprType& xpr, Index i)
326 : Base(internal::const_cast_ptr(&xpr.coeffRef(
327 (BlockRows==1) && (BlockCols==XprType::ColsAtCompileTime) ? i : 0,
328 (BlockRows==XprType::RowsAtCompileTime) && (BlockCols==1) ? i : 0)),
329 BlockRows==1 ? 1 : xpr.rows(),
330 BlockCols==1 ? 1 : xpr.cols()),
331 m_xpr(xpr)
332 {
333 init();
334 }
335
338 inline BlockImpl_dense(XprType& xpr, Index startRow, Index startCol)
339 : Base(internal::const_cast_ptr(&xpr.coeffRef(startRow,startCol))), m_xpr(xpr)
340 {
341 init();
342 }
343
346 inline BlockImpl_dense(XprType& xpr,
347 Index startRow, Index startCol,
348 Index blockRows, Index blockCols)
349 : Base(internal::const_cast_ptr(&xpr.coeffRef(startRow,startCol)), blockRows, blockCols),
350 m_xpr(xpr)
351 {
352 init();
353 }
354
355 const typename internal::remove_all<typename XprType::Nested>::type& nestedExpression() const
356 {
357 return m_xpr;
358 }
359
361 inline Index innerStride() const
362 {
363 return internal::traits<BlockType>::HasSameStorageOrderAsXprType
364 ? m_xpr.innerStride()
365 : m_xpr.outerStride();
366 }
367
369 inline Index outerStride() const
370 {
371 return m_outerStride;
372 }
373
374 #ifndef __SUNPRO_CC
375 // FIXME sunstudio is not friendly with the above friend...
376 // META-FIXME there is no 'friend' keyword around here. Is this obsolete?
377 protected:
378 #endif
379
380 #ifndef EIGEN_PARSED_BY_DOXYGEN
382 inline BlockImpl_dense(XprType& xpr, const Scalar* data, Index blockRows, Index blockCols)
383 : Base(data, blockRows, blockCols), m_xpr(xpr)
384 {
385 init();
386 }
387 #endif
388
389 protected:
390 void init()
391 {
392 m_outerStride = internal::traits<BlockType>::HasSameStorageOrderAsXprType
393 ? m_xpr.outerStride()
394 : m_xpr.innerStride();
395 }
396
397 typename XprType::Nested m_xpr;
398 Index m_outerStride;
399};
400
401} // end namespace internal
402
403} // end namespace Eigen
404
405#endif // EIGEN_BLOCK_H
Block(XprType &xpr, Index a_startRow, Index a_startCol)
Definition Block.h:123
Block(XprType &xpr, Index i)
Definition Block.h:114
Block(XprType &xpr, Index a_startRow, Index a_startCol, Index blockRows, Index blockCols)
Definition Block.h:133
const unsigned int DirectAccessBit
Definition Constants.h:142
const unsigned int LvalueBit
Definition Constants.h:131
const unsigned int RowMajorBit
Definition Constants.h:53
const unsigned int AlignedBit
Definition Constants.h:147
const unsigned int PacketAccessBit
Definition Constants.h:81
const unsigned int LinearAccessBit
Definition Constants.h:117
Definition LDLT.h:18