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
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 */
48 template<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. */
diagonalSize()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 #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>
167 Derived& lazyAssign(const MatrixPowerProduct<MatrixPower, Lhs,Rhs>& other);
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>
176 const typename ProductReturnType<Derived,OtherDerived>::Type
177 operator*(const MatrixBase<OtherDerived> &other) const;
178
179 template<typename OtherDerived>
180 const typename LazyProductReturnType<Derived,OtherDerived>::Type
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>
193 const DiagonalProduct<Derived, DiagonalDerived, OnTheRight>
194 operator*(const DiagonalBase<DiagonalDerived> &diagonal) const;
195
196 template<typename OtherDerived>
197 typename internal::scalar_product_traits<typename internal::traits<Derived>::Scalar,typename internal::traits<OtherDerived>::Scalar>::ReturnType
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
216 typedef Diagonal<Derived> DiagonalReturnType;
217 DiagonalReturnType diagonal();
218 typedef typename internal::add_const<Diagonal<const Derived> >::type ConstDiagonalReturnType;
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 typedef Diagonal<Derived,DynamicIndex> DiagonalDynamicIndexReturnType;
228 typedef typename internal::add_const<Diagonal<const Derived,DynamicIndex> >::type ConstDiagonalDynamicIndexReturnType;
229
230 DiagonalDynamicIndexReturnType diagonal(Index index);
231 ConstDiagonalDynamicIndexReturnType diagonal(Index index) const;
232
233 #ifdef EIGEN2_SUPPORT
234 template<unsigned int Mode> typename internal::eigen2_part_return_type<Derived, Mode>::type part();
235 template<unsigned int Mode> const typename internal::eigen2_part_return_type<Derived, Mode>::type part() const;
236
237 // huuuge hack. make Eigen2's matrix.part<Diagonal>() work in eigen3. Problem: Diagonal is now a class template instead
238 // of an integer constant. Solution: overload the part() method template wrt template parameters list.
239 template<template<typename T, int N> class U>
part()240 const DiagonalWrapper<ConstDiagonalReturnType> part() const
241 { return diagonal().asDiagonal(); }
242 #endif // EIGEN2_SUPPORT
243
244 template<unsigned int Mode> struct TriangularViewReturnType { typedef TriangularView<Derived, Mode> Type; };
245 template<unsigned int Mode> struct ConstTriangularViewReturnType { typedef const TriangularView<const Derived, Mode> Type; };
246
247 template<unsigned int Mode> typename TriangularViewReturnType<Mode>::Type triangularView();
248 template<unsigned int Mode> typename ConstTriangularViewReturnType<Mode>::Type triangularView() const;
249
250 template<unsigned int UpLo> struct SelfAdjointViewReturnType { typedef SelfAdjointView<Derived, UpLo> Type; };
251 template<unsigned int UpLo> struct ConstSelfAdjointViewReturnType { typedef const SelfAdjointView<const Derived, UpLo> Type; };
252
253 template<unsigned int UpLo> typename SelfAdjointViewReturnType<UpLo>::Type selfadjointView();
254 template<unsigned int UpLo> typename ConstSelfAdjointViewReturnType<UpLo>::Type selfadjointView() const;
255
256 const SparseView<Derived> sparseView(const Scalar& m_reference = Scalar(0),
257 const typename NumTraits<Scalar>::Real& m_epsilon = NumTraits<Scalar>::dummy_precision()) const;
258 static const IdentityReturnType Identity();
259 static const IdentityReturnType Identity(Index rows, Index cols);
260 static const BasisReturnType Unit(Index size, Index i);
261 static const BasisReturnType Unit(Index i);
262 static const BasisReturnType UnitX();
263 static const BasisReturnType UnitY();
264 static const BasisReturnType UnitZ();
265 static const BasisReturnType UnitW();
266
267 const DiagonalWrapper<const Derived> asDiagonal() const;
268 const PermutationWrapper<const Derived> asPermutation() const;
269
270 Derived& setIdentity();
271 Derived& setIdentity(Index rows, Index cols);
272
273 bool isIdentity(const RealScalar& prec = NumTraits<Scalar>::dummy_precision()) const;
274 bool isDiagonal(const RealScalar& prec = NumTraits<Scalar>::dummy_precision()) const;
275
276 bool isUpperTriangular(const RealScalar& prec = NumTraits<Scalar>::dummy_precision()) const;
277 bool isLowerTriangular(const RealScalar& prec = NumTraits<Scalar>::dummy_precision()) const;
278
279 template<typename OtherDerived>
280 bool isOrthogonal(const MatrixBase<OtherDerived>& other,
281 const RealScalar& prec = NumTraits<Scalar>::dummy_precision()) const;
282 bool isUnitary(const RealScalar& prec = NumTraits<Scalar>::dummy_precision()) const;
283
284 /** \returns true if each coefficients of \c *this and \a other are all exactly equal.
285 * \warning When using floating point scalar values you probably should rather use a
286 * fuzzy comparison such as isApprox()
287 * \sa isApprox(), operator!= */
288 template<typename OtherDerived>
289 inline bool operator==(const MatrixBase<OtherDerived>& other) const
290 { return cwiseEqual(other).all(); }
291
292 /** \returns true if at least one pair of coefficients of \c *this and \a other are not exactly equal to each other.
293 * \warning When using floating point scalar values you probably should rather use a
294 * fuzzy comparison such as isApprox()
295 * \sa isApprox(), operator== */
296 template<typename OtherDerived>
297 inline bool operator!=(const MatrixBase<OtherDerived>& other) const
298 { return cwiseNotEqual(other).any(); }
299
300 NoAlias<Derived,Eigen::MatrixBase > noalias();
301
302 inline const ForceAlignedAccess<Derived> forceAlignedAccess() const;
303 inline ForceAlignedAccess<Derived> forceAlignedAccess();
304 template<bool Enable> inline typename internal::add_const_on_value_type<typename internal::conditional<Enable,ForceAlignedAccess<Derived>,Derived&>::type>::type forceAlignedAccessIf() const;
305 template<bool Enable> inline typename internal::conditional<Enable,ForceAlignedAccess<Derived>,Derived&>::type forceAlignedAccessIf();
306
307 Scalar trace() const;
308
309 /////////// Array module ///////////
310
311 template<int p> RealScalar lpNorm() const;
312
matrix()313 MatrixBase<Derived>& matrix() { return *this; }
matrix()314 const MatrixBase<Derived>& matrix() const { return *this; }
315
316 /** \returns an \link Eigen::ArrayBase Array \endlink expression of this matrix
317 * \sa ArrayBase::matrix() */
array()318 ArrayWrapper<Derived> array() { return derived(); }
array()319 const ArrayWrapper<const Derived> array() const { return derived(); }
320
321 /////////// LU module ///////////
322
323 const FullPivLU<PlainObject> fullPivLu() const;
324 const PartialPivLU<PlainObject> partialPivLu() const;
325
326 #if EIGEN2_SUPPORT_STAGE < STAGE20_RESOLVE_API_CONFLICTS
327 const LU<PlainObject> lu() const;
328 #endif
329
330 #ifdef EIGEN2_SUPPORT
331 const LU<PlainObject> eigen2_lu() const;
332 #endif
333
334 #if EIGEN2_SUPPORT_STAGE > STAGE20_RESOLVE_API_CONFLICTS
335 const PartialPivLU<PlainObject> lu() const;
336 #endif
337
338 #ifdef EIGEN2_SUPPORT
339 template<typename ResultType>
computeInverse(MatrixBase<ResultType> * result)340 void computeInverse(MatrixBase<ResultType> *result) const {
341 *result = this->inverse();
342 }
343 #endif
344
345 const internal::inverse_impl<Derived> inverse() const;
346 template<typename ResultType>
347 void computeInverseAndDetWithCheck(
348 ResultType& inverse,
349 typename ResultType::Scalar& determinant,
350 bool& invertible,
351 const RealScalar& absDeterminantThreshold = NumTraits<Scalar>::dummy_precision()
352 ) const;
353 template<typename ResultType>
354 void computeInverseWithCheck(
355 ResultType& inverse,
356 bool& invertible,
357 const RealScalar& absDeterminantThreshold = NumTraits<Scalar>::dummy_precision()
358 ) const;
359 Scalar determinant() const;
360
361 /////////// Cholesky module ///////////
362
363 const LLT<PlainObject> llt() const;
364 const LDLT<PlainObject> ldlt() const;
365
366 /////////// QR module ///////////
367
368 const HouseholderQR<PlainObject> householderQr() const;
369 const ColPivHouseholderQR<PlainObject> colPivHouseholderQr() const;
370 const FullPivHouseholderQR<PlainObject> fullPivHouseholderQr() const;
371
372 #ifdef EIGEN2_SUPPORT
373 const QR<PlainObject> qr() const;
374 #endif
375
376 EigenvaluesReturnType eigenvalues() const;
377 RealScalar operatorNorm() const;
378
379 /////////// SVD module ///////////
380
381 JacobiSVD<PlainObject> jacobiSvd(unsigned int computationOptions = 0) const;
382
383 #ifdef EIGEN2_SUPPORT
384 SVD<PlainObject> svd() const;
385 #endif
386
387 /////////// Geometry module ///////////
388
389 #ifndef EIGEN_PARSED_BY_DOXYGEN
390 /// \internal helper struct to form the return type of the cross product
391 template<typename OtherDerived> struct cross_product_return_type {
392 typedef typename internal::scalar_product_traits<typename internal::traits<Derived>::Scalar,typename internal::traits<OtherDerived>::Scalar>::ReturnType Scalar;
393 typedef Matrix<Scalar,MatrixBase::RowsAtCompileTime,MatrixBase::ColsAtCompileTime> type;
394 };
395 #endif // EIGEN_PARSED_BY_DOXYGEN
396 template<typename OtherDerived>
397 typename cross_product_return_type<OtherDerived>::type
398 cross(const MatrixBase<OtherDerived>& other) const;
399 template<typename OtherDerived>
400 PlainObject cross3(const MatrixBase<OtherDerived>& other) const;
401 PlainObject unitOrthogonal(void) const;
402 Matrix<Scalar,3,1> eulerAngles(Index a0, Index a1, Index a2) const;
403
404 #if EIGEN2_SUPPORT_STAGE > STAGE20_RESOLVE_API_CONFLICTS
405 ScalarMultipleReturnType operator*(const UniformScaling<Scalar>& s) const;
406 // put this as separate enum value to work around possible GCC 4.3 bug (?)
407 enum { HomogeneousReturnTypeDirection = ColsAtCompileTime==1?Vertical:Horizontal };
408 typedef Homogeneous<Derived, HomogeneousReturnTypeDirection> HomogeneousReturnType;
409 HomogeneousReturnType homogeneous() const;
410 #endif
411
412 enum {
413 SizeMinusOne = SizeAtCompileTime==Dynamic ? Dynamic : SizeAtCompileTime-1
414 };
415 typedef Block<const Derived,
416 internal::traits<Derived>::ColsAtCompileTime==1 ? SizeMinusOne : 1,
417 internal::traits<Derived>::ColsAtCompileTime==1 ? 1 : SizeMinusOne> ConstStartMinusOne;
418 typedef CwiseUnaryOp<internal::scalar_quotient1_op<typename internal::traits<Derived>::Scalar>,
419 const ConstStartMinusOne > HNormalizedReturnType;
420
421 const HNormalizedReturnType hnormalized() const;
422
423 ////////// Householder module ///////////
424
425 void makeHouseholderInPlace(Scalar& tau, RealScalar& beta);
426 template<typename EssentialPart>
427 void makeHouseholder(EssentialPart& essential,
428 Scalar& tau, RealScalar& beta) const;
429 template<typename EssentialPart>
430 void applyHouseholderOnTheLeft(const EssentialPart& essential,
431 const Scalar& tau,
432 Scalar* workspace);
433 template<typename EssentialPart>
434 void applyHouseholderOnTheRight(const EssentialPart& essential,
435 const Scalar& tau,
436 Scalar* workspace);
437
438 ///////// Jacobi module /////////
439
440 template<typename OtherScalar>
441 void applyOnTheLeft(Index p, Index q, const JacobiRotation<OtherScalar>& j);
442 template<typename OtherScalar>
443 void applyOnTheRight(Index p, Index q, const JacobiRotation<OtherScalar>& j);
444
445 ///////// MatrixFunctions module /////////
446
447 typedef typename internal::stem_function<Scalar>::type StemFunction;
448 const MatrixExponentialReturnValue<Derived> exp() const;
449 const MatrixFunctionReturnValue<Derived> matrixFunction(StemFunction f) const;
450 const MatrixFunctionReturnValue<Derived> cosh() const;
451 const MatrixFunctionReturnValue<Derived> sinh() const;
452 const MatrixFunctionReturnValue<Derived> cos() const;
453 const MatrixFunctionReturnValue<Derived> sin() const;
454 const MatrixSquareRootReturnValue<Derived> sqrt() const;
455 const MatrixLogarithmReturnValue<Derived> log() const;
456 const MatrixPowerReturnValue<Derived> pow(const RealScalar& p) const;
457
458 #ifdef EIGEN2_SUPPORT
459 template<typename ProductDerived, typename Lhs, typename Rhs>
460 Derived& operator+=(const Flagged<ProductBase<ProductDerived, Lhs,Rhs>, 0,
461 EvalBeforeAssigningBit>& other);
462
463 template<typename ProductDerived, typename Lhs, typename Rhs>
464 Derived& operator-=(const Flagged<ProductBase<ProductDerived, Lhs,Rhs>, 0,
465 EvalBeforeAssigningBit>& other);
466
467 /** \deprecated because .lazy() is deprecated
468 * Overloaded for cache friendly product evaluation */
469 template<typename OtherDerived>
lazyAssign(const Flagged<OtherDerived,0,EvalBeforeAssigningBit> & other)470 Derived& lazyAssign(const Flagged<OtherDerived, 0, EvalBeforeAssigningBit>& other)
471 { return lazyAssign(other._expression()); }
472
473 template<unsigned int Added>
474 const Flagged<Derived, Added, 0> marked() const;
475 const Flagged<Derived, 0, EvalBeforeAssigningBit> lazy() const;
476
477 inline const Cwise<Derived> cwise() const;
478 inline Cwise<Derived> cwise();
479
480 VectorBlock<Derived> start(Index size);
481 const VectorBlock<const Derived> start(Index size) const;
482 VectorBlock<Derived> end(Index size);
483 const VectorBlock<const Derived> end(Index size) const;
484 template<int Size> VectorBlock<Derived,Size> start();
485 template<int Size> const VectorBlock<const Derived,Size> start() const;
486 template<int Size> VectorBlock<Derived,Size> end();
487 template<int Size> const VectorBlock<const Derived,Size> end() const;
488
489 Minor<Derived> minor(Index row, Index col);
490 const Minor<Derived> minor(Index row, Index col) const;
491 #endif
492
493 protected:
MatrixBase()494 MatrixBase() : Base() {}
495
496 private:
497 explicit MatrixBase(int);
498 MatrixBase(int,int);
499 template<typename OtherDerived> explicit MatrixBase(const MatrixBase<OtherDerived>&);
500 protected:
501 // mixing arrays and matrices is not legal
502 template<typename OtherDerived> Derived& operator+=(const ArrayBase<OtherDerived>& )
503 {EIGEN_STATIC_ASSERT(std::ptrdiff_t(sizeof(typename OtherDerived::Scalar))==-1,YOU_CANNOT_MIX_ARRAYS_AND_MATRICES); return *this;}
504 // mixing arrays and matrices is not legal
505 template<typename OtherDerived> Derived& operator-=(const ArrayBase<OtherDerived>& )
506 {EIGEN_STATIC_ASSERT(std::ptrdiff_t(sizeof(typename OtherDerived::Scalar))==-1,YOU_CANNOT_MIX_ARRAYS_AND_MATRICES); return *this;}
507 };
508
509
510 /***************************************************************************
511 * Implementation of matrix base methods
512 ***************************************************************************/
513
514 /** replaces \c *this by \c *this * \a other.
515 *
516 * \returns a reference to \c *this
517 *
518 * Example: \include MatrixBase_applyOnTheRight.cpp
519 * Output: \verbinclude MatrixBase_applyOnTheRight.out
520 */
521 template<typename Derived>
522 template<typename OtherDerived>
523 inline Derived&
524 MatrixBase<Derived>::operator*=(const EigenBase<OtherDerived> &other)
525 {
526 other.derived().applyThisOnTheRight(derived());
527 return derived();
528 }
529
530 /** replaces \c *this by \c *this * \a other. It is equivalent to MatrixBase::operator*=().
531 *
532 * Example: \include MatrixBase_applyOnTheRight.cpp
533 * Output: \verbinclude MatrixBase_applyOnTheRight.out
534 */
535 template<typename Derived>
536 template<typename OtherDerived>
applyOnTheRight(const EigenBase<OtherDerived> & other)537 inline void MatrixBase<Derived>::applyOnTheRight(const EigenBase<OtherDerived> &other)
538 {
539 other.derived().applyThisOnTheRight(derived());
540 }
541
542 /** replaces \c *this by \a other * \c *this.
543 *
544 * Example: \include MatrixBase_applyOnTheLeft.cpp
545 * Output: \verbinclude MatrixBase_applyOnTheLeft.out
546 */
547 template<typename Derived>
548 template<typename OtherDerived>
applyOnTheLeft(const EigenBase<OtherDerived> & other)549 inline void MatrixBase<Derived>::applyOnTheLeft(const EigenBase<OtherDerived> &other)
550 {
551 other.derived().applyThisOnTheLeft(derived());
552 }
553
554 } // end namespace Eigen
555
556 #endif // EIGEN_MATRIXBASE_H
557