• Home
  • Line#
  • Scopes#
  • Navigate#
  • Raw
  • Download
1namespace Eigen {
2
3/** \page TopicCustomizing_CustomScalar Using custom scalar types
4\anchor user_defined_scalars
5
6By 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.
7On x86-64 systems, \c long \c double permits to locally enforces the use of x87 registers with extended accuracy (in comparison to SSE).
8
9In order to add support for a custom type \c T you need:
10-# make sure the common operator (+,-,*,/,etc.) are supported by the type \c T
11-# add a specialization of struct Eigen::NumTraits<T> (see \ref NumTraits)
12-# 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.
13     (see the file Eigen/src/Core/MathFunctions.h)
14
15The 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.
16
17Here 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.
18
19\code
20#ifndef ADOLCSUPPORT_H
21#define ADOLCSUPPORT_H
22
23#define ADOLC_TAPELESS
24#include <adolc/adouble.h>
25#include <Eigen/Core>
26
27namespace Eigen {
28
29template<> struct NumTraits<adtl::adouble>
30 : NumTraits<double> // permits to get the epsilon, dummy_precision, lowest, highest functions
31{
32  typedef adtl::adouble Real;
33  typedef adtl::adouble NonInteger;
34  typedef adtl::adouble Nested;
35
36  enum {
37    IsComplex = 0,
38    IsInteger = 0,
39    IsSigned = 1,
40    RequireInitialization = 1,
41    ReadCost = 1,
42    AddCost = 3,
43    MulCost = 3
44  };
45};
46
47}
48
49namespace adtl {
50
51inline const adouble& conj(const adouble& x)  { return x; }
52inline const adouble& real(const adouble& x)  { return x; }
53inline adouble imag(const adouble&)    { return 0.; }
54inline adouble abs(const adouble&  x)  { return fabs(x); }
55inline adouble abs2(const adouble& x)  { return x*x; }
56
57}
58
59#endif // ADOLCSUPPORT_H
60\endcode
61
62This 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.
63
64\code
65#include <gmpxx.h>
66#include <Eigen/Core>
67#include <boost/operators.hpp>
68
69namespace Eigen {
70  template<> struct NumTraits<mpq_class> : GenericNumTraits<mpq_class>
71  {
72    typedef mpq_class Real;
73    typedef mpq_class NonInteger;
74    typedef mpq_class Nested;
75
76    static inline Real epsilon() { return 0; }
77    static inline Real dummy_precision() { return 0; }
78    static inline Real digits10() { return 0; }
79
80    enum {
81      IsInteger = 0,
82      IsSigned = 1,
83      IsComplex = 0,
84      RequireInitialization = 1,
85      ReadCost = 6,
86      AddCost = 150,
87      MulCost = 100
88    };
89  };
90
91  namespace internal {
92
93    template<> struct scalar_score_coeff_op<mpq_class> {
94      struct result_type : boost::totally_ordered1<result_type> {
95        std::size_t len;
96        result_type(int i = 0) : len(i) {} // Eigen uses Score(0) and Score()
97        result_type(mpq_class const& q) :
98          len(mpz_size(q.get_num_mpz_t())+
99              mpz_size(q.get_den_mpz_t())-1) {}
100        friend bool operator<(result_type x, result_type y) {
101          // 0 is the worst possible pivot
102          if (x.len == 0) return y.len > 0;
103          if (y.len == 0) return false;
104          // Prefer a pivot with a small representation
105          return x.len > y.len;
106        }
107        friend bool operator==(result_type x, result_type y) {
108          // Only used to test if the score is 0
109          return x.len == y.len;
110        }
111      };
112      result_type operator()(mpq_class const& x) const { return x; }
113    };
114  }
115}
116\endcode
117
118*/
119
120}
121