// This file is part of Eigen, a lightweight C++ template library // for linear algebra. // // Copyright (C) 2006-2009 Benoit Jacob // Copyright (C) 2008 Gael Guennebaud // // This Source Code Form is subject to the terms of the Mozilla // Public License v. 2.0. If a copy of the MPL was not distributed // with this file, You can obtain one at http://mozilla.org/MPL/2.0/. #ifndef EIGEN_MATRIXBASE_H #define EIGEN_MATRIXBASE_H namespace Eigen { /** \class MatrixBase * \ingroup Core_Module * * \brief Base class for all dense matrices, vectors, and expressions * * This class is the base that is inherited by all matrix, vector, and related expression * types. Most of the Eigen API is contained in this class, and its base classes. Other important * classes for the Eigen API are Matrix, and VectorwiseOp. * * Note that some methods are defined in other modules such as the \ref LU_Module LU module * for all functions related to matrix inversions. * * \tparam Derived is the derived type, e.g. a matrix type, or an expression, etc. * * When writing a function taking Eigen objects as argument, if you want your function * to take as argument any matrix, vector, or expression, just let it take a * MatrixBase argument. As an example, here is a function printFirstRow which, given * a matrix, vector, or expression \a x, prints the first row of \a x. * * \code template void printFirstRow(const Eigen::MatrixBase& x) { cout << x.row(0) << endl; } * \endcode * * This class can be extended with the help of the plugin mechanism described on the page * \ref TopicCustomizingEigen by defining the preprocessor symbol \c EIGEN_MATRIXBASE_PLUGIN. * * \sa \ref TopicClassHierarchy */ template class MatrixBase : public DenseBase { public: #ifndef EIGEN_PARSED_BY_DOXYGEN typedef MatrixBase StorageBaseType; typedef typename internal::traits::StorageKind StorageKind; typedef typename internal::traits::Index Index; typedef typename internal::traits::Scalar Scalar; typedef typename internal::packet_traits::type PacketScalar; typedef typename NumTraits::Real RealScalar; typedef DenseBase Base; using Base::RowsAtCompileTime; using Base::ColsAtCompileTime; using Base::SizeAtCompileTime; using Base::MaxRowsAtCompileTime; using Base::MaxColsAtCompileTime; using Base::MaxSizeAtCompileTime; using Base::IsVectorAtCompileTime; using Base::Flags; using Base::CoeffReadCost; using Base::derived; using Base::const_cast_derived; using Base::rows; using Base::cols; using Base::size; using Base::coeff; using Base::coeffRef; using Base::lazyAssign; using Base::eval; using Base::operator+=; using Base::operator-=; using Base::operator*=; using Base::operator/=; typedef typename Base::CoeffReturnType CoeffReturnType; typedef typename Base::ConstTransposeReturnType ConstTransposeReturnType; typedef typename Base::RowXpr RowXpr; typedef typename Base::ColXpr ColXpr; #endif // not EIGEN_PARSED_BY_DOXYGEN #ifndef EIGEN_PARSED_BY_DOXYGEN /** type of the equivalent square matrix */ typedef Matrix SquareMatrixType; #endif // not EIGEN_PARSED_BY_DOXYGEN /** \returns the size of the main diagonal, which is min(rows(),cols()). * \sa rows(), cols(), SizeAtCompileTime. */ inline Index diagonalSize() const { return (std::min)(rows(),cols()); } /** \brief The plain matrix type corresponding to this expression. * * This is not necessarily exactly the return type of eval(). In the case of plain matrices, * the return type of eval() is a const reference to a matrix, not a matrix! It is however guaranteed * that the return type of eval() is either PlainObject or const PlainObject&. */ typedef Matrix::Scalar, internal::traits::RowsAtCompileTime, internal::traits::ColsAtCompileTime, AutoAlign | (internal::traits::Flags&RowMajorBit ? RowMajor : ColMajor), internal::traits::MaxRowsAtCompileTime, internal::traits::MaxColsAtCompileTime > PlainObject; #ifndef EIGEN_PARSED_BY_DOXYGEN /** \internal Represents a matrix with all coefficients equal to one another*/ typedef CwiseNullaryOp,Derived> ConstantReturnType; /** \internal the return type of MatrixBase::adjoint() */ typedef typename internal::conditional::IsComplex, CwiseUnaryOp, ConstTransposeReturnType>, ConstTransposeReturnType >::type AdjointReturnType; /** \internal Return type of eigenvalues() */ typedef Matrix, internal::traits::ColsAtCompileTime, 1, ColMajor> EigenvaluesReturnType; /** \internal the return type of identity */ typedef CwiseNullaryOp,Derived> IdentityReturnType; /** \internal the return type of unit vectors */ typedef Block, SquareMatrixType>, internal::traits::RowsAtCompileTime, internal::traits::ColsAtCompileTime> BasisReturnType; #endif // not EIGEN_PARSED_BY_DOXYGEN #define EIGEN_CURRENT_STORAGE_BASE_CLASS Eigen::MatrixBase # include "../plugins/CommonCwiseUnaryOps.h" # include "../plugins/CommonCwiseBinaryOps.h" # include "../plugins/MatrixCwiseUnaryOps.h" # include "../plugins/MatrixCwiseBinaryOps.h" # ifdef EIGEN_MATRIXBASE_PLUGIN # include EIGEN_MATRIXBASE_PLUGIN # endif #undef EIGEN_CURRENT_STORAGE_BASE_CLASS /** Special case of the template operator=, in order to prevent the compiler * from generating a default operator= (issue hit with g++ 4.1) */ Derived& operator=(const MatrixBase& other); // We cannot inherit here via Base::operator= since it is causing // trouble with MSVC. template Derived& operator=(const DenseBase& other); template Derived& operator=(const EigenBase& other); template Derived& operator=(const ReturnByValue& other); #ifndef EIGEN_PARSED_BY_DOXYGEN template Derived& lazyAssign(const ProductBase& other); template Derived& lazyAssign(const MatrixPowerProduct& other); #endif // not EIGEN_PARSED_BY_DOXYGEN template Derived& operator+=(const MatrixBase& other); template Derived& operator-=(const MatrixBase& other); template const typename ProductReturnType::Type operator*(const MatrixBase &other) const; template const typename LazyProductReturnType::Type lazyProduct(const MatrixBase &other) const; template Derived& operator*=(const EigenBase& other); template void applyOnTheLeft(const EigenBase& other); template void applyOnTheRight(const EigenBase& other); template const DiagonalProduct operator*(const DiagonalBase &diagonal) const; template typename internal::scalar_product_traits::Scalar,typename internal::traits::Scalar>::ReturnType dot(const MatrixBase& other) const; #ifdef EIGEN2_SUPPORT template Scalar eigen2_dot(const MatrixBase& other) const; #endif RealScalar squaredNorm() const; RealScalar norm() const; RealScalar stableNorm() const; RealScalar blueNorm() const; RealScalar hypotNorm() const; const PlainObject normalized() const; void normalize(); const AdjointReturnType adjoint() const; void adjointInPlace(); typedef Diagonal DiagonalReturnType; DiagonalReturnType diagonal(); typedef typename internal::add_const >::type ConstDiagonalReturnType; ConstDiagonalReturnType diagonal() const; template struct DiagonalIndexReturnType { typedef Diagonal Type; }; template struct ConstDiagonalIndexReturnType { typedef const Diagonal Type; }; template typename DiagonalIndexReturnType::Type diagonal(); template typename ConstDiagonalIndexReturnType::Type diagonal() const; typedef Diagonal DiagonalDynamicIndexReturnType; typedef typename internal::add_const >::type ConstDiagonalDynamicIndexReturnType; DiagonalDynamicIndexReturnType diagonal(Index index); ConstDiagonalDynamicIndexReturnType diagonal(Index index) const; #ifdef EIGEN2_SUPPORT template typename internal::eigen2_part_return_type::type part(); template const typename internal::eigen2_part_return_type::type part() const; // huuuge hack. make Eigen2's matrix.part() work in eigen3. Problem: Diagonal is now a class template instead // of an integer constant. Solution: overload the part() method template wrt template parameters list. template class U> const DiagonalWrapper part() const { return diagonal().asDiagonal(); } #endif // EIGEN2_SUPPORT template struct TriangularViewReturnType { typedef TriangularView Type; }; template struct ConstTriangularViewReturnType { typedef const TriangularView Type; }; template typename TriangularViewReturnType::Type triangularView(); template typename ConstTriangularViewReturnType::Type triangularView() const; template struct SelfAdjointViewReturnType { typedef SelfAdjointView Type; }; template struct ConstSelfAdjointViewReturnType { typedef const SelfAdjointView Type; }; template typename SelfAdjointViewReturnType::Type selfadjointView(); template typename ConstSelfAdjointViewReturnType::Type selfadjointView() const; const SparseView sparseView(const Scalar& m_reference = Scalar(0), const typename NumTraits::Real& m_epsilon = NumTraits::dummy_precision()) const; static const IdentityReturnType Identity(); static const IdentityReturnType Identity(Index rows, Index cols); static const BasisReturnType Unit(Index size, Index i); static const BasisReturnType Unit(Index i); static const BasisReturnType UnitX(); static const BasisReturnType UnitY(); static const BasisReturnType UnitZ(); static const BasisReturnType UnitW(); const DiagonalWrapper asDiagonal() const; const PermutationWrapper asPermutation() const; Derived& setIdentity(); Derived& setIdentity(Index rows, Index cols); bool isIdentity(const RealScalar& prec = NumTraits::dummy_precision()) const; bool isDiagonal(const RealScalar& prec = NumTraits::dummy_precision()) const; bool isUpperTriangular(const RealScalar& prec = NumTraits::dummy_precision()) const; bool isLowerTriangular(const RealScalar& prec = NumTraits::dummy_precision()) const; template bool isOrthogonal(const MatrixBase& other, const RealScalar& prec = NumTraits::dummy_precision()) const; bool isUnitary(const RealScalar& prec = NumTraits::dummy_precision()) const; /** \returns true if each coefficients of \c *this and \a other are all exactly equal. * \warning When using floating point scalar values you probably should rather use a * fuzzy comparison such as isApprox() * \sa isApprox(), operator!= */ template inline bool operator==(const MatrixBase& other) const { return cwiseEqual(other).all(); } /** \returns true if at least one pair of coefficients of \c *this and \a other are not exactly equal to each other. * \warning When using floating point scalar values you probably should rather use a * fuzzy comparison such as isApprox() * \sa isApprox(), operator== */ template inline bool operator!=(const MatrixBase& other) const { return cwiseNotEqual(other).any(); } NoAlias noalias(); inline const ForceAlignedAccess forceAlignedAccess() const; inline ForceAlignedAccess forceAlignedAccess(); template inline typename internal::add_const_on_value_type,Derived&>::type>::type forceAlignedAccessIf() const; template inline typename internal::conditional,Derived&>::type forceAlignedAccessIf(); Scalar trace() const; /////////// Array module /////////// template RealScalar lpNorm() const; MatrixBase& matrix() { return *this; } const MatrixBase& matrix() const { return *this; } /** \returns an \link Eigen::ArrayBase Array \endlink expression of this matrix * \sa ArrayBase::matrix() */ ArrayWrapper array() { return derived(); } const ArrayWrapper array() const { return derived(); } /////////// LU module /////////// const FullPivLU fullPivLu() const; const PartialPivLU partialPivLu() const; #if EIGEN2_SUPPORT_STAGE < STAGE20_RESOLVE_API_CONFLICTS const LU lu() const; #endif #ifdef EIGEN2_SUPPORT const LU eigen2_lu() const; #endif #if EIGEN2_SUPPORT_STAGE > STAGE20_RESOLVE_API_CONFLICTS const PartialPivLU lu() const; #endif #ifdef EIGEN2_SUPPORT template void computeInverse(MatrixBase *result) const { *result = this->inverse(); } #endif const internal::inverse_impl inverse() const; template void computeInverseAndDetWithCheck( ResultType& inverse, typename ResultType::Scalar& determinant, bool& invertible, const RealScalar& absDeterminantThreshold = NumTraits::dummy_precision() ) const; template void computeInverseWithCheck( ResultType& inverse, bool& invertible, const RealScalar& absDeterminantThreshold = NumTraits::dummy_precision() ) const; Scalar determinant() const; /////////// Cholesky module /////////// const LLT llt() const; const LDLT ldlt() const; /////////// QR module /////////// const HouseholderQR householderQr() const; const ColPivHouseholderQR colPivHouseholderQr() const; const FullPivHouseholderQR fullPivHouseholderQr() const; #ifdef EIGEN2_SUPPORT const QR qr() const; #endif EigenvaluesReturnType eigenvalues() const; RealScalar operatorNorm() const; /////////// SVD module /////////// JacobiSVD jacobiSvd(unsigned int computationOptions = 0) const; #ifdef EIGEN2_SUPPORT SVD svd() const; #endif /////////// Geometry module /////////// #ifndef EIGEN_PARSED_BY_DOXYGEN /// \internal helper struct to form the return type of the cross product template struct cross_product_return_type { typedef typename internal::scalar_product_traits::Scalar,typename internal::traits::Scalar>::ReturnType Scalar; typedef Matrix type; }; #endif // EIGEN_PARSED_BY_DOXYGEN template typename cross_product_return_type::type cross(const MatrixBase& other) const; template PlainObject cross3(const MatrixBase& other) const; PlainObject unitOrthogonal(void) const; Matrix eulerAngles(Index a0, Index a1, Index a2) const; #if EIGEN2_SUPPORT_STAGE > STAGE20_RESOLVE_API_CONFLICTS ScalarMultipleReturnType operator*(const UniformScaling& s) const; // put this as separate enum value to work around possible GCC 4.3 bug (?) enum { HomogeneousReturnTypeDirection = ColsAtCompileTime==1?Vertical:Horizontal }; typedef Homogeneous HomogeneousReturnType; HomogeneousReturnType homogeneous() const; #endif enum { SizeMinusOne = SizeAtCompileTime==Dynamic ? Dynamic : SizeAtCompileTime-1 }; typedef Block::ColsAtCompileTime==1 ? SizeMinusOne : 1, internal::traits::ColsAtCompileTime==1 ? 1 : SizeMinusOne> ConstStartMinusOne; typedef CwiseUnaryOp::Scalar>, const ConstStartMinusOne > HNormalizedReturnType; const HNormalizedReturnType hnormalized() const; ////////// Householder module /////////// void makeHouseholderInPlace(Scalar& tau, RealScalar& beta); template void makeHouseholder(EssentialPart& essential, Scalar& tau, RealScalar& beta) const; template void applyHouseholderOnTheLeft(const EssentialPart& essential, const Scalar& tau, Scalar* workspace); template void applyHouseholderOnTheRight(const EssentialPart& essential, const Scalar& tau, Scalar* workspace); ///////// Jacobi module ///////// template void applyOnTheLeft(Index p, Index q, const JacobiRotation& j); template void applyOnTheRight(Index p, Index q, const JacobiRotation& j); ///////// MatrixFunctions module ///////// typedef typename internal::stem_function::type StemFunction; const MatrixExponentialReturnValue exp() const; const MatrixFunctionReturnValue matrixFunction(StemFunction f) const; const MatrixFunctionReturnValue cosh() const; const MatrixFunctionReturnValue sinh() const; const MatrixFunctionReturnValue cos() const; const MatrixFunctionReturnValue sin() const; const MatrixSquareRootReturnValue sqrt() const; const MatrixLogarithmReturnValue log() const; const MatrixPowerReturnValue pow(const RealScalar& p) const; #ifdef EIGEN2_SUPPORT template Derived& operator+=(const Flagged, 0, EvalBeforeAssigningBit>& other); template Derived& operator-=(const Flagged, 0, EvalBeforeAssigningBit>& other); /** \deprecated because .lazy() is deprecated * Overloaded for cache friendly product evaluation */ template Derived& lazyAssign(const Flagged& other) { return lazyAssign(other._expression()); } template const Flagged marked() const; const Flagged lazy() const; inline const Cwise cwise() const; inline Cwise cwise(); VectorBlock start(Index size); const VectorBlock start(Index size) const; VectorBlock end(Index size); const VectorBlock end(Index size) const; template VectorBlock start(); template const VectorBlock start() const; template VectorBlock end(); template const VectorBlock end() const; Minor minor(Index row, Index col); const Minor minor(Index row, Index col) const; #endif protected: MatrixBase() : Base() {} private: explicit MatrixBase(int); MatrixBase(int,int); template explicit MatrixBase(const MatrixBase&); protected: // mixing arrays and matrices is not legal template Derived& operator+=(const ArrayBase& ) {EIGEN_STATIC_ASSERT(std::ptrdiff_t(sizeof(typename OtherDerived::Scalar))==-1,YOU_CANNOT_MIX_ARRAYS_AND_MATRICES); return *this;} // mixing arrays and matrices is not legal template Derived& operator-=(const ArrayBase& ) {EIGEN_STATIC_ASSERT(std::ptrdiff_t(sizeof(typename OtherDerived::Scalar))==-1,YOU_CANNOT_MIX_ARRAYS_AND_MATRICES); return *this;} }; /*************************************************************************** * Implementation of matrix base methods ***************************************************************************/ /** replaces \c *this by \c *this * \a other. * * \returns a reference to \c *this * * Example: \include MatrixBase_applyOnTheRight.cpp * Output: \verbinclude MatrixBase_applyOnTheRight.out */ template template inline Derived& MatrixBase::operator*=(const EigenBase &other) { other.derived().applyThisOnTheRight(derived()); return derived(); } /** replaces \c *this by \c *this * \a other. It is equivalent to MatrixBase::operator*=(). * * Example: \include MatrixBase_applyOnTheRight.cpp * Output: \verbinclude MatrixBase_applyOnTheRight.out */ template template inline void MatrixBase::applyOnTheRight(const EigenBase &other) { other.derived().applyThisOnTheRight(derived()); } /** replaces \c *this by \a other * \c *this. * * Example: \include MatrixBase_applyOnTheLeft.cpp * Output: \verbinclude MatrixBase_applyOnTheLeft.out */ template template inline void MatrixBase::applyOnTheLeft(const EigenBase &other) { other.derived().applyThisOnTheLeft(derived()); } } // end namespace Eigen #endif // EIGEN_MATRIXBASE_H