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