1 //  Copyright John Maddock 2009.
2 //  Use, modification and distribution are subject to the
3 //  Boost Software License, Version 1.0. (See accompanying file
4 //  LICENSE_1_0.txt or copy at http://www.boost.org/LICENSE_1_0.txt)
5 
6 #include "required_defines.hpp"
7 
8 #include "performance_measure.hpp"
9 
10 #include <boost/math/special_functions/bessel.hpp>
11 #include <boost/array.hpp>
12 
13 #define T double
14 #include "../test/bessel_j_data.ipp"
15 #include "../test/bessel_j_int_data.ipp"
16 #include "../test/bessel_j_large_data.ipp"
17 #include "../test/bessel_y01_data.ipp"
18 #include "../test/bessel_yn_data.ipp"
19 #include "../test/bessel_yv_data.ipp"
20 #include "../test/bessel_k_int_data.ipp"
21 #include "../test/bessel_k_data.ipp"
22 #include "../test/bessel_i_int_data.ipp"
23 #include "../test/bessel_i_data.ipp"
24 #include "../test/sph_bessel_data.ipp"
25 #include "../test/sph_neumann_data.ipp"
26 
27 
28 template <std::size_t N>
bessel_evaluate2(const boost::array<boost::array<T,3>,N> & data)29 double bessel_evaluate2(const boost::array<boost::array<T, 3>, N>& data)
30 {
31    double result = 0;
32    for(unsigned i = 0; i < N; ++i)
33       result += boost::math::cyl_bessel_j(data[i][0], data[i][1]);
34    return result;
35 }
36 
37 BOOST_MATH_PERFORMANCE_TEST(bessel_test, "cyl_bessel_j")
38 {
39    double result= bessel_evaluate2(bessel_j_data);
40    result += bessel_evaluate2(bessel_j_int_data);
41    result += bessel_evaluate2(bessel_j_large_data);
42 
43    consume_result(result);
44    set_call_count(
45       (sizeof(bessel_j_data)
46       + sizeof(bessel_j_int_data)
47       + sizeof(bessel_j_large_data)) / sizeof(bessel_j_data[0]));
48 }
49 
50 template <std::size_t N>
bessel_y_evaluate2(const boost::array<boost::array<T,3>,N> & data)51 double bessel_y_evaluate2(const boost::array<boost::array<T, 3>, N>& data)
52 {
53    double result = 0;
54    for(unsigned i = 0; i < N; ++i)
55       result += boost::math::cyl_neumann(data[i][0], data[i][1]);
56    return result;
57 }
58 
59 BOOST_MATH_PERFORMANCE_TEST(bessel_test, "cyl_neumann")
60 {
61    double result= bessel_y_evaluate2(bessel_y01_data);
62    result += bessel_y_evaluate2(bessel_yn_data);
63    result += bessel_y_evaluate2(bessel_yv_data);
64 
65    consume_result(result);
66    set_call_count(
67       (sizeof(bessel_y01_data)
68       + sizeof(bessel_yn_data)
69       + sizeof(bessel_yv_data)) / sizeof(bessel_j_data[0]));
70 }
71 
72 template <std::size_t N>
bessel_i_evaluate(const boost::array<boost::array<T,3>,N> & data)73 double bessel_i_evaluate(const boost::array<boost::array<T, 3>, N>& data)
74 {
75    double result = 0;
76    for(unsigned i = 0; i < N; ++i)
77       result += boost::math::cyl_bessel_i(data[i][0], data[i][1]);
78    return result;
79 }
80 
81 BOOST_MATH_PERFORMANCE_TEST(bessel_test, "cyl_bessel_i")
82 {
83    double result= bessel_i_evaluate(bessel_i_int_data);
84    result += bessel_i_evaluate(bessel_i_data);
85 
86    consume_result(result);
87    set_call_count(
88       (sizeof(bessel_i_data)
89       + sizeof(bessel_i_int_data)) / sizeof(bessel_i_data[0]));
90 }
91 
92 template <std::size_t N>
bessel_k_evaluate(const boost::array<boost::array<T,3>,N> & data)93 double bessel_k_evaluate(const boost::array<boost::array<T, 3>, N>& data)
94 {
95    double result = 0;
96    for(unsigned i = 0; i < N; ++i)
97       result += boost::math::cyl_bessel_k(data[i][0], data[i][1]);
98    return result;
99 }
100 
101 BOOST_MATH_PERFORMANCE_TEST(bessel_test, "cyl_bessel_k")
102 {
103    double result= bessel_k_evaluate(bessel_k_data);
104    result += bessel_k_evaluate(bessel_k_int_data);
105 
106    consume_result(result);
107    set_call_count(
108       (sizeof(bessel_k_data)
109       + sizeof(bessel_k_int_data)) / sizeof(bessel_k_data[0]));
110 }
111 
112 template <std::size_t N>
bessel_sph_j_evaluate(const boost::array<boost::array<T,3>,N> & data)113 double bessel_sph_j_evaluate(const boost::array<boost::array<T, 3>, N>& data)
114 {
115    double result = 0;
116    for(unsigned i = 0; i < N; ++i)
117       result += boost::math::sph_bessel(static_cast<unsigned>(data[i][0]), data[i][1]);
118    return result;
119 }
120 
121 BOOST_MATH_PERFORMANCE_TEST(bessel_test, "sph_bessel")
122 {
123    double result= bessel_sph_j_evaluate(sph_bessel_data);
124 
125    consume_result(result);
126    set_call_count(
127       sizeof(sph_bessel_data) / sizeof(sph_bessel_data[0]));
128 }
129 
130 template <std::size_t N>
bessel_sph_y_evaluate(const boost::array<boost::array<T,3>,N> & data)131 double bessel_sph_y_evaluate(const boost::array<boost::array<T, 3>, N>& data)
132 {
133    double result = 0;
134    for(unsigned i = 0; i < N; ++i)
135       result += boost::math::sph_neumann(static_cast<unsigned>(data[i][0]), data[i][1]);
136    return result;
137 }
138 
139 BOOST_MATH_PERFORMANCE_TEST(bessel_test, "sph_neumann")
140 {
141    double result= bessel_sph_y_evaluate(sph_neumann_data);
142 
143    consume_result(result);
144    set_call_count(
145       sizeof(sph_neumann_data) / sizeof(sph_neumann_data[0]));
146 }
147 
148 #ifdef TEST_DCDFLIB
149 
150 #endif
151 
152 #ifdef TEST_CEPHES
153 
154 extern "C" double jv(double, double);
155 extern "C" double yv(double, double);
156 extern "C" double iv(double, double);
157 
158 template <std::size_t N>
bessel_evaluate2_cephes(const boost::array<boost::array<T,3>,N> & data)159 double bessel_evaluate2_cephes(const boost::array<boost::array<T, 3>, N>& data)
160 {
161    double result = 0;
162    for(unsigned i = 0; i < N; ++i)
163       result += jv(data[i][0], data[i][1]);
164    return result;
165 }
166 
167 BOOST_MATH_PERFORMANCE_TEST(bessel_test, "cyl_bessel_j-cephes")
168 {
169    double result = bessel_evaluate2_cephes(bessel_j_data);
170    result += bessel_evaluate2_cephes(bessel_j_int_data);
171    result += bessel_evaluate2_cephes(bessel_j_large_data);
172 
173    consume_result(result);
174    set_call_count(
175       (sizeof(bessel_j_data)
176       + sizeof(bessel_j_int_data)
177       + sizeof(bessel_j_large_data)) / sizeof(bessel_j_data[0]));
178 }
179 
180 template <std::size_t N>
bessel_y_evaluate2_cephes(const boost::array<boost::array<T,3>,N> & data)181 double bessel_y_evaluate2_cephes(const boost::array<boost::array<T, 3>, N>& data)
182 {
183    double result = 0;
184    for(unsigned i = 0; i < N; ++i)
185       result += yv(data[i][0], data[i][1]);
186    return result;
187 }
188 
189 BOOST_MATH_PERFORMANCE_TEST(bessel_test, "cyl_neumann-cephes")
190 {
191    double result= bessel_y_evaluate2_cephes(bessel_y01_data);
192    result += bessel_y_evaluate2_cephes(bessel_yn_data);
193    result += bessel_y_evaluate2_cephes(bessel_yv_data);
194 
195    consume_result(result);
196    set_call_count(
197       (sizeof(bessel_y01_data)
198       + sizeof(bessel_yn_data)
199       + sizeof(bessel_yv_data)) / sizeof(bessel_j_data[0]));
200 }
201 
202 template <std::size_t N>
bessel_i_evaluate_cephes(const boost::array<boost::array<T,3>,N> & data)203 double bessel_i_evaluate_cephes(const boost::array<boost::array<T, 3>, N>& data)
204 {
205    double result = 0;
206    for(unsigned i = 0; i < N; ++i)
207       result += iv(data[i][0], data[i][1]);
208    return result;
209 }
210 
211 BOOST_MATH_PERFORMANCE_TEST(bessel_test, "cyl_bessel_i-cephes")
212 {
213    double result= bessel_i_evaluate_cephes(bessel_i_int_data);
214    result += bessel_i_evaluate_cephes(bessel_i_data);
215 
216    consume_result(result);
217    set_call_count(
218       (sizeof(bessel_i_data)
219       + sizeof(bessel_i_int_data)) / sizeof(bessel_i_data[0]));
220 }
221 
222 #endif
223 
224 #ifdef TEST_GSL
225 
226 #include <gsl/gsl_sf_bessel.h>
227 
228 template <std::size_t N>
bessel_evaluate2_gsl(const boost::array<boost::array<T,3>,N> & data)229 double bessel_evaluate2_gsl(const boost::array<boost::array<T, 3>, N>& data)
230 {
231    double result = 0;
232    for(unsigned i = 0; i < N; ++i)
233       result += gsl_sf_bessel_Jnu(data[i][0], data[i][1]);
234    return result;
235 }
236 
237 BOOST_MATH_PERFORMANCE_TEST(bessel_test, "cyl_bessel_j-gsl")
238 {
239    double result = bessel_evaluate2_gsl(bessel_j_data);
240    result += bessel_evaluate2_gsl(bessel_j_int_data);
241    result += bessel_evaluate2_gsl(bessel_j_large_data);
242 
243    consume_result(result);
244    set_call_count(
245       (sizeof(bessel_j_data)
246       + sizeof(bessel_j_int_data)
247       + sizeof(bessel_j_large_data)) / sizeof(bessel_j_data[0]));
248 }
249 
250 template <std::size_t N>
bessel_y_evaluate2_gsl(const boost::array<boost::array<T,3>,N> & data)251 double bessel_y_evaluate2_gsl(const boost::array<boost::array<T, 3>, N>& data)
252 {
253    double result = 0;
254    for(unsigned i = 0; i < N; ++i)
255       result += gsl_sf_bessel_Ynu(data[i][0], data[i][1]);
256    return result;
257 }
258 
259 BOOST_MATH_PERFORMANCE_TEST(bessel_test, "cyl_neumann-gsl")
260 {
261    double result= bessel_y_evaluate2_gsl(bessel_y01_data);
262    result += bessel_y_evaluate2_gsl(bessel_yn_data);
263    result += bessel_y_evaluate2_gsl(bessel_yv_data);
264 
265    consume_result(result);
266    set_call_count(
267       (sizeof(bessel_y01_data)
268       + sizeof(bessel_yn_data)
269       + sizeof(bessel_yv_data)) / sizeof(bessel_j_data[0]));
270 }
271 
272 template <std::size_t N>
bessel_i_evaluate_gsl(const boost::array<boost::array<T,3>,N> & data)273 double bessel_i_evaluate_gsl(const boost::array<boost::array<T, 3>, N>& data)
274 {
275    double result = 0;
276    for(unsigned i = 0; i < N; ++i)
277       result += gsl_sf_bessel_Inu(data[i][0], data[i][1]);
278    return result;
279 }
280 
281 BOOST_MATH_PERFORMANCE_TEST(bessel_test, "cyl_bessel_i-gsl")
282 {
283    double result= bessel_i_evaluate_gsl(bessel_i_int_data);
284    result += bessel_i_evaluate_gsl(bessel_i_data);
285 
286    consume_result(result);
287    set_call_count(
288       (sizeof(bessel_i_data)
289       + sizeof(bessel_i_int_data)) / sizeof(bessel_i_data[0]));
290 }
291 
292 template <std::size_t N>
bessel_k_evaluate_gsl(const boost::array<boost::array<T,3>,N> & data)293 double bessel_k_evaluate_gsl(const boost::array<boost::array<T, 3>, N>& data)
294 {
295    double result = 0;
296    for(unsigned i = 0; i < N; ++i)
297       result += gsl_sf_bessel_Knu(data[i][0], data[i][1]);
298    return result;
299 }
300 
301 BOOST_MATH_PERFORMANCE_TEST(bessel_test, "cyl_bessel_k-gsl")
302 {
303    double result= bessel_k_evaluate_gsl(bessel_k_data);
304    result += bessel_k_evaluate_gsl(bessel_k_int_data);
305 
306    consume_result(result);
307    set_call_count(
308       (sizeof(bessel_k_data)
309       + sizeof(bessel_k_int_data)) / sizeof(bessel_k_data[0]));
310 }
311 
312 #endif
313 
314