MatrixBase.h
Go to the documentation of this file.
1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2006-2009 Benoit Jacob <jacob.benoit.1@gmail.com>
5 // Copyright (C) 2008 Gael Guennebaud <gael.guennebaud@inria.fr>
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_MATRIXBASE_H
12 #define EIGEN_MATRIXBASE_H
13 
14 namespace Eigen {
15 
48 template<typename Derived> class MatrixBase
49  : public DenseBase<Derived>
50 {
51  public:
52 #ifndef EIGEN_PARSED_BY_DOXYGEN
59 
68  using Base::Flags;
69  using Base::CoeffReadCost;
70 
71  using Base::derived;
72  using Base::const_cast_derived;
73  using Base::rows;
74  using Base::cols;
75  using Base::size;
76  using Base::coeff;
77  using Base::coeffRef;
78  using Base::lazyAssign;
79  using Base::eval;
80  using Base::operator+=;
81  using Base::operator-=;
82  using Base::operator*=;
83  using Base::operator/=;
84 
87  typedef typename Base::RowXpr RowXpr;
88  typedef typename Base::ColXpr ColXpr;
89 #endif // not EIGEN_PARSED_BY_DOXYGEN
90 
91 
92 
93 #ifndef EIGEN_PARSED_BY_DOXYGEN
94 
97 #endif // not EIGEN_PARSED_BY_DOXYGEN
98 
101  inline Index diagonalSize() const { return (std::min)(rows(),cols()); }
116 
117 #ifndef EIGEN_PARSED_BY_DOXYGEN
118 
122  CwiseUnaryOp<internal::scalar_conjugate_op<Scalar>, ConstTransposeReturnType>,
123  ConstTransposeReturnType
131  internal::traits<Derived>::RowsAtCompileTime,
132  internal::traits<Derived>::ColsAtCompileTime> BasisReturnType;
133 #endif // not EIGEN_PARSED_BY_DOXYGEN
134 
135 #define EIGEN_CURRENT_STORAGE_BASE_CLASS Eigen::MatrixBase
136 # include "../plugins/CommonCwiseUnaryOps.h"
137 # include "../plugins/CommonCwiseBinaryOps.h"
138 # include "../plugins/MatrixCwiseUnaryOps.h"
139 # include "../plugins/MatrixCwiseBinaryOps.h"
140 # ifdef EIGEN_MATRIXBASE_PLUGIN
141 # include EIGEN_MATRIXBASE_PLUGIN
142 # endif
143 #undef EIGEN_CURRENT_STORAGE_BASE_CLASS
144 
148  Derived& operator=(const MatrixBase& other);
149 
150  // We cannot inherit here via Base::operator= since it is causing
151  // trouble with MSVC.
152 
153  template <typename OtherDerived>
154  Derived& operator=(const DenseBase<OtherDerived>& other);
155 
156  template <typename OtherDerived>
157  Derived& operator=(const EigenBase<OtherDerived>& other);
159  template<typename OtherDerived>
160  Derived& operator=(const ReturnByValue<OtherDerived>& other);
161 
162 #ifndef EIGEN_PARSED_BY_DOXYGEN
163  template<typename ProductDerived, typename Lhs, typename Rhs>
164  Derived& lazyAssign(const ProductBase<ProductDerived, Lhs,Rhs>& other);
165 
166  template<typename MatrixPower, typename Lhs, typename Rhs>
168 #endif // not EIGEN_PARSED_BY_DOXYGEN
169 
170  template<typename OtherDerived>
171  Derived& operator+=(const MatrixBase<OtherDerived>& other);
172  template<typename OtherDerived>
173  Derived& operator-=(const MatrixBase<OtherDerived>& other);
174 
175  template<typename OtherDerived>
177  operator*(const MatrixBase<OtherDerived> &other) const;
178 
179  template<typename OtherDerived>
181  lazyProduct(const MatrixBase<OtherDerived> &other) const;
182 
183  template<typename OtherDerived>
184  Derived& operator*=(const EigenBase<OtherDerived>& other);
185 
186  template<typename OtherDerived>
187  void applyOnTheLeft(const EigenBase<OtherDerived>& other);
188 
189  template<typename OtherDerived>
190  void applyOnTheRight(const EigenBase<OtherDerived>& other);
191 
192  template<typename DiagonalDerived>
195 
196  template<typename OtherDerived>
198  dot(const MatrixBase<OtherDerived>& other) const;
199 
200  #ifdef EIGEN2_SUPPORT
201  template<typename OtherDerived>
202  Scalar eigen2_dot(const MatrixBase<OtherDerived>& other) const;
203  #endif
204 
205  RealScalar squaredNorm() const;
206  RealScalar norm() const;
207  RealScalar stableNorm() const;
208  RealScalar blueNorm() const;
209  RealScalar hypotNorm() const;
210  const PlainObject normalized() const;
211  void normalize();
212 
213  const AdjointReturnType adjoint() const;
214  void adjointInPlace();
215 
217  DiagonalReturnType diagonal();
219  ConstDiagonalReturnType diagonal() const;
220 
221  template<int Index> struct DiagonalIndexReturnType { typedef Diagonal<Derived,Index> Type; };
222  template<int Index> struct ConstDiagonalIndexReturnType { typedef const Diagonal<const Derived,Index> Type; };
223 
224  template<int Index> typename DiagonalIndexReturnType<Index>::Type diagonal();
225  template<int Index> typename ConstDiagonalIndexReturnType<Index>::Type diagonal() const;
226 
227  // Note: The "MatrixBase::" prefixes are added to help MSVC9 to match these declarations with the later implementations.
228  // On the other hand they confuse MSVC8...
229  #if (defined _MSC_VER) && (_MSC_VER >= 1500) // 2008 or later
230  typename MatrixBase::template DiagonalIndexReturnType<DynamicIndex>::Type diagonal(Index index);
231  typename MatrixBase::template ConstDiagonalIndexReturnType<DynamicIndex>::Type diagonal(Index index) const;
232  #else
235  #endif
236 
237  #ifdef EIGEN2_SUPPORT
238  template<unsigned int Mode> typename internal::eigen2_part_return_type<Derived, Mode>::type part();
239  template<unsigned int Mode> const typename internal::eigen2_part_return_type<Derived, Mode>::type part() const;
240 
241  // huuuge hack. make Eigen2's matrix.part<Diagonal>() work in eigen3. Problem: Diagonal is now a class template instead
242  // of an integer constant. Solution: overload the part() method template wrt template parameters list.
243  template<template<typename T, int N> class U>
245  { return diagonal().asDiagonal(); }
246  #endif // EIGEN2_SUPPORT
247 
248  template<unsigned int Mode> struct TriangularViewReturnType { typedef TriangularView<Derived, Mode> Type; };
249  template<unsigned int Mode> struct ConstTriangularViewReturnType { typedef const TriangularView<const Derived, Mode> Type; };
250 
251  template<unsigned int Mode> typename TriangularViewReturnType<Mode>::Type triangularView();
252  template<unsigned int Mode> typename ConstTriangularViewReturnType<Mode>::Type triangularView() const;
253 
254  template<unsigned int UpLo> struct SelfAdjointViewReturnType { typedef SelfAdjointView<Derived, UpLo> Type; };
255  template<unsigned int UpLo> struct ConstSelfAdjointViewReturnType { typedef const SelfAdjointView<const Derived, UpLo> Type; };
256 
257  template<unsigned int UpLo> typename SelfAdjointViewReturnType<UpLo>::Type selfadjointView();
258  template<unsigned int UpLo> typename ConstSelfAdjointViewReturnType<UpLo>::Type selfadjointView() const;
259 
260  const SparseView<Derived> sparseView(const Scalar& m_reference = Scalar(0),
261  const typename NumTraits<Scalar>::Real& m_epsilon = NumTraits<Scalar>::dummy_precision()) const;
262  static const IdentityReturnType Identity();
263  static const IdentityReturnType Identity(Index rows, Index cols);
264  static const BasisReturnType Unit(Index size, Index i);
265  static const BasisReturnType Unit(Index i);
266  static const BasisReturnType UnitX();
267  static const BasisReturnType UnitY();
268  static const BasisReturnType UnitZ();
269  static const BasisReturnType UnitW();
270 
273 
274  Derived& setIdentity();
275  Derived& setIdentity(Index rows, Index cols);
276 
277  bool isIdentity(const RealScalar& prec = NumTraits<Scalar>::dummy_precision()) const;
278  bool isDiagonal(const RealScalar& prec = NumTraits<Scalar>::dummy_precision()) const;
279 
280  bool isUpperTriangular(const RealScalar& prec = NumTraits<Scalar>::dummy_precision()) const;
281  bool isLowerTriangular(const RealScalar& prec = NumTraits<Scalar>::dummy_precision()) const;
282 
283  template<typename OtherDerived>
284  bool isOrthogonal(const MatrixBase<OtherDerived>& other,
285  const RealScalar& prec = NumTraits<Scalar>::dummy_precision()) const;
286  bool isUnitary(const RealScalar& prec = NumTraits<Scalar>::dummy_precision()) const;
287 
292  template<typename OtherDerived>
293  inline bool operator==(const MatrixBase<OtherDerived>& other) const
294  { return cwiseEqual(other).all(); }
295 
300  template<typename OtherDerived>
301  inline bool operator!=(const MatrixBase<OtherDerived>& other) const
302  { return cwiseNotEqual(other).any(); }
303 
305 
306  inline const ForceAlignedAccess<Derived> forceAlignedAccess() const;
309  template<bool Enable> inline typename internal::conditional<Enable,ForceAlignedAccess<Derived>,Derived&>::type forceAlignedAccessIf();
310 
311  Scalar trace() const;
312 
314 
315  template<int p> RealScalar lpNorm() const;
316 
317  MatrixBase<Derived>& matrix() { return *this; }
318  const MatrixBase<Derived>& matrix() const { return *this; }
319 
322  ArrayWrapper<Derived> array() { return derived(); }
323  const ArrayWrapper<const Derived> array() const { return derived(); }
324 
326 
327  const FullPivLU<PlainObject> fullPivLu() const;
329 
330  #if EIGEN2_SUPPORT_STAGE < STAGE20_RESOLVE_API_CONFLICTS
331  const LU<PlainObject> lu() const;
332  #endif
333 
334  #ifdef EIGEN2_SUPPORT
335  const LU<PlainObject> eigen2_lu() const;
336  #endif
337 
338  #if EIGEN2_SUPPORT_STAGE > STAGE20_RESOLVE_API_CONFLICTS
339  const PartialPivLU<PlainObject> lu() const;
340  #endif
341 
342  #ifdef EIGEN2_SUPPORT
343  template<typename ResultType>
344  void computeInverse(MatrixBase<ResultType> *result) const {
345  *result = this->inverse();
346  }
347  #endif
348 
350  template<typename ResultType>
352  ResultType& inverse,
353  typename ResultType::Scalar& determinant,
354  bool& invertible,
355  const RealScalar& absDeterminantThreshold = NumTraits<Scalar>::dummy_precision()
356  ) const;
357  template<typename ResultType>
359  ResultType& inverse,
360  bool& invertible,
361  const RealScalar& absDeterminantThreshold = NumTraits<Scalar>::dummy_precision()
362  ) const;
363  Scalar determinant() const;
364 
366 
367  const LLT<PlainObject> llt() const;
368  const LDLT<PlainObject> ldlt() const;
369 
371 
375 
376  #ifdef EIGEN2_SUPPORT
377  const QR<PlainObject> qr() const;
378  #endif
379 
380  EigenvaluesReturnType eigenvalues() const;
381  RealScalar operatorNorm() const;
382 
384 
385  JacobiSVD<PlainObject> jacobiSvd(unsigned int computationOptions = 0) const;
386 
387  #ifdef EIGEN2_SUPPORT
388  SVD<PlainObject> svd() const;
389  #endif
390 
392 
393  #ifndef EIGEN_PARSED_BY_DOXYGEN
394  template<typename OtherDerived> struct cross_product_return_type {
398  };
399  #endif // EIGEN_PARSED_BY_DOXYGEN
400  template<typename OtherDerived>
402  cross(const MatrixBase<OtherDerived>& other) const;
403  template<typename OtherDerived>
404  PlainObject cross3(const MatrixBase<OtherDerived>& other) const;
405  PlainObject unitOrthogonal(void) const;
406  Matrix<Scalar,3,1> eulerAngles(Index a0, Index a1, Index a2) const;
407 
408  #if EIGEN2_SUPPORT_STAGE > STAGE20_RESOLVE_API_CONFLICTS
410  // put this as separate enum value to work around possible GCC 4.3 bug (?)
411  enum { HomogeneousReturnTypeDirection = ColsAtCompileTime==1?Vertical:Horizontal };
412  typedef Homogeneous<Derived, HomogeneousReturnTypeDirection> HomogeneousReturnType;
413  HomogeneousReturnType homogeneous() const;
414  #endif
415 
416  enum {
418  };
419  typedef Block<const Derived,
420  internal::traits<Derived>::ColsAtCompileTime==1 ? SizeMinusOne : 1,
421  internal::traits<Derived>::ColsAtCompileTime==1 ? 1 : SizeMinusOne> ConstStartMinusOne;
423  const ConstStartMinusOne > HNormalizedReturnType;
424 
425  const HNormalizedReturnType hnormalized() const;
426 
428 
429  void makeHouseholderInPlace(Scalar& tau, RealScalar& beta);
430  template<typename EssentialPart>
431  void makeHouseholder(EssentialPart& essential,
432  Scalar& tau, RealScalar& beta) const;
433  template<typename EssentialPart>
434  void applyHouseholderOnTheLeft(const EssentialPart& essential,
435  const Scalar& tau,
436  Scalar* workspace);
437  template<typename EssentialPart>
438  void applyHouseholderOnTheRight(const EssentialPart& essential,
439  const Scalar& tau,
440  Scalar* workspace);
441 
443 
444  template<typename OtherScalar>
445  void applyOnTheLeft(Index p, Index q, const JacobiRotation<OtherScalar>& j);
446  template<typename OtherScalar>
447  void applyOnTheRight(Index p, Index q, const JacobiRotation<OtherScalar>& j);
448 
450 
453  const MatrixFunctionReturnValue<Derived> matrixFunction(StemFunction f) const;
460  const MatrixPowerReturnValue<Derived> pow(const RealScalar& p) const;
461 
462 #ifdef EIGEN2_SUPPORT
463  template<typename ProductDerived, typename Lhs, typename Rhs>
465  EvalBeforeAssigningBit>& other);
466 
467  template<typename ProductDerived, typename Lhs, typename Rhs>
469  EvalBeforeAssigningBit>& other);
470 
473  template<typename OtherDerived>
475  { return lazyAssign(other._expression()); }
476 
477  template<unsigned int Added>
478  const Flagged<Derived, Added, 0> marked() const;
480 
481  inline const Cwise<Derived> cwise() const;
482  inline Cwise<Derived> cwise();
483 
484  VectorBlock<Derived> start(Index size);
485  const VectorBlock<const Derived> start(Index size) const;
486  VectorBlock<Derived> end(Index size);
487  const VectorBlock<const Derived> end(Index size) const;
488  template<int Size> VectorBlock<Derived,Size> start();
489  template<int Size> const VectorBlock<const Derived,Size> start() const;
490  template<int Size> VectorBlock<Derived,Size> end();
491  template<int Size> const VectorBlock<const Derived,Size> end() const;
492 
493  Minor<Derived> minor(Index row, Index col);
494  const Minor<Derived> minor(Index row, Index col) const;
495 #endif
496 
497  protected:
498  MatrixBase() : Base() {}
499 
500  private:
501  explicit MatrixBase(int);
502  MatrixBase(int,int);
503  template<typename OtherDerived> explicit MatrixBase(const MatrixBase<OtherDerived>&);
504  protected:
505  // mixing arrays and matrices is not legal
506  template<typename OtherDerived> Derived& operator+=(const ArrayBase<OtherDerived>& )
507  {EIGEN_STATIC_ASSERT(std::ptrdiff_t(sizeof(typename OtherDerived::Scalar))==-1,YOU_CANNOT_MIX_ARRAYS_AND_MATRICES); return *this;}
508  // mixing arrays and matrices is not legal
509  template<typename OtherDerived> Derived& operator-=(const ArrayBase<OtherDerived>& )
510  {EIGEN_STATIC_ASSERT(std::ptrdiff_t(sizeof(typename OtherDerived::Scalar))==-1,YOU_CANNOT_MIX_ARRAYS_AND_MATRICES); return *this;}
511 };
512 
513 } // end namespace Eigen
514 
515 #endif // EIGEN_MATRIXBASE_H
Expression of the product of two general matrices or vectors.
Generic expression of a matrix where all coefficients are defined by a functor.
Robust Cholesky decomposition of a matrix with pivoting.
Definition: LDLT.h:45
return ReturnType(abs2(x.value()), x.derivatives()*(Scalar(2)*x.value()))
CwiseUnaryOp< internal::scalar_quotient1_op< typename internal::traits< Derived >::Scalar >, const ConstStartMinusOne > HNormalizedReturnType
Definition: MatrixBase.h:423
ColXpr col(Index i)
Definition: DenseBase.h:709
Scalar determinant() const
Definition: Determinant.h:92
PlainObject unitOrthogonal(void) const
Definition: OrthoMethods.h:210
Enforce aligned packet loads and stores regardless of what is requested.
const LDLT< PlainObject > ldlt() const
Definition: LDLT.h:593
Matrix< Scalar, MatrixBase::RowsAtCompileTime, MatrixBase::ColsAtCompileTime > type
Definition: MatrixBase.h:397
static const BasisReturnType UnitW()
bool operator!=(const MatrixBase< OtherDerived > &other) const
Definition: MatrixBase.h:301
Expression of a mathematical vector or matrix as an array object.
Definition: ArrayWrapper.h:36
Householder rank-revealing QR decomposition of a matrix with full pivoting.
internal::traits< Derived >::Index Index
Definition: MatrixBase.h:55
Base::CoeffReturnType CoeffReturnType
Definition: DenseBase.h:98
bool operator==(const MatrixBase< OtherDerived > &other) const
Definition: MatrixBase.h:293
Expression with modified flags.
Definition: Flagged.h:39
void applyOnTheLeft(const EigenBase< OtherDerived > &other)
Definition: EigenBase.h:154
Definition: LU.h:16
Derived & setIdentity()
Pseudo expression providing an operator = assuming no aliasing.
Definition: NoAlias.h:31
Scalar trace() const
Definition: Redux.h:401
internal::stem_function< Scalar >::type StemFunction
Definition: MatrixBase.h:451
Diagonal< Derived, Index > Type
Definition: MatrixBase.h:221
const LazyProductReturnType< Derived, OtherDerived >::Type lazyProduct(const MatrixBase< OtherDerived > &other) const
SelfAdjointView< Derived, UpLo > Type
Definition: MatrixBase.h:254
Proxy for the matrix square root of some matrix (expression).
DiagonalReturnType diagonal()
Definition: Diagonal.h:168
MatrixBase StorageBaseType
Definition: MatrixBase.h:53
internal::add_const< Transpose< const Derived > >::type ConstTransposeReturnType
Definition: DenseBase.h:284
RealScalar operatorNorm() const
Computes the L2 operator norm.
bool isLowerTriangular(const RealScalar &prec=NumTraits< Scalar >::dummy_precision()) const
internal::scalar_product_traits< typename internal::traits< Derived >::Scalar, typename internal::traits< OtherDerived >::Scalar >::ReturnType dot(const MatrixBase< OtherDerived > &other) const
Definition: Dot.h:63
void applyOnTheRight(const EigenBase< OtherDerived > &other)
Definition: EigenBase.h:146
LU decomposition of a matrix with partial pivoting, and related features.
CwiseNullaryOp< internal::scalar_identity_op< Scalar >, Derived > IdentityReturnType
Definition: MatrixBase.h:128
const MatrixLogarithmReturnValue< Derived > log() const
Derived & operator+=(const MatrixBase< OtherDerived > &other)
iterative scaling algorithm to equilibrate rows and column norms in matrices
Definition: matrix.hpp:471
Block< const Derived, internal::traits< Derived >::ColsAtCompileTime==1?SizeMinusOne:1, internal::traits< Derived >::ColsAtCompileTime==1?1:SizeMinusOne > ConstStartMinusOne
Definition: MatrixBase.h:421
Rotation given by a cosine-sine pair.
const ForceAlignedAccess< Derived > forceAlignedAccess() const
Holds information about the various numeric (i.e. scalar) types allowed by Eigen. ...
Definition: NumTraits.h:88
RowXpr row(Index i)
Definition: DenseBase.h:726
void makeHouseholder(EssentialPart &essential, Scalar &tau, RealScalar &beta) const
Definition: Householder.h:65
#define EIGEN_STATIC_ASSERT(CONDITION, MSG)
Definition: StaticAssert.h:111
EigenvaluesReturnType eigenvalues() const
Computes the eigenvalues of a matrix.
Base::ConstTransposeReturnType ConstTransposeReturnType
Definition: MatrixBase.h:86
internal::traits< Derived >::Scalar Scalar
Definition: MatrixBase.h:56
static const BasisReturnType UnitY()
const unsigned int RowMajorBit
Base class for all dense matrices, vectors, and arrays.
Definition: DenseBase.h:41
const MatrixFunctionReturnValue< Derived > sinh() const
const CwiseUnaryOp< std::binder1st< std::equal_to< Scalar > >, const Derived > cwiseEqual(const Scalar &s) const
Definition: MatrixBase.h:64
const SparseView< Derived > sparseView(const Scalar &m_reference=Scalar(0), const typename NumTraits< Scalar >::Real &m_epsilon=NumTraits< Scalar >::dummy_precision()) const
Definition: SparseView.h:91
const TriangularView< const Derived, Mode > Type
Definition: MatrixBase.h:249
const ArrayWrapper< const Derived > array() const
Definition: MatrixBase.h:323
Proxy for the matrix function of some matrix (expression).
Block< Derived, internal::traits< Derived >::RowsAtCompileTime, 1,!IsRowMajor > ColXpr
Definition: DenseBase.h:16
internal::packet_traits< Scalar >::type PacketScalar
Definition: MatrixBase.h:57
Base::ColXpr ColXpr
Definition: MatrixBase.h:88
Matrix< Scalar, 3, 1 > eulerAngles(Index a0, Index a1, Index a2) const
Definition: EulerAngles.h:37
const PartialPivLU< PlainObject > partialPivLu() const
Definition: PartialPivLU.h:477
MatrixBase< Derived > & matrix()
Definition: MatrixBase.h:317
const HNormalizedReturnType hnormalized() const
Definition: Homogeneous.h:158
const MatrixFunctionReturnValue< Derived > cos() const
void applyHouseholderOnTheLeft(const EssentialPart &essential, const Scalar &tau, Scalar *workspace)
Definition: Householder.h:112
void normalize()
Definition: Dot.h:154
const MatrixExponentialReturnValue< Derived > exp() const
Base::RowXpr RowXpr
Definition: MatrixBase.h:87
Expression of a fixed-size or dynamic-size sub-vector.
internal::scalar_product_traits< typename internal::traits< Derived >::Scalar, typename internal::traits< OtherDerived >::Scalar >::ReturnType Scalar
Definition: MatrixBase.h:396
const MatrixBase< Derived > & matrix() const
Definition: MatrixBase.h:318
Derived & lazyAssign(const DenseBase< OtherDerived > &other)
Pseudo expression providing additional coefficient-wise operations.
Definition: Cwise.h:50
const MatrixPowerReturnValue< Derived > pow(const RealScalar &p) const
Definition: MatrixPower.h:504
RealScalar blueNorm() const
Definition: StableNorm.h:171
ArrayWrapper< Derived > array()
Definition: MatrixBase.h:322
Householder rank-revealing QR decomposition of a matrix with column-pivoting.
const Diagonal< const Derived, Index > Type
Definition: MatrixBase.h:222
Standard Cholesky decomposition (LL^T) of a matrix and associated features.
Definition: LLT.h:50
const MatrixFunctionReturnValue< Derived > cosh() const
Block< const CwiseNullaryOp< internal::scalar_identity_op< Scalar >, SquareMatrixType >, internal::traits< Derived >::RowsAtCompileTime, internal::traits< Derived >::ColsAtCompileTime > BasisReturnType
Definition: MatrixBase.h:132
const SelfAdjointView< const Derived, UpLo > Type
Definition: MatrixBase.h:255
DenseBase< Derived > Base
Definition: MatrixBase.h:60
NumTraits< Scalar >::Real RealScalar
Definition: MatrixBase.h:58
const PlainObject normalized() const
Definition: Dot.h:139
Standard SVD decomposition of a matrix and associated features.
Definition: SVD.h:30
internal::traits< Derived >::StorageKind StorageKind
Definition: MatrixBase.h:54
Expression of a minor.
Definition: Minor.h:53
void computeInverseAndDetWithCheck(ResultType &inverse, typename ResultType::Scalar &determinant, bool &invertible, const RealScalar &absDeterminantThreshold=NumTraits< Scalar >::dummy_precision()) const
Definition: Inverse.h:347
Expression of a selfadjoint matrix from a triangular part of a dense matrix.
Block< Derived, 1, internal::traits< Derived >::ColsAtCompileTime, IsRowMajor > RowXpr
Definition: DenseBase.h:19
internal::conditional< NumTraits< Scalar >::IsComplex, CwiseUnaryOp< internal::scalar_conjugate_op< Scalar >, ConstTransposeReturnType >, ConstTransposeReturnType >::type AdjointReturnType
Definition: MatrixBase.h:124
PlainObject cross3(const MatrixBase< OtherDerived > &other) const
Definition: OrthoMethods.h:74
Base class for all 1D and 2D array, and related expressions.
Definition: ArrayBase.h:39
Class to view a vector of integers as a permutation matrix.
static const BasisReturnType UnitZ()
const unsigned int EvalBeforeAssigningBit
void adjointInPlace()
Definition: Transpose.h:321
const MatrixSquareRootReturnValue< Derived > sqrt() const
RealScalar stableNorm() const
Definition: StableNorm.h:140
const ColPivHouseholderQR< PlainObject > colPivHouseholderQr() const
Index diagonalSize() const
Definition: MatrixBase.h:101
bool isIdentity(const RealScalar &prec=NumTraits< Scalar >::dummy_precision()) const
Derived & operator+=(const ArrayBase< OtherDerived > &)
Definition: MatrixBase.h:506
RealScalar squaredNorm() const
Definition: Dot.h:113
static const IdentityReturnType Identity()
Derived & operator-=(const ArrayBase< OtherDerived > &)
Definition: MatrixBase.h:509
Definition: QR.h:17
Derived & operator-=(const MatrixBase< OtherDerived > &other)
const PermutationWrapper< const Derived > asPermutation() const
internal::add_const< Diagonal< const Derived > >::type ConstDiagonalReturnType
Definition: MatrixBase.h:218
SelfAdjointViewReturnType< UpLo >::Type selfadjointView()
void computeInverseWithCheck(ResultType &inverse, bool &invertible, const RealScalar &absDeterminantThreshold=NumTraits< Scalar >::dummy_precision()) const
Definition: Inverse.h:386
bool isUnitary(const RealScalar &prec=NumTraits< Scalar >::dummy_precision()) const
Definition: Dot.h:247
#define EIGEN_SIZE_MAX(a, b)
cross_product_return_type< OtherDerived >::type cross(const MatrixBase< OtherDerived > &other) const
Expression of a fixed-size or dynamic-size block.
Definition: Core/Block.h:102
Matrix< std::complex< RealScalar >, internal::traits< Derived >::ColsAtCompileTime, 1, ColMajor > EigenvaluesReturnType
Definition: MatrixBase.h:126
RealScalar hypotNorm() const
Definition: StableNorm.h:183
LU decomposition of a matrix with complete pivoting, and related features.
EIGEN_STRONG_INLINE EvalReturnType eval() const
Definition: DenseBase.h:362
void applyHouseholderOnTheRight(const EssentialPart &essential, const Scalar &tau, Scalar *workspace)
Definition: Householder.h:149
const FullPivLU< PlainObject > fullPivLu() const
Definition: FullPivLU.h:735
Householder QR decomposition of a matrix.
bool isUpperTriangular(const RealScalar &prec=NumTraits< Scalar >::dummy_precision()) const
Two-sided Jacobi SVD decomposition of a rectangular matrix.
Base class for triangular part in a matrix.
const CwiseBinaryOp< std::not_equal_to< Scalar >, const Derived, const OtherDerived > cwiseNotEqual(const EIGEN_CURRENT_STORAGE_BASE_CLASS< OtherDerived > &other) const
Definition: MatrixBase.h:61
Expression of a diagonal matrix.
const ExpressionType & _expression() const
Definition: Flagged.h:112
Proxy for the matrix exponential of some matrix (expression).
Derived & operator*=(const EigenBase< OtherDerived > &other)
Definition: EigenBase.h:136
NoAlias< Derived, Eigen::MatrixBase > noalias()
Definition: NoAlias.h:127
bool isOrthogonal(const MatrixBase< OtherDerived > &other, const RealScalar &prec=NumTraits< Scalar >::dummy_precision()) const
Definition: Dot.h:228
Expression of a diagonal/subdiagonal/superdiagonal in a matrix.
Definition: Diagonal.h:64
const MatrixFunctionReturnValue< Derived > matrixFunction(StemFunction f) const
const LLT< PlainObject > llt() const
Definition: LLT.h:473
CwiseNullaryOp< internal::scalar_constant_op< Scalar >, Derived > ConstantReturnType
Definition: MatrixBase.h:119
Derived & lazyAssign(const ProductBase< ProductDerived, Lhs, Rhs > &other)
Definition: ProductBase.h:270
void makeHouseholderInPlace(Scalar &tau, RealScalar &beta)
Definition: Householder.h:42
internal::add_const_on_value_type< typename internal::conditional< Enable, ForceAlignedAccess< Derived >, Derived & >::type >::type forceAlignedAccessIf() const
static const BasisReturnType UnitX()
Generic expression where a coefficient-wise unary operator is applied to an expression.
Definition: CwiseUnaryOp.h:59
bool isDiagonal(const RealScalar &prec=NumTraits< Scalar >::dummy_precision()) const
const DiagonalWrapper< const Derived > asDiagonal() const
The matrix class, also used for vectors and row-vectors.
Definition: Matrix.h:127
Diagonal< Derived > DiagonalReturnType
Definition: MatrixBase.h:216
TriangularViewReturnType< Mode >::Type triangularView()
const AdjointReturnType adjoint() const
Definition: Transpose.h:237
RealScalar lpNorm() const
TriangularView< Derived, Mode > Type
Definition: MatrixBase.h:248
Derived & operator=(const MatrixBase &other)
Definition: Assign.h:555
const internal::inverse_impl< Derived > inverse() const
Definition: Inverse.h:320
Proxy for the matrix power of some matrix (expression).
static const BasisReturnType Unit(Index size, Index i)
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:48
Base::CoeffReturnType CoeffReturnType
Definition: MatrixBase.h:85
Proxy for the matrix logarithm of some matrix (expression).
RealScalar norm() const
Definition: Dot.h:125
const HouseholderQR< PlainObject > householderQr() const
JacobiSVD< PlainObject > jacobiSvd(unsigned int computationOptions=0) const
Definition: JacobiSVD.h:872
const MatrixFunctionReturnValue< Derived > sin() const
const FullPivHouseholderQR< PlainObject > fullPivHouseholderQr() const
Matrix< Scalar, EIGEN_SIZE_MAX(RowsAtCompileTime, ColsAtCompileTime), EIGEN_SIZE_MAX(RowsAtCompileTime, ColsAtCompileTime)> SquareMatrixType
Definition: MatrixBase.h:96
Expression of one (or a set of) homogeneous vector(s)
Definition: Homogeneous.h:61


acado
Author(s): Milan Vukov, Rien Quirynen
autogenerated on Mon Jun 10 2019 12:34:52