summaryrefslogtreecommitdiff
path: root/libs/math/performance/test_igamma.cpp
blob: 91150443fdd7910dbb859927428c462ff8660993 (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
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
//  Copyright John Maddock 2007.
//  Use, modification and distribution are subject to the
//  Boost Software License, Version 1.0. (See accompanying file
//  LICENSE_1_0.txt or copy at http://www.boost.org/LICENSE_1_0.txt)

#include "required_defines.hpp"

#include "performance_measure.hpp"

#include <boost/math/special_functions/gamma.hpp>
#include <boost/array.hpp>

#define T double
#include "../test/igamma_big_data.ipp"
#include "../test/igamma_int_data.ipp"
#include "../test/igamma_med_data.ipp"
#include "../test/igamma_small_data.ipp"

template <std::size_t N>
double igamma_evaluate2(const boost::array<boost::array<T, 6>, N>& data)
{
   double result = 0;
   for(unsigned i = 0; i < N; ++i)
   {
      result += boost::math::gamma_p(data[i][0], data[i][1]);
      result += boost::math::gamma_q(data[i][0], data[i][1]);
   }
   return result;
}

BOOST_MATH_PERFORMANCE_TEST(igamma_test, "igamma")
{
   double result = igamma_evaluate2(igamma_big_data);
   result += igamma_evaluate2(igamma_int_data);
   result += igamma_evaluate2(igamma_med_data);
   result += igamma_evaluate2(igamma_small_data);

   consume_result(result);
   set_call_count(
      2 * (sizeof(igamma_big_data) 
      + sizeof(igamma_int_data) 
      + sizeof(igamma_med_data)
      + sizeof(igamma_small_data)) / sizeof(igamma_big_data[0]));
}

template <std::size_t N>
double igamma_inv_evaluate2(const boost::array<boost::array<T, 6>, N>& data)
{
   double result = 0;
   for(unsigned i = 0; i < N; ++i)
   {
      result += boost::math::gamma_p_inv(data[i][0], data[i][5]);
      result += boost::math::gamma_q_inv(data[i][0], data[i][3]);
   }
   return result;
}

BOOST_MATH_PERFORMANCE_TEST(igamma_inv_test, "igamma_inv")
{
   double result = igamma_inv_evaluate2(igamma_big_data);
   result += igamma_inv_evaluate2(igamma_int_data);
   result += igamma_inv_evaluate2(igamma_med_data);
   result += igamma_inv_evaluate2(igamma_small_data);

   consume_result(result);
   set_call_count(
      2 * (sizeof(igamma_big_data) 
      + sizeof(igamma_int_data) 
      + sizeof(igamma_med_data)
      + sizeof(igamma_small_data)) / sizeof(igamma_big_data[0]));
}

template <std::size_t N>
double igamma_inva_evaluate2(const boost::array<boost::array<T, 6>, N>& data)
{
   double result = 0;
   for(unsigned i = 0; i < N; ++i)
   {
      result += boost::math::gamma_p_inva(data[i][1], data[i][5]);
      result += boost::math::gamma_q_inva(data[i][1], data[i][3]);
   }
   return result;
}

BOOST_MATH_PERFORMANCE_TEST(igamma_inva_test, "igamma_inva")
{
   double result = igamma_inva_evaluate2(igamma_big_data);
   result += igamma_inva_evaluate2(igamma_int_data);
   result += igamma_inva_evaluate2(igamma_med_data);
   result += igamma_inva_evaluate2(igamma_small_data);

   consume_result(result);
   set_call_count(
      2 * (sizeof(igamma_big_data) 
      + sizeof(igamma_int_data) 
      + sizeof(igamma_med_data)
      + sizeof(igamma_small_data)) / sizeof(igamma_big_data[0]));
}

#ifdef TEST_CEPHES

extern "C" {

double igam(double, double);
double igami(double, double);

}

template <std::size_t N>
double igamma_evaluate_cephes(const boost::array<boost::array<T, 6>, N>& data)
{
   double result = 0;
   for(unsigned i = 0; i < N; ++i)
      result += igam(data[i][0], data[i][1]);
   return result;
}

BOOST_MATH_PERFORMANCE_TEST(igamma_test, "igamma-cephes")
{
   double result = igamma_evaluate_cephes(igamma_big_data);
   result += igamma_evaluate_cephes(igamma_int_data);
   result += igamma_evaluate_cephes(igamma_med_data);
   result += igamma_evaluate_cephes(igamma_small_data);

   consume_result(result);
   set_call_count(
      (sizeof(igamma_big_data) 
      + sizeof(igamma_int_data) 
      + sizeof(igamma_med_data)
      + sizeof(igamma_small_data)) / sizeof(igamma_big_data[0]));
}

template <std::size_t N>
double igamma_inv_evaluate_cephes(const boost::array<boost::array<T, 6>, N>& data)
{
   double result = 0;
   for(unsigned i = 0; i < N; ++i)
      result += igami(data[i][0], data[i][3]); // note needs complement of probability!!
   return result;
}

//
// This test does not run to completion, gets stuck
// in infinite loop inside cephes....
//
#if 0
BOOST_MATH_PERFORMANCE_TEST(igamma_inv_test, "igamma_inv-cephes")
{
   double result = igamma_inv_evaluate_cephes(igamma_big_data);
   result += igamma_inv_evaluate_cephes(igamma_int_data);
   result += igamma_inv_evaluate_cephes(igamma_med_data);
   result += igamma_inv_evaluate_cephes(igamma_small_data);

   consume_result(result);
   set_call_count(
      (sizeof(igamma_big_data) 
      + sizeof(igamma_int_data) 
      + sizeof(igamma_med_data)
      + sizeof(igamma_small_data)) / sizeof(igamma_big_data[0]));
}
#endif

#endif

#ifdef TEST_GSL

#include <gsl/gsl_sf_gamma.h>

template <std::size_t N>
double igamma_evaluate_gsl(const boost::array<boost::array<T, 6>, N>& data)
{
   double result = 0;
   for(unsigned i = 0; i < N; ++i)
      result += gsl_sf_gamma_inc_P(data[i][0], data[i][1]);
   return result;
}

BOOST_MATH_PERFORMANCE_TEST(igamma_test, "igamma-gsl")
{
   double result = igamma_evaluate_gsl(igamma_big_data);
   result += igamma_evaluate_gsl(igamma_int_data);
   result += igamma_evaluate_gsl(igamma_med_data);
   result += igamma_evaluate_gsl(igamma_small_data);

   consume_result(result);
   set_call_count(
      (sizeof(igamma_big_data) 
      + sizeof(igamma_int_data) 
      + sizeof(igamma_med_data)
      + sizeof(igamma_small_data)) / sizeof(igamma_big_data[0]));
}

#endif

#ifdef TEST_DCDFLIB
#include <dcdflib.h>
namespace dcd{

inline double gamma_q(double x, double y)
{ 
   double ans, qans;
   int i = 0;
   gamma_inc (&x, &y, &ans, &qans, &i); 
   return qans;
}

inline double gamma_p(double x, double y)
{ 
   double ans, qans;
   int i = 0;
   gamma_inc (&x, &y, &ans, &qans, &i); 
   return ans;
}

inline double gamma_q_inv(double x, double y)
{ 
   double ans, p, nul;
   int i = 0;
   p = 1 - y;
   nul = 0;
   gamma_inc_inv (&x, &ans, &nul, &p, &y, &i); 
   return ans;
}

inline double gamma_p_inv(double x, double y)
{ 
   double ans, p, nul;
   int i = 0;
   p = 1 - y;
   nul = 0;
   gamma_inc_inv (&x, &ans, &nul, &y, &p, &i); 
   return ans;
}

}

template <std::size_t N>
double igamma_evaluate2_dcd(const boost::array<boost::array<T, 6>, N>& data)
{
   double result = 0;
   for(unsigned i = 0; i < N; ++i)
   {
      result += dcd::gamma_p(data[i][0], data[i][1]);
      result += dcd::gamma_q(data[i][0], data[i][1]);
   }
   return result;
}

BOOST_MATH_PERFORMANCE_TEST(igamma_test, "igamma-dcd")
{
   double result = igamma_evaluate2_dcd(igamma_big_data);
   result += igamma_evaluate2_dcd(igamma_int_data);
   result += igamma_evaluate2_dcd(igamma_med_data);
   result += igamma_evaluate2_dcd(igamma_small_data);

   consume_result(result);
   set_call_count(
      2 * (sizeof(igamma_big_data) 
      + sizeof(igamma_int_data) 
      + sizeof(igamma_med_data)
      + sizeof(igamma_small_data)) / sizeof(igamma_big_data[0]));
}

template <std::size_t N>
double igamma_inv_evaluate2_dcd(const boost::array<boost::array<T, 6>, N>& data)
{
   double result = 0;
   for(unsigned i = 0; i < N; ++i)
   {
      result += dcd::gamma_p_inv(data[i][0], data[i][5]);
      result += dcd::gamma_q_inv(data[i][0], data[i][3]);
   }
   return result;
}

BOOST_MATH_PERFORMANCE_TEST(igamma_inv_test, "igamma_inv-dcd")
{
   double result = igamma_inv_evaluate2_dcd(igamma_big_data);
   result += igamma_inv_evaluate2_dcd(igamma_int_data);
   result += igamma_inv_evaluate2_dcd(igamma_med_data);
   result += igamma_inv_evaluate2_dcd(igamma_small_data);

   consume_result(result);
   set_call_count(
      2 * (sizeof(igamma_big_data) 
      + sizeof(igamma_int_data) 
      + sizeof(igamma_med_data)
      + sizeof(igamma_small_data)) / sizeof(igamma_big_data[0]));
}

#endif