Go to the documentation of this file.00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026 #ifndef EIGEN_VECTORBLOCK_H
00027 #define EIGEN_VECTORBLOCK_H
00028
00060 namespace internal {
00061 template<typename VectorType, int Size>
00062 struct traits<VectorBlock<VectorType, Size> >
00063 : public traits<Block<VectorType,
00064 traits<VectorType>::Flags & RowMajorBit ? 1 : Size,
00065 traits<VectorType>::Flags & RowMajorBit ? Size : 1> >
00066 {
00067 };
00068 }
00069
00070 template<typename VectorType, int Size> class VectorBlock
00071 : public Block<VectorType,
00072 internal::traits<VectorType>::Flags & RowMajorBit ? 1 : Size,
00073 internal::traits<VectorType>::Flags & RowMajorBit ? Size : 1>
00074 {
00075 typedef Block<VectorType,
00076 internal::traits<VectorType>::Flags & RowMajorBit ? 1 : Size,
00077 internal::traits<VectorType>::Flags & RowMajorBit ? Size : 1> Base;
00078 enum {
00079 IsColVector = !(internal::traits<VectorType>::Flags & RowMajorBit)
00080 };
00081 public:
00082 EIGEN_DENSE_PUBLIC_INTERFACE(VectorBlock)
00083
00084 using Base::operator=;
00085
00088 inline VectorBlock(VectorType& vector, Index start, Index size)
00089 : Base(vector,
00090 IsColVector ? start : 0, IsColVector ? 0 : start,
00091 IsColVector ? size : 1, IsColVector ? 1 : size)
00092 {
00093 EIGEN_STATIC_ASSERT_VECTOR_ONLY(VectorBlock);
00094 }
00095
00098 inline VectorBlock(VectorType& vector, Index start)
00099 : Base(vector, IsColVector ? start : 0, IsColVector ? 0 : start)
00100 {
00101 EIGEN_STATIC_ASSERT_VECTOR_ONLY(VectorBlock);
00102 }
00103 };
00104
00105
00122 template<typename Derived>
00123 inline typename DenseBase<Derived>::SegmentReturnType
00124 DenseBase<Derived>::segment(Index start, Index size)
00125 {
00126 EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived)
00127 return SegmentReturnType(derived(), start, size);
00128 }
00129
00131 template<typename Derived>
00132 inline typename DenseBase<Derived>::ConstSegmentReturnType
00133 DenseBase<Derived>::segment(Index start, Index size) const
00134 {
00135 EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived)
00136 return ConstSegmentReturnType(derived(), start, size);
00137 }
00138
00154 template<typename Derived>
00155 inline typename DenseBase<Derived>::SegmentReturnType
00156 DenseBase<Derived>::head(Index size)
00157 {
00158 EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived)
00159 return SegmentReturnType(derived(), 0, size);
00160 }
00161
00163 template<typename Derived>
00164 inline typename DenseBase<Derived>::ConstSegmentReturnType
00165 DenseBase<Derived>::head(Index size) const
00166 {
00167 EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived)
00168 return ConstSegmentReturnType(derived(), 0, size);
00169 }
00170
00186 template<typename Derived>
00187 inline typename DenseBase<Derived>::SegmentReturnType
00188 DenseBase<Derived>::tail(Index size)
00189 {
00190 EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived)
00191 return SegmentReturnType(derived(), this->size() - size, size);
00192 }
00193
00195 template<typename Derived>
00196 inline typename DenseBase<Derived>::ConstSegmentReturnType
00197 DenseBase<Derived>::tail(Index size) const
00198 {
00199 EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived)
00200 return ConstSegmentReturnType(derived(), this->size() - size, size);
00201 }
00202
00216 template<typename Derived>
00217 template<int Size>
00218 inline typename DenseBase<Derived>::template FixedSegmentReturnType<Size>::Type
00219 DenseBase<Derived>::segment(Index start)
00220 {
00221 EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived)
00222 return typename FixedSegmentReturnType<Size>::Type(derived(), start);
00223 }
00224
00226 template<typename Derived>
00227 template<int Size>
00228 inline typename DenseBase<Derived>::template ConstFixedSegmentReturnType<Size>::Type
00229 DenseBase<Derived>::segment(Index start) const
00230 {
00231 EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived)
00232 return typename ConstFixedSegmentReturnType<Size>::Type(derived(), start);
00233 }
00234
00246 template<typename Derived>
00247 template<int Size>
00248 inline typename DenseBase<Derived>::template FixedSegmentReturnType<Size>::Type
00249 DenseBase<Derived>::head()
00250 {
00251 EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived)
00252 return typename FixedSegmentReturnType<Size>::Type(derived(), 0);
00253 }
00254
00256 template<typename Derived>
00257 template<int Size>
00258 inline typename DenseBase<Derived>::template ConstFixedSegmentReturnType<Size>::Type
00259 DenseBase<Derived>::head() const
00260 {
00261 EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived)
00262 return typename ConstFixedSegmentReturnType<Size>::Type(derived(), 0);
00263 }
00264
00276 template<typename Derived>
00277 template<int Size>
00278 inline typename DenseBase<Derived>::template FixedSegmentReturnType<Size>::Type
00279 DenseBase<Derived>::tail()
00280 {
00281 EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived)
00282 return typename FixedSegmentReturnType<Size>::Type(derived(), size() - Size);
00283 }
00284
00286 template<typename Derived>
00287 template<int Size>
00288 inline typename DenseBase<Derived>::template ConstFixedSegmentReturnType<Size>::Type
00289 DenseBase<Derived>::tail() const
00290 {
00291 EIGEN_STATIC_ASSERT_VECTOR_ONLY(Derived)
00292 return typename ConstFixedSegmentReturnType<Size>::Type(derived(), size() - Size);
00293 }
00294
00295
00296 #endif // EIGEN_VECTORBLOCK_H