Home | History | Annotate | Download | only in Core
      1 // This file is part of Eigen, a lightweight C++ template library
      2 // for linear algebra.
      3 //
      4 // Copyright (C) 2008-2009 Gael Guennebaud <gael.guennebaud (at) inria.fr>
      5 // Copyright (C) 2006-2008 Benoit Jacob <jacob.benoit.1 (at) 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_DENSESTORAGEBASE_H
     12 #define EIGEN_DENSESTORAGEBASE_H
     13 
     14 #if defined(EIGEN_INITIALIZE_MATRICES_BY_ZERO)
     15 # define EIGEN_INITIALIZE_COEFFS
     16 # define EIGEN_INITIALIZE_COEFFS_IF_THAT_OPTION_IS_ENABLED for(int i=0;i<base().size();++i) coeffRef(i)=Scalar(0);
     17 #elif defined(EIGEN_INITIALIZE_MATRICES_BY_NAN)
     18 # define EIGEN_INITIALIZE_COEFFS
     19 # define EIGEN_INITIALIZE_COEFFS_IF_THAT_OPTION_IS_ENABLED for(int i=0;i<base().size();++i) coeffRef(i)=std::numeric_limits<Scalar>::quiet_NaN();
     20 #else
     21 # undef EIGEN_INITIALIZE_COEFFS
     22 # define EIGEN_INITIALIZE_COEFFS_IF_THAT_OPTION_IS_ENABLED
     23 #endif
     24 
     25 namespace Eigen {
     26 
     27 namespace internal {
     28 
     29 template<int MaxSizeAtCompileTime> struct check_rows_cols_for_overflow {
     30   template<typename Index>
     31   static EIGEN_ALWAYS_INLINE void run(Index, Index)
     32   {
     33   }
     34 };
     35 
     36 template<> struct check_rows_cols_for_overflow<Dynamic> {
     37   template<typename Index>
     38   static EIGEN_ALWAYS_INLINE void run(Index rows, Index cols)
     39   {
     40     // http://hg.mozilla.org/mozilla-central/file/6c8a909977d3/xpcom/ds/CheckedInt.h#l242
     41     // we assume Index is signed
     42     Index max_index = (size_t(1) << (8 * sizeof(Index) - 1)) - 1; // assume Index is signed
     43     bool error = (rows == 0 || cols == 0) ? false
     44                : (rows > max_index / cols);
     45     if (error)
     46       throw_std_bad_alloc();
     47   }
     48 };
     49 
     50 template <typename Derived,
     51           typename OtherDerived = Derived,
     52           bool IsVector = bool(Derived::IsVectorAtCompileTime) && bool(OtherDerived::IsVectorAtCompileTime)>
     53 struct conservative_resize_like_impl;
     54 
     55 template<typename MatrixTypeA, typename MatrixTypeB, bool SwapPointers> struct matrix_swap_impl;
     56 
     57 } // end namespace internal
     58 
     59 /** \class PlainObjectBase
     60   * \brief %Dense storage base class for matrices and arrays.
     61   *
     62   * This class can be extended with the help of the plugin mechanism described on the page
     63   * \ref TopicCustomizingEigen by defining the preprocessor symbol \c EIGEN_PLAINOBJECTBASE_PLUGIN.
     64   *
     65   * \sa \ref TopicClassHierarchy
     66   */
     67 #ifdef EIGEN_PARSED_BY_DOXYGEN
     68 namespace internal {
     69 
     70 // this is a warkaround to doxygen not being able to understand the inheritence logic
     71 // when it is hidden by the dense_xpr_base helper struct.
     72 template<typename Derived> struct dense_xpr_base_dispatcher_for_doxygen;// : public MatrixBase<Derived> {};
     73 /** This class is just a workaround for Doxygen and it does not not actually exist. */
     74 template<typename _Scalar, int _Rows, int _Cols, int _Options, int _MaxRows, int _MaxCols>
     75 struct dense_xpr_base_dispatcher_for_doxygen<Matrix<_Scalar, _Rows, _Cols, _Options, _MaxRows, _MaxCols> >
     76     : public MatrixBase<Matrix<_Scalar, _Rows, _Cols, _Options, _MaxRows, _MaxCols> > {};
     77 /** This class is just a workaround for Doxygen and it does not not actually exist. */
     78 template<typename _Scalar, int _Rows, int _Cols, int _Options, int _MaxRows, int _MaxCols>
     79 struct dense_xpr_base_dispatcher_for_doxygen<Array<_Scalar, _Rows, _Cols, _Options, _MaxRows, _MaxCols> >
     80     : public ArrayBase<Array<_Scalar, _Rows, _Cols, _Options, _MaxRows, _MaxCols> > {};
     81 
     82 } // namespace internal
     83 
     84 template<typename Derived>
     85 class PlainObjectBase : public internal::dense_xpr_base_dispatcher_for_doxygen<Derived>
     86 #else
     87 template<typename Derived>
     88 class PlainObjectBase : public internal::dense_xpr_base<Derived>::type
     89 #endif
     90 {
     91   public:
     92     enum { Options = internal::traits<Derived>::Options };
     93     typedef typename internal::dense_xpr_base<Derived>::type Base;
     94 
     95     typedef typename internal::traits<Derived>::StorageKind StorageKind;
     96     typedef typename internal::traits<Derived>::Index Index;
     97     typedef typename internal::traits<Derived>::Scalar Scalar;
     98     typedef typename internal::packet_traits<Scalar>::type PacketScalar;
     99     typedef typename NumTraits<Scalar>::Real RealScalar;
    100     typedef Derived DenseType;
    101 
    102     using Base::RowsAtCompileTime;
    103     using Base::ColsAtCompileTime;
    104     using Base::SizeAtCompileTime;
    105     using Base::MaxRowsAtCompileTime;
    106     using Base::MaxColsAtCompileTime;
    107     using Base::MaxSizeAtCompileTime;
    108     using Base::IsVectorAtCompileTime;
    109     using Base::Flags;
    110 
    111     template<typename PlainObjectType, int MapOptions, typename StrideType> friend class Eigen::Map;
    112     friend  class Eigen::Map<Derived, Unaligned>;
    113     typedef Eigen::Map<Derived, Unaligned>  MapType;
    114     friend  class Eigen::Map<const Derived, Unaligned>;
    115     typedef const Eigen::Map<const Derived, Unaligned> ConstMapType;
    116     friend  class Eigen::Map<Derived, Aligned>;
    117     typedef Eigen::Map<Derived, Aligned> AlignedMapType;
    118     friend  class Eigen::Map<const Derived, Aligned>;
    119     typedef const Eigen::Map<const Derived, Aligned> ConstAlignedMapType;
    120     template<typename StrideType> struct StridedMapType { typedef Eigen::Map<Derived, Unaligned, StrideType> type; };
    121     template<typename StrideType> struct StridedConstMapType { typedef Eigen::Map<const Derived, Unaligned, StrideType> type; };
    122     template<typename StrideType> struct StridedAlignedMapType { typedef Eigen::Map<Derived, Aligned, StrideType> type; };
    123     template<typename StrideType> struct StridedConstAlignedMapType { typedef Eigen::Map<const Derived, Aligned, StrideType> type; };
    124 
    125   protected:
    126     DenseStorage<Scalar, Base::MaxSizeAtCompileTime, Base::RowsAtCompileTime, Base::ColsAtCompileTime, Options> m_storage;
    127 
    128   public:
    129     enum { NeedsToAlign = SizeAtCompileTime != Dynamic && (internal::traits<Derived>::Flags & AlignedBit) != 0 };
    130     EIGEN_MAKE_ALIGNED_OPERATOR_NEW_IF(NeedsToAlign)
    131 
    132     Base& base() { return *static_cast<Base*>(this); }
    133     const Base& base() const { return *static_cast<const Base*>(this); }
    134 
    135     EIGEN_STRONG_INLINE Index rows() const { return m_storage.rows(); }
    136     EIGEN_STRONG_INLINE Index cols() const { return m_storage.cols(); }
    137 
    138     EIGEN_STRONG_INLINE const Scalar& coeff(Index rowId, Index colId) const
    139     {
    140       if(Flags & RowMajorBit)
    141         return m_storage.data()[colId + rowId * m_storage.cols()];
    142       else // column-major
    143         return m_storage.data()[rowId + colId * m_storage.rows()];
    144     }
    145 
    146     EIGEN_STRONG_INLINE const Scalar& coeff(Index index) const
    147     {
    148       return m_storage.data()[index];
    149     }
    150 
    151     EIGEN_STRONG_INLINE Scalar& coeffRef(Index rowId, Index colId)
    152     {
    153       if(Flags & RowMajorBit)
    154         return m_storage.data()[colId + rowId * m_storage.cols()];
    155       else // column-major
    156         return m_storage.data()[rowId + colId * m_storage.rows()];
    157     }
    158 
    159     EIGEN_STRONG_INLINE Scalar& coeffRef(Index index)
    160     {
    161       return m_storage.data()[index];
    162     }
    163 
    164     EIGEN_STRONG_INLINE const Scalar& coeffRef(Index rowId, Index colId) const
    165     {
    166       if(Flags & RowMajorBit)
    167         return m_storage.data()[colId + rowId * m_storage.cols()];
    168       else // column-major
    169         return m_storage.data()[rowId + colId * m_storage.rows()];
    170     }
    171 
    172     EIGEN_STRONG_INLINE const Scalar& coeffRef(Index index) const
    173     {
    174       return m_storage.data()[index];
    175     }
    176 
    177     /** \internal */
    178     template<int LoadMode>
    179     EIGEN_STRONG_INLINE PacketScalar packet(Index rowId, Index colId) const
    180     {
    181       return internal::ploadt<PacketScalar, LoadMode>
    182                (m_storage.data() + (Flags & RowMajorBit
    183                                    ? colId + rowId * m_storage.cols()
    184                                    : rowId + colId * m_storage.rows()));
    185     }
    186 
    187     /** \internal */
    188     template<int LoadMode>
    189     EIGEN_STRONG_INLINE PacketScalar packet(Index index) const
    190     {
    191       return internal::ploadt<PacketScalar, LoadMode>(m_storage.data() + index);
    192     }
    193 
    194     /** \internal */
    195     template<int StoreMode>
    196     EIGEN_STRONG_INLINE void writePacket(Index rowId, Index colId, const PacketScalar& val)
    197     {
    198       internal::pstoret<Scalar, PacketScalar, StoreMode>
    199               (m_storage.data() + (Flags & RowMajorBit
    200                                    ? colId + rowId * m_storage.cols()
    201                                    : rowId + colId * m_storage.rows()), val);
    202     }
    203 
    204     /** \internal */
    205     template<int StoreMode>
    206     EIGEN_STRONG_INLINE void writePacket(Index index, const PacketScalar& val)
    207     {
    208       internal::pstoret<Scalar, PacketScalar, StoreMode>(m_storage.data() + index, val);
    209     }
    210 
    211     /** \returns a const pointer to the data array of this matrix */
    212     EIGEN_STRONG_INLINE const Scalar *data() const
    213     { return m_storage.data(); }
    214 
    215     /** \returns a pointer to the data array of this matrix */
    216     EIGEN_STRONG_INLINE Scalar *data()
    217     { return m_storage.data(); }
    218 
    219     /** Resizes \c *this to a \a rows x \a cols matrix.
    220       *
    221       * This method is intended for dynamic-size matrices, although it is legal to call it on any
    222       * matrix as long as fixed dimensions are left unchanged. If you only want to change the number
    223       * of rows and/or of columns, you can use resize(NoChange_t, Index), resize(Index, NoChange_t).
    224       *
    225       * If the current number of coefficients of \c *this exactly matches the
    226       * product \a rows * \a cols, then no memory allocation is performed and
    227       * the current values are left unchanged. In all other cases, including
    228       * shrinking, the data is reallocated and all previous values are lost.
    229       *
    230       * Example: \include Matrix_resize_int_int.cpp
    231       * Output: \verbinclude Matrix_resize_int_int.out
    232       *
    233       * \sa resize(Index) for vectors, resize(NoChange_t, Index), resize(Index, NoChange_t)
    234       */
    235     EIGEN_STRONG_INLINE void resize(Index nbRows, Index nbCols)
    236     {
    237       eigen_assert(   EIGEN_IMPLIES(RowsAtCompileTime!=Dynamic,nbRows==RowsAtCompileTime)
    238                    && EIGEN_IMPLIES(ColsAtCompileTime!=Dynamic,nbCols==ColsAtCompileTime)
    239                    && EIGEN_IMPLIES(RowsAtCompileTime==Dynamic && MaxRowsAtCompileTime!=Dynamic,nbRows<=MaxRowsAtCompileTime)
    240                    && EIGEN_IMPLIES(ColsAtCompileTime==Dynamic && MaxColsAtCompileTime!=Dynamic,nbCols<=MaxColsAtCompileTime)
    241                    && nbRows>=0 && nbCols>=0 && "Invalid sizes when resizing a matrix or array.");
    242       internal::check_rows_cols_for_overflow<MaxSizeAtCompileTime>::run(nbRows, nbCols);
    243       #ifdef EIGEN_INITIALIZE_COEFFS
    244         Index size = nbRows*nbCols;
    245         bool size_changed = size != this->size();
    246         m_storage.resize(size, nbRows, nbCols);
    247         if(size_changed) EIGEN_INITIALIZE_COEFFS_IF_THAT_OPTION_IS_ENABLED
    248       #else
    249         internal::check_rows_cols_for_overflow<MaxSizeAtCompileTime>::run(nbRows, nbCols);
    250         m_storage.resize(nbRows*nbCols, nbRows, nbCols);
    251       #endif
    252     }
    253 
    254     /** Resizes \c *this to a vector of length \a size
    255       *
    256       * \only_for_vectors. This method does not work for
    257       * partially dynamic matrices when the static dimension is anything other
    258       * than 1. For example it will not work with Matrix<double, 2, Dynamic>.
    259       *
    260       * Example: \include Matrix_resize_int.cpp
    261       * Output: \verbinclude Matrix_resize_int.out
    262       *
    263       * \sa resize(Index,Index), resize(NoChange_t, Index), resize(Index, NoChange_t)
    264       */
    265     inline void resize(Index size)
    266     {
    267       EIGEN_STATIC_ASSERT_VECTOR_ONLY(PlainObjectBase)
    268       eigen_assert(((SizeAtCompileTime == Dynamic && (MaxSizeAtCompileTime==Dynamic || size<=MaxSizeAtCompileTime)) || SizeAtCompileTime == size) && size>=0);
    269       #ifdef EIGEN_INITIALIZE_COEFFS
    270         bool size_changed = size != this->size();
    271       #endif
    272       if(RowsAtCompileTime == 1)
    273         m_storage.resize(size, 1, size);
    274       else
    275         m_storage.resize(size, size, 1);
    276       #ifdef EIGEN_INITIALIZE_COEFFS
    277         if(size_changed) EIGEN_INITIALIZE_COEFFS_IF_THAT_OPTION_IS_ENABLED
    278       #endif
    279     }
    280 
    281     /** Resizes the matrix, changing only the number of columns. For the parameter of type NoChange_t, just pass the special value \c NoChange
    282       * as in the example below.
    283       *
    284       * Example: \include Matrix_resize_NoChange_int.cpp
    285       * Output: \verbinclude Matrix_resize_NoChange_int.out
    286       *
    287       * \sa resize(Index,Index)
    288       */
    289     inline void resize(NoChange_t, Index nbCols)
    290     {
    291       resize(rows(), nbCols);
    292     }
    293 
    294     /** Resizes the matrix, changing only the number of rows. For the parameter of type NoChange_t, just pass the special value \c NoChange
    295       * as in the example below.
    296       *
    297       * Example: \include Matrix_resize_int_NoChange.cpp
    298       * Output: \verbinclude Matrix_resize_int_NoChange.out
    299       *
    300       * \sa resize(Index,Index)
    301       */
    302     inline void resize(Index nbRows, NoChange_t)
    303     {
    304       resize(nbRows, cols());
    305     }
    306 
    307     /** Resizes \c *this to have the same dimensions as \a other.
    308       * Takes care of doing all the checking that's needed.
    309       *
    310       * Note that copying a row-vector into a vector (and conversely) is allowed.
    311       * The resizing, if any, is then done in the appropriate way so that row-vectors
    312       * remain row-vectors and vectors remain vectors.
    313       */
    314     template<typename OtherDerived>
    315     EIGEN_STRONG_INLINE void resizeLike(const EigenBase<OtherDerived>& _other)
    316     {
    317       const OtherDerived& other = _other.derived();
    318       internal::check_rows_cols_for_overflow<MaxSizeAtCompileTime>::run(other.rows(), other.cols());
    319       const Index othersize = other.rows()*other.cols();
    320       if(RowsAtCompileTime == 1)
    321       {
    322         eigen_assert(other.rows() == 1 || other.cols() == 1);
    323         resize(1, othersize);
    324       }
    325       else if(ColsAtCompileTime == 1)
    326       {
    327         eigen_assert(other.rows() == 1 || other.cols() == 1);
    328         resize(othersize, 1);
    329       }
    330       else resize(other.rows(), other.cols());
    331     }
    332 
    333     /** Resizes the matrix to \a rows x \a cols while leaving old values untouched.
    334       *
    335       * The method is intended for matrices of dynamic size. If you only want to change the number
    336       * of rows and/or of columns, you can use conservativeResize(NoChange_t, Index) or
    337       * conservativeResize(Index, NoChange_t).
    338       *
    339       * Matrices are resized relative to the top-left element. In case values need to be
    340       * appended to the matrix they will be uninitialized.
    341       */
    342     EIGEN_STRONG_INLINE void conservativeResize(Index nbRows, Index nbCols)
    343     {
    344       internal::conservative_resize_like_impl<Derived>::run(*this, nbRows, nbCols);
    345     }
    346 
    347     /** Resizes the matrix to \a rows x \a cols while leaving old values untouched.
    348       *
    349       * As opposed to conservativeResize(Index rows, Index cols), this version leaves
    350       * the number of columns unchanged.
    351       *
    352       * In case the matrix is growing, new rows will be uninitialized.
    353       */
    354     EIGEN_STRONG_INLINE void conservativeResize(Index nbRows, NoChange_t)
    355     {
    356       // Note: see the comment in conservativeResize(Index,Index)
    357       conservativeResize(nbRows, cols());
    358     }
    359 
    360     /** Resizes the matrix to \a rows x \a cols while leaving old values untouched.
    361       *
    362       * As opposed to conservativeResize(Index rows, Index cols), this version leaves
    363       * the number of rows unchanged.
    364       *
    365       * In case the matrix is growing, new columns will be uninitialized.
    366       */
    367     EIGEN_STRONG_INLINE void conservativeResize(NoChange_t, Index nbCols)
    368     {
    369       // Note: see the comment in conservativeResize(Index,Index)
    370       conservativeResize(rows(), nbCols);
    371     }
    372 
    373     /** Resizes the vector to \a size while retaining old values.
    374       *
    375       * \only_for_vectors. This method does not work for
    376       * partially dynamic matrices when the static dimension is anything other
    377       * than 1. For example it will not work with Matrix<double, 2, Dynamic>.
    378       *
    379       * When values are appended, they will be uninitialized.
    380       */
    381     EIGEN_STRONG_INLINE void conservativeResize(Index size)
    382     {
    383       internal::conservative_resize_like_impl<Derived>::run(*this, size);
    384     }
    385 
    386     /** Resizes the matrix to \a rows x \a cols of \c other, while leaving old values untouched.
    387       *
    388       * The method is intended for matrices of dynamic size. If you only want to change the number
    389       * of rows and/or of columns, you can use conservativeResize(NoChange_t, Index) or
    390       * conservativeResize(Index, NoChange_t).
    391       *
    392       * Matrices are resized relative to the top-left element. In case values need to be
    393       * appended to the matrix they will copied from \c other.
    394       */
    395     template<typename OtherDerived>
    396     EIGEN_STRONG_INLINE void conservativeResizeLike(const DenseBase<OtherDerived>& other)
    397     {
    398       internal::conservative_resize_like_impl<Derived,OtherDerived>::run(*this, other);
    399     }
    400 
    401     /** This is a special case of the templated operator=. Its purpose is to
    402       * prevent a default operator= from hiding the templated operator=.
    403       */
    404     EIGEN_STRONG_INLINE Derived& operator=(const PlainObjectBase& other)
    405     {
    406       return _set(other);
    407     }
    408 
    409     /** \sa MatrixBase::lazyAssign() */
    410     template<typename OtherDerived>
    411     EIGEN_STRONG_INLINE Derived& lazyAssign(const DenseBase<OtherDerived>& other)
    412     {
    413       _resize_to_match(other);
    414       return Base::lazyAssign(other.derived());
    415     }
    416 
    417     template<typename OtherDerived>
    418     EIGEN_STRONG_INLINE Derived& operator=(const ReturnByValue<OtherDerived>& func)
    419     {
    420       resize(func.rows(), func.cols());
    421       return Base::operator=(func);
    422     }
    423 
    424     EIGEN_STRONG_INLINE PlainObjectBase() : m_storage()
    425     {
    426 //       _check_template_params();
    427 //       EIGEN_INITIALIZE_COEFFS_IF_THAT_OPTION_IS_ENABLED
    428     }
    429 
    430 #ifndef EIGEN_PARSED_BY_DOXYGEN
    431     // FIXME is it still needed ?
    432     /** \internal */
    433     PlainObjectBase(internal::constructor_without_unaligned_array_assert)
    434       : m_storage(internal::constructor_without_unaligned_array_assert())
    435     {
    436 //       _check_template_params(); EIGEN_INITIALIZE_COEFFS_IF_THAT_OPTION_IS_ENABLED
    437     }
    438 #endif
    439 
    440     EIGEN_STRONG_INLINE PlainObjectBase(Index a_size, Index nbRows, Index nbCols)
    441       : m_storage(a_size, nbRows, nbCols)
    442     {
    443 //       _check_template_params();
    444 //       EIGEN_INITIALIZE_COEFFS_IF_THAT_OPTION_IS_ENABLED
    445     }
    446 
    447     /** \copydoc MatrixBase::operator=(const EigenBase<OtherDerived>&)
    448       */
    449     template<typename OtherDerived>
    450     EIGEN_STRONG_INLINE Derived& operator=(const EigenBase<OtherDerived> &other)
    451     {
    452       _resize_to_match(other);
    453       Base::operator=(other.derived());
    454       return this->derived();
    455     }
    456 
    457     /** \sa MatrixBase::operator=(const EigenBase<OtherDerived>&) */
    458     template<typename OtherDerived>
    459     EIGEN_STRONG_INLINE PlainObjectBase(const EigenBase<OtherDerived> &other)
    460       : m_storage(other.derived().rows() * other.derived().cols(), other.derived().rows(), other.derived().cols())
    461     {
    462       _check_template_params();
    463       internal::check_rows_cols_for_overflow<MaxSizeAtCompileTime>::run(other.derived().rows(), other.derived().cols());
    464       Base::operator=(other.derived());
    465     }
    466 
    467     /** \name Map
    468       * These are convenience functions returning Map objects. The Map() static functions return unaligned Map objects,
    469       * while the AlignedMap() functions return aligned Map objects and thus should be called only with 16-byte-aligned
    470       * \a data pointers.
    471       *
    472       * \see class Map
    473       */
    474     //@{
    475     static inline ConstMapType Map(const Scalar* data)
    476     { return ConstMapType(data); }
    477     static inline MapType Map(Scalar* data)
    478     { return MapType(data); }
    479     static inline ConstMapType Map(const Scalar* data, Index size)
    480     { return ConstMapType(data, size); }
    481     static inline MapType Map(Scalar* data, Index size)
    482     { return MapType(data, size); }
    483     static inline ConstMapType Map(const Scalar* data, Index rows, Index cols)
    484     { return ConstMapType(data, rows, cols); }
    485     static inline MapType Map(Scalar* data, Index rows, Index cols)
    486     { return MapType(data, rows, cols); }
    487 
    488     static inline ConstAlignedMapType MapAligned(const Scalar* data)
    489     { return ConstAlignedMapType(data); }
    490     static inline AlignedMapType MapAligned(Scalar* data)
    491     { return AlignedMapType(data); }
    492     static inline ConstAlignedMapType MapAligned(const Scalar* data, Index size)
    493     { return ConstAlignedMapType(data, size); }
    494     static inline AlignedMapType MapAligned(Scalar* data, Index size)
    495     { return AlignedMapType(data, size); }
    496     static inline ConstAlignedMapType MapAligned(const Scalar* data, Index rows, Index cols)
    497     { return ConstAlignedMapType(data, rows, cols); }
    498     static inline AlignedMapType MapAligned(Scalar* data, Index rows, Index cols)
    499     { return AlignedMapType(data, rows, cols); }
    500 
    501     template<int Outer, int Inner>
    502     static inline typename StridedConstMapType<Stride<Outer, Inner> >::type Map(const Scalar* data, const Stride<Outer, Inner>& stride)
    503     { return typename StridedConstMapType<Stride<Outer, Inner> >::type(data, stride); }
    504     template<int Outer, int Inner>
    505     static inline typename StridedMapType<Stride<Outer, Inner> >::type Map(Scalar* data, const Stride<Outer, Inner>& stride)
    506     { return typename StridedMapType<Stride<Outer, Inner> >::type(data, stride); }
    507     template<int Outer, int Inner>
    508     static inline typename StridedConstMapType<Stride<Outer, Inner> >::type Map(const Scalar* data, Index size, const Stride<Outer, Inner>& stride)
    509     { return typename StridedConstMapType<Stride<Outer, Inner> >::type(data, size, stride); }
    510     template<int Outer, int Inner>
    511     static inline typename StridedMapType<Stride<Outer, Inner> >::type Map(Scalar* data, Index size, const Stride<Outer, Inner>& stride)
    512     { return typename StridedMapType<Stride<Outer, Inner> >::type(data, size, stride); }
    513     template<int Outer, int Inner>
    514     static inline typename StridedConstMapType<Stride<Outer, Inner> >::type Map(const Scalar* data, Index rows, Index cols, const Stride<Outer, Inner>& stride)
    515     { return typename StridedConstMapType<Stride<Outer, Inner> >::type(data, rows, cols, stride); }
    516     template<int Outer, int Inner>
    517     static inline typename StridedMapType<Stride<Outer, Inner> >::type Map(Scalar* data, Index rows, Index cols, const Stride<Outer, Inner>& stride)
    518     { return typename StridedMapType<Stride<Outer, Inner> >::type(data, rows, cols, stride); }
    519 
    520     template<int Outer, int Inner>
    521     static inline typename StridedConstAlignedMapType<Stride<Outer, Inner> >::type MapAligned(const Scalar* data, const Stride<Outer, Inner>& stride)
    522     { return typename StridedConstAlignedMapType<Stride<Outer, Inner> >::type(data, stride); }
    523     template<int Outer, int Inner>
    524     static inline typename StridedAlignedMapType<Stride<Outer, Inner> >::type MapAligned(Scalar* data, const Stride<Outer, Inner>& stride)
    525     { return typename StridedAlignedMapType<Stride<Outer, Inner> >::type(data, stride); }
    526     template<int Outer, int Inner>
    527     static inline typename StridedConstAlignedMapType<Stride<Outer, Inner> >::type MapAligned(const Scalar* data, Index size, const Stride<Outer, Inner>& stride)
    528     { return typename StridedConstAlignedMapType<Stride<Outer, Inner> >::type(data, size, stride); }
    529     template<int Outer, int Inner>
    530     static inline typename StridedAlignedMapType<Stride<Outer, Inner> >::type MapAligned(Scalar* data, Index size, const Stride<Outer, Inner>& stride)
    531     { return typename StridedAlignedMapType<Stride<Outer, Inner> >::type(data, size, stride); }
    532     template<int Outer, int Inner>
    533     static inline typename StridedConstAlignedMapType<Stride<Outer, Inner> >::type MapAligned(const Scalar* data, Index rows, Index cols, const Stride<Outer, Inner>& stride)
    534     { return typename StridedConstAlignedMapType<Stride<Outer, Inner> >::type(data, rows, cols, stride); }
    535     template<int Outer, int Inner>
    536     static inline typename StridedAlignedMapType<Stride<Outer, Inner> >::type MapAligned(Scalar* data, Index rows, Index cols, const Stride<Outer, Inner>& stride)
    537     { return typename StridedAlignedMapType<Stride<Outer, Inner> >::type(data, rows, cols, stride); }
    538     //@}
    539 
    540     using Base::setConstant;
    541     Derived& setConstant(Index size, const Scalar& value);
    542     Derived& setConstant(Index rows, Index cols, const Scalar& value);
    543 
    544     using Base::setZero;
    545     Derived& setZero(Index size);
    546     Derived& setZero(Index rows, Index cols);
    547 
    548     using Base::setOnes;
    549     Derived& setOnes(Index size);
    550     Derived& setOnes(Index rows, Index cols);
    551 
    552     using Base::setRandom;
    553     Derived& setRandom(Index size);
    554     Derived& setRandom(Index rows, Index cols);
    555 
    556     #ifdef EIGEN_PLAINOBJECTBASE_PLUGIN
    557     #include EIGEN_PLAINOBJECTBASE_PLUGIN
    558     #endif
    559 
    560   protected:
    561     /** \internal Resizes *this in preparation for assigning \a other to it.
    562       * Takes care of doing all the checking that's needed.
    563       *
    564       * Note that copying a row-vector into a vector (and conversely) is allowed.
    565       * The resizing, if any, is then done in the appropriate way so that row-vectors
    566       * remain row-vectors and vectors remain vectors.
    567       */
    568     template<typename OtherDerived>
    569     EIGEN_STRONG_INLINE void _resize_to_match(const EigenBase<OtherDerived>& other)
    570     {
    571       #ifdef EIGEN_NO_AUTOMATIC_RESIZING
    572       eigen_assert((this->size()==0 || (IsVectorAtCompileTime ? (this->size() == other.size())
    573                  : (rows() == other.rows() && cols() == other.cols())))
    574         && "Size mismatch. Automatic resizing is disabled because EIGEN_NO_AUTOMATIC_RESIZING is defined");
    575       EIGEN_ONLY_USED_FOR_DEBUG(other);
    576       if(this->size()==0)
    577         resizeLike(other);
    578       #else
    579       resizeLike(other);
    580       #endif
    581     }
    582 
    583     /**
    584       * \brief Copies the value of the expression \a other into \c *this with automatic resizing.
    585       *
    586       * *this might be resized to match the dimensions of \a other. If *this was a null matrix (not already initialized),
    587       * it will be initialized.
    588       *
    589       * Note that copying a row-vector into a vector (and conversely) is allowed.
    590       * The resizing, if any, is then done in the appropriate way so that row-vectors
    591       * remain row-vectors and vectors remain vectors.
    592       *
    593       * \sa operator=(const MatrixBase<OtherDerived>&), _set_noalias()
    594       *
    595       * \internal
    596       */
    597     template<typename OtherDerived>
    598     EIGEN_STRONG_INLINE Derived& _set(const DenseBase<OtherDerived>& other)
    599     {
    600       _set_selector(other.derived(), typename internal::conditional<static_cast<bool>(int(OtherDerived::Flags) & EvalBeforeAssigningBit), internal::true_type, internal::false_type>::type());
    601       return this->derived();
    602     }
    603 
    604     template<typename OtherDerived>
    605     EIGEN_STRONG_INLINE void _set_selector(const OtherDerived& other, const internal::true_type&) { _set_noalias(other.eval()); }
    606 
    607     template<typename OtherDerived>
    608     EIGEN_STRONG_INLINE void _set_selector(const OtherDerived& other, const internal::false_type&) { _set_noalias(other); }
    609 
    610     /** \internal Like _set() but additionally makes the assumption that no aliasing effect can happen (which
    611       * is the case when creating a new matrix) so one can enforce lazy evaluation.
    612       *
    613       * \sa operator=(const MatrixBase<OtherDerived>&), _set()
    614       */
    615     template<typename OtherDerived>
    616     EIGEN_STRONG_INLINE Derived& _set_noalias(const DenseBase<OtherDerived>& other)
    617     {
    618       // I don't think we need this resize call since the lazyAssign will anyways resize
    619       // and lazyAssign will be called by the assign selector.
    620       //_resize_to_match(other);
    621       // the 'false' below means to enforce lazy evaluation. We don't use lazyAssign() because
    622       // it wouldn't allow to copy a row-vector into a column-vector.
    623       return internal::assign_selector<Derived,OtherDerived,false>::run(this->derived(), other.derived());
    624     }
    625 
    626     template<typename T0, typename T1>
    627     EIGEN_STRONG_INLINE void _init2(Index nbRows, Index nbCols, typename internal::enable_if<Base::SizeAtCompileTime!=2,T0>::type* = 0)
    628     {
    629       EIGEN_STATIC_ASSERT(bool(NumTraits<T0>::IsInteger) &&
    630                           bool(NumTraits<T1>::IsInteger),
    631                           FLOATING_POINT_ARGUMENT_PASSED__INTEGER_WAS_EXPECTED)
    632       resize(nbRows,nbCols);
    633     }
    634     template<typename T0, typename T1>
    635     EIGEN_STRONG_INLINE void _init2(const Scalar& val0, const Scalar& val1, typename internal::enable_if<Base::SizeAtCompileTime==2,T0>::type* = 0)
    636     {
    637       EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(PlainObjectBase, 2)
    638       m_storage.data()[0] = val0;
    639       m_storage.data()[1] = val1;
    640     }
    641 
    642     template<typename MatrixTypeA, typename MatrixTypeB, bool SwapPointers>
    643     friend struct internal::matrix_swap_impl;
    644 
    645     /** \internal generic implementation of swap for dense storage since for dynamic-sized matrices of same type it is enough to swap the
    646       * data pointers.
    647       */
    648     template<typename OtherDerived>
    649     void _swap(DenseBase<OtherDerived> const & other)
    650     {
    651       enum { SwapPointers = internal::is_same<Derived, OtherDerived>::value && Base::SizeAtCompileTime==Dynamic };
    652       internal::matrix_swap_impl<Derived, OtherDerived, bool(SwapPointers)>::run(this->derived(), other.const_cast_derived());
    653     }
    654 
    655   public:
    656 #ifndef EIGEN_PARSED_BY_DOXYGEN
    657     static EIGEN_STRONG_INLINE void _check_template_params()
    658     {
    659       EIGEN_STATIC_ASSERT((EIGEN_IMPLIES(MaxRowsAtCompileTime==1 && MaxColsAtCompileTime!=1, (Options&RowMajor)==RowMajor)
    660                         && EIGEN_IMPLIES(MaxColsAtCompileTime==1 && MaxRowsAtCompileTime!=1, (Options&RowMajor)==0)
    661                         && ((RowsAtCompileTime == Dynamic) || (RowsAtCompileTime >= 0))
    662                         && ((ColsAtCompileTime == Dynamic) || (ColsAtCompileTime >= 0))
    663                         && ((MaxRowsAtCompileTime == Dynamic) || (MaxRowsAtCompileTime >= 0))
    664                         && ((MaxColsAtCompileTime == Dynamic) || (MaxColsAtCompileTime >= 0))
    665                         && (MaxRowsAtCompileTime == RowsAtCompileTime || RowsAtCompileTime==Dynamic)
    666                         && (MaxColsAtCompileTime == ColsAtCompileTime || ColsAtCompileTime==Dynamic)
    667                         && (Options & (DontAlign|RowMajor)) == Options),
    668         INVALID_MATRIX_TEMPLATE_PARAMETERS)
    669     }
    670 #endif
    671 
    672 private:
    673     enum { ThisConstantIsPrivateInPlainObjectBase };
    674 };
    675 
    676 namespace internal {
    677 
    678 template <typename Derived, typename OtherDerived, bool IsVector>
    679 struct conservative_resize_like_impl
    680 {
    681   typedef typename Derived::Index Index;
    682   static void run(DenseBase<Derived>& _this, Index rows, Index cols)
    683   {
    684     if (_this.rows() == rows && _this.cols() == cols) return;
    685     EIGEN_STATIC_ASSERT_DYNAMIC_SIZE(Derived)
    686 
    687     if ( ( Derived::IsRowMajor && _this.cols() == cols) || // row-major and we change only the number of rows
    688          (!Derived::IsRowMajor && _this.rows() == rows) )  // column-major and we change only the number of columns
    689     {
    690       internal::check_rows_cols_for_overflow<Derived::MaxSizeAtCompileTime>::run(rows, cols);
    691       _this.derived().m_storage.conservativeResize(rows*cols,rows,cols);
    692     }
    693     else
    694     {
    695       // The storage order does not allow us to use reallocation.
    696       typename Derived::PlainObject tmp(rows,cols);
    697       const Index common_rows = (std::min)(rows, _this.rows());
    698       const Index common_cols = (std::min)(cols, _this.cols());
    699       tmp.block(0,0,common_rows,common_cols) = _this.block(0,0,common_rows,common_cols);
    700       _this.derived().swap(tmp);
    701     }
    702   }
    703 
    704   static void run(DenseBase<Derived>& _this, const DenseBase<OtherDerived>& other)
    705   {
    706     if (_this.rows() == other.rows() && _this.cols() == other.cols()) return;
    707 
    708     // Note: Here is space for improvement. Basically, for conservativeResize(Index,Index),
    709     // neither RowsAtCompileTime or ColsAtCompileTime must be Dynamic. If only one of the
    710     // dimensions is dynamic, one could use either conservativeResize(Index rows, NoChange_t) or
    711     // conservativeResize(NoChange_t, Index cols). For these methods new static asserts like
    712     // EIGEN_STATIC_ASSERT_DYNAMIC_ROWS and EIGEN_STATIC_ASSERT_DYNAMIC_COLS would be good.
    713     EIGEN_STATIC_ASSERT_DYNAMIC_SIZE(Derived)
    714     EIGEN_STATIC_ASSERT_DYNAMIC_SIZE(OtherDerived)
    715 
    716     if ( ( Derived::IsRowMajor && _this.cols() == other.cols()) || // row-major and we change only the number of rows
    717          (!Derived::IsRowMajor && _this.rows() == other.rows()) )  // column-major and we change only the number of columns
    718     {
    719       const Index new_rows = other.rows() - _this.rows();
    720       const Index new_cols = other.cols() - _this.cols();
    721       _this.derived().m_storage.conservativeResize(other.size(),other.rows(),other.cols());
    722       if (new_rows>0)
    723         _this.bottomRightCorner(new_rows, other.cols()) = other.bottomRows(new_rows);
    724       else if (new_cols>0)
    725         _this.bottomRightCorner(other.rows(), new_cols) = other.rightCols(new_cols);
    726     }
    727     else
    728     {
    729       // The storage order does not allow us to use reallocation.
    730       typename Derived::PlainObject tmp(other);
    731       const Index common_rows = (std::min)(tmp.rows(), _this.rows());
    732       const Index common_cols = (std::min)(tmp.cols(), _this.cols());
    733       tmp.block(0,0,common_rows,common_cols) = _this.block(0,0,common_rows,common_cols);
    734       _this.derived().swap(tmp);
    735     }
    736   }
    737 };
    738 
    739 // Here, the specialization for vectors inherits from the general matrix case
    740 // to allow calling .conservativeResize(rows,cols) on vectors.
    741 template <typename Derived, typename OtherDerived>
    742 struct conservative_resize_like_impl<Derived,OtherDerived,true>
    743   : conservative_resize_like_impl<Derived,OtherDerived,false>
    744 {
    745   using conservative_resize_like_impl<Derived,OtherDerived,false>::run;
    746 
    747   typedef typename Derived::Index Index;
    748   static void run(DenseBase<Derived>& _this, Index size)
    749   {
    750     const Index new_rows = Derived::RowsAtCompileTime==1 ? 1 : size;
    751     const Index new_cols = Derived::RowsAtCompileTime==1 ? size : 1;
    752     _this.derived().m_storage.conservativeResize(size,new_rows,new_cols);
    753   }
    754 
    755   static void run(DenseBase<Derived>& _this, const DenseBase<OtherDerived>& other)
    756   {
    757     if (_this.rows() == other.rows() && _this.cols() == other.cols()) return;
    758 
    759     const Index num_new_elements = other.size() - _this.size();
    760 
    761     const Index new_rows = Derived::RowsAtCompileTime==1 ? 1 : other.rows();
    762     const Index new_cols = Derived::RowsAtCompileTime==1 ? other.cols() : 1;
    763     _this.derived().m_storage.conservativeResize(other.size(),new_rows,new_cols);
    764 
    765     if (num_new_elements > 0)
    766       _this.tail(num_new_elements) = other.tail(num_new_elements);
    767   }
    768 };
    769 
    770 template<typename MatrixTypeA, typename MatrixTypeB, bool SwapPointers>
    771 struct matrix_swap_impl
    772 {
    773   static inline void run(MatrixTypeA& a, MatrixTypeB& b)
    774   {
    775     a.base().swap(b);
    776   }
    777 };
    778 
    779 template<typename MatrixTypeA, typename MatrixTypeB>
    780 struct matrix_swap_impl<MatrixTypeA, MatrixTypeB, true>
    781 {
    782   static inline void run(MatrixTypeA& a, MatrixTypeB& b)
    783   {
    784     static_cast<typename MatrixTypeA::Base&>(a).m_storage.swap(static_cast<typename MatrixTypeB::Base&>(b).m_storage);
    785   }
    786 };
    787 
    788 } // end namespace internal
    789 
    790 } // end namespace Eigen
    791 
    792 #endif // EIGEN_DENSESTORAGEBASE_H
    793