00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011 #ifndef EIGEN_DIAGONALMATRIX_H
00012 #define EIGEN_DIAGONALMATRIX_H
00013
00014 namespace Eigen {
00015
00016 #ifndef EIGEN_PARSED_BY_DOXYGEN
00017 template<typename Derived>
00018 class DiagonalBase : public EigenBase<Derived>
00019 {
00020 public:
00021 typedef typename internal::traits<Derived>::DiagonalVectorType DiagonalVectorType;
00022 typedef typename DiagonalVectorType::Scalar Scalar;
00023 typedef typename DiagonalVectorType::RealScalar RealScalar;
00024 typedef typename internal::traits<Derived>::StorageKind StorageKind;
00025 typedef typename internal::traits<Derived>::Index Index;
00026
00027 enum {
00028 RowsAtCompileTime = DiagonalVectorType::SizeAtCompileTime,
00029 ColsAtCompileTime = DiagonalVectorType::SizeAtCompileTime,
00030 MaxRowsAtCompileTime = DiagonalVectorType::MaxSizeAtCompileTime,
00031 MaxColsAtCompileTime = DiagonalVectorType::MaxSizeAtCompileTime,
00032 IsVectorAtCompileTime = 0,
00033 Flags = 0
00034 };
00035
00036 typedef Matrix<Scalar, RowsAtCompileTime, ColsAtCompileTime, 0, MaxRowsAtCompileTime, MaxColsAtCompileTime> DenseMatrixType;
00037 typedef DenseMatrixType DenseType;
00038 typedef DiagonalMatrix<Scalar,DiagonalVectorType::SizeAtCompileTime,DiagonalVectorType::MaxSizeAtCompileTime> PlainObject;
00039
00040 inline const Derived& derived() const { return *static_cast<const Derived*>(this); }
00041 inline Derived& derived() { return *static_cast<Derived*>(this); }
00042
00043 DenseMatrixType toDenseMatrix() const { return derived(); }
00044 template<typename DenseDerived>
00045 void evalTo(MatrixBase<DenseDerived> &other) const;
00046 template<typename DenseDerived>
00047 void addTo(MatrixBase<DenseDerived> &other) const
00048 { other.diagonal() += diagonal(); }
00049 template<typename DenseDerived>
00050 void subTo(MatrixBase<DenseDerived> &other) const
00051 { other.diagonal() -= diagonal(); }
00052
00053 inline const DiagonalVectorType& diagonal() const { return derived().diagonal(); }
00054 inline DiagonalVectorType& diagonal() { return derived().diagonal(); }
00055
00056 inline Index rows() const { return diagonal().size(); }
00057 inline Index cols() const { return diagonal().size(); }
00058
00059 template<typename MatrixDerived>
00060 const DiagonalProduct<MatrixDerived, Derived, OnTheLeft>
00061 operator*(const MatrixBase<MatrixDerived> &matrix) const;
00062
00063 inline const DiagonalWrapper<const CwiseUnaryOp<internal::scalar_inverse_op<Scalar>, const DiagonalVectorType> >
00064 inverse() const
00065 {
00066 return diagonal().cwiseInverse();
00067 }
00068
00069 inline const DiagonalWrapper<const CwiseUnaryOp<internal::scalar_multiple_op<Scalar>, const DiagonalVectorType> >
00070 operator*(const Scalar& scalar) const
00071 {
00072 return diagonal() * scalar;
00073 }
00074 friend inline const DiagonalWrapper<const CwiseUnaryOp<internal::scalar_multiple_op<Scalar>, const DiagonalVectorType> >
00075 operator*(const Scalar& scalar, const DiagonalBase& other)
00076 {
00077 return other.diagonal() * scalar;
00078 }
00079
00080 #ifdef EIGEN2_SUPPORT
00081 template<typename OtherDerived>
00082 bool isApprox(const DiagonalBase<OtherDerived>& other, typename NumTraits<Scalar>::Real precision = NumTraits<Scalar>::dummy_precision()) const
00083 {
00084 return diagonal().isApprox(other.diagonal(), precision);
00085 }
00086 template<typename OtherDerived>
00087 bool isApprox(const MatrixBase<OtherDerived>& other, typename NumTraits<Scalar>::Real precision = NumTraits<Scalar>::dummy_precision()) const
00088 {
00089 return toDenseMatrix().isApprox(other, precision);
00090 }
00091 #endif
00092 };
00093
00094 template<typename Derived>
00095 template<typename DenseDerived>
00096 void DiagonalBase<Derived>::evalTo(MatrixBase<DenseDerived> &other) const
00097 {
00098 other.setZero();
00099 other.diagonal() = diagonal();
00100 }
00101 #endif
00102
00116 namespace internal {
00117 template<typename _Scalar, int SizeAtCompileTime, int MaxSizeAtCompileTime>
00118 struct traits<DiagonalMatrix<_Scalar,SizeAtCompileTime,MaxSizeAtCompileTime> >
00119 : traits<Matrix<_Scalar,SizeAtCompileTime,SizeAtCompileTime,0,MaxSizeAtCompileTime,MaxSizeAtCompileTime> >
00120 {
00121 typedef Matrix<_Scalar,SizeAtCompileTime,1,0,MaxSizeAtCompileTime,1> DiagonalVectorType;
00122 typedef Dense StorageKind;
00123 typedef DenseIndex Index;
00124 enum {
00125 Flags = LvalueBit
00126 };
00127 };
00128 }
00129 template<typename _Scalar, int SizeAtCompileTime, int MaxSizeAtCompileTime>
00130 class DiagonalMatrix
00131 : public DiagonalBase<DiagonalMatrix<_Scalar,SizeAtCompileTime,MaxSizeAtCompileTime> >
00132 {
00133 public:
00134 #ifndef EIGEN_PARSED_BY_DOXYGEN
00135 typedef typename internal::traits<DiagonalMatrix>::DiagonalVectorType DiagonalVectorType;
00136 typedef const DiagonalMatrix& Nested;
00137 typedef _Scalar Scalar;
00138 typedef typename internal::traits<DiagonalMatrix>::StorageKind StorageKind;
00139 typedef typename internal::traits<DiagonalMatrix>::Index Index;
00140 #endif
00141
00142 protected:
00143
00144 DiagonalVectorType m_diagonal;
00145
00146 public:
00147
00149 inline const DiagonalVectorType& diagonal() const { return m_diagonal; }
00151 inline DiagonalVectorType& diagonal() { return m_diagonal; }
00152
00154 inline DiagonalMatrix() {}
00155
00157 inline DiagonalMatrix(Index dim) : m_diagonal(dim) {}
00158
00160 inline DiagonalMatrix(const Scalar& x, const Scalar& y) : m_diagonal(x,y) {}
00161
00163 inline DiagonalMatrix(const Scalar& x, const Scalar& y, const Scalar& z) : m_diagonal(x,y,z) {}
00164
00166 template<typename OtherDerived>
00167 inline DiagonalMatrix(const DiagonalBase<OtherDerived>& other) : m_diagonal(other.diagonal()) {}
00168
00169 #ifndef EIGEN_PARSED_BY_DOXYGEN
00170
00171 inline DiagonalMatrix(const DiagonalMatrix& other) : m_diagonal(other.diagonal()) {}
00172 #endif
00173
00175 template<typename OtherDerived>
00176 explicit inline DiagonalMatrix(const MatrixBase<OtherDerived>& other) : m_diagonal(other)
00177 {}
00178
00180 template<typename OtherDerived>
00181 DiagonalMatrix& operator=(const DiagonalBase<OtherDerived>& other)
00182 {
00183 m_diagonal = other.diagonal();
00184 return *this;
00185 }
00186
00187 #ifndef EIGEN_PARSED_BY_DOXYGEN
00188
00191 DiagonalMatrix& operator=(const DiagonalMatrix& other)
00192 {
00193 m_diagonal = other.diagonal();
00194 return *this;
00195 }
00196 #endif
00197
00199 inline void resize(Index size) { m_diagonal.resize(size); }
00201 inline void setZero() { m_diagonal.setZero(); }
00203 inline void setZero(Index size) { m_diagonal.setZero(size); }
00205 inline void setIdentity() { m_diagonal.setOnes(); }
00207 inline void setIdentity(Index size) { m_diagonal.setOnes(size); }
00208 };
00209
00224 namespace internal {
00225 template<typename _DiagonalVectorType>
00226 struct traits<DiagonalWrapper<_DiagonalVectorType> >
00227 {
00228 typedef _DiagonalVectorType DiagonalVectorType;
00229 typedef typename DiagonalVectorType::Scalar Scalar;
00230 typedef typename DiagonalVectorType::Index Index;
00231 typedef typename DiagonalVectorType::StorageKind StorageKind;
00232 enum {
00233 RowsAtCompileTime = DiagonalVectorType::SizeAtCompileTime,
00234 ColsAtCompileTime = DiagonalVectorType::SizeAtCompileTime,
00235 MaxRowsAtCompileTime = DiagonalVectorType::SizeAtCompileTime,
00236 MaxColsAtCompileTime = DiagonalVectorType::SizeAtCompileTime,
00237 Flags = traits<DiagonalVectorType>::Flags & LvalueBit
00238 };
00239 };
00240 }
00241
00242 template<typename _DiagonalVectorType>
00243 class DiagonalWrapper
00244 : public DiagonalBase<DiagonalWrapper<_DiagonalVectorType> >, internal::no_assignment_operator
00245 {
00246 public:
00247 #ifndef EIGEN_PARSED_BY_DOXYGEN
00248 typedef _DiagonalVectorType DiagonalVectorType;
00249 typedef DiagonalWrapper Nested;
00250 #endif
00251
00253 inline DiagonalWrapper(DiagonalVectorType& diagonal) : m_diagonal(diagonal) {}
00254
00256 const DiagonalVectorType& diagonal() const { return m_diagonal; }
00257
00258 protected:
00259 typename DiagonalVectorType::Nested m_diagonal;
00260 };
00261
00271 template<typename Derived>
00272 inline const DiagonalWrapper<const Derived>
00273 MatrixBase<Derived>::asDiagonal() const
00274 {
00275 return derived();
00276 }
00277
00286 template<typename Derived>
00287 bool MatrixBase<Derived>::isDiagonal(RealScalar prec) const
00288 {
00289 if(cols() != rows()) return false;
00290 RealScalar maxAbsOnDiagonal = static_cast<RealScalar>(-1);
00291 for(Index j = 0; j < cols(); ++j)
00292 {
00293 RealScalar absOnDiagonal = internal::abs(coeff(j,j));
00294 if(absOnDiagonal > maxAbsOnDiagonal) maxAbsOnDiagonal = absOnDiagonal;
00295 }
00296 for(Index j = 0; j < cols(); ++j)
00297 for(Index i = 0; i < j; ++i)
00298 {
00299 if(!internal::isMuchSmallerThan(coeff(i, j), maxAbsOnDiagonal, prec)) return false;
00300 if(!internal::isMuchSmallerThan(coeff(j, i), maxAbsOnDiagonal, prec)) return false;
00301 }
00302 return true;
00303 }
00304
00305 }
00306
00307 #endif // EIGEN_DIAGONALMATRIX_H