Eigen-unsupported  3.4.90 (git rev 9589cc4e7fd8e4538bedef80dd36c7738977a8be)
 
Loading...
Searching...
No Matches
MPRealSupport
1// This file is part of a joint effort between Eigen, a lightweight C++ template library
2// for linear algebra, and MPFR C++, a C++ interface to MPFR library (http://www.holoborodko.com/pavel/)
3//
4// Copyright (C) 2010-2012 Pavel Holoborodko <pavel@holoborodko.com>
5// Copyright (C) 2010 Konstantin Holoborodko <konstantin@holoborodko.com>
6// Copyright (C) 2010 Gael Guennebaud <gael.guennebaud@inria.fr>
7//
8// This Source Code Form is subject to the terms of the Mozilla
9// Public License v. 2.0. If a copy of the MPL was not distributed
10// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
11
12#ifndef EIGEN_MPREALSUPPORT_MODULE_H
13#define EIGEN_MPREALSUPPORT_MODULE_H
14
15#include "../../Eigen/Core"
16#include <mpreal.h>
17
18namespace Eigen {
19
61
62template <>
63struct NumTraits<mpfr::mpreal> : GenericNumTraits<mpfr::mpreal> {
64 enum {
65 IsInteger = 0,
66 IsSigned = 1,
67 IsComplex = 0,
68 RequireInitialization = 1,
69 ReadCost = HugeCost,
70 AddCost = HugeCost,
71 MulCost = HugeCost
72 };
73
74 typedef mpfr::mpreal Real;
75 typedef mpfr::mpreal NonInteger;
76
77 static inline Real highest(long Precision = mpfr::mpreal::get_default_prec()) { return mpfr::maxval(Precision); }
78 static inline Real lowest(long Precision = mpfr::mpreal::get_default_prec()) { return -mpfr::maxval(Precision); }
79
80 // Constants
81 static inline Real Pi(long Precision = mpfr::mpreal::get_default_prec()) { return mpfr::const_pi(Precision); }
82 static inline Real Euler(long Precision = mpfr::mpreal::get_default_prec()) { return mpfr::const_euler(Precision); }
83 static inline Real Log2(long Precision = mpfr::mpreal::get_default_prec()) { return mpfr::const_log2(Precision); }
84 static inline Real Catalan(long Precision = mpfr::mpreal::get_default_prec()) {
85 return mpfr::const_catalan(Precision);
86 }
87
88 static inline Real epsilon(long Precision = mpfr::mpreal::get_default_prec()) {
89 return mpfr::machine_epsilon(Precision);
90 }
91 static inline Real epsilon(const Real& x) { return mpfr::machine_epsilon(x); }
92
93#ifdef MPREAL_HAVE_DYNAMIC_STD_NUMERIC_LIMITS
94 static inline int digits10(long Precision = mpfr::mpreal::get_default_prec()) {
95 return std::numeric_limits<Real>::digits10(Precision);
96 }
97 static inline int digits10(const Real& x) { return std::numeric_limits<Real>::digits10(x); }
98
99 static inline int max_digits10(long Precision = mpfr::mpreal::get_default_prec()) {
100 return std::numeric_limits<Real>::max_digits10(Precision);
101 }
102
103 static inline int digits() { return std::numeric_limits<Real>::digits(); }
104 static inline int digits(const Real& x) { return std::numeric_limits<Real>::digits(x); }
105#endif
106
107 static inline Real dummy_precision() {
108 mpfr_prec_t weak_prec = ((mpfr::mpreal::get_default_prec() - 1) * 90) / 100;
109 return mpfr::machine_epsilon(weak_prec);
110 }
111};
112
113namespace internal {
114
115template <>
116inline mpfr::mpreal random<mpfr::mpreal>() {
117 return mpfr::random();
118}
119
120template <>
121inline mpfr::mpreal random<mpfr::mpreal>(const mpfr::mpreal& a, const mpfr::mpreal& b) {
122 return a + (b - a) * random<mpfr::mpreal>();
123}
124
125inline bool isMuchSmallerThan(const mpfr::mpreal& a, const mpfr::mpreal& b, const mpfr::mpreal& eps) {
126 return mpfr::abs(a) <= mpfr::abs(b) * eps;
127}
128
129inline bool isApprox(const mpfr::mpreal& a, const mpfr::mpreal& b, const mpfr::mpreal& eps) {
130 return mpfr::isEqualFuzzy(a, b, eps);
131}
132
133inline bool isApproxOrLessThan(const mpfr::mpreal& a, const mpfr::mpreal& b, const mpfr::mpreal& eps) {
134 return a <= b || mpfr::isEqualFuzzy(a, b, eps);
135}
136
137template <>
138inline long double cast<mpfr::mpreal, long double>(const mpfr::mpreal& x) {
139 return x.toLDouble();
140}
141
142template <>
143inline double cast<mpfr::mpreal, double>(const mpfr::mpreal& x) {
144 return x.toDouble();
145}
146
147template <>
148inline long cast<mpfr::mpreal, long>(const mpfr::mpreal& x) {
149 return x.toLong();
150}
151
152template <>
153inline int cast<mpfr::mpreal, int>(const mpfr::mpreal& x) {
154 return int(x.toLong());
155}
156
157// Specialize GEBP kernel and traits for mpreal (no need for peeling, nor complicated stuff)
158// This also permits to directly call mpfr's routines and avoid many temporaries produced by mpreal
159template <>
160class gebp_traits<mpfr::mpreal, mpfr::mpreal, false, false> {
161 public:
162 typedef mpfr::mpreal ResScalar;
163 enum {
164 Vectorizable = false,
165 LhsPacketSize = 1,
166 RhsPacketSize = 1,
167 ResPacketSize = 1,
168 NumberOfRegisters = 1,
169 nr = 1,
170 mr = 1,
171 LhsProgress = 1,
172 RhsProgress = 1
173 };
174 typedef ResScalar LhsPacket;
175 typedef ResScalar RhsPacket;
176 typedef ResScalar ResPacket;
177 typedef LhsPacket LhsPacket4Packing;
178};
179
180template <typename Index, typename DataMapper, bool ConjugateLhs, bool ConjugateRhs>
181struct gebp_kernel<mpfr::mpreal, mpfr::mpreal, Index, DataMapper, 1, 1, ConjugateLhs, ConjugateRhs> {
182 typedef mpfr::mpreal mpreal;
183
184 EIGEN_DONT_INLINE void operator()(const DataMapper& res, const mpreal* blockA, const mpreal* blockB, Index rows,
185 Index depth, Index cols, const mpreal& alpha, Index strideA = -1,
186 Index strideB = -1, Index offsetA = 0, Index offsetB = 0) {
187 if (rows == 0 || cols == 0 || depth == 0) return;
188
189 mpreal acc1(0, mpfr_get_prec(blockA[0].mpfr_srcptr())), tmp(0, mpfr_get_prec(blockA[0].mpfr_srcptr()));
190
191 if (strideA == -1) strideA = depth;
192 if (strideB == -1) strideB = depth;
193
194 for (Index i = 0; i < rows; ++i) {
195 for (Index j = 0; j < cols; ++j) {
196 const mpreal* A = blockA + i * strideA + offsetA;
197 const mpreal* B = blockB + j * strideB + offsetB;
198
199 acc1 = 0;
200 for (Index k = 0; k < depth; k++) {
201 mpfr_mul(tmp.mpfr_ptr(), A[k].mpfr_srcptr(), B[k].mpfr_srcptr(), mpreal::get_default_rnd());
202 mpfr_add(acc1.mpfr_ptr(), acc1.mpfr_ptr(), tmp.mpfr_ptr(), mpreal::get_default_rnd());
203 }
204
205 mpfr_mul(acc1.mpfr_ptr(), acc1.mpfr_srcptr(), alpha.mpfr_srcptr(), mpreal::get_default_rnd());
206 mpfr_add(res(i, j).mpfr_ptr(), res(i, j).mpfr_srcptr(), acc1.mpfr_srcptr(), mpreal::get_default_rnd());
207 }
208 }
209 }
210};
211} // end namespace internal
212} // namespace Eigen
213
214#endif // EIGEN_MPREALSUPPORT_MODULE_H
Namespace containing all symbols from the Eigen library.
const int HugeCost
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index