source: pacpussensors/trunk/Vislab/lib3dv/eigen/Eigen/src/Core/MatrixBase.h@ 136

Last change on this file since 136 was 136, checked in by ldecherf, 7 years ago

Doc

File size: 23.7 KB
Line 
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
14namespace Eigen {
15
16/** \class MatrixBase
17 * \ingroup Core_Module
18 *
19 * \brief Base class for all dense matrices, vectors, and expressions
20 *
21 * This class is the base that is inherited by all matrix, vector, and related expression
22 * types. Most of the Eigen API is contained in this class, and its base classes. Other important
23 * classes for the Eigen API are Matrix, and VectorwiseOp.
24 *
25 * Note that some methods are defined in other modules such as the \ref LU_Module LU module
26 * for all functions related to matrix inversions.
27 *
28 * \tparam Derived is the derived type, e.g. a matrix type, or an expression, etc.
29 *
30 * When writing a function taking Eigen objects as argument, if you want your function
31 * to take as argument any matrix, vector, or expression, just let it take a
32 * MatrixBase argument. As an example, here is a function printFirstRow which, given
33 * a matrix, vector, or expression \a x, prints the first row of \a x.
34 *
35 * \code
36 template<typename Derived>
37 void printFirstRow(const Eigen::MatrixBase<Derived>& x)
38 {
39 cout << x.row(0) << endl;
40 }
41 * \endcode
42 *
43 * This class can be extended with the help of the plugin mechanism described on the page
44 * \ref TopicCustomizingEigen by defining the preprocessor symbol \c EIGEN_MATRIXBASE_PLUGIN.
45 *
46 * \sa \ref TopicClassHierarchy
47 */
48template<typename Derived> class MatrixBase
49 : public DenseBase<Derived>
50{
51 public:
52#ifndef EIGEN_PARSED_BY_DOXYGEN
53 typedef MatrixBase StorageBaseType;
54 typedef typename internal::traits<Derived>::StorageKind StorageKind;
55 typedef typename internal::traits<Derived>::Index Index;
56 typedef typename internal::traits<Derived>::Scalar Scalar;
57 typedef typename internal::packet_traits<Scalar>::type PacketScalar;
58 typedef typename NumTraits<Scalar>::Real RealScalar;
59
60 typedef DenseBase<Derived> Base;
61 using Base::RowsAtCompileTime;
62 using Base::ColsAtCompileTime;
63 using Base::SizeAtCompileTime;
64 using Base::MaxRowsAtCompileTime;
65 using Base::MaxColsAtCompileTime;
66 using Base::MaxSizeAtCompileTime;
67 using Base::IsVectorAtCompileTime;
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
85 typedef typename Base::CoeffReturnType CoeffReturnType;
86 typedef typename Base::ConstTransposeReturnType ConstTransposeReturnType;
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 /** type of the equivalent square matrix */
95 typedef Matrix<Scalar,EIGEN_SIZE_MAX(RowsAtCompileTime,ColsAtCompileTime),
96 EIGEN_SIZE_MAX(RowsAtCompileTime,ColsAtCompileTime)> SquareMatrixType;
97#endif // not EIGEN_PARSED_BY_DOXYGEN
98
99 /** \returns the size of the main diagonal, which is min(rows(),cols()).
100 * \sa rows(), cols(), SizeAtCompileTime. */
101 inline Index diagonalSize() const { return (std::min)(rows(),cols()); }
102
103 /** \brief The plain matrix type corresponding to this expression.
104 *
105 * This is not necessarily exactly the return type of eval(). In the case of plain matrices,
106 * the return type of eval() is a const reference to a matrix, not a matrix! It is however guaranteed
107 * that the return type of eval() is either PlainObject or const PlainObject&.
108 */
109 typedef Matrix<typename internal::traits<Derived>::Scalar,
110 internal::traits<Derived>::RowsAtCompileTime,
111 internal::traits<Derived>::ColsAtCompileTime,
112 AutoAlign | (internal::traits<Derived>::Flags&RowMajorBit ? RowMajor : ColMajor),
113 internal::traits<Derived>::MaxRowsAtCompileTime,
114 internal::traits<Derived>::MaxColsAtCompileTime
115 > PlainObject;
116
117#ifndef EIGEN_PARSED_BY_DOXYGEN
118 /** \internal Represents a matrix with all coefficients equal to one another*/
119 typedef CwiseNullaryOp<internal::scalar_constant_op<Scalar>,Derived> ConstantReturnType;
120 /** \internal the return type of MatrixBase::adjoint() */
121 typedef typename internal::conditional<NumTraits<Scalar>::IsComplex,
122 CwiseUnaryOp<internal::scalar_conjugate_op<Scalar>, ConstTransposeReturnType>,
123 ConstTransposeReturnType
124 >::type AdjointReturnType;
125 /** \internal Return type of eigenvalues() */
126 typedef Matrix<std::complex<RealScalar>, internal::traits<Derived>::ColsAtCompileTime, 1, ColMajor> EigenvaluesReturnType;
127 /** \internal the return type of identity */
128 typedef CwiseNullaryOp<internal::scalar_identity_op<Scalar>,Derived> IdentityReturnType;
129 /** \internal the return type of unit vectors */
130 typedef Block<const CwiseNullaryOp<internal::scalar_identity_op<Scalar>, SquareMatrixType>,
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
145 /** Special case of the template operator=, in order to prevent the compiler
146 * from generating a default operator= (issue hit with g++ 4.1)
147 */
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);
158
159 template<typename OtherDerived>
160 Derived& operator=(const ReturnByValue<OtherDerived>& other);
161
162 template<typename ProductDerived, typename Lhs, typename Rhs>
163 Derived& lazyAssign(const ProductBase<ProductDerived, Lhs,Rhs>& other);
164
165 template<typename MatrixPower, typename Lhs, typename Rhs>
166 Derived& lazyAssign(const MatrixPowerProduct<MatrixPower, Lhs,Rhs>& other);
167
168 template<typename OtherDerived>
169 Derived& operator+=(const MatrixBase<OtherDerived>& other);
170 template<typename OtherDerived>
171 Derived& operator-=(const MatrixBase<OtherDerived>& other);
172
173 template<typename OtherDerived>
174 const typename ProductReturnType<Derived,OtherDerived>::Type
175 operator*(const MatrixBase<OtherDerived> &other) const;
176
177 template<typename OtherDerived>
178 const typename LazyProductReturnType<Derived,OtherDerived>::Type
179 lazyProduct(const MatrixBase<OtherDerived> &other) const;
180
181 template<typename OtherDerived>
182 Derived& operator*=(const EigenBase<OtherDerived>& other);
183
184 template<typename OtherDerived>
185 void applyOnTheLeft(const EigenBase<OtherDerived>& other);
186
187 template<typename OtherDerived>
188 void applyOnTheRight(const EigenBase<OtherDerived>& other);
189
190 template<typename DiagonalDerived>
191 const DiagonalProduct<Derived, DiagonalDerived, OnTheRight>
192 operator*(const DiagonalBase<DiagonalDerived> &diagonal) const;
193
194 template<typename OtherDerived>
195 typename internal::scalar_product_traits<typename internal::traits<Derived>::Scalar,typename internal::traits<OtherDerived>::Scalar>::ReturnType
196 dot(const MatrixBase<OtherDerived>& other) const;
197
198 #ifdef EIGEN2_SUPPORT
199 template<typename OtherDerived>
200 Scalar eigen2_dot(const MatrixBase<OtherDerived>& other) const;
201 #endif
202
203 RealScalar squaredNorm() const;
204 RealScalar norm() const;
205 RealScalar stableNorm() const;
206 RealScalar blueNorm() const;
207 RealScalar hypotNorm() const;
208 const PlainObject normalized() const;
209 void normalize();
210
211 const AdjointReturnType adjoint() const;
212 void adjointInPlace();
213
214 typedef Diagonal<Derived> DiagonalReturnType;
215 DiagonalReturnType diagonal();
216 typedef typename internal::add_const<Diagonal<const Derived> >::type ConstDiagonalReturnType;
217 ConstDiagonalReturnType diagonal() const;
218
219 template<int Index> struct DiagonalIndexReturnType { typedef Diagonal<Derived,Index> Type; };
220 template<int Index> struct ConstDiagonalIndexReturnType { typedef const Diagonal<const Derived,Index> Type; };
221
222 template<int Index> typename DiagonalIndexReturnType<Index>::Type diagonal();
223 template<int Index> typename ConstDiagonalIndexReturnType<Index>::Type diagonal() const;
224
225 typedef Diagonal<Derived,DynamicIndex> DiagonalDynamicIndexReturnType;
226 typedef typename internal::add_const<Diagonal<const Derived,DynamicIndex> >::type ConstDiagonalDynamicIndexReturnType;
227
228 DiagonalDynamicIndexReturnType diagonal(Index index);
229 ConstDiagonalDynamicIndexReturnType diagonal(Index index) const;
230
231 #ifdef EIGEN2_SUPPORT
232 template<unsigned int Mode> typename internal::eigen2_part_return_type<Derived, Mode>::type part();
233 template<unsigned int Mode> const typename internal::eigen2_part_return_type<Derived, Mode>::type part() const;
234
235 // huuuge hack. make Eigen2's matrix.part<Diagonal>() work in eigen3. Problem: Diagonal is now a class template instead
236 // of an integer constant. Solution: overload the part() method template wrt template parameters list.
237 template<template<typename T, int N> class U>
238 const DiagonalWrapper<ConstDiagonalReturnType> part() const
239 { return diagonal().asDiagonal(); }
240 #endif // EIGEN2_SUPPORT
241
242 template<unsigned int Mode> struct TriangularViewReturnType { typedef TriangularView<Derived, Mode> Type; };
243 template<unsigned int Mode> struct ConstTriangularViewReturnType { typedef const TriangularView<const Derived, Mode> Type; };
244
245 template<unsigned int Mode> typename TriangularViewReturnType<Mode>::Type triangularView();
246 template<unsigned int Mode> typename ConstTriangularViewReturnType<Mode>::Type triangularView() const;
247
248 template<unsigned int UpLo> struct SelfAdjointViewReturnType { typedef SelfAdjointView<Derived, UpLo> Type; };
249 template<unsigned int UpLo> struct ConstSelfAdjointViewReturnType { typedef const SelfAdjointView<const Derived, UpLo> Type; };
250
251 template<unsigned int UpLo> typename SelfAdjointViewReturnType<UpLo>::Type selfadjointView();
252 template<unsigned int UpLo> typename ConstSelfAdjointViewReturnType<UpLo>::Type selfadjointView() const;
253
254 const SparseView<Derived> sparseView(const Scalar& m_reference = Scalar(0),
255 const typename NumTraits<Scalar>::Real& m_epsilon = NumTraits<Scalar>::dummy_precision()) const;
256 static const IdentityReturnType Identity();
257 static const IdentityReturnType Identity(Index rows, Index cols);
258 static const BasisReturnType Unit(Index size, Index i);
259 static const BasisReturnType Unit(Index i);
260 static const BasisReturnType UnitX();
261 static const BasisReturnType UnitY();
262 static const BasisReturnType UnitZ();
263 static const BasisReturnType UnitW();
264
265 const DiagonalWrapper<const Derived> asDiagonal() const;
266 const PermutationWrapper<const Derived> asPermutation() const;
267
268 Derived& setIdentity();
269 Derived& setIdentity(Index rows, Index cols);
270
271 bool isIdentity(const RealScalar& prec = NumTraits<Scalar>::dummy_precision()) const;
272 bool isDiagonal(const RealScalar& prec = NumTraits<Scalar>::dummy_precision()) const;
273
274 bool isUpperTriangular(const RealScalar& prec = NumTraits<Scalar>::dummy_precision()) const;
275 bool isLowerTriangular(const RealScalar& prec = NumTraits<Scalar>::dummy_precision()) const;
276
277 template<typename OtherDerived>
278 bool isOrthogonal(const MatrixBase<OtherDerived>& other,
279 const RealScalar& prec = NumTraits<Scalar>::dummy_precision()) const;
280 bool isUnitary(const RealScalar& prec = NumTraits<Scalar>::dummy_precision()) const;
281
282 /** \returns true if each coefficients of \c *this and \a other are all exactly equal.
283 * \warning When using floating point scalar values you probably should rather use a
284 * fuzzy comparison such as isApprox()
285 * \sa isApprox(), operator!= */
286 template<typename OtherDerived>
287 inline bool operator==(const MatrixBase<OtherDerived>& other) const
288 { return cwiseEqual(other).all(); }
289
290 /** \returns true if at least one pair of coefficients of \c *this and \a other are not exactly equal to each other.
291 * \warning When using floating point scalar values you probably should rather use a
292 * fuzzy comparison such as isApprox()
293 * \sa isApprox(), operator== */
294 template<typename OtherDerived>
295 inline bool operator!=(const MatrixBase<OtherDerived>& other) const
296 { return cwiseNotEqual(other).any(); }
297
298 NoAlias<Derived,Eigen::MatrixBase > noalias();
299
300 inline const ForceAlignedAccess<Derived> forceAlignedAccess() const;
301 inline ForceAlignedAccess<Derived> forceAlignedAccess();
302 template<bool Enable> inline typename internal::add_const_on_value_type<typename internal::conditional<Enable,ForceAlignedAccess<Derived>,Derived&>::type>::type forceAlignedAccessIf() const;
303 template<bool Enable> inline typename internal::conditional<Enable,ForceAlignedAccess<Derived>,Derived&>::type forceAlignedAccessIf();
304
305 Scalar trace() const;
306
307/////////// Array module ///////////
308
309 template<int p> RealScalar lpNorm() const;
310
311 MatrixBase<Derived>& matrix() { return *this; }
312 const MatrixBase<Derived>& matrix() const { return *this; }
313
314 /** \returns an \link Eigen::ArrayBase Array \endlink expression of this matrix
315 * \sa ArrayBase::matrix() */
316 ArrayWrapper<Derived> array() { return derived(); }
317 const ArrayWrapper<const Derived> array() const { return derived(); }
318
319/////////// LU module ///////////
320
321 const FullPivLU<PlainObject> fullPivLu() const;
322 const PartialPivLU<PlainObject> partialPivLu() const;
323
324 #if EIGEN2_SUPPORT_STAGE < STAGE20_RESOLVE_API_CONFLICTS
325 const LU<PlainObject> lu() const;
326 #endif
327
328 #ifdef EIGEN2_SUPPORT
329 const LU<PlainObject> eigen2_lu() const;
330 #endif
331
332 #if EIGEN2_SUPPORT_STAGE > STAGE20_RESOLVE_API_CONFLICTS
333 const PartialPivLU<PlainObject> lu() const;
334 #endif
335
336 #ifdef EIGEN2_SUPPORT
337 template<typename ResultType>
338 void computeInverse(MatrixBase<ResultType> *result) const {
339 *result = this->inverse();
340 }
341 #endif
342
343 const internal::inverse_impl<Derived> inverse() const;
344 template<typename ResultType>
345 void computeInverseAndDetWithCheck(
346 ResultType& inverse,
347 typename ResultType::Scalar& determinant,
348 bool& invertible,
349 const RealScalar& absDeterminantThreshold = NumTraits<Scalar>::dummy_precision()
350 ) const;
351 template<typename ResultType>
352 void computeInverseWithCheck(
353 ResultType& inverse,
354 bool& invertible,
355 const RealScalar& absDeterminantThreshold = NumTraits<Scalar>::dummy_precision()
356 ) const;
357 Scalar determinant() const;
358
359/////////// Cholesky module ///////////
360
361 const LLT<PlainObject> llt() const;
362 const LDLT<PlainObject> ldlt() const;
363
364/////////// QR module ///////////
365
366 const HouseholderQR<PlainObject> householderQr() const;
367 const ColPivHouseholderQR<PlainObject> colPivHouseholderQr() const;
368 const FullPivHouseholderQR<PlainObject> fullPivHouseholderQr() const;
369
370 #ifdef EIGEN2_SUPPORT
371 const QR<PlainObject> qr() const;
372 #endif
373
374 EigenvaluesReturnType eigenvalues() const;
375 RealScalar operatorNorm() const;
376
377/////////// SVD module ///////////
378
379 JacobiSVD<PlainObject> jacobiSvd(unsigned int computationOptions = 0) const;
380
381 #ifdef EIGEN2_SUPPORT
382 SVD<PlainObject> svd() const;
383 #endif
384
385/////////// Geometry module ///////////
386
387 #ifndef EIGEN_PARSED_BY_DOXYGEN
388 /// \internal helper struct to form the return type of the cross product
389 template<typename OtherDerived> struct cross_product_return_type {
390 typedef typename internal::scalar_product_traits<typename internal::traits<Derived>::Scalar,typename internal::traits<OtherDerived>::Scalar>::ReturnType Scalar;
391 typedef Matrix<Scalar,MatrixBase::RowsAtCompileTime,MatrixBase::ColsAtCompileTime> type;
392 };
393 #endif // EIGEN_PARSED_BY_DOXYGEN
394 template<typename OtherDerived>
395 typename cross_product_return_type<OtherDerived>::type
396 cross(const MatrixBase<OtherDerived>& other) const;
397 template<typename OtherDerived>
398 PlainObject cross3(const MatrixBase<OtherDerived>& other) const;
399 PlainObject unitOrthogonal(void) const;
400 Matrix<Scalar,3,1> eulerAngles(Index a0, Index a1, Index a2) const;
401
402 #if EIGEN2_SUPPORT_STAGE > STAGE20_RESOLVE_API_CONFLICTS
403 ScalarMultipleReturnType operator*(const UniformScaling<Scalar>& s) const;
404 // put this as separate enum value to work around possible GCC 4.3 bug (?)
405 enum { HomogeneousReturnTypeDirection = ColsAtCompileTime==1?Vertical:Horizontal };
406 typedef Homogeneous<Derived, HomogeneousReturnTypeDirection> HomogeneousReturnType;
407 HomogeneousReturnType homogeneous() const;
408 #endif
409
410 enum {
411 SizeMinusOne = SizeAtCompileTime==Dynamic ? Dynamic : SizeAtCompileTime-1
412 };
413 typedef Block<const Derived,
414 internal::traits<Derived>::ColsAtCompileTime==1 ? SizeMinusOne : 1,
415 internal::traits<Derived>::ColsAtCompileTime==1 ? 1 : SizeMinusOne> ConstStartMinusOne;
416 typedef CwiseUnaryOp<internal::scalar_quotient1_op<typename internal::traits<Derived>::Scalar>,
417 const ConstStartMinusOne > HNormalizedReturnType;
418
419 const HNormalizedReturnType hnormalized() const;
420
421////////// Householder module ///////////
422
423 void makeHouseholderInPlace(Scalar& tau, RealScalar& beta);
424 template<typename EssentialPart>
425 void makeHouseholder(EssentialPart& essential,
426 Scalar& tau, RealScalar& beta) const;
427 template<typename EssentialPart>
428 void applyHouseholderOnTheLeft(const EssentialPart& essential,
429 const Scalar& tau,
430 Scalar* workspace);
431 template<typename EssentialPart>
432 void applyHouseholderOnTheRight(const EssentialPart& essential,
433 const Scalar& tau,
434 Scalar* workspace);
435
436///////// Jacobi module /////////
437
438 template<typename OtherScalar>
439 void applyOnTheLeft(Index p, Index q, const JacobiRotation<OtherScalar>& j);
440 template<typename OtherScalar>
441 void applyOnTheRight(Index p, Index q, const JacobiRotation<OtherScalar>& j);
442
443///////// SparseCore module /////////
444
445 template<typename OtherDerived>
446 EIGEN_STRONG_INLINE const typename SparseMatrixBase<OtherDerived>::template CwiseProductDenseReturnType<Derived>::Type
447 cwiseProduct(const SparseMatrixBase<OtherDerived> &other) const
448 {
449 return other.cwiseProduct(derived());
450 }
451
452///////// MatrixFunctions module /////////
453
454 typedef typename internal::stem_function<Scalar>::type StemFunction;
455 const MatrixExponentialReturnValue<Derived> exp() const;
456 const MatrixFunctionReturnValue<Derived> matrixFunction(StemFunction f) const;
457 const MatrixFunctionReturnValue<Derived> cosh() const;
458 const MatrixFunctionReturnValue<Derived> sinh() const;
459 const MatrixFunctionReturnValue<Derived> cos() const;
460 const MatrixFunctionReturnValue<Derived> sin() const;
461 const MatrixSquareRootReturnValue<Derived> sqrt() const;
462 const MatrixLogarithmReturnValue<Derived> log() const;
463 const MatrixPowerReturnValue<Derived> pow(const RealScalar& p) const;
464
465#ifdef EIGEN2_SUPPORT
466 template<typename ProductDerived, typename Lhs, typename Rhs>
467 Derived& operator+=(const Flagged<ProductBase<ProductDerived, Lhs,Rhs>, 0,
468 EvalBeforeAssigningBit>& other);
469
470 template<typename ProductDerived, typename Lhs, typename Rhs>
471 Derived& operator-=(const Flagged<ProductBase<ProductDerived, Lhs,Rhs>, 0,
472 EvalBeforeAssigningBit>& other);
473
474 /** \deprecated because .lazy() is deprecated
475 * Overloaded for cache friendly product evaluation */
476 template<typename OtherDerived>
477 Derived& lazyAssign(const Flagged<OtherDerived, 0, EvalBeforeAssigningBit>& other)
478 { return lazyAssign(other._expression()); }
479
480 template<unsigned int Added>
481 const Flagged<Derived, Added, 0> marked() const;
482 const Flagged<Derived, 0, EvalBeforeAssigningBit> lazy() const;
483
484 inline const Cwise<Derived> cwise() const;
485 inline Cwise<Derived> cwise();
486
487 VectorBlock<Derived> start(Index size);
488 const VectorBlock<const Derived> start(Index size) const;
489 VectorBlock<Derived> end(Index size);
490 const VectorBlock<const Derived> end(Index size) const;
491 template<int Size> VectorBlock<Derived,Size> start();
492 template<int Size> const VectorBlock<const Derived,Size> start() const;
493 template<int Size> VectorBlock<Derived,Size> end();
494 template<int Size> const VectorBlock<const Derived,Size> end() const;
495
496 Minor<Derived> minor(Index row, Index col);
497 const Minor<Derived> minor(Index row, Index col) const;
498#endif
499
500 protected:
501 MatrixBase() : Base() {}
502
503 private:
504 explicit MatrixBase(int);
505 MatrixBase(int,int);
506 template<typename OtherDerived> explicit MatrixBase(const MatrixBase<OtherDerived>&);
507 protected:
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 // mixing arrays and matrices is not legal
512 template<typename OtherDerived> Derived& operator-=(const ArrayBase<OtherDerived>& )
513 {EIGEN_STATIC_ASSERT(std::ptrdiff_t(sizeof(typename OtherDerived::Scalar))==-1,YOU_CANNOT_MIX_ARRAYS_AND_MATRICES); return *this;}
514};
515
516
517/***************************************************************************
518* Implementation of matrix base methods
519***************************************************************************/
520
521/** replaces \c *this by \c *this * \a other.
522 *
523 * \returns a reference to \c *this
524 *
525 * Example: \include MatrixBase_applyOnTheRight.cpp
526 * Output: \verbinclude MatrixBase_applyOnTheRight.out
527 */
528template<typename Derived>
529template<typename OtherDerived>
530inline Derived&
531MatrixBase<Derived>::operator*=(const EigenBase<OtherDerived> &other)
532{
533 other.derived().applyThisOnTheRight(derived());
534 return derived();
535}
536
537/** replaces \c *this by \c *this * \a other. It is equivalent to MatrixBase::operator*=().
538 *
539 * Example: \include MatrixBase_applyOnTheRight.cpp
540 * Output: \verbinclude MatrixBase_applyOnTheRight.out
541 */
542template<typename Derived>
543template<typename OtherDerived>
544inline void MatrixBase<Derived>::applyOnTheRight(const EigenBase<OtherDerived> &other)
545{
546 other.derived().applyThisOnTheRight(derived());
547}
548
549/** replaces \c *this by \a other * \c *this.
550 *
551 * Example: \include MatrixBase_applyOnTheLeft.cpp
552 * Output: \verbinclude MatrixBase_applyOnTheLeft.out
553 */
554template<typename Derived>
555template<typename OtherDerived>
556inline void MatrixBase<Derived>::applyOnTheLeft(const EigenBase<OtherDerived> &other)
557{
558 other.derived().applyThisOnTheLeft(derived());
559}
560
561} // end namespace Eigen
562
563#endif // EIGEN_MATRIXBASE_H
Note: See TracBrowser for help on using the repository browser.