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 TopicCustomizing_Plugins by defining the preprocessor symbol \c EIGEN_MATRIXBASE_PLUGIN.
45 *
46 * \sa \blank \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>::StorageIndex StorageIndex;
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
70 using Base::derived;
71 using Base::const_cast_derived;
72 using Base::rows;
73 using Base::cols;
74 using Base::size;
75 using Base::coeff;
76 using Base::coeffRef;
77 using Base::lazyAssign;
78 using Base::eval;
79 using Base::operator+=;
80 using Base::operator-=;
81 using Base::operator*=;
82 using Base::operator/=;
83
84 typedef typename Base::CoeffReturnType CoeffReturnType;
85 typedef typename Base::ConstTransposeReturnType ConstTransposeReturnType;
86 typedef typename Base::RowXpr RowXpr;
87 typedef typename Base::ColXpr ColXpr;
88#endif // not EIGEN_PARSED_BY_DOXYGEN
89
90
91
92#ifndef EIGEN_PARSED_BY_DOXYGEN
93 /** type of the equivalent square matrix */
94 typedef Matrix<Scalar,EIGEN_SIZE_MAX(RowsAtCompileTime,ColsAtCompileTime),
95 EIGEN_SIZE_MAX(RowsAtCompileTime,ColsAtCompileTime)> SquareMatrixType;
96#endif // not EIGEN_PARSED_BY_DOXYGEN
97
98 /** \returns the size of the main diagonal, which is min(rows(),cols()).
99 * \sa rows(), cols(), SizeAtCompileTime. */
100 EIGEN_DEVICE_FUNC
101 inline Index diagonalSize() const { return (numext::mini)(rows(),cols()); }
102
103 typedef typename Base::PlainObject PlainObject;
104
105#ifndef EIGEN_PARSED_BY_DOXYGEN
106 /** \internal Represents a matrix with all coefficients equal to one another*/
107 typedef CwiseNullaryOp<internal::scalar_constant_op<Scalar>,PlainObject> ConstantReturnType;
108 /** \internal the return type of MatrixBase::adjoint() */
109 typedef typename internal::conditional<NumTraits<Scalar>::IsComplex,
110 CwiseUnaryOp<internal::scalar_conjugate_op<Scalar>, ConstTransposeReturnType>,
111 ConstTransposeReturnType
112 >::type AdjointReturnType;
113 /** \internal Return type of eigenvalues() */
114 typedef Matrix<std::complex<RealScalar>, internal::traits<Derived>::ColsAtCompileTime, 1, ColMajor> EigenvaluesReturnType;
115 /** \internal the return type of identity */
116 typedef CwiseNullaryOp<internal::scalar_identity_op<Scalar>,PlainObject> IdentityReturnType;
117 /** \internal the return type of unit vectors */
118 typedef Block<const CwiseNullaryOp<internal::scalar_identity_op<Scalar>, SquareMatrixType>,
119 internal::traits<Derived>::RowsAtCompileTime,
120 internal::traits<Derived>::ColsAtCompileTime> BasisReturnType;
121#endif // not EIGEN_PARSED_BY_DOXYGEN
122
123#define EIGEN_CURRENT_STORAGE_BASE_CLASS Eigen::MatrixBase
124#define EIGEN_DOC_UNARY_ADDONS(X,Y)
125# include "../plugins/CommonCwiseUnaryOps.h"
126# include "../plugins/CommonCwiseBinaryOps.h"
127# include "../plugins/MatrixCwiseUnaryOps.h"
128# include "../plugins/MatrixCwiseBinaryOps.h"
129# ifdef EIGEN_MATRIXBASE_PLUGIN
130# include EIGEN_MATRIXBASE_PLUGIN
131# endif
132#undef EIGEN_CURRENT_STORAGE_BASE_CLASS
133#undef EIGEN_DOC_UNARY_ADDONS
134
135 /** Special case of the template operator=, in order to prevent the compiler
136 * from generating a default operator= (issue hit with g++ 4.1)
137 */
138 EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE
139 Derived& operator=(const MatrixBase& other);
140
141 // We cannot inherit here via Base::operator= since it is causing
142 // trouble with MSVC.
143
144 template <typename OtherDerived>
145 EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE
146 Derived& operator=(const DenseBase<OtherDerived>& other);
147
148 template <typename OtherDerived>
149 EIGEN_DEVICE_FUNC
150 Derived& operator=(const EigenBase<OtherDerived>& other);
151
152 template<typename OtherDerived>
153 EIGEN_DEVICE_FUNC
154 Derived& operator=(const ReturnByValue<OtherDerived>& other);
155
156 template<typename OtherDerived>
157 EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE
158 Derived& operator+=(const MatrixBase<OtherDerived>& other);
159 template<typename OtherDerived>
160 EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE
161 Derived& operator-=(const MatrixBase<OtherDerived>& other);
162
163 template<typename OtherDerived>
164 EIGEN_DEVICE_FUNC
165 const Product<Derived,OtherDerived>
166 operator*(const MatrixBase<OtherDerived> &other) const;
167
168 template<typename OtherDerived>
169 EIGEN_DEVICE_FUNC
170 const Product<Derived,OtherDerived,LazyProduct>
171 lazyProduct(const MatrixBase<OtherDerived> &other) const;
172
173 template<typename OtherDerived>
174 Derived& operator*=(const EigenBase<OtherDerived>& other);
175
176 template<typename OtherDerived>
177 void applyOnTheLeft(const EigenBase<OtherDerived>& other);
178
179 template<typename OtherDerived>
180 void applyOnTheRight(const EigenBase<OtherDerived>& other);
181
182 template<typename DiagonalDerived>
183 EIGEN_DEVICE_FUNC
184 const Product<Derived, DiagonalDerived, LazyProduct>
185 operator*(const DiagonalBase<DiagonalDerived> &diagonal) const;
186
187 template<typename OtherDerived>
188 EIGEN_DEVICE_FUNC
189 typename ScalarBinaryOpTraits<typename internal::traits<Derived>::Scalar,typename internal::traits<OtherDerived>::Scalar>::ReturnType
190 dot(const MatrixBase<OtherDerived>& other) const;
191
192 EIGEN_DEVICE_FUNC RealScalar squaredNorm() const;
193 EIGEN_DEVICE_FUNC RealScalar norm() const;
194 RealScalar stableNorm() const;
195 RealScalar blueNorm() const;
196 RealScalar hypotNorm() const;
197 EIGEN_DEVICE_FUNC const PlainObject normalized() const;
198 EIGEN_DEVICE_FUNC const PlainObject stableNormalized() const;
199 EIGEN_DEVICE_FUNC void normalize();
200 EIGEN_DEVICE_FUNC void stableNormalize();
201
202 EIGEN_DEVICE_FUNC const AdjointReturnType adjoint() const;
203 EIGEN_DEVICE_FUNC void adjointInPlace();
204
205 typedef Diagonal<Derived> DiagonalReturnType;
206 EIGEN_DEVICE_FUNC
207 DiagonalReturnType diagonal();
208
209 typedef typename internal::add_const<Diagonal<const Derived> >::type ConstDiagonalReturnType;
210 EIGEN_DEVICE_FUNC
211 ConstDiagonalReturnType diagonal() const;
212
213 template<int Index> struct DiagonalIndexReturnType { typedef Diagonal<Derived,Index> Type; };
214 template<int Index> struct ConstDiagonalIndexReturnType { typedef const Diagonal<const Derived,Index> Type; };
215
216 template<int Index>
217 EIGEN_DEVICE_FUNC
218 typename DiagonalIndexReturnType<Index>::Type diagonal();
219
220 template<int Index>
221 EIGEN_DEVICE_FUNC
222 typename ConstDiagonalIndexReturnType<Index>::Type diagonal() const;
223
224 typedef Diagonal<Derived,DynamicIndex> DiagonalDynamicIndexReturnType;
225 typedef typename internal::add_const<Diagonal<const Derived,DynamicIndex> >::type ConstDiagonalDynamicIndexReturnType;
226
227 EIGEN_DEVICE_FUNC
228 DiagonalDynamicIndexReturnType diagonal(Index index);
229 EIGEN_DEVICE_FUNC
230 ConstDiagonalDynamicIndexReturnType diagonal(Index index) const;
231
232 template<unsigned int Mode> struct TriangularViewReturnType { typedef TriangularView<Derived, Mode> Type; };
233 template<unsigned int Mode> struct ConstTriangularViewReturnType { typedef const TriangularView<const Derived, Mode> Type; };
234
235 template<unsigned int Mode>
236 EIGEN_DEVICE_FUNC
237 typename TriangularViewReturnType<Mode>::Type triangularView();
238 template<unsigned int Mode>
239 EIGEN_DEVICE_FUNC
240 typename ConstTriangularViewReturnType<Mode>::Type triangularView() const;
241
242 template<unsigned int UpLo> struct SelfAdjointViewReturnType { typedef SelfAdjointView<Derived, UpLo> Type; };
243 template<unsigned int UpLo> struct ConstSelfAdjointViewReturnType { typedef const SelfAdjointView<const Derived, UpLo> Type; };
244
245 template<unsigned int UpLo>
246 EIGEN_DEVICE_FUNC
247 typename SelfAdjointViewReturnType<UpLo>::Type selfadjointView();
248 template<unsigned int UpLo>
249 EIGEN_DEVICE_FUNC
250 typename ConstSelfAdjointViewReturnType<UpLo>::Type selfadjointView() const;
251
252 const SparseView<Derived> sparseView(const Scalar& m_reference = Scalar(0),
253 const typename NumTraits<Scalar>::Real& m_epsilon = NumTraits<Scalar>::dummy_precision()) const;
254 EIGEN_DEVICE_FUNC static const IdentityReturnType Identity();
255 EIGEN_DEVICE_FUNC static const IdentityReturnType Identity(Index rows, Index cols);
256 EIGEN_DEVICE_FUNC static const BasisReturnType Unit(Index size, Index i);
257 EIGEN_DEVICE_FUNC static const BasisReturnType Unit(Index i);
258 EIGEN_DEVICE_FUNC static const BasisReturnType UnitX();
259 EIGEN_DEVICE_FUNC static const BasisReturnType UnitY();
260 EIGEN_DEVICE_FUNC static const BasisReturnType UnitZ();
261 EIGEN_DEVICE_FUNC static const BasisReturnType UnitW();
262
263 EIGEN_DEVICE_FUNC
264 const DiagonalWrapper<const Derived> asDiagonal() const;
265 const PermutationWrapper<const Derived> asPermutation() const;
266
267 EIGEN_DEVICE_FUNC
268 Derived& setIdentity();
269 EIGEN_DEVICE_FUNC
270 Derived& setIdentity(Index rows, Index cols);
271
272 bool isIdentity(const RealScalar& prec = NumTraits<Scalar>::dummy_precision()) const;
273 bool isDiagonal(const RealScalar& prec = NumTraits<Scalar>::dummy_precision()) const;
274
275 bool isUpperTriangular(const RealScalar& prec = NumTraits<Scalar>::dummy_precision()) const;
276 bool isLowerTriangular(const RealScalar& prec = NumTraits<Scalar>::dummy_precision()) const;
277
278 template<typename OtherDerived>
279 bool isOrthogonal(const MatrixBase<OtherDerived>& other,
280 const RealScalar& prec = NumTraits<Scalar>::dummy_precision()) const;
281 bool isUnitary(const RealScalar& prec = NumTraits<Scalar>::dummy_precision()) const;
282
283 /** \returns true if each coefficients of \c *this and \a other are all exactly equal.
284 * \warning When using floating point scalar values you probably should rather use a
285 * fuzzy comparison such as isApprox()
286 * \sa isApprox(), operator!= */
287 template<typename OtherDerived>
288 EIGEN_DEVICE_FUNC inline bool operator==(const MatrixBase<OtherDerived>& other) const
289 { return cwiseEqual(other).all(); }
290
291 /** \returns true if at least one pair of coefficients of \c *this and \a other are not exactly equal to each other.
292 * \warning When using floating point scalar values you probably should rather use a
293 * fuzzy comparison such as isApprox()
294 * \sa isApprox(), operator== */
295 template<typename OtherDerived>
296 EIGEN_DEVICE_FUNC inline bool operator!=(const MatrixBase<OtherDerived>& other) const
297 { return cwiseNotEqual(other).any(); }
298
299 NoAlias<Derived,Eigen::MatrixBase > noalias();
300
301 // TODO forceAlignedAccess is temporarily disabled
302 // Need to find a nicer workaround.
303 inline const Derived& forceAlignedAccess() const { return derived(); }
304 inline Derived& forceAlignedAccess() { return derived(); }
305 template<bool Enable> inline const Derived& forceAlignedAccessIf() const { return derived(); }
306 template<bool Enable> inline Derived& forceAlignedAccessIf() { return derived(); }
307
308 EIGEN_DEVICE_FUNC Scalar trace() const;
309
310 template<int p> EIGEN_DEVICE_FUNC RealScalar lpNorm() const;
311
312 EIGEN_DEVICE_FUNC MatrixBase<Derived>& matrix() { return *this; }
313 EIGEN_DEVICE_FUNC const MatrixBase<Derived>& matrix() const { return *this; }
314
315 /** \returns an \link Eigen::ArrayBase Array \endlink expression of this matrix
316 * \sa ArrayBase::matrix() */
317 EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE ArrayWrapper<Derived> array() { return ArrayWrapper<Derived>(derived()); }
318 /** \returns a const \link Eigen::ArrayBase Array \endlink expression of this matrix
319 * \sa ArrayBase::matrix() */
320 EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE const ArrayWrapper<const Derived> array() const { return ArrayWrapper<const Derived>(derived()); }
321
322/////////// LU module ///////////
323
324 inline const FullPivLU<PlainObject> fullPivLu() const;
325 inline const PartialPivLU<PlainObject> partialPivLu() const;
326
327 inline const PartialPivLU<PlainObject> lu() const;
328
329 inline const Inverse<Derived> inverse() const;
330
331 template<typename ResultType>
332 inline void computeInverseAndDetWithCheck(
333 ResultType& inverse,
334 typename ResultType::Scalar& determinant,
335 bool& invertible,
336 const RealScalar& absDeterminantThreshold = NumTraits<Scalar>::dummy_precision()
337 ) const;
338 template<typename ResultType>
339 inline void computeInverseWithCheck(
340 ResultType& inverse,
341 bool& invertible,
342 const RealScalar& absDeterminantThreshold = NumTraits<Scalar>::dummy_precision()
343 ) const;
344 Scalar determinant() const;
345
346/////////// Cholesky module ///////////
347
348 inline const LLT<PlainObject> llt() const;
349 inline const LDLT<PlainObject> ldlt() const;
350
351/////////// QR module ///////////
352
353 inline const HouseholderQR<PlainObject> householderQr() const;
354 inline const ColPivHouseholderQR<PlainObject> colPivHouseholderQr() const;
355 inline const FullPivHouseholderQR<PlainObject> fullPivHouseholderQr() const;
356 inline const CompleteOrthogonalDecomposition<PlainObject> completeOrthogonalDecomposition() const;
357
358/////////// Eigenvalues module ///////////
359
360 inline EigenvaluesReturnType eigenvalues() const;
361 inline RealScalar operatorNorm() const;
362
363/////////// SVD module ///////////
364
365 inline JacobiSVD<PlainObject> jacobiSvd(unsigned int computationOptions = 0) const;
366 inline BDCSVD<PlainObject> bdcSvd(unsigned int computationOptions = 0) const;
367
368/////////// Geometry module ///////////
369
370 #ifndef EIGEN_PARSED_BY_DOXYGEN
371 /// \internal helper struct to form the return type of the cross product
372 template<typename OtherDerived> struct cross_product_return_type {
373 typedef typename ScalarBinaryOpTraits<typename internal::traits<Derived>::Scalar,typename internal::traits<OtherDerived>::Scalar>::ReturnType Scalar;
374 typedef Matrix<Scalar,MatrixBase::RowsAtCompileTime,MatrixBase::ColsAtCompileTime> type;
375 };
376 #endif // EIGEN_PARSED_BY_DOXYGEN
377 template<typename OtherDerived>
378 EIGEN_DEVICE_FUNC
379#ifndef EIGEN_PARSED_BY_DOXYGEN
380 inline typename cross_product_return_type<OtherDerived>::type
381#else
382 inline PlainObject
383#endif
384 cross(const MatrixBase<OtherDerived>& other) const;
385
386 template<typename OtherDerived>
387 EIGEN_DEVICE_FUNC
388 inline PlainObject cross3(const MatrixBase<OtherDerived>& other) const;
389
390 EIGEN_DEVICE_FUNC
391 inline PlainObject unitOrthogonal(void) const;
392
393 EIGEN_DEVICE_FUNC
394 inline Matrix<Scalar,3,1> eulerAngles(Index a0, Index a1, Index a2) const;
395
396 // put this as separate enum value to work around possible GCC 4.3 bug (?)
397 enum { HomogeneousReturnTypeDirection = ColsAtCompileTime==1&&RowsAtCompileTime==1 ? ((internal::traits<Derived>::Flags&RowMajorBit)==RowMajorBit ? Horizontal : Vertical)
398 : ColsAtCompileTime==1 ? Vertical : Horizontal };
399 typedef Homogeneous<Derived, HomogeneousReturnTypeDirection> HomogeneousReturnType;
400 EIGEN_DEVICE_FUNC
401 inline HomogeneousReturnType homogeneous() const;
402
403 enum {
404 SizeMinusOne = SizeAtCompileTime==Dynamic ? Dynamic : SizeAtCompileTime-1
405 };
406 typedef Block<const Derived,
407 internal::traits<Derived>::ColsAtCompileTime==1 ? SizeMinusOne : 1,
408 internal::traits<Derived>::ColsAtCompileTime==1 ? 1 : SizeMinusOne> ConstStartMinusOne;
409 typedef EIGEN_EXPR_BINARYOP_SCALAR_RETURN_TYPE(ConstStartMinusOne,Scalar,quotient) HNormalizedReturnType;
410 EIGEN_DEVICE_FUNC
411 inline const HNormalizedReturnType hnormalized() const;
412
413////////// Householder module ///////////
414
415 void makeHouseholderInPlace(Scalar& tau, RealScalar& beta);
416 template<typename EssentialPart>
417 void makeHouseholder(EssentialPart& essential,
418 Scalar& tau, RealScalar& beta) const;
419 template<typename EssentialPart>
420 void applyHouseholderOnTheLeft(const EssentialPart& essential,
421 const Scalar& tau,
422 Scalar* workspace);
423 template<typename EssentialPart>
424 void applyHouseholderOnTheRight(const EssentialPart& essential,
425 const Scalar& tau,
426 Scalar* workspace);
427
428///////// Jacobi module /////////
429
430 template<typename OtherScalar>
431 void applyOnTheLeft(Index p, Index q, const JacobiRotation<OtherScalar>& j);
432 template<typename OtherScalar>
433 void applyOnTheRight(Index p, Index q, const JacobiRotation<OtherScalar>& j);
434
435///////// SparseCore module /////////
436
437 template<typename OtherDerived>
438 EIGEN_STRONG_INLINE const typename SparseMatrixBase<OtherDerived>::template CwiseProductDenseReturnType<Derived>::Type
439 cwiseProduct(const SparseMatrixBase<OtherDerived> &other) const
440 {
441 return other.cwiseProduct(derived());
442 }
443
444///////// MatrixFunctions module /////////
445
446 typedef typename internal::stem_function<Scalar>::type StemFunction;
447#define EIGEN_MATRIX_FUNCTION(ReturnType, Name, Description) \
448 /** \returns an expression of the matrix Description of \c *this. \brief This function requires the <a href="unsupported/group__MatrixFunctions__Module.html"> unsupported MatrixFunctions module</a>. To compute the coefficient-wise Description use ArrayBase::##Name . */ \
449 const ReturnType<Derived> Name() const;
450#define EIGEN_MATRIX_FUNCTION_1(ReturnType, Name, Description, Argument) \
451 /** \returns an expression of the matrix Description of \c *this. \brief This function requires the <a href="unsupported/group__MatrixFunctions__Module.html"> unsupported MatrixFunctions module</a>. To compute the coefficient-wise Description use ArrayBase::##Name . */ \
452 const ReturnType<Derived> Name(Argument) const;
453
454 EIGEN_MATRIX_FUNCTION(MatrixExponentialReturnValue, exp, exponential)
455 /** \brief Helper function for the <a href="unsupported/group__MatrixFunctions__Module.html"> unsupported MatrixFunctions module</a>.*/
456 const MatrixFunctionReturnValue<Derived> matrixFunction(StemFunction f) const;
457 EIGEN_MATRIX_FUNCTION(MatrixFunctionReturnValue, cosh, hyperbolic cosine)
458 EIGEN_MATRIX_FUNCTION(MatrixFunctionReturnValue, sinh, hyperbolic sine)
459 EIGEN_MATRIX_FUNCTION(MatrixFunctionReturnValue, cos, cosine)
460 EIGEN_MATRIX_FUNCTION(MatrixFunctionReturnValue, sin, sine)
461 EIGEN_MATRIX_FUNCTION(MatrixSquareRootReturnValue, sqrt, square root)
462 EIGEN_MATRIX_FUNCTION(MatrixLogarithmReturnValue, log, logarithm)
463 EIGEN_MATRIX_FUNCTION_1(MatrixPowerReturnValue, pow, power to \c p, const RealScalar& p)
464 EIGEN_MATRIX_FUNCTION_1(MatrixComplexPowerReturnValue, pow, power to \c p, const std::complex<RealScalar>& p)
465
466 protected:
467 EIGEN_DEVICE_FUNC MatrixBase() : Base() {}
468
469 private:
470 EIGEN_DEVICE_FUNC explicit MatrixBase(int);
471 EIGEN_DEVICE_FUNC MatrixBase(int,int);
472 template<typename OtherDerived> EIGEN_DEVICE_FUNC explicit MatrixBase(const MatrixBase<OtherDerived>&);
473 protected:
474 // mixing arrays and matrices is not legal
475 template<typename OtherDerived> Derived& operator+=(const ArrayBase<OtherDerived>& )
476 {EIGEN_STATIC_ASSERT(std::ptrdiff_t(sizeof(typename OtherDerived::Scalar))==-1,YOU_CANNOT_MIX_ARRAYS_AND_MATRICES); return *this;}
477 // mixing arrays and matrices is not legal
478 template<typename OtherDerived> Derived& operator-=(const ArrayBase<OtherDerived>& )
479 {EIGEN_STATIC_ASSERT(std::ptrdiff_t(sizeof(typename OtherDerived::Scalar))==-1,YOU_CANNOT_MIX_ARRAYS_AND_MATRICES); return *this;}
480};
481
482
483/***************************************************************************
484* Implementation of matrix base methods
485***************************************************************************/
486
487/** replaces \c *this by \c *this * \a other.
488 *
489 * \returns a reference to \c *this
490 *
491 * Example: \include MatrixBase_applyOnTheRight.cpp
492 * Output: \verbinclude MatrixBase_applyOnTheRight.out
493 */
494template<typename Derived>
495template<typename OtherDerived>
496inline Derived&
497MatrixBase<Derived>::operator*=(const EigenBase<OtherDerived> &other)
498{
499 other.derived().applyThisOnTheRight(derived());
500 return derived();
501}
502
503/** replaces \c *this by \c *this * \a other. It is equivalent to MatrixBase::operator*=().
504 *
505 * Example: \include MatrixBase_applyOnTheRight.cpp
506 * Output: \verbinclude MatrixBase_applyOnTheRight.out
507 */
508template<typename Derived>
509template<typename OtherDerived>
510inline void MatrixBase<Derived>::applyOnTheRight(const EigenBase<OtherDerived> &other)
511{
512 other.derived().applyThisOnTheRight(derived());
513}
514
515/** replaces \c *this by \a other * \c *this.
516 *
517 * Example: \include MatrixBase_applyOnTheLeft.cpp
518 * Output: \verbinclude MatrixBase_applyOnTheLeft.out
519 */
520template<typename Derived>
521template<typename OtherDerived>
522inline void MatrixBase<Derived>::applyOnTheLeft(const EigenBase<OtherDerived> &other)
523{
524 other.derived().applyThisOnTheLeft(derived());
525}
526
527} // end namespace Eigen
528
529#endif // EIGEN_MATRIXBASE_H
530