aboutsummaryrefslogtreecommitdiffhomepage
path: root/test/boostmultiprec.cpp
blob: 3e16aeabd245079639cf91ce3c35367725428dc4 (plain)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
// This file is part of Eigen, a lightweight C++ template library
// for linear algebra.
//
// Copyright (C) 2016 Gael Guennebaud <gael.guennebaud@inria.fr>
//
// This Source Code Form is subject to the terms of the Mozilla
// Public License v. 2.0. If a copy of the MPL was not distributed
// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.

#include <sstream>

#ifdef EIGEN_TEST_MAX_SIZE
#undef EIGEN_TEST_MAX_SIZE
#endif

#define EIGEN_TEST_MAX_SIZE 50

#ifdef EIGEN_TEST_PART_1
#include "cholesky.cpp"
#endif

#ifdef EIGEN_TEST_PART_2
#include "lu.cpp"
#endif

#ifdef EIGEN_TEST_PART_3
#include "qr.cpp"
#endif

#ifdef EIGEN_TEST_PART_4
#include "qr_colpivoting.cpp"
#endif

#ifdef EIGEN_TEST_PART_5
#include "qr_fullpivoting.cpp"
#endif

#ifdef EIGEN_TEST_PART_6
#include "eigensolver_selfadjoint.cpp"
#endif

#ifdef EIGEN_TEST_PART_7
#include "jacobisvd.cpp"
#endif

#ifdef EIGEN_TEST_PART_8
#include "bdcsvd.cpp"
#endif

#include <Eigen/Dense>

#undef min
#undef max
#undef isnan
#undef isinf
#undef isfinite

#include <boost/multiprecision/cpp_dec_float.hpp>
#include <boost/multiprecision/number.hpp>

namespace mp = boost::multiprecision;
typedef mp::number<mp::cpp_dec_float<100>, mp::et_off> Real; // swith to et_on for testing with expression templates

namespace Eigen {
  template<> struct NumTraits<Real> : GenericNumTraits<Real> {
    static inline Real dummy_precision() { return 1e-50; }
  };

  template<typename T1,typename T2,typename T3,typename T4,typename T5>
  struct NumTraits<boost::multiprecision::detail::expression<T1,T2,T3,T4,T5> > : NumTraits<Real> {};

  template<>
  Real test_precision<Real>() { return 1e-50; }

  namespace internal {
    template<typename NewType>
    struct cast_impl<Real,NewType>
    {
      static inline NewType run(const Real& x)
      {
        return x.template convert_to<NewType>();
      }
    };
  }
}

namespace boost {
namespace multiprecision {
  // to make ADL works as expected:
  using boost::math::isfinite;

  // some specialization for the unit tests:
  inline bool test_isMuchSmallerThan(const Real& a, const Real& b) {
    return internal::isMuchSmallerThan(a, b, test_precision<Real>());
  }

  inline bool test_isApprox(const Real& a, const Real& b) {
    return internal::isApprox(a, b, test_precision<Real>());
  }

  inline bool test_isApproxOrLessThan(const Real& a, const Real& b) {
    return internal::isApproxOrLessThan(a, b, test_precision<Real>());
  }

  Real get_test_precision(const Real&) {
    return test_precision<Real>();
  }

  Real test_relative_error(const Real &a, const Real &b) {
    return Eigen::numext::sqrt(Real(Eigen::numext::abs2(a-b))/Real((Eigen::numext::mini)(Eigen::numext::abs2(a),Eigen::numext::abs2(b))));
  }
}
}

namespace Eigen {

}

void test_boostmultiprec()
{
  typedef Matrix<Real,Dynamic,Dynamic> Mat;

  std::cout << "NumTraits<Real>::epsilon()         = " << NumTraits<Real>::epsilon() << std::endl;
  std::cout << "NumTraits<Real>::dummy_precision() = " << NumTraits<Real>::dummy_precision() << std::endl;
  std::cout << "NumTraits<Real>::lowest()          = " << NumTraits<Real>::lowest() << std::endl;
  std::cout << "NumTraits<Real>::highest()         = " << NumTraits<Real>::highest() << std::endl;

  // chekc stream output
  {
    Mat A(10,10);
    A.setRandom();
    std::stringstream ss;
    ss << A;
  }

  for(int i = 0; i < g_repeat; i++) {
    int s = internal::random<int>(1,EIGEN_TEST_MAX_SIZE);

    CALL_SUBTEST_1( cholesky(Mat(s,s)) );

    CALL_SUBTEST_2( lu_non_invertible<Mat>() );
    CALL_SUBTEST_2( lu_invertible<Mat>() );

    CALL_SUBTEST_3( qr(Mat(internal::random<int>(1,EIGEN_TEST_MAX_SIZE),internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
    CALL_SUBTEST_3( qr_invertible<Mat>() );

    CALL_SUBTEST_4( qr<Mat>() );
    CALL_SUBTEST_4( cod<Mat>() );
    CALL_SUBTEST_4( qr_invertible<Mat>() );

    CALL_SUBTEST_5( qr<Mat>() );
    CALL_SUBTEST_5( qr_invertible<Mat>() );

    CALL_SUBTEST_6( selfadjointeigensolver(Mat(s,s)) );

    TEST_SET_BUT_UNUSED_VARIABLE(s)
  }

  CALL_SUBTEST_7(( jacobisvd(Mat(internal::random<int>(EIGEN_TEST_MAX_SIZE/4, EIGEN_TEST_MAX_SIZE), internal::random<int>(EIGEN_TEST_MAX_SIZE/4, EIGEN_TEST_MAX_SIZE/2))) ));
  CALL_SUBTEST_8(( bdcsvd(Mat(internal::random<int>(EIGEN_TEST_MAX_SIZE/4, EIGEN_TEST_MAX_SIZE), internal::random<int>(EIGEN_TEST_MAX_SIZE/4, EIGEN_TEST_MAX_SIZE/2))) ));
}