By default, Eigen currently supports standard floating-point types (float
, double
, std::complex<float>
, std::complex<double>
, long
double
), as well as all native integer types (e.g., int
, unsigned
int
, short
, etc.), and bool
. On x86-64 systems, long
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 T
you need:
- make sure the common operator (+,-,*,/,etc.) are supported by the type
T
- add a specialization of struct Eigen::NumTraits<T> (see 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 T
, or in the std
namespace though that second approach is not recommended.
Here is a concrete example adding support for the Adolc's adouble
type. Adolc is an automatic differentiation library. The type adouble
is basically a real value tracking the values of any number of partial derivatives.
#ifndef ADOLCSUPPORT_H
#define ADOLCSUPPORT_H
#define ADOLC_TAPELESS
#include <adolc/adouble.h>
template<> struct NumTraits<adtl::adouble>
: NumTraits<double>
{
typedef adtl::adouble
Real;
typedef adtl::adouble NonInteger;
typedef adtl::adouble Nested;
enum {
};
};
}
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
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_imag_op< typename Derived::Scalar >, const Derived > imag(const Eigen::ArrayBase< Derived > &x)
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_abs2_op< typename Derived::Scalar >, const Derived > abs2(const Eigen::ArrayBase< Derived > &x)
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_conjugate_op< typename Derived::Scalar >, const Derived > conj(const Eigen::ArrayBase< Derived > &x)
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_real_op< typename Derived::Scalar >, const Derived > real(const Eigen::ArrayBase< Derived > &x)
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_abs_op< typename Derived::Scalar >, const Derived > abs(const Eigen::ArrayBase< Derived > &x)
This other example adds support for the mpq_class
type from GMP. 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.
#include <gmpxx.h>
#include <boost/operators.hpp>
template<> struct NumTraits<mpq_class> : GenericNumTraits<mpq_class>
{
typedef mpq_class NonInteger;
typedef mpq_class Nested;
static inline Real epsilon() {
return 0; }
static inline Real dummy_precision() {
return 0; }
static inline int digits10() { return 0; }
enum {
};
};
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) {}
result_type(mpq_class const& q) :
len(mpz_size(q.get_num_mpz_t())+
mpz_size(q.get_den_mpz_t())-1) {}
if (
x.len == 0)
return y.len > 0;
if (
y.len == 0)
return false;
}
}
};
result_type
operator()(mpq_class
const&
x)
const {
return x; }
};
}
}
IndexedView_or_Block operator()(const RowIndices &rowIndices, const ColIndices &colIndices)
bool operator==(const bfloat16 &a, const bfloat16 &b)
bool operator<(const bfloat16 &a, const bfloat16 &b)