source: pacpussensors/trunk/Vislab/lib3dv-1.2.0/lib3dv/eigen/doc/CustomizingEigen.dox

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

Doc

File size: 8.0 KB
Line 
1namespace Eigen {
2
3/** \page TopicCustomizingEigen Customizing/Extending Eigen
4
5Eigen 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\eigenAutoToc
8
9\section ExtendingMatrixBase Extending MatrixBase (and other classes)
10
11In 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.
12
13You 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:
14\code
15class MatrixBase {
16 // ...
17 #ifdef EIGEN_MATRIXBASE_PLUGIN
18 #include EIGEN_MATRIXBASE_PLUGIN
19 #endif
20};
21\endcode
22Therefore 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.
23
24You 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.
25
26Here is an example of an extension file for adding methods to MatrixBase: \n
27\b MatrixBaseAddons.h
28\code
29inline Scalar at(uint i, uint j) const { return this->operator()(i,j); }
30inline Scalar& at(uint i, uint j) { return this->operator()(i,j); }
31inline Scalar at(uint i) const { return this->operator[](i); }
32inline Scalar& at(uint i) { return this->operator[](i); }
33
34inline RealScalar squaredLength() const { return squaredNorm(); }
35inline RealScalar length() const { return norm(); }
36inline RealScalar invLength(void) const { return fast_inv_sqrt(squaredNorm()); }
37
38template<typename OtherDerived>
39inline Scalar squaredDistanceTo(const MatrixBase<OtherDerived>& other) const
40{ return (derived() - other.derived()).squaredNorm(); }
41
42template<typename OtherDerived>
43inline RealScalar distanceTo(const MatrixBase<OtherDerived>& other) const
44{ return internal::sqrt(derived().squaredDistanceTo(other)); }
45
46inline void scaleTo(RealScalar l) { RealScalar vl = norm(); if (vl>1e-9) derived() *= (l/vl); }
47
48inline Transpose<Derived> transposed() {return this->transpose();}
49inline const Transpose<Derived> transposed() const {return this->transpose();}
50
51inline uint minComponentId(void) const { int i; this->minCoeff(&i); return i; }
52inline uint maxComponentId(void) const { int i; this->maxCoeff(&i); return i; }
53
54template<typename OtherDerived>
55void makeFloor(const MatrixBase<OtherDerived>& other) { derived() = derived().cwiseMin(other.derived()); }
56template<typename OtherDerived>
57void makeCeil(const MatrixBase<OtherDerived>& other) { derived() = derived().cwiseMax(other.derived()); }
58
59const CwiseUnaryOp<internal::scalar_add_op<Scalar>, Derived>
60operator+(const Scalar& scalar) const
61{ return CwiseUnaryOp<internal::scalar_add_op<Scalar>, Derived>(derived(), internal::scalar_add_op<Scalar>(scalar)); }
62
63friend const CwiseUnaryOp<internal::scalar_add_op<Scalar>, Derived>
64operator+(const Scalar& scalar, const MatrixBase<Derived>& mat)
65{ return CwiseUnaryOp<internal::scalar_add_op<Scalar>, Derived>(mat.derived(), internal::scalar_add_op<Scalar>(scalar)); }
66\endcode
67
68Then one can the following declaration in the config.h or whatever prerequisites header file of his project:
69\code
70#define EIGEN_MATRIXBASE_PLUGIN "MatrixBaseAddons.h"
71\endcode
72
73\section InheritingFromMatrix Inheriting from Matrix
74
75Before inheriting from Matrix, be really, i mean REALLY sure that using
76EIGEN_MATRIX_PLUGIN is not what you really want (see previous section).
77If you just need to add few members to Matrix, this is the way to go.
78
79An example of when you actually need to inherit Matrix, is when you have
80several layers of heritage such as MyVerySpecificVector1,MyVerySpecificVector1 -> MyVector1 -> Matrix and.
81MyVerySpecificVector3,MyVerySpecificVector4 -> MyVector2 -> Matrix.
82
83In order for your object to work within the %Eigen framework, you need to
84define a few members in your inherited class.
85
86Here is a minimalistic example:\n
87\code
88class MyVectorType : public Eigen::VectorXd
89{
90public:
91 MyVectorType(void):Eigen::VectorXd() {}
92
93 typedef Eigen::VectorXd Base;
94
95 // This constructor allows you to construct MyVectorType from Eigen expressions
96 template<typename OtherDerived>
97 MyVectorType(const Eigen::MatrixBase<OtherDerived>& other)
98 : Eigen::Vector3d(other)
99 { }
100
101 // This method allows you to assign Eigen expressions to MyVectorType
102 template<typename OtherDerived>
103 MyVectorType & operator= (const Eigen::MatrixBase <OtherDerived>& other)
104 {
105 this->Base::operator=(other);
106 return *this;
107 }
108};
109\endcode
110
111This is the kind of error you can get if you don't provide those methods
112\code
113error: no match for ‘operator=’ in ‘delta =
114(((Eigen::MatrixBase<Eigen::Matrix<std::complex<float>, 10000, 1, 2, 10000,
1151> >*)(& delta)) + 8u)->Eigen::MatrixBase<Derived>::cwise [with Derived =
116Eigen::Matrix<std::complex<float>, 10000, 1, 2, 10000,
1171>]().Eigen::Cwise<ExpressionType>::operator* [with OtherDerived =
118Eigen::Matrix<std::complex<float>, 10000, 1, 2, 10000, 1>, ExpressionType =
119Eigen::Matrix<std::complex<float>, 10000, 1, 2, 10000, 1>](((const
120Eigen::MatrixBase<Eigen::Matrix<std::complex<float>, 10000, 1, 2, 10000, 1>
121>&)(((const Eigen::MatrixBase<Eigen::Matrix<std::complex<float>, 10000, 1,
122>2, 10000, 1> >*)((const spectral1d*)where)) + 8u)))’
123\endcode
124
125\anchor user_defined_scalars \section CustomScalarType Using custom scalar types
126
127By 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.
128On x86-64 systems, \c long \c double permits to locally enforces the use of x87 registers with extended accuracy (in comparison to SSE).
129
130In order to add support for a custom type \c T you need:
131-# make sure the common operator (+,-,*,/,etc.) are supported by the type \c T
132-# add a specialization of struct Eigen::NumTraits<T> (see \ref NumTraits)
133-# 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.
134 (see the file Eigen/src/Core/MathFunctions.h)
135
136The 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.
137
138Here 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.
139
140\code
141#ifndef ADOLCSUPPORT_H
142#define ADOLCSUPPORT_H
143
144#define ADOLC_TAPELESS
145#include <adolc/adouble.h>
146#include <Eigen/Core>
147
148namespace Eigen {
149
150template<> struct NumTraits<adtl::adouble>
151 : NumTraits<double> // permits to get the epsilon, dummy_precision, lowest, highest functions
152{
153 typedef adtl::adouble Real;
154 typedef adtl::adouble NonInteger;
155 typedef adtl::adouble Nested;
156
157 enum {
158 IsComplex = 0,
159 IsInteger = 0,
160 IsSigned = 1,
161 RequireInitialization = 1,
162 ReadCost = 1,
163 AddCost = 3,
164 MulCost = 3
165 };
166};
167
168}
169
170namespace adtl {
171
172inline const adouble& conj(const adouble& x) { return x; }
173inline const adouble& real(const adouble& x) { return x; }
174inline adouble imag(const adouble&) { return 0.; }
175inline adouble abs(const adouble& x) { return fabs(x); }
176inline adouble abs2(const adouble& x) { return x*x; }
177
178}
179
180#endif // ADOLCSUPPORT_H
181\endcode
182
183
184\sa \ref TopicPreprocessorDirectives
185
186*/
187
188}
Note: See TracBrowser for help on using the repository browser.