1 | // This file is part of Eigen, a lightweight C++ template library |
---|
2 | // for linear algebra. |
---|
3 | // |
---|
4 | // Copyright (C) 2008-2010 Gael Guennebaud <gael.guennebaud@inria.fr> |
---|
5 | // Copyright (C) 2006-2008 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_VECTORBLOCK_H |
---|
12 | #define EIGEN_VECTORBLOCK_H |
---|
13 | |
---|
14 | namespace Eigen { |
---|
15 | |
---|
16 | /** \class VectorBlock |
---|
17 | * \ingroup Core_Module |
---|
18 | * |
---|
19 | * \brief Expression of a fixed-size or dynamic-size sub-vector |
---|
20 | * |
---|
21 | * \param VectorType the type of the object in which we are taking a sub-vector |
---|
22 | * \param Size size of the sub-vector we are taking at compile time (optional) |
---|
23 | * |
---|
24 | * This class represents an expression of either a fixed-size or dynamic-size sub-vector. |
---|
25 | * It is the return type of DenseBase::segment(Index,Index) and DenseBase::segment<int>(Index) and |
---|
26 | * most of the time this is the only way it is used. |
---|
27 | * |
---|
28 | * However, if you want to directly maniputate sub-vector expressions, |
---|
29 | * for instance if you want to write a function returning such an expression, you |
---|
30 | * will need to use this class. |
---|
31 | * |
---|
32 | * Here is an example illustrating the dynamic case: |
---|
33 | * \include class_VectorBlock.cpp |
---|
34 | * Output: \verbinclude class_VectorBlock.out |
---|
35 | * |
---|
36 | * \note Even though this expression has dynamic size, in the case where \a VectorType |
---|
37 | * has fixed size, this expression inherits a fixed maximal size which means that evaluating |
---|
38 | * it does not cause a dynamic memory allocation. |
---|
39 | * |
---|
40 | * Here is an example illustrating the fixed-size case: |
---|
41 | * \include class_FixedVectorBlock.cpp |
---|
42 | * Output: \verbinclude class_FixedVectorBlock.out |
---|
43 | * |
---|
44 | * \sa class Block, DenseBase::segment(Index,Index,Index,Index), DenseBase::segment(Index,Index) |
---|
45 | */ |
---|
46 | |
---|
47 | namespace internal { |
---|
48 | template<typename VectorType, int Size> |
---|
49 | struct traits<VectorBlock<VectorType, Size> > |
---|
50 | : public traits<Block<VectorType, |
---|
51 | traits<VectorType>::Flags & RowMajorBit ? 1 : Size, |
---|
52 | traits<VectorType>::Flags & RowMajorBit ? Size : 1> > |
---|
53 | { |
---|
54 | }; |
---|
55 | } |
---|
56 | |
---|
57 | template<typename VectorType, int Size> class VectorBlock |
---|
58 | : public Block<VectorType, |
---|
59 | internal::traits<VectorType>::Flags & RowMajorBit ? 1 : Size, |
---|
60 | internal::traits<VectorType>::Flags & RowMajorBit ? Size : 1> |
---|
61 | { |
---|
62 | typedef Block<VectorType, |
---|
63 | internal::traits<VectorType>::Flags & RowMajorBit ? 1 : Size, |
---|
64 | internal::traits<VectorType>::Flags & RowMajorBit ? Size : 1> Base; |
---|
65 | enum { |
---|
66 | IsColVector = !(internal::traits<VectorType>::Flags & RowMajorBit) |
---|
67 | }; |
---|
68 | public: |
---|
69 | EIGEN_DENSE_PUBLIC_INTERFACE(VectorBlock) |
---|
70 | |
---|
71 | using Base::operator=; |
---|
72 | |
---|
73 | /** Dynamic-size constructor |
---|
74 | */ |
---|
75 | inline VectorBlock(VectorType& vector, Index start, Index size) |
---|
76 | : Base(vector, |
---|
77 | IsColVector ? start : 0, IsColVector ? 0 : start, |
---|
78 | IsColVector ? size : 1, IsColVector ? 1 : size) |
---|
79 | { |
---|
80 | EIGEN_STATIC_ASSERT_VECTOR_ONLY(VectorBlock); |
---|
81 | } |
---|
82 | |
---|
83 | /** Fixed-size constructor |
---|
84 | */ |
---|
85 | inline VectorBlock(VectorType& vector, Index start) |
---|
86 | : Base(vector, IsColVector ? start : 0, IsColVector ? 0 : start) |
---|
87 | { |
---|
88 | EIGEN_STATIC_ASSERT_VECTOR_ONLY(VectorBlock); |
---|
89 | } |
---|
90 | }; |
---|
91 | |
---|
92 | |
---|
93 | /** \returns a dynamic-size expression of a segment (i.e. a vector block) in *this. |
---|
94 | * |
---|
95 | * \only_for_vectors |
---|
96 | * |
---|
97 | * \param start the first coefficient in the segment |
---|
98 | * \param size the number of coefficients in the segment |
---|
99 | * |
---|
100 | * Example: \include MatrixBase_segment_int_int.cpp |
---|
101 | * Output: \verbinclude MatrixBase_segment_int_int.out |
---|
102 | * |
---|
103 | * \note Even though the returned expression has dynamic size, in the case |
---|
104 | * when it is applied to a fixed-size vector, it inherits a fixed maximal size, |
---|
105 | * which means that evaluating it does not cause a dynamic memory allocation. |
---|
106 | * |
---|
107 | * \sa class Block, segment(Index) |
---|
108 | */ |
---|
109 | template<typename Derived> |
---|
110 | inline typename DenseBase<Derived>::SegmentReturnType |
---|
111 | DenseBase<Derived>::segment(Index start, Index size) |
---|
112 | { |
---|
113 | EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived) |
---|
114 | return SegmentReturnType(derived(), start, size); |
---|
115 | } |
---|
116 | |
---|
117 | /** This is the const version of segment(Index,Index).*/ |
---|
118 | template<typename Derived> |
---|
119 | inline typename DenseBase<Derived>::ConstSegmentReturnType |
---|
120 | DenseBase<Derived>::segment(Index start, Index size) const |
---|
121 | { |
---|
122 | EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived) |
---|
123 | return ConstSegmentReturnType(derived(), start, size); |
---|
124 | } |
---|
125 | |
---|
126 | /** \returns a dynamic-size expression of the first coefficients of *this. |
---|
127 | * |
---|
128 | * \only_for_vectors |
---|
129 | * |
---|
130 | * \param size the number of coefficients in the block |
---|
131 | * |
---|
132 | * Example: \include MatrixBase_start_int.cpp |
---|
133 | * Output: \verbinclude MatrixBase_start_int.out |
---|
134 | * |
---|
135 | * \note Even though the returned expression has dynamic size, in the case |
---|
136 | * when it is applied to a fixed-size vector, it inherits a fixed maximal size, |
---|
137 | * which means that evaluating it does not cause a dynamic memory allocation. |
---|
138 | * |
---|
139 | * \sa class Block, block(Index,Index) |
---|
140 | */ |
---|
141 | template<typename Derived> |
---|
142 | inline typename DenseBase<Derived>::SegmentReturnType |
---|
143 | DenseBase<Derived>::head(Index size) |
---|
144 | { |
---|
145 | EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived) |
---|
146 | return SegmentReturnType(derived(), 0, size); |
---|
147 | } |
---|
148 | |
---|
149 | /** This is the const version of head(Index).*/ |
---|
150 | template<typename Derived> |
---|
151 | inline typename DenseBase<Derived>::ConstSegmentReturnType |
---|
152 | DenseBase<Derived>::head(Index size) const |
---|
153 | { |
---|
154 | EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived) |
---|
155 | return ConstSegmentReturnType(derived(), 0, size); |
---|
156 | } |
---|
157 | |
---|
158 | /** \returns a dynamic-size expression of the last coefficients of *this. |
---|
159 | * |
---|
160 | * \only_for_vectors |
---|
161 | * |
---|
162 | * \param size the number of coefficients in the block |
---|
163 | * |
---|
164 | * Example: \include MatrixBase_end_int.cpp |
---|
165 | * Output: \verbinclude MatrixBase_end_int.out |
---|
166 | * |
---|
167 | * \note Even though the returned expression has dynamic size, in the case |
---|
168 | * when it is applied to a fixed-size vector, it inherits a fixed maximal size, |
---|
169 | * which means that evaluating it does not cause a dynamic memory allocation. |
---|
170 | * |
---|
171 | * \sa class Block, block(Index,Index) |
---|
172 | */ |
---|
173 | template<typename Derived> |
---|
174 | inline typename DenseBase<Derived>::SegmentReturnType |
---|
175 | DenseBase<Derived>::tail(Index size) |
---|
176 | { |
---|
177 | EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived) |
---|
178 | return SegmentReturnType(derived(), this->size() - size, size); |
---|
179 | } |
---|
180 | |
---|
181 | /** This is the const version of tail(Index).*/ |
---|
182 | template<typename Derived> |
---|
183 | inline typename DenseBase<Derived>::ConstSegmentReturnType |
---|
184 | DenseBase<Derived>::tail(Index size) const |
---|
185 | { |
---|
186 | EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived) |
---|
187 | return ConstSegmentReturnType(derived(), this->size() - size, size); |
---|
188 | } |
---|
189 | |
---|
190 | /** \returns a fixed-size expression of a segment (i.e. a vector block) in \c *this |
---|
191 | * |
---|
192 | * \only_for_vectors |
---|
193 | * |
---|
194 | * The template parameter \a Size is the number of coefficients in the block |
---|
195 | * |
---|
196 | * \param start the index of the first element of the sub-vector |
---|
197 | * |
---|
198 | * Example: \include MatrixBase_template_int_segment.cpp |
---|
199 | * Output: \verbinclude MatrixBase_template_int_segment.out |
---|
200 | * |
---|
201 | * \sa class Block |
---|
202 | */ |
---|
203 | template<typename Derived> |
---|
204 | template<int Size> |
---|
205 | inline typename DenseBase<Derived>::template FixedSegmentReturnType<Size>::Type |
---|
206 | DenseBase<Derived>::segment(Index start) |
---|
207 | { |
---|
208 | EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived) |
---|
209 | return typename FixedSegmentReturnType<Size>::Type(derived(), start); |
---|
210 | } |
---|
211 | |
---|
212 | /** This is the const version of segment<int>(Index).*/ |
---|
213 | template<typename Derived> |
---|
214 | template<int Size> |
---|
215 | inline typename DenseBase<Derived>::template ConstFixedSegmentReturnType<Size>::Type |
---|
216 | DenseBase<Derived>::segment(Index start) const |
---|
217 | { |
---|
218 | EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived) |
---|
219 | return typename ConstFixedSegmentReturnType<Size>::Type(derived(), start); |
---|
220 | } |
---|
221 | |
---|
222 | /** \returns a fixed-size expression of the first coefficients of *this. |
---|
223 | * |
---|
224 | * \only_for_vectors |
---|
225 | * |
---|
226 | * The template parameter \a Size is the number of coefficients in the block |
---|
227 | * |
---|
228 | * Example: \include MatrixBase_template_int_start.cpp |
---|
229 | * Output: \verbinclude MatrixBase_template_int_start.out |
---|
230 | * |
---|
231 | * \sa class Block |
---|
232 | */ |
---|
233 | template<typename Derived> |
---|
234 | template<int Size> |
---|
235 | inline typename DenseBase<Derived>::template FixedSegmentReturnType<Size>::Type |
---|
236 | DenseBase<Derived>::head() |
---|
237 | { |
---|
238 | EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived) |
---|
239 | return typename FixedSegmentReturnType<Size>::Type(derived(), 0); |
---|
240 | } |
---|
241 | |
---|
242 | /** This is the const version of head<int>().*/ |
---|
243 | template<typename Derived> |
---|
244 | template<int Size> |
---|
245 | inline typename DenseBase<Derived>::template ConstFixedSegmentReturnType<Size>::Type |
---|
246 | DenseBase<Derived>::head() const |
---|
247 | { |
---|
248 | EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived) |
---|
249 | return typename ConstFixedSegmentReturnType<Size>::Type(derived(), 0); |
---|
250 | } |
---|
251 | |
---|
252 | /** \returns a fixed-size expression of the last coefficients of *this. |
---|
253 | * |
---|
254 | * \only_for_vectors |
---|
255 | * |
---|
256 | * The template parameter \a Size is the number of coefficients in the block |
---|
257 | * |
---|
258 | * Example: \include MatrixBase_template_int_end.cpp |
---|
259 | * Output: \verbinclude MatrixBase_template_int_end.out |
---|
260 | * |
---|
261 | * \sa class Block |
---|
262 | */ |
---|
263 | template<typename Derived> |
---|
264 | template<int Size> |
---|
265 | inline typename DenseBase<Derived>::template FixedSegmentReturnType<Size>::Type |
---|
266 | DenseBase<Derived>::tail() |
---|
267 | { |
---|
268 | EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived) |
---|
269 | return typename FixedSegmentReturnType<Size>::Type(derived(), size() - Size); |
---|
270 | } |
---|
271 | |
---|
272 | /** This is the const version of tail<int>.*/ |
---|
273 | template<typename Derived> |
---|
274 | template<int Size> |
---|
275 | inline typename DenseBase<Derived>::template ConstFixedSegmentReturnType<Size>::Type |
---|
276 | DenseBase<Derived>::tail() const |
---|
277 | { |
---|
278 | EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived) |
---|
279 | return typename ConstFixedSegmentReturnType<Size>::Type(derived(), size() - Size); |
---|
280 | } |
---|
281 | |
---|
282 | } // end namespace Eigen |
---|
283 | |
---|
284 | #endif // EIGEN_VECTORBLOCK_H |
---|