1 namespace Eigen { 2 3 /** \page TopicCustomizingEigen Customizing/Extending Eigen 4 5 Eigen can be extended in several ways, for instance, by defining global methods, \ref ExtendingMatrixBase "by adding custom methods to MatrixBase", adding support to \ref CustomScalarType "custom types" etc. 6 7 \b Table \b of \b contents 8 - \ref ExtendingMatrixBase 9 - \ref InheritingFromMatrix 10 - \ref CustomScalarType 11 12 \section ExtendingMatrixBase Extending MatrixBase (and other classes) 13 14 In this section we will see how to add custom methods to MatrixBase. Since all expressions and matrix types inherit MatrixBase, adding a method to MatrixBase make it immediately available to all expressions ! A typical use case is, for instance, to make Eigen compatible with another API. 15 16 You certainly know that in C++ it is not possible to add methods to an existing class. So how that's possible ? Here the trick is to include in the declaration of MatrixBase a file defined by the preprocessor token \c EIGEN_MATRIXBASE_PLUGIN: 17 \code 18 class MatrixBase { 19 // ... 20 #ifdef EIGEN_MATRIXBASE_PLUGIN 21 #include EIGEN_MATRIXBASE_PLUGIN 22 #endif 23 }; 24 \endcode 25 Therefore to extend MatrixBase with your own methods you just have to create a file with your method declaration and define EIGEN_MATRIXBASE_PLUGIN before you include any Eigen's header file. 26 27 You can extend many of the other classes used in Eigen by defining similarly named preprocessor symbols. For instance, define \c EIGEN_ARRAYBASE_PLUGIN if you want to extend the ArrayBase class. A full list of classes that can be extended in this way and the corresponding preprocessor symbols can be found on our page \ref TopicPreprocessorDirectives. 28 29 Here is an example of an extension file for adding methods to MatrixBase: \n 30 \b MatrixBaseAddons.h 31 \code 32 inline Scalar at(uint i, uint j) const { return this->operator()(i,j); } 33 inline Scalar& at(uint i, uint j) { return this->operator()(i,j); } 34 inline Scalar at(uint i) const { return this->operator[](i); } 35 inline Scalar& at(uint i) { return this->operator[](i); } 36 37 inline RealScalar squaredLength() const { return squaredNorm(); } 38 inline RealScalar length() const { return norm(); } 39 inline RealScalar invLength(void) const { return fast_inv_sqrt(squaredNorm()); } 40 41 template<typename OtherDerived> 42 inline Scalar squaredDistanceTo(const MatrixBase<OtherDerived>& other) const 43 { return (derived() - other.derived()).squaredNorm(); } 44 45 template<typename OtherDerived> 46 inline RealScalar distanceTo(const MatrixBase<OtherDerived>& other) const 47 { return internal::sqrt(derived().squaredDistanceTo(other)); } 48 49 inline void scaleTo(RealScalar l) { RealScalar vl = norm(); if (vl>1e-9) derived() *= (l/vl); } 50 51 inline Transpose<Derived> transposed() {return this->transpose();} 52 inline const Transpose<Derived> transposed() const {return this->transpose();} 53 54 inline uint minComponentId(void) const { int i; this->minCoeff(&i); return i; } 55 inline uint maxComponentId(void) const { int i; this->maxCoeff(&i); return i; } 56 57 template<typename OtherDerived> 58 void makeFloor(const MatrixBase<OtherDerived>& other) { derived() = derived().cwiseMin(other.derived()); } 59 template<typename OtherDerived> 60 void makeCeil(const MatrixBase<OtherDerived>& other) { derived() = derived().cwiseMax(other.derived()); } 61 62 const CwiseUnaryOp<internal::scalar_add_op<Scalar>, Derived> 63 operator+(const Scalar& scalar) const 64 { return CwiseUnaryOp<internal::scalar_add_op<Scalar>, Derived>(derived(), internal::scalar_add_op<Scalar>(scalar)); } 65 66 friend const CwiseUnaryOp<internal::scalar_add_op<Scalar>, Derived> 67 operator+(const Scalar& scalar, const MatrixBase<Derived>& mat) 68 { return CwiseUnaryOp<internal::scalar_add_op<Scalar>, Derived>(mat.derived(), internal::scalar_add_op<Scalar>(scalar)); } 69 \endcode 70 71 Then one can the following declaration in the config.h or whatever prerequisites header file of his project: 72 \code 73 #define EIGEN_MATRIXBASE_PLUGIN "MatrixBaseAddons.h" 74 \endcode 75 76 \section InheritingFromMatrix Inheriting from Matrix 77 78 Before inheriting from Matrix, be really, i mean REALLY sure that using 79 EIGEN_MATRIX_PLUGIN is not what you really want (see previous section). 80 If you just need to add few members to Matrix, this is the way to go. 81 82 An example of when you actually need to inherit Matrix, is when you have 83 several layers of heritage such as MyVerySpecificVector1,MyVerySpecificVector1 -> MyVector1 -> Matrix and. 84 MyVerySpecificVector3,MyVerySpecificVector4 -> MyVector2 -> Matrix. 85 86 In order for your object to work within the %Eigen framework, you need to 87 define a few members in your inherited class. 88 89 Here is a minimalistic example:\n 90 \code 91 class MyVectorType : public Eigen::VectorXd 92 { 93 public: 94 MyVectorType(void):Eigen::VectorXd() {} 95 96 typedef Eigen::VectorXd Base; 97 98 // This constructor allows you to construct MyVectorType from Eigen expressions 99 template<typename OtherDerived> 100 MyVectorType(const Eigen::MatrixBase<OtherDerived>& other) 101 : Eigen::Vector3d(other) 102 { } 103 104 // This method allows you to assign Eigen expressions to MyVectorType 105 template<typename OtherDerived> 106 MyVectorType & operator= (const Eigen::MatrixBase <OtherDerived>& other) 107 { 108 this->Base::operator=(other); 109 return *this; 110 } 111 }; 112 \endcode 113 114 This is the kind of error you can get if you don't provide those methods 115 \code 116 error: no match for operator= in delta = 117 (((Eigen::MatrixBase<Eigen::Matrix<std::complex<float>, 10000, 1, 2, 10000, 118 1> >*)(& delta)) + 8u)->Eigen::MatrixBase<Derived>::cwise [with Derived = 119 Eigen::Matrix<std::complex<float>, 10000, 1, 2, 10000, 120 1>]().Eigen::Cwise<ExpressionType>::operator* [with OtherDerived = 121 Eigen::Matrix<std::complex<float>, 10000, 1, 2, 10000, 1>, ExpressionType = 122 Eigen::Matrix<std::complex<float>, 10000, 1, 2, 10000, 1>](((const 123 Eigen::MatrixBase<Eigen::Matrix<std::complex<float>, 10000, 1, 2, 10000, 1> 124 >&)(((const Eigen::MatrixBase<Eigen::Matrix<std::complex<float>, 10000, 1, 125 >2, 10000, 1> >*)((const spectral1d*)where)) + 8u))) 126 \endcode 127 128 \anchor user_defined_scalars \section CustomScalarType Using custom scalar types 129 130 By default, Eigen currently supports standard floating-point types (\c float, \c double, \c std::complex<float>, \c std::complex<double>, \c long \c double), as well as all integrale types (e.g., \c int, \c unsigned \c int, \c short, etc.), and \c bool. 131 On x86-64 systems, \c long \c double permits to locally enforces the use of x87 registers with extended accuracy (in comparison to SSE). 132 133 In order to add support for a custom type \c T you need: 134 -# make sure the common operator (+,-,*,/,etc.) are supported by the type \c T 135 -# add a specialization of struct Eigen::NumTraits<T> (see \ref NumTraits) 136 -# define the math functions that makes sense for your type. This includes standard ones like sqrt, pow, sin, tan, conj, real, imag, etc, as well as abs2 which is Eigen specific. 137 (see the file Eigen/src/Core/MathFunctions.h) 138 139 The math function should be defined in the same namespace than \c T, or in the \c std namespace though that second appraoch is not recommended. 140 141 Here is a concrete example adding support for the Adolc's \c adouble type. <a href="https://projects.coin-or.org/ADOL-C">Adolc</a> is an automatic differentiation library. The type \c adouble is basically a real value tracking the values of any number of partial derivatives. 142 143 \code 144 #ifndef ADOLCSUPPORT_H 145 #define ADOLCSUPPORT_H 146 147 #define ADOLC_TAPELESS 148 #include <adolc/adouble.h> 149 #include <Eigen/Core> 150 151 namespace Eigen { 152 153 template<> struct NumTraits<adtl::adouble> 154 : NumTraits<double> // permits to get the epsilon, dummy_precision, lowest, highest functions 155 { 156 typedef adtl::adouble Real; 157 typedef adtl::adouble NonInteger; 158 typedef adtl::adouble Nested; 159 160 enum { 161 IsComplex = 0, 162 IsInteger = 0, 163 IsSigned = 1, 164 RequireInitialization = 1, 165 ReadCost = 1, 166 AddCost = 3, 167 MulCost = 3 168 }; 169 }; 170 171 } 172 173 namespace adtl { 174 175 inline const adouble& conj(const adouble& x) { return x; } 176 inline const adouble& real(const adouble& x) { return x; } 177 inline adouble imag(const adouble&) { return 0.; } 178 inline adouble abs(const adouble& x) { return fabs(x); } 179 inline adouble abs2(const adouble& x) { return x*x; } 180 181 } 182 183 #endif // ADOLCSUPPORT_H 184 \endcode 185 186 187 \sa \ref TopicPreprocessorDirectives 188 189 */ 190 191 } 192