|
|
namespace Eigen {
/** \page TopicCustomizingEigen Customizing/Extending Eigen
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.
\eigenAutoToc
\section ExtendingMatrixBase Extending MatrixBase (and other classes)
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.
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: \code class MatrixBase { // ... #ifdef EIGEN_MATRIXBASE_PLUGIN #include EIGEN_MATRIXBASE_PLUGIN #endif }; \endcode 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.
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.
Here is an example of an extension file for adding methods to MatrixBase: \n \b MatrixBaseAddons.h \code inline Scalar at(uint i, uint j) const { return this->operator()(i,j); } inline Scalar& at(uint i, uint j) { return this->operator()(i,j); } inline Scalar at(uint i) const { return this->operator[](i); } inline Scalar& at(uint i) { return this->operator[](i); }
inline RealScalar squaredLength() const { return squaredNorm(); } inline RealScalar length() const { return norm(); } inline RealScalar invLength(void) const { return fast_inv_sqrt(squaredNorm()); }
template<typename OtherDerived> inline Scalar squaredDistanceTo(const MatrixBase<OtherDerived>& other) const { return (derived() - other.derived()).squaredNorm(); }
template<typename OtherDerived> inline RealScalar distanceTo(const MatrixBase<OtherDerived>& other) const { return internal::sqrt(derived().squaredDistanceTo(other)); }
inline void scaleTo(RealScalar l) { RealScalar vl = norm(); if (vl>1e-9) derived() *= (l/vl); }
inline Transpose<Derived> transposed() {return this->transpose();} inline const Transpose<Derived> transposed() const {return this->transpose();}
inline uint minComponentId(void) const { int i; this->minCoeff(&i); return i; } inline uint maxComponentId(void) const { int i; this->maxCoeff(&i); return i; }
template<typename OtherDerived> void makeFloor(const MatrixBase<OtherDerived>& other) { derived() = derived().cwiseMin(other.derived()); } template<typename OtherDerived> void makeCeil(const MatrixBase<OtherDerived>& other) { derived() = derived().cwiseMax(other.derived()); }
const CwiseUnaryOp<internal::scalar_add_op<Scalar>, Derived> operator+(const Scalar& scalar) const { return CwiseUnaryOp<internal::scalar_add_op<Scalar>, Derived>(derived(), internal::scalar_add_op<Scalar>(scalar)); }
friend const CwiseUnaryOp<internal::scalar_add_op<Scalar>, Derived> operator+(const Scalar& scalar, const MatrixBase<Derived>& mat) { return CwiseUnaryOp<internal::scalar_add_op<Scalar>, Derived>(mat.derived(), internal::scalar_add_op<Scalar>(scalar)); } \endcode
Then one can the following declaration in the config.h or whatever prerequisites header file of his project: \code #define EIGEN_MATRIXBASE_PLUGIN "MatrixBaseAddons.h" \endcode
\section InheritingFromMatrix Inheriting from Matrix
Before inheriting from Matrix, be really, I mean REALLY, sure that using EIGEN_MATRIX_PLUGIN is not what you really want (see previous section). If you just need to add few members to Matrix, this is the way to go.
An example of when you actually need to inherit Matrix, is when you have several layers of heritage such as MyVerySpecificVector1, MyVerySpecificVector2 -> MyVector1 -> Matrix and MyVerySpecificVector3, MyVerySpecificVector4 -> MyVector2 -> Matrix.
In order for your object to work within the %Eigen framework, you need to define a few members in your inherited class.
Here is a minimalistic example:
\include CustomizingEigen_Inheritance.cpp
Output: \verbinclude CustomizingEigen_Inheritance.out
This is the kind of error you can get if you don't provide those methods \verbatim error: no match for ‘operator=’ in ‘v = Eigen::operator*( const Eigen::MatrixBase<Eigen::Matrix<double, -0x000000001, 1, 0, -0x000000001, 1> >::Scalar&, const Eigen::MatrixBase<Eigen::Matrix<double, -0x000000001, 1> >::StorageBaseType&) (((const Eigen::MatrixBase<Eigen::Matrix<double, -0x000000001, 1> >::StorageBaseType&) ((const Eigen::MatrixBase<Eigen::Matrix<double, -0x000000001, 1> >::StorageBaseType*)(& v))))’ \endverbatim
\anchor user_defined_scalars \section CustomScalarType Using custom scalar types
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 native integer types (e.g., \c int, \c unsigned \c int, \c short, etc.), and \c bool. On x86-64 systems, \c long \c double permits to locally enforces the use of x87 registers with extended accuracy (in comparison to SSE).
In order to add support for a custom type \c T you need: -# make sure the common operator (+,-,*,/,etc.) are supported by the type \c T -# add a specialization of struct Eigen::NumTraits<T> (see \ref NumTraits) -# 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. (see the file Eigen/src/Core/MathFunctions.h)
The math function should be defined in the same namespace than \c T, or in the \c std namespace though that second approach is not recommended.
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.
\code #ifndef ADOLCSUPPORT_H #define ADOLCSUPPORT_H
#define ADOLC_TAPELESS #include <adolc/adouble.h> #include <Eigen/Core>
namespace Eigen {
template<> struct NumTraits<adtl::adouble> : NumTraits<double> // permits to get the epsilon, dummy_precision, lowest, highest functions { typedef adtl::adouble Real; typedef adtl::adouble NonInteger; typedef adtl::adouble Nested;
enum { IsComplex = 0, IsInteger = 0, IsSigned = 1, RequireInitialization = 1, ReadCost = 1, AddCost = 3, MulCost = 3 }; };
}
namespace adtl {
inline const adouble& conj(const adouble& x) { return x; } inline const adouble& real(const adouble& x) { return x; } inline adouble imag(const adouble&) { return 0.; } inline adouble abs(const adouble& x) { return fabs(x); } inline adouble abs2(const adouble& x) { return x*x; }
}
#endif // ADOLCSUPPORT_H \endcode
This other example adds support for the \c mpq_class type from <a href="https://gmplib.org/">GMP</a>. It shows in particular how to change the way Eigen picks the best pivot during LU factorization. It selects the coefficient with the highest score, where the score is by default the absolute value of a number, but we can define a different score, for instance to prefer pivots with a more compact representation (this is an example, not a recommendation). Note that the scores should always be non-negative and only zero is allowed to have a score of zero. Also, this can interact badly with thresholds for inexact scalar types.
\code #include <gmpxx.h> #include <Eigen/Core> #include <boost/operators.hpp>
namespace Eigen { template<class> struct NumTraits; template<> struct NumTraits<mpq_class> { typedef mpq_class Real; typedef mpq_class NonInteger; typedef mpq_class Nested;
static inline Real epsilon() { return 0; } static inline Real dummy_precision() { return 0; }
enum { IsInteger = 0, IsSigned = 1, IsComplex = 0, RequireInitialization = 1, ReadCost = 6, AddCost = 150, MulCost = 100 }; };
namespace internal { template<> struct significant_decimals_impl<mpq_class> { // Infinite precision when printing static inline int run() { return 0; } };
template<> struct scalar_score_coeff_op<mpq_class> { struct result_type : boost::totally_ordered1<result_type> { std::size_t len; result_type(int i = 0) : len(i) {} // Eigen uses Score(0) and Score() result_type(mpq_class const& q) : len(mpz_size(q.get_num_mpz_t())+ mpz_size(q.get_den_mpz_t())-1) {} friend bool operator<(result_type x, result_type y) { // 0 is the worst possible pivot if (x.len == 0) return y.len > 0; if (y.len == 0) return false; // Prefer a pivot with a small representation return x.len > y.len; } friend bool operator==(result_type x, result_type y) { // Only used to test if the score is 0 return x.len == y.len; } }; result_type operator()(mpq_class const& x) const { return x; } }; } } \endcode
\sa \ref TopicPreprocessorDirectives
*/
}
|