| 1 | // Copyright John Maddock 2007, 2014. |
| 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 | #ifndef BOOST_MATH_ZETA_HPP |
| 7 | #define BOOST_MATH_ZETA_HPP |
| 8 | |
| 9 | #ifdef _MSC_VER |
| 10 | #pragma once |
| 11 | #endif |
| 12 | |
| 13 | #include <boost/math/special_functions/math_fwd.hpp> |
| 14 | #include <boost/math/tools/precision.hpp> |
| 15 | #include <boost/math/tools/series.hpp> |
| 16 | #include <boost/math/tools/big_constant.hpp> |
| 17 | #include <boost/math/policies/error_handling.hpp> |
| 18 | #include <boost/math/special_functions/gamma.hpp> |
| 19 | #include <boost/math/special_functions/factorials.hpp> |
| 20 | #include <boost/math/special_functions/sin_pi.hpp> |
| 21 | |
| 22 | #if defined(__GNUC__) && defined(BOOST_MATH_USE_FLOAT128) |
| 23 | // |
| 24 | // This is the only way we can avoid |
| 25 | // warning: non-standard suffix on floating constant [-Wpedantic] |
| 26 | // when building with -Wall -pedantic. Neither __extension__ |
| 27 | // nor #pragma diagnostic ignored work :( |
| 28 | // |
| 29 | #pragma GCC system_header |
| 30 | #endif |
| 31 | |
| 32 | namespace boost{ namespace math{ namespace detail{ |
| 33 | |
| 34 | #if 0 |
| 35 | // |
| 36 | // This code is commented out because we have a better more rapidly converging series |
| 37 | // now. Retained for future reference and in case the new code causes any issues down the line.... |
| 38 | // |
| 39 | |
| 40 | template <class T, class Policy> |
| 41 | struct zeta_series_cache_size |
| 42 | { |
| 43 | // |
| 44 | // Work how large to make our cache size when evaluating the series |
| 45 | // evaluation: normally this is just large enough for the series |
| 46 | // to have converged, but for arbitrary precision types we need a |
| 47 | // really large cache to achieve reasonable precision in a reasonable |
| 48 | // time. This is important when constructing rational approximations |
| 49 | // to zeta for example. |
| 50 | // |
| 51 | typedef typename boost::math::policies::precision<T,Policy>::type precision_type; |
| 52 | typedef typename mpl::if_< |
| 53 | mpl::less_equal<precision_type, boost::integral_constant<int, 0> >, |
| 54 | boost::integral_constant<int, 5000>, |
| 55 | typename mpl::if_< |
| 56 | mpl::less_equal<precision_type, boost::integral_constant<int, 64> >, |
| 57 | boost::integral_constant<int, 70>, |
| 58 | typename mpl::if_< |
| 59 | mpl::less_equal<precision_type, boost::integral_constant<int, 113> >, |
| 60 | boost::integral_constant<int, 100>, |
| 61 | boost::integral_constant<int, 5000> |
| 62 | >::type |
| 63 | >::type |
| 64 | >::type type; |
| 65 | }; |
| 66 | |
| 67 | template <class T, class Policy> |
| 68 | T zeta_series_imp(T s, T sc, const Policy&) |
| 69 | { |
| 70 | // |
| 71 | // Series evaluation from: |
| 72 | // Havil, J. Gamma: Exploring Euler's Constant. |
| 73 | // Princeton, NJ: Princeton University Press, 2003. |
| 74 | // |
| 75 | // See also http://mathworld.wolfram.com/RiemannZetaFunction.html |
| 76 | // |
| 77 | BOOST_MATH_STD_USING |
| 78 | T sum = 0; |
| 79 | T mult = 0.5; |
| 80 | T change; |
| 81 | typedef typename zeta_series_cache_size<T,Policy>::type cache_size; |
| 82 | T powers[cache_size::value] = { 0, }; |
| 83 | unsigned n = 0; |
| 84 | do{ |
| 85 | T binom = -static_cast<T>(n); |
| 86 | T nested_sum = 1; |
| 87 | if(n < sizeof(powers) / sizeof(powers[0])) |
| 88 | powers[n] = pow(static_cast<T>(n + 1), -s); |
| 89 | for(unsigned k = 1; k <= n; ++k) |
| 90 | { |
| 91 | T p; |
| 92 | if(k < sizeof(powers) / sizeof(powers[0])) |
| 93 | { |
| 94 | p = powers[k]; |
| 95 | //p = pow(k + 1, -s); |
| 96 | } |
| 97 | else |
| 98 | p = pow(static_cast<T>(k + 1), -s); |
| 99 | nested_sum += binom * p; |
| 100 | binom *= (k - static_cast<T>(n)) / (k + 1); |
| 101 | } |
| 102 | change = mult * nested_sum; |
| 103 | sum += change; |
| 104 | mult /= 2; |
| 105 | ++n; |
| 106 | }while(fabs(change / sum) > tools::epsilon<T>()); |
| 107 | |
| 108 | return sum * 1 / -boost::math::powm1(T(2), sc); |
| 109 | } |
| 110 | |
| 111 | // |
| 112 | // Classical p-series: |
| 113 | // |
| 114 | template <class T> |
| 115 | struct zeta_series2 |
| 116 | { |
| 117 | typedef T result_type; |
| 118 | zeta_series2(T _s) : s(-_s), k(1){} |
| 119 | T operator()() |
| 120 | { |
| 121 | BOOST_MATH_STD_USING |
| 122 | return pow(static_cast<T>(k++), s); |
| 123 | } |
| 124 | private: |
| 125 | T s; |
| 126 | unsigned k; |
| 127 | }; |
| 128 | |
| 129 | template <class T, class Policy> |
| 130 | inline T zeta_series2_imp(T s, const Policy& pol) |
| 131 | { |
| 132 | boost::uintmax_t max_iter = policies::get_max_series_iterations<Policy>();; |
| 133 | zeta_series2<T> f(s); |
| 134 | T result = tools::sum_series( |
| 135 | f, |
| 136 | policies::get_epsilon<T, Policy>(), |
| 137 | max_iter); |
| 138 | policies::check_series_iterations<T>("boost::math::zeta_series2<%1%>(%1%)" , max_iter, pol); |
| 139 | return result; |
| 140 | } |
| 141 | #endif |
| 142 | |
| 143 | template <class T, class Policy> |
| 144 | T zeta_polynomial_series(T s, T sc, Policy const &) |
| 145 | { |
| 146 | // |
| 147 | // This is algorithm 3 from: |
| 148 | // |
| 149 | // "An Efficient Algorithm for the Riemann Zeta Function", P. Borwein, |
| 150 | // Canadian Mathematical Society, Conference Proceedings. |
| 151 | // See: http://www.cecm.sfu.ca/personal/pborwein/PAPERS/P155.pdf |
| 152 | // |
| 153 | BOOST_MATH_STD_USING |
| 154 | int n = itrunc(T(log(boost::math::tools::epsilon<T>()) / -2)); |
| 155 | T sum = 0; |
| 156 | T two_n = ldexp(T(1), n); |
| 157 | int ej_sign = 1; |
| 158 | for(int j = 0; j < n; ++j) |
| 159 | { |
| 160 | sum += ej_sign * -two_n / pow(T(j + 1), s); |
| 161 | ej_sign = -ej_sign; |
| 162 | } |
| 163 | T ej_sum = 1; |
| 164 | T ej_term = 1; |
| 165 | for(int j = n; j <= 2 * n - 1; ++j) |
| 166 | { |
| 167 | sum += ej_sign * (ej_sum - two_n) / pow(T(j + 1), s); |
| 168 | ej_sign = -ej_sign; |
| 169 | ej_term *= 2 * n - j; |
| 170 | ej_term /= j - n + 1; |
| 171 | ej_sum += ej_term; |
| 172 | } |
| 173 | return -sum / (two_n * (-powm1(T(2), sc))); |
| 174 | } |
| 175 | |
| 176 | template <class T, class Policy> |
| 177 | T zeta_imp_prec(T s, T sc, const Policy& pol, const boost::integral_constant<int, 0>&) |
| 178 | { |
| 179 | BOOST_MATH_STD_USING |
| 180 | T result; |
| 181 | if(s >= policies::digits<T, Policy>()) |
| 182 | return 1; |
| 183 | result = zeta_polynomial_series(s, sc, pol); |
| 184 | #if 0 |
| 185 | // Old code archived for future reference: |
| 186 | |
| 187 | // |
| 188 | // Only use power series if it will converge in 100 |
| 189 | // iterations or less: the more iterations it consumes |
| 190 | // the slower convergence becomes so we have to be very |
| 191 | // careful in it's usage. |
| 192 | // |
| 193 | if (s > -log(tools::epsilon<T>()) / 4.5) |
| 194 | result = detail::zeta_series2_imp(s, pol); |
| 195 | else |
| 196 | result = detail::zeta_series_imp(s, sc, pol); |
| 197 | #endif |
| 198 | return result; |
| 199 | } |
| 200 | |
| 201 | template <class T, class Policy> |
| 202 | inline T zeta_imp_prec(T s, T sc, const Policy&, const boost::integral_constant<int, 53>&) |
| 203 | { |
| 204 | BOOST_MATH_STD_USING |
| 205 | T result; |
| 206 | if(s < 1) |
| 207 | { |
| 208 | // Rational Approximation |
| 209 | // Maximum Deviation Found: 2.020e-18 |
| 210 | // Expected Error Term: -2.020e-18 |
| 211 | // Max error found at double precision: 3.994987e-17 |
| 212 | static const T P[6] = { |
| 213 | static_cast<T>(0.24339294433593750202L), |
| 214 | static_cast<T>(-0.49092470516353571651L), |
| 215 | static_cast<T>(0.0557616214776046784287L), |
| 216 | static_cast<T>(-0.00320912498879085894856L), |
| 217 | static_cast<T>(0.000451534528645796438704L), |
| 218 | static_cast<T>(-0.933241270357061460782e-5L), |
| 219 | }; |
| 220 | static const T Q[6] = { |
| 221 | static_cast<T>(1L), |
| 222 | static_cast<T>(-0.279960334310344432495L), |
| 223 | static_cast<T>(0.0419676223309986037706L), |
| 224 | static_cast<T>(-0.00413421406552171059003L), |
| 225 | static_cast<T>(0.00024978985622317935355L), |
| 226 | static_cast<T>(-0.101855788418564031874e-4L), |
| 227 | }; |
| 228 | result = tools::evaluate_polynomial(P, sc) / tools::evaluate_polynomial(Q, sc); |
| 229 | result -= 1.2433929443359375F; |
| 230 | result += (sc); |
| 231 | result /= (sc); |
| 232 | } |
| 233 | else if(s <= 2) |
| 234 | { |
| 235 | // Maximum Deviation Found: 9.007e-20 |
| 236 | // Expected Error Term: 9.007e-20 |
| 237 | static const T P[6] = { |
| 238 | static_cast<T>(0.577215664901532860516L), |
| 239 | static_cast<T>(0.243210646940107164097L), |
| 240 | static_cast<T>(0.0417364673988216497593L), |
| 241 | static_cast<T>(0.00390252087072843288378L), |
| 242 | static_cast<T>(0.000249606367151877175456L), |
| 243 | static_cast<T>(0.110108440976732897969e-4L), |
| 244 | }; |
| 245 | static const T Q[6] = { |
| 246 | static_cast<T>(1.0), |
| 247 | static_cast<T>(0.295201277126631761737L), |
| 248 | static_cast<T>(0.043460910607305495864L), |
| 249 | static_cast<T>(0.00434930582085826330659L), |
| 250 | static_cast<T>(0.000255784226140488490982L), |
| 251 | static_cast<T>(0.10991819782396112081e-4L), |
| 252 | }; |
| 253 | result = tools::evaluate_polynomial(P, T(-sc)) / tools::evaluate_polynomial(Q, T(-sc)); |
| 254 | result += 1 / (-sc); |
| 255 | } |
| 256 | else if(s <= 4) |
| 257 | { |
| 258 | // Maximum Deviation Found: 5.946e-22 |
| 259 | // Expected Error Term: -5.946e-22 |
| 260 | static const float Y = 0.6986598968505859375; |
| 261 | static const T P[6] = { |
| 262 | static_cast<T>(-0.0537258300023595030676L), |
| 263 | static_cast<T>(0.0445163473292365591906L), |
| 264 | static_cast<T>(0.0128677673534519952905L), |
| 265 | static_cast<T>(0.00097541770457391752726L), |
| 266 | static_cast<T>(0.769875101573654070925e-4L), |
| 267 | static_cast<T>(0.328032510000383084155e-5L), |
| 268 | }; |
| 269 | static const T Q[7] = { |
| 270 | 1.0f, |
| 271 | static_cast<T>(0.33383194553034051422L), |
| 272 | static_cast<T>(0.0487798431291407621462L), |
| 273 | static_cast<T>(0.00479039708573558490716L), |
| 274 | static_cast<T>(0.000270776703956336357707L), |
| 275 | static_cast<T>(0.106951867532057341359e-4L), |
| 276 | static_cast<T>(0.236276623974978646399e-7L), |
| 277 | }; |
| 278 | result = tools::evaluate_polynomial(P, T(s - 2)) / tools::evaluate_polynomial(Q, T(s - 2)); |
| 279 | result += Y + 1 / (-sc); |
| 280 | } |
| 281 | else if(s <= 7) |
| 282 | { |
| 283 | // Maximum Deviation Found: 2.955e-17 |
| 284 | // Expected Error Term: 2.955e-17 |
| 285 | // Max error found at double precision: 2.009135e-16 |
| 286 | |
| 287 | static const T P[6] = { |
| 288 | static_cast<T>(-2.49710190602259410021L), |
| 289 | static_cast<T>(-2.60013301809475665334L), |
| 290 | static_cast<T>(-0.939260435377109939261L), |
| 291 | static_cast<T>(-0.138448617995741530935L), |
| 292 | static_cast<T>(-0.00701721240549802377623L), |
| 293 | static_cast<T>(-0.229257310594893932383e-4L), |
| 294 | }; |
| 295 | static const T Q[9] = { |
| 296 | 1.0f, |
| 297 | static_cast<T>(0.706039025937745133628L), |
| 298 | static_cast<T>(0.15739599649558626358L), |
| 299 | static_cast<T>(0.0106117950976845084417L), |
| 300 | static_cast<T>(-0.36910273311764618902e-4L), |
| 301 | static_cast<T>(0.493409563927590008943e-5L), |
| 302 | static_cast<T>(-0.234055487025287216506e-6L), |
| 303 | static_cast<T>(0.718833729365459760664e-8L), |
| 304 | static_cast<T>(-0.1129200113474947419e-9L), |
| 305 | }; |
| 306 | result = tools::evaluate_polynomial(P, T(s - 4)) / tools::evaluate_polynomial(Q, T(s - 4)); |
| 307 | result = 1 + exp(result); |
| 308 | } |
| 309 | else if(s < 15) |
| 310 | { |
| 311 | // Maximum Deviation Found: 7.117e-16 |
| 312 | // Expected Error Term: 7.117e-16 |
| 313 | // Max error found at double precision: 9.387771e-16 |
| 314 | static const T P[7] = { |
| 315 | static_cast<T>(-4.78558028495135619286L), |
| 316 | static_cast<T>(-1.89197364881972536382L), |
| 317 | static_cast<T>(-0.211407134874412820099L), |
| 318 | static_cast<T>(-0.000189204758260076688518L), |
| 319 | static_cast<T>(0.00115140923889178742086L), |
| 320 | static_cast<T>(0.639949204213164496988e-4L), |
| 321 | static_cast<T>(0.139348932445324888343e-5L), |
| 322 | }; |
| 323 | static const T Q[9] = { |
| 324 | 1.0f, |
| 325 | static_cast<T>(0.244345337378188557777L), |
| 326 | static_cast<T>(0.00873370754492288653669L), |
| 327 | static_cast<T>(-0.00117592765334434471562L), |
| 328 | static_cast<T>(-0.743743682899933180415e-4L), |
| 329 | static_cast<T>(-0.21750464515767984778e-5L), |
| 330 | static_cast<T>(0.471001264003076486547e-8L), |
| 331 | static_cast<T>(-0.833378440625385520576e-10L), |
| 332 | static_cast<T>(0.699841545204845636531e-12L), |
| 333 | }; |
| 334 | result = tools::evaluate_polynomial(P, T(s - 7)) / tools::evaluate_polynomial(Q, T(s - 7)); |
| 335 | result = 1 + exp(result); |
| 336 | } |
| 337 | else if(s < 36) |
| 338 | { |
| 339 | // Max error in interpolated form: 1.668e-17 |
| 340 | // Max error found at long double precision: 1.669714e-17 |
| 341 | static const T P[8] = { |
| 342 | static_cast<T>(-10.3948950573308896825L), |
| 343 | static_cast<T>(-2.85827219671106697179L), |
| 344 | static_cast<T>(-0.347728266539245787271L), |
| 345 | static_cast<T>(-0.0251156064655346341766L), |
| 346 | static_cast<T>(-0.00119459173416968685689L), |
| 347 | static_cast<T>(-0.382529323507967522614e-4L), |
| 348 | static_cast<T>(-0.785523633796723466968e-6L), |
| 349 | static_cast<T>(-0.821465709095465524192e-8L), |
| 350 | }; |
| 351 | static const T Q[10] = { |
| 352 | 1.0f, |
| 353 | static_cast<T>(0.208196333572671890965L), |
| 354 | static_cast<T>(0.0195687657317205033485L), |
| 355 | static_cast<T>(0.00111079638102485921877L), |
| 356 | static_cast<T>(0.408507746266039256231e-4L), |
| 357 | static_cast<T>(0.955561123065693483991e-6L), |
| 358 | static_cast<T>(0.118507153474022900583e-7L), |
| 359 | static_cast<T>(0.222609483627352615142e-14L), |
| 360 | }; |
| 361 | result = tools::evaluate_polynomial(P, T(s - 15)) / tools::evaluate_polynomial(Q, T(s - 15)); |
| 362 | result = 1 + exp(result); |
| 363 | } |
| 364 | else if(s < 56) |
| 365 | { |
| 366 | result = 1 + pow(T(2), -s); |
| 367 | } |
| 368 | else |
| 369 | { |
| 370 | result = 1; |
| 371 | } |
| 372 | return result; |
| 373 | } |
| 374 | |
| 375 | template <class T, class Policy> |
| 376 | T zeta_imp_prec(T s, T sc, const Policy&, const boost::integral_constant<int, 64>&) |
| 377 | { |
| 378 | BOOST_MATH_STD_USING |
| 379 | T result; |
| 380 | if(s < 1) |
| 381 | { |
| 382 | // Rational Approximation |
| 383 | // Maximum Deviation Found: 3.099e-20 |
| 384 | // Expected Error Term: 3.099e-20 |
| 385 | // Max error found at long double precision: 5.890498e-20 |
| 386 | static const T P[6] = { |
| 387 | BOOST_MATH_BIG_CONSTANT(T, 64, 0.243392944335937499969), |
| 388 | BOOST_MATH_BIG_CONSTANT(T, 64, -0.496837806864865688082), |
| 389 | BOOST_MATH_BIG_CONSTANT(T, 64, 0.0680008039723709987107), |
| 390 | BOOST_MATH_BIG_CONSTANT(T, 64, -0.00511620413006619942112), |
| 391 | BOOST_MATH_BIG_CONSTANT(T, 64, 0.000455369899250053003335), |
| 392 | BOOST_MATH_BIG_CONSTANT(T, 64, -0.279496685273033761927e-4), |
| 393 | }; |
| 394 | static const T Q[7] = { |
| 395 | BOOST_MATH_BIG_CONSTANT(T, 64, 1.0), |
| 396 | BOOST_MATH_BIG_CONSTANT(T, 64, -0.30425480068225790522), |
| 397 | BOOST_MATH_BIG_CONSTANT(T, 64, 0.050052748580371598736), |
| 398 | BOOST_MATH_BIG_CONSTANT(T, 64, -0.00519355671064700627862), |
| 399 | BOOST_MATH_BIG_CONSTANT(T, 64, 0.000360623385771198350257), |
| 400 | BOOST_MATH_BIG_CONSTANT(T, 64, -0.159600883054550987633e-4), |
| 401 | BOOST_MATH_BIG_CONSTANT(T, 64, 0.339770279812410586032e-6), |
| 402 | }; |
| 403 | result = tools::evaluate_polynomial(P, sc) / tools::evaluate_polynomial(Q, sc); |
| 404 | result -= 1.2433929443359375F; |
| 405 | result += (sc); |
| 406 | result /= (sc); |
| 407 | } |
| 408 | else if(s <= 2) |
| 409 | { |
| 410 | // Maximum Deviation Found: 1.059e-21 |
| 411 | // Expected Error Term: 1.059e-21 |
| 412 | // Max error found at long double precision: 1.626303e-19 |
| 413 | |
| 414 | static const T P[6] = { |
| 415 | BOOST_MATH_BIG_CONSTANT(T, 64, 0.577215664901532860605), |
| 416 | BOOST_MATH_BIG_CONSTANT(T, 64, 0.222537368917162139445), |
| 417 | BOOST_MATH_BIG_CONSTANT(T, 64, 0.0356286324033215682729), |
| 418 | BOOST_MATH_BIG_CONSTANT(T, 64, 0.00304465292366350081446), |
| 419 | BOOST_MATH_BIG_CONSTANT(T, 64, 0.000178102511649069421904), |
| 420 | BOOST_MATH_BIG_CONSTANT(T, 64, 0.700867470265983665042e-5), |
| 421 | }; |
| 422 | static const T Q[7] = { |
| 423 | BOOST_MATH_BIG_CONSTANT(T, 64, 1.0), |
| 424 | BOOST_MATH_BIG_CONSTANT(T, 64, 0.259385759149531030085), |
| 425 | BOOST_MATH_BIG_CONSTANT(T, 64, 0.0373974962106091316854), |
| 426 | BOOST_MATH_BIG_CONSTANT(T, 64, 0.00332735159183332820617), |
| 427 | BOOST_MATH_BIG_CONSTANT(T, 64, 0.000188690420706998606469), |
| 428 | BOOST_MATH_BIG_CONSTANT(T, 64, 0.635994377921861930071e-5), |
| 429 | BOOST_MATH_BIG_CONSTANT(T, 64, 0.226583954978371199405e-7), |
| 430 | }; |
| 431 | result = tools::evaluate_polynomial(P, T(-sc)) / tools::evaluate_polynomial(Q, T(-sc)); |
| 432 | result += 1 / (-sc); |
| 433 | } |
| 434 | else if(s <= 4) |
| 435 | { |
| 436 | // Maximum Deviation Found: 5.946e-22 |
| 437 | // Expected Error Term: -5.946e-22 |
| 438 | static const float Y = 0.6986598968505859375; |
| 439 | static const T P[7] = { |
| 440 | BOOST_MATH_BIG_CONSTANT(T, 64, -0.053725830002359501027), |
| 441 | BOOST_MATH_BIG_CONSTANT(T, 64, 0.0470551187571475844778), |
| 442 | BOOST_MATH_BIG_CONSTANT(T, 64, 0.0101339410415759517471), |
| 443 | BOOST_MATH_BIG_CONSTANT(T, 64, 0.00100240326666092854528), |
| 444 | BOOST_MATH_BIG_CONSTANT(T, 64, 0.685027119098122814867e-4), |
| 445 | BOOST_MATH_BIG_CONSTANT(T, 64, 0.390972820219765942117e-5), |
| 446 | BOOST_MATH_BIG_CONSTANT(T, 64, 0.540319769113543934483e-7), |
| 447 | }; |
| 448 | static const T Q[8] = { |
| 449 | BOOST_MATH_BIG_CONSTANT(T, 64, 1.0), |
| 450 | BOOST_MATH_BIG_CONSTANT(T, 64, 0.286577739726542730421), |
| 451 | BOOST_MATH_BIG_CONSTANT(T, 64, 0.0447355811517733225843), |
| 452 | BOOST_MATH_BIG_CONSTANT(T, 64, 0.00430125107610252363302), |
| 453 | BOOST_MATH_BIG_CONSTANT(T, 64, 0.000284956969089786662045), |
| 454 | BOOST_MATH_BIG_CONSTANT(T, 64, 0.116188101609848411329e-4), |
| 455 | BOOST_MATH_BIG_CONSTANT(T, 64, 0.278090318191657278204e-6), |
| 456 | BOOST_MATH_BIG_CONSTANT(T, 64, -0.19683620233222028478e-8), |
| 457 | }; |
| 458 | result = tools::evaluate_polynomial(P, T(s - 2)) / tools::evaluate_polynomial(Q, T(s - 2)); |
| 459 | result += Y + 1 / (-sc); |
| 460 | } |
| 461 | else if(s <= 7) |
| 462 | { |
| 463 | // Max error found at long double precision: 8.132216e-19 |
| 464 | static const T P[8] = { |
| 465 | BOOST_MATH_BIG_CONSTANT(T, 64, -2.49710190602259407065), |
| 466 | BOOST_MATH_BIG_CONSTANT(T, 64, -3.36664913245960625334), |
| 467 | BOOST_MATH_BIG_CONSTANT(T, 64, -1.77180020623777595452), |
| 468 | BOOST_MATH_BIG_CONSTANT(T, 64, -0.464717885249654313933), |
| 469 | BOOST_MATH_BIG_CONSTANT(T, 64, -0.0643694921293579472583), |
| 470 | BOOST_MATH_BIG_CONSTANT(T, 64, -0.00464265386202805715487), |
| 471 | BOOST_MATH_BIG_CONSTANT(T, 64, -0.000165556579779704340166), |
| 472 | BOOST_MATH_BIG_CONSTANT(T, 64, -0.252884970740994069582e-5), |
| 473 | }; |
| 474 | static const T Q[9] = { |
| 475 | BOOST_MATH_BIG_CONSTANT(T, 64, 1.0), |
| 476 | BOOST_MATH_BIG_CONSTANT(T, 64, 1.01300131390690459085), |
| 477 | BOOST_MATH_BIG_CONSTANT(T, 64, 0.387898115758643503827), |
| 478 | BOOST_MATH_BIG_CONSTANT(T, 64, 0.0695071490045701135188), |
| 479 | BOOST_MATH_BIG_CONSTANT(T, 64, 0.00586908595251442839291), |
| 480 | BOOST_MATH_BIG_CONSTANT(T, 64, 0.000217752974064612188616), |
| 481 | BOOST_MATH_BIG_CONSTANT(T, 64, 0.397626583349419011731e-5), |
| 482 | BOOST_MATH_BIG_CONSTANT(T, 64, -0.927884739284359700764e-8), |
| 483 | BOOST_MATH_BIG_CONSTANT(T, 64, 0.119810501805618894381e-9), |
| 484 | }; |
| 485 | result = tools::evaluate_polynomial(P, T(s - 4)) / tools::evaluate_polynomial(Q, T(s - 4)); |
| 486 | result = 1 + exp(result); |
| 487 | } |
| 488 | else if(s < 15) |
| 489 | { |
| 490 | // Max error in interpolated form: 1.133e-18 |
| 491 | // Max error found at long double precision: 2.183198e-18 |
| 492 | static const T P[9] = { |
| 493 | BOOST_MATH_BIG_CONSTANT(T, 64, -4.78558028495135548083), |
| 494 | BOOST_MATH_BIG_CONSTANT(T, 64, -3.23873322238609358947), |
| 495 | BOOST_MATH_BIG_CONSTANT(T, 64, -0.892338582881021799922), |
| 496 | BOOST_MATH_BIG_CONSTANT(T, 64, -0.131326296217965913809), |
| 497 | BOOST_MATH_BIG_CONSTANT(T, 64, -0.0115651591773783712996), |
| 498 | BOOST_MATH_BIG_CONSTANT(T, 64, -0.000657728968362695775205), |
| 499 | BOOST_MATH_BIG_CONSTANT(T, 64, -0.252051328129449973047e-4), |
| 500 | BOOST_MATH_BIG_CONSTANT(T, 64, -0.626503445372641798925e-6), |
| 501 | BOOST_MATH_BIG_CONSTANT(T, 64, -0.815696314790853893484e-8), |
| 502 | }; |
| 503 | static const T Q[9] = { |
| 504 | BOOST_MATH_BIG_CONSTANT(T, 64, 1.0), |
| 505 | BOOST_MATH_BIG_CONSTANT(T, 64, 0.525765665400123515036), |
| 506 | BOOST_MATH_BIG_CONSTANT(T, 64, 0.10852641753657122787), |
| 507 | BOOST_MATH_BIG_CONSTANT(T, 64, 0.0115669945375362045249), |
| 508 | BOOST_MATH_BIG_CONSTANT(T, 64, 0.000732896513858274091966), |
| 509 | BOOST_MATH_BIG_CONSTANT(T, 64, 0.30683952282420248448e-4), |
| 510 | BOOST_MATH_BIG_CONSTANT(T, 64, 0.819649214609633126119e-6), |
| 511 | BOOST_MATH_BIG_CONSTANT(T, 64, 0.117957556472335968146e-7), |
| 512 | BOOST_MATH_BIG_CONSTANT(T, 64, -0.193432300973017671137e-12), |
| 513 | }; |
| 514 | result = tools::evaluate_polynomial(P, T(s - 7)) / tools::evaluate_polynomial(Q, T(s - 7)); |
| 515 | result = 1 + exp(result); |
| 516 | } |
| 517 | else if(s < 42) |
| 518 | { |
| 519 | // Max error in interpolated form: 1.668e-17 |
| 520 | // Max error found at long double precision: 1.669714e-17 |
| 521 | static const T P[9] = { |
| 522 | BOOST_MATH_BIG_CONSTANT(T, 64, -10.3948950573308861781), |
| 523 | BOOST_MATH_BIG_CONSTANT(T, 64, -2.82646012777913950108), |
| 524 | BOOST_MATH_BIG_CONSTANT(T, 64, -0.342144362739570333665), |
| 525 | BOOST_MATH_BIG_CONSTANT(T, 64, -0.0249285145498722647472), |
| 526 | BOOST_MATH_BIG_CONSTANT(T, 64, -0.00122493108848097114118), |
| 527 | BOOST_MATH_BIG_CONSTANT(T, 64, -0.423055371192592850196e-4), |
| 528 | BOOST_MATH_BIG_CONSTANT(T, 64, -0.1025215577185967488e-5), |
| 529 | BOOST_MATH_BIG_CONSTANT(T, 64, -0.165096762663509467061e-7), |
| 530 | BOOST_MATH_BIG_CONSTANT(T, 64, -0.145392555873022044329e-9), |
| 531 | }; |
| 532 | static const T Q[10] = { |
| 533 | BOOST_MATH_BIG_CONSTANT(T, 64, 1.0), |
| 534 | BOOST_MATH_BIG_CONSTANT(T, 64, 0.205135978585281988052), |
| 535 | BOOST_MATH_BIG_CONSTANT(T, 64, 0.0192359357875879453602), |
| 536 | BOOST_MATH_BIG_CONSTANT(T, 64, 0.00111496452029715514119), |
| 537 | BOOST_MATH_BIG_CONSTANT(T, 64, 0.434928449016693986857e-4), |
| 538 | BOOST_MATH_BIG_CONSTANT(T, 64, 0.116911068726610725891e-5), |
| 539 | BOOST_MATH_BIG_CONSTANT(T, 64, 0.206704342290235237475e-7), |
| 540 | BOOST_MATH_BIG_CONSTANT(T, 64, 0.209772836100827647474e-9), |
| 541 | BOOST_MATH_BIG_CONSTANT(T, 64, -0.939798249922234703384e-16), |
| 542 | BOOST_MATH_BIG_CONSTANT(T, 64, 0.264584017421245080294e-18), |
| 543 | }; |
| 544 | result = tools::evaluate_polynomial(P, T(s - 15)) / tools::evaluate_polynomial(Q, T(s - 15)); |
| 545 | result = 1 + exp(result); |
| 546 | } |
| 547 | else if(s < 63) |
| 548 | { |
| 549 | result = 1 + pow(T(2), -s); |
| 550 | } |
| 551 | else |
| 552 | { |
| 553 | result = 1; |
| 554 | } |
| 555 | return result; |
| 556 | } |
| 557 | |
| 558 | template <class T, class Policy> |
| 559 | T zeta_imp_prec(T s, T sc, const Policy&, const boost::integral_constant<int, 113>&) |
| 560 | { |
| 561 | BOOST_MATH_STD_USING |
| 562 | T result; |
| 563 | if(s < 1) |
| 564 | { |
| 565 | // Rational Approximation |
| 566 | // Maximum Deviation Found: 9.493e-37 |
| 567 | // Expected Error Term: 9.492e-37 |
| 568 | // Max error found at long double precision: 7.281332e-31 |
| 569 | |
| 570 | static const T P[10] = { |
| 571 | BOOST_MATH_BIG_CONSTANT(T, 113, -1.0), |
| 572 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.0353008629988648122808504280990313668), |
| 573 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.0107795651204927743049369868548706909), |
| 574 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.000523961870530500751114866884685172975), |
| 575 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.661805838304910731947595897966487515e-4), |
| 576 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.658932670403818558510656304189164638e-5), |
| 577 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.103437265642266106533814021041010453e-6), |
| 578 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.116818787212666457105375746642927737e-7), |
| 579 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.660690993901506912123512551294239036e-9), |
| 580 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.113103113698388531428914333768142527e-10), |
| 581 | }; |
| 582 | static const T Q[11] = { |
| 583 | BOOST_MATH_BIG_CONSTANT(T, 113, 1.0), |
| 584 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.387483472099602327112637481818565459), |
| 585 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.0802265315091063135271497708694776875), |
| 586 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.0110727276164171919280036408995078164), |
| 587 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.00112552716946286252000434849173787243), |
| 588 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.874554160748626916455655180296834352e-4), |
| 589 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.530097847491828379568636739662278322e-5), |
| 590 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.248461553590496154705565904497247452e-6), |
| 591 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.881834921354014787309644951507523899e-8), |
| 592 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.217062446168217797598596496310953025e-9), |
| 593 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.315823200002384492377987848307151168e-11), |
| 594 | }; |
| 595 | result = tools::evaluate_polynomial(P, sc) / tools::evaluate_polynomial(Q, sc); |
| 596 | result += (sc); |
| 597 | result /= (sc); |
| 598 | } |
| 599 | else if(s <= 2) |
| 600 | { |
| 601 | // Maximum Deviation Found: 1.616e-37 |
| 602 | // Expected Error Term: -1.615e-37 |
| 603 | |
| 604 | static const T P[10] = { |
| 605 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.577215664901532860606512090082402431), |
| 606 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.255597968739771510415479842335906308), |
| 607 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.0494056503552807274142218876983542205), |
| 608 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.00551372778611700965268920983472292325), |
| 609 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.00043667616723970574871427830895192731), |
| 610 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.268562259154821957743669387915239528e-4), |
| 611 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.109249633923016310141743084480436612e-5), |
| 612 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.273895554345300227466534378753023924e-7), |
| 613 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.583103205551702720149237384027795038e-9), |
| 614 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.835774625259919268768735944711219256e-11), |
| 615 | }; |
| 616 | static const T Q[11] = { |
| 617 | BOOST_MATH_BIG_CONSTANT(T, 113, 1.0), |
| 618 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.316661751179735502065583176348292881), |
| 619 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.0540401806533507064453851182728635272), |
| 620 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.00598621274107420237785899476374043797), |
| 621 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.000474907812321704156213038740142079615), |
| 622 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.272125421722314389581695715835862418e-4), |
| 623 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.112649552156479800925522445229212933e-5), |
| 624 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.301838975502992622733000078063330461e-7), |
| 625 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.422960728687211282539769943184270106e-9), |
| 626 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.377105263588822468076813329270698909e-11), |
| 627 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.581926559304525152432462127383600681e-13), |
| 628 | }; |
| 629 | result = tools::evaluate_polynomial(P, T(-sc)) / tools::evaluate_polynomial(Q, T(-sc)); |
| 630 | result += 1 / (-sc); |
| 631 | } |
| 632 | else if(s <= 4) |
| 633 | { |
| 634 | // Maximum Deviation Found: 1.891e-36 |
| 635 | // Expected Error Term: -1.891e-36 |
| 636 | // Max error found: 2.171527e-35 |
| 637 | |
| 638 | static const float Y = 0.6986598968505859375; |
| 639 | static const T P[11] = { |
| 640 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.0537258300023595010275848333539748089), |
| 641 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.0429086930802630159457448174466342553), |
| 642 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.0136148228754303412510213395034056857), |
| 643 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.00190231601036042925183751238033763915), |
| 644 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.000186880390916311438818302549192456581), |
| 645 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.145347370745893262394287982691323657e-4), |
| 646 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.805843276446813106414036600485884885e-6), |
| 647 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.340818159286739137503297172091882574e-7), |
| 648 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.115762357488748996526167305116837246e-8), |
| 649 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.231904754577648077579913403645767214e-10), |
| 650 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.340169592866058506675897646629036044e-12), |
| 651 | }; |
| 652 | static const T Q[12] = { |
| 653 | BOOST_MATH_BIG_CONSTANT(T, 113, 1.0), |
| 654 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.363755247765087100018556983050520554), |
| 655 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.0696581979014242539385695131258321598), |
| 656 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.00882208914484611029571547753782014817), |
| 657 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.000815405623261946661762236085660996718), |
| 658 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.571366167062457197282642344940445452e-4), |
| 659 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.309278269271853502353954062051797838e-5), |
| 660 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.12822982083479010834070516053794262e-6), |
| 661 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.397876357325018976733953479182110033e-8), |
| 662 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.8484432107648683277598472295289279e-10), |
| 663 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.105677416606909614301995218444080615e-11), |
| 664 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.547223964564003701979951154093005354e-15), |
| 665 | }; |
| 666 | result = tools::evaluate_polynomial(P, T(s - 2)) / tools::evaluate_polynomial(Q, T(s - 2)); |
| 667 | result += Y + 1 / (-sc); |
| 668 | } |
| 669 | else if(s <= 6) |
| 670 | { |
| 671 | // Max error in interpolated form: 1.510e-37 |
| 672 | // Max error found at long double precision: 2.769266e-34 |
| 673 | |
| 674 | static const T Y = 3.28348541259765625F; |
| 675 | |
| 676 | static const T P[13] = { |
| 677 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.786383506575062179339611614117697622), |
| 678 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.495766593395271370974685959652073976), |
| 679 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.409116737851754766422360889037532228), |
| 680 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.57340744006238263817895456842655987), |
| 681 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.280479899797421910694892949057963111), |
| 682 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.0753148409447590257157585696212649869), |
| 683 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.0122934003684672788499099362823748632), |
| 684 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.00126148398446193639247961370266962927), |
| 685 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.828465038179772939844657040917364896e-4), |
| 686 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.361008916706050977143208468690645684e-5), |
| 687 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.109879825497910544424797771195928112e-6), |
| 688 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.214539416789686920918063075528797059e-8), |
| 689 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.15090220092460596872172844424267351e-10), |
| 690 | }; |
| 691 | static const T Q[14] = { |
| 692 | BOOST_MATH_BIG_CONSTANT(T, 113, 1.0), |
| 693 | BOOST_MATH_BIG_CONSTANT(T, 113, 1.69490865837142338462982225731926485), |
| 694 | BOOST_MATH_BIG_CONSTANT(T, 113, 1.22697696630994080733321401255942464), |
| 695 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.495409420862526540074366618006341533), |
| 696 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.122368084916843823462872905024259633), |
| 697 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.0191412993625268971656513890888208623), |
| 698 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.00191401538628980617753082598351559642), |
| 699 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.000123318142456272424148930280876444459), |
| 700 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.531945488232526067889835342277595709e-5), |
| 701 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.161843184071894368337068779669116236e-6), |
| 702 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.305796079600152506743828859577462778e-8), |
| 703 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.233582592298450202680170811044408894e-10), |
| 704 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.275363878344548055574209713637734269e-13), |
| 705 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.221564186807357535475441900517843892e-15), |
| 706 | }; |
| 707 | result = tools::evaluate_polynomial(P, T(s - 4)) / tools::evaluate_polynomial(Q, T(s - 4)); |
| 708 | result -= Y; |
| 709 | result = 1 + exp(result); |
| 710 | } |
| 711 | else if(s < 10) |
| 712 | { |
| 713 | // Max error in interpolated form: 1.999e-34 |
| 714 | // Max error found at long double precision: 2.156186e-33 |
| 715 | |
| 716 | static const T P[13] = { |
| 717 | BOOST_MATH_BIG_CONSTANT(T, 113, -4.0545627381873738086704293881227365), |
| 718 | BOOST_MATH_BIG_CONSTANT(T, 113, -4.70088348734699134347906176097717782), |
| 719 | BOOST_MATH_BIG_CONSTANT(T, 113, -2.36921550900925512951976617607678789), |
| 720 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.684322583796369508367726293719322866), |
| 721 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.126026534540165129870721937592996324), |
| 722 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.015636903921778316147260572008619549), |
| 723 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.00135442294754728549644376325814460807), |
| 724 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.842793965853572134365031384646117061e-4), |
| 725 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.385602133791111663372015460784978351e-5), |
| 726 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.130458500394692067189883214401478539e-6), |
| 727 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.315861074947230418778143153383660035e-8), |
| 728 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.500334720512030826996373077844707164e-10), |
| 729 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.420204769185233365849253969097184005e-12), |
| 730 | }; |
| 731 | static const T Q[14] = { |
| 732 | BOOST_MATH_BIG_CONSTANT(T, 113, 1.0), |
| 733 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.97663511666410096104783358493318814), |
| 734 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.40878780231201806504987368939673249), |
| 735 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.0963890666609396058945084107597727252), |
| 736 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.0142207619090854604824116070866614505), |
| 737 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.00139010220902667918476773423995750877), |
| 738 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.940669540194694997889636696089994734e-4), |
| 739 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.458220848507517004399292480807026602e-5), |
| 740 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.16345521617741789012782420625435495e-6), |
| 741 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.414007452533083304371566316901024114e-8), |
| 742 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.68701473543366328016953742622661377e-10), |
| 743 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.603461891080716585087883971886075863e-12), |
| 744 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.294670713571839023181857795866134957e-16), |
| 745 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.147003914536437243143096875069813451e-18), |
| 746 | }; |
| 747 | result = tools::evaluate_polynomial(P, T(s - 6)) / tools::evaluate_polynomial(Q, T(s - 6)); |
| 748 | result = 1 + exp(result); |
| 749 | } |
| 750 | else if(s < 17) |
| 751 | { |
| 752 | // Max error in interpolated form: 1.641e-32 |
| 753 | // Max error found at long double precision: 1.696121e-32 |
| 754 | static const T P[13] = { |
| 755 | BOOST_MATH_BIG_CONSTANT(T, 113, -6.91319491921722925920883787894829678), |
| 756 | BOOST_MATH_BIG_CONSTANT(T, 113, -3.65491257639481960248690596951049048), |
| 757 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.813557553449954526442644544105257881), |
| 758 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.0994317301685870959473658713841138083), |
| 759 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.00726896610245676520248617014211734906), |
| 760 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.000317253318715075854811266230916762929), |
| 761 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.66851422826636750855184211580127133e-5), |
| 762 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.879464154730985406003332577806849971e-7), |
| 763 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.113838903158254250631678791998294628e-7), |
| 764 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.379184410304927316385211327537817583e-9), |
| 765 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.612992858643904887150527613446403867e-11), |
| 766 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.347873737198164757035457841688594788e-13), |
| 767 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.289187187441625868404494665572279364e-15), |
| 768 | }; |
| 769 | static const T Q[14] = { |
| 770 | BOOST_MATH_BIG_CONSTANT(T, 113, 1.0), |
| 771 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.427310044448071818775721584949868806), |
| 772 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.074602514873055756201435421385243062), |
| 773 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.00688651562174480772901425121653945942), |
| 774 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.000360174847635115036351323894321880445), |
| 775 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.973556847713307543918865405758248777e-5), |
| 776 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.853455848314516117964634714780874197e-8), |
| 777 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.118203513654855112421673192194622826e-7), |
| 778 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.462521662511754117095006543363328159e-9), |
| 779 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.834212591919475633107355719369463143e-11), |
| 780 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.5354594751002702935740220218582929e-13), |
| 781 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.406451690742991192964889603000756203e-15), |
| 782 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.887948682401000153828241615760146728e-19), |
| 783 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.34980761098820347103967203948619072e-21), |
| 784 | }; |
| 785 | result = tools::evaluate_polynomial(P, T(s - 10)) / tools::evaluate_polynomial(Q, T(s - 10)); |
| 786 | result = 1 + exp(result); |
| 787 | } |
| 788 | else if(s < 30) |
| 789 | { |
| 790 | // Max error in interpolated form: 1.563e-31 |
| 791 | // Max error found at long double precision: 1.562725e-31 |
| 792 | |
| 793 | static const T P[13] = { |
| 794 | BOOST_MATH_BIG_CONSTANT(T, 113, -11.7824798233959252791987402769438322), |
| 795 | BOOST_MATH_BIG_CONSTANT(T, 113, -4.36131215284987731928174218354118102), |
| 796 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.732260980060982349410898496846972204), |
| 797 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.0744985185694913074484248803015717388), |
| 798 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.00517228281320594683022294996292250527), |
| 799 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.000260897206152101522569969046299309939), |
| 800 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.989553462123121764865178453128769948e-5), |
| 801 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.286916799741891410827712096608826167e-6), |
| 802 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.637262477796046963617949532211619729e-8), |
| 803 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.106796831465628373325491288787760494e-9), |
| 804 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.129343095511091870860498356205376823e-11), |
| 805 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.102397936697965977221267881716672084e-13), |
| 806 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.402663128248642002351627980255756363e-16), |
| 807 | }; |
| 808 | static const T Q[14] = { |
| 809 | BOOST_MATH_BIG_CONSTANT(T, 113, 1.0), |
| 810 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.311288325355705609096155335186466508), |
| 811 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.0438318468940415543546769437752132748), |
| 812 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.00374396349183199548610264222242269536), |
| 813 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.000218707451200585197339671707189281302), |
| 814 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.927578767487930747532953583797351219e-5), |
| 815 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.294145760625753561951137473484889639e-6), |
| 816 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.704618586690874460082739479535985395e-8), |
| 817 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.126333332872897336219649130062221257e-9), |
| 818 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.16317315713773503718315435769352765e-11), |
| 819 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.137846712823719515148344938160275695e-13), |
| 820 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.580975420554224366450994232723910583e-16), |
| 821 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.291354445847552426900293580511392459e-22), |
| 822 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.73614324724785855925025452085443636e-25), |
| 823 | }; |
| 824 | result = tools::evaluate_polynomial(P, T(s - 17)) / tools::evaluate_polynomial(Q, T(s - 17)); |
| 825 | result = 1 + exp(result); |
| 826 | } |
| 827 | else if(s < 74) |
| 828 | { |
| 829 | // Max error in interpolated form: 2.311e-27 |
| 830 | // Max error found at long double precision: 2.297544e-27 |
| 831 | static const T P[14] = { |
| 832 | BOOST_MATH_BIG_CONSTANT(T, 113, -20.7944102007844314586649688802236072), |
| 833 | BOOST_MATH_BIG_CONSTANT(T, 113, -4.95759941987499442499908748130192187), |
| 834 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.563290752832461751889194629200298688), |
| 835 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.0406197001137935911912457120706122877), |
| 836 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.0020846534789473022216888863613422293), |
| 837 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.808095978462109173749395599401375667e-4), |
| 838 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.244706022206249301640890603610060959e-5), |
| 839 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.589477682919645930544382616501666572e-7), |
| 840 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.113699573675553496343617442433027672e-8), |
| 841 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.174767860183598149649901223128011828e-10), |
| 842 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.210051620306761367764549971980026474e-12), |
| 843 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.189187969537370950337212675466400599e-14), |
| 844 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.116313253429564048145641663778121898e-16), |
| 845 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.376708747782400769427057630528578187e-19), |
| 846 | }; |
| 847 | static const T Q[16] = { |
| 848 | BOOST_MATH_BIG_CONSTANT(T, 113, 1.0), |
| 849 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.205076752981410805177554569784219717), |
| 850 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.0202526722696670378999575738524540269), |
| 851 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.001278305290005994980069466658219057), |
| 852 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.576404779858501791742255670403304787e-4), |
| 853 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.196477049872253010859712483984252067e-5), |
| 854 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.521863830500876189501054079974475762e-7), |
| 855 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.109524209196868135198775445228552059e-8), |
| 856 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.181698713448644481083966260949267825e-10), |
| 857 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.234793316975091282090312036524695562e-12), |
| 858 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.227490441461460571047545264251399048e-14), |
| 859 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.151500292036937400913870642638520668e-16), |
| 860 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.543475775154780935815530649335936121e-19), |
| 861 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.241647013434111434636554455083309352e-28), |
| 862 | BOOST_MATH_BIG_CONSTANT(T, 113, -0.557103423021951053707162364713587374e-31), |
| 863 | BOOST_MATH_BIG_CONSTANT(T, 113, 0.618708773442584843384712258199645166e-34), |
| 864 | }; |
| 865 | result = tools::evaluate_polynomial(P, T(s - 30)) / tools::evaluate_polynomial(Q, T(s - 30)); |
| 866 | result = 1 + exp(result); |
| 867 | } |
| 868 | else if(s < 117) |
| 869 | { |
| 870 | result = 1 + pow(T(2), -s); |
| 871 | } |
| 872 | else |
| 873 | { |
| 874 | result = 1; |
| 875 | } |
| 876 | return result; |
| 877 | } |
| 878 | |
| 879 | template <class T, class Policy> |
| 880 | T zeta_imp_odd_integer(int s, const T&, const Policy&, const true_type&) |
| 881 | { |
| 882 | static const T results[] = { |
| 883 | BOOST_MATH_BIG_CONSTANT(T, 113, 1.2020569031595942853997381615114500), BOOST_MATH_BIG_CONSTANT(T, 113, 1.0369277551433699263313654864570342), BOOST_MATH_BIG_CONSTANT(T, 113, 1.0083492773819228268397975498497968), BOOST_MATH_BIG_CONSTANT(T, 113, 1.0020083928260822144178527692324121), BOOST_MATH_BIG_CONSTANT(T, 113, 1.0004941886041194645587022825264699), BOOST_MATH_BIG_CONSTANT(T, 113, 1.0001227133475784891467518365263574), BOOST_MATH_BIG_CONSTANT(T, 113, 1.0000305882363070204935517285106451), BOOST_MATH_BIG_CONSTANT(T, 113, 1.0000076371976378997622736002935630), BOOST_MATH_BIG_CONSTANT(T, 113, 1.0000019082127165539389256569577951), BOOST_MATH_BIG_CONSTANT(T, 113, 1.0000004769329867878064631167196044), BOOST_MATH_BIG_CONSTANT(T, 113, 1.0000001192199259653110730677887189), BOOST_MATH_BIG_CONSTANT(T, 113, 1.0000000298035035146522801860637051), BOOST_MATH_BIG_CONSTANT(T, 113, 1.0000000074507117898354294919810042), BOOST_MATH_BIG_CONSTANT(T, 113, 1.0000000018626597235130490064039099), BOOST_MATH_BIG_CONSTANT(T, 113, 1.0000000004656629065033784072989233), BOOST_MATH_BIG_CONSTANT(T, 113, 1.0000000001164155017270051977592974), BOOST_MATH_BIG_CONSTANT(T, 113, 1.0000000000291038504449709968692943), BOOST_MATH_BIG_CONSTANT(T, 113, 1.0000000000072759598350574810145209), BOOST_MATH_BIG_CONSTANT(T, 113, 1.0000000000018189896503070659475848), BOOST_MATH_BIG_CONSTANT(T, 113, 1.0000000000004547473783042154026799), BOOST_MATH_BIG_CONSTANT(T, 113, 1.0000000000001136868407680227849349), BOOST_MATH_BIG_CONSTANT(T, 113, 1.0000000000000284217097688930185546), BOOST_MATH_BIG_CONSTANT(T, 113, 1.0000000000000071054273952108527129), BOOST_MATH_BIG_CONSTANT(T, 113, 1.0000000000000017763568435791203275), BOOST_MATH_BIG_CONSTANT(T, 113, 1.0000000000000004440892103143813364), BOOST_MATH_BIG_CONSTANT(T, 113, 1.0000000000000001110223025141066134), BOOST_MATH_BIG_CONSTANT(T, 113, 1.0000000000000000277555756213612417), BOOST_MATH_BIG_CONSTANT(T, 113, 1.0000000000000000069388939045441537), BOOST_MATH_BIG_CONSTANT(T, 113, 1.0000000000000000017347234760475766), BOOST_MATH_BIG_CONSTANT(T, 113, 1.0000000000000000004336808690020650), BOOST_MATH_BIG_CONSTANT(T, 113, 1.0000000000000000001084202172494241), BOOST_MATH_BIG_CONSTANT(T, 113, 1.0000000000000000000271050543122347), BOOST_MATH_BIG_CONSTANT(T, 113, 1.0000000000000000000067762635780452), BOOST_MATH_BIG_CONSTANT(T, 113, 1.0000000000000000000016940658945098), BOOST_MATH_BIG_CONSTANT(T, 113, 1.0000000000000000000004235164736273), BOOST_MATH_BIG_CONSTANT(T, 113, 1.0000000000000000000001058791184068), BOOST_MATH_BIG_CONSTANT(T, 113, 1.0000000000000000000000264697796017), BOOST_MATH_BIG_CONSTANT(T, 113, 1.0000000000000000000000066174449004), BOOST_MATH_BIG_CONSTANT(T, 113, 1.0000000000000000000000016543612251), BOOST_MATH_BIG_CONSTANT(T, 113, 1.0000000000000000000000004135903063), BOOST_MATH_BIG_CONSTANT(T, 113, 1.0000000000000000000000001033975766), BOOST_MATH_BIG_CONSTANT(T, 113, 1.0000000000000000000000000258493941), BOOST_MATH_BIG_CONSTANT(T, 113, 1.0000000000000000000000000064623485), BOOST_MATH_BIG_CONSTANT(T, 113, 1.0000000000000000000000000016155871), BOOST_MATH_BIG_CONSTANT(T, 113, 1.0000000000000000000000000004038968), BOOST_MATH_BIG_CONSTANT(T, 113, 1.0000000000000000000000000001009742), BOOST_MATH_BIG_CONSTANT(T, 113, 1.0000000000000000000000000000252435), BOOST_MATH_BIG_CONSTANT(T, 113, 1.0000000000000000000000000000063109), BOOST_MATH_BIG_CONSTANT(T, 113, 1.0000000000000000000000000000015777), BOOST_MATH_BIG_CONSTANT(T, 113, 1.0000000000000000000000000000003944), BOOST_MATH_BIG_CONSTANT(T, 113, 1.0000000000000000000000000000000986), BOOST_MATH_BIG_CONSTANT(T, 113, 1.0000000000000000000000000000000247), BOOST_MATH_BIG_CONSTANT(T, 113, 1.0000000000000000000000000000000062), BOOST_MATH_BIG_CONSTANT(T, 113, 1.0000000000000000000000000000000015), BOOST_MATH_BIG_CONSTANT(T, 113, 1.0000000000000000000000000000000004), BOOST_MATH_BIG_CONSTANT(T, 113, 1.0000000000000000000000000000000001), |
| 884 | }; |
| 885 | return s > 113 ? 1 : results[(s - 3) / 2]; |
| 886 | } |
| 887 | |
| 888 | template <class T, class Policy> |
| 889 | T zeta_imp_odd_integer(int s, const T& sc, const Policy& pol, const false_type&) |
| 890 | { |
| 891 | static BOOST_MATH_THREAD_LOCAL bool is_init = false; |
| 892 | static BOOST_MATH_THREAD_LOCAL T results[50] = {}; |
| 893 | static BOOST_MATH_THREAD_LOCAL int digits = tools::digits<T>(); |
| 894 | int current_digits = tools::digits<T>(); |
| 895 | if(digits != current_digits) |
| 896 | { |
| 897 | // Oh my precision has changed... |
| 898 | is_init = false; |
| 899 | } |
| 900 | if(!is_init) |
| 901 | { |
| 902 | is_init = true; |
| 903 | digits = current_digits; |
| 904 | for(unsigned k = 0; k < sizeof(results) / sizeof(results[0]); ++k) |
| 905 | { |
| 906 | T arg = k * 2 + 3; |
| 907 | T c_arg = 1 - arg; |
| 908 | results[k] = zeta_polynomial_series(arg, c_arg, pol); |
| 909 | } |
| 910 | } |
| 911 | unsigned index = (s - 3) / 2; |
| 912 | return index >= sizeof(results) / sizeof(results[0]) ? zeta_polynomial_series(T(s), sc, pol): results[index]; |
| 913 | } |
| 914 | |
| 915 | template <class T, class Policy, class Tag> |
| 916 | T zeta_imp(T s, T sc, const Policy& pol, const Tag& tag) |
| 917 | { |
| 918 | BOOST_MATH_STD_USING |
| 919 | static const char* function = "boost::math::zeta<%1%>" ; |
| 920 | if(sc == 0) |
| 921 | return policies::raise_pole_error<T>( |
| 922 | function, |
| 923 | "Evaluation of zeta function at pole %1%" , |
| 924 | s, pol); |
| 925 | T result; |
| 926 | // |
| 927 | // Trivial case: |
| 928 | // |
| 929 | if(s > policies::digits<T, Policy>()) |
| 930 | return 1; |
| 931 | // |
| 932 | // Start by seeing if we have a simple closed form: |
| 933 | // |
| 934 | if(floor(s) == s) |
| 935 | { |
| 936 | #ifndef BOOST_NO_EXCEPTIONS |
| 937 | // Without exceptions we expect itrunc to return INT_MAX on overflow |
| 938 | // and we fall through anyway. |
| 939 | try |
| 940 | { |
| 941 | #endif |
| 942 | int v = itrunc(s); |
| 943 | if(v == s) |
| 944 | { |
| 945 | if(v < 0) |
| 946 | { |
| 947 | if(((-v) & 1) == 0) |
| 948 | return 0; |
| 949 | int n = (-v + 1) / 2; |
| 950 | if(n <= (int)boost::math::max_bernoulli_b2n<T>::value) |
| 951 | return T((-v & 1) ? -1 : 1) * boost::math::unchecked_bernoulli_b2n<T>(n) / (1 - v); |
| 952 | } |
| 953 | else if((v & 1) == 0) |
| 954 | { |
| 955 | if(((v / 2) <= (int)boost::math::max_bernoulli_b2n<T>::value) && (v <= (int)boost::math::max_factorial<T>::value)) |
| 956 | return T(((v / 2 - 1) & 1) ? -1 : 1) * ldexp(T(1), v - 1) * pow(constants::pi<T, Policy>(), v) * |
| 957 | boost::math::unchecked_bernoulli_b2n<T>(v / 2) / boost::math::unchecked_factorial<T>(v); |
| 958 | return T(((v / 2 - 1) & 1) ? -1 : 1) * ldexp(T(1), v - 1) * pow(constants::pi<T, Policy>(), v) * |
| 959 | boost::math::bernoulli_b2n<T>(v / 2) / boost::math::factorial<T>(v); |
| 960 | } |
| 961 | else |
| 962 | return zeta_imp_odd_integer(v, sc, pol, boost::integral_constant<bool, (Tag::value <= 113) && Tag::value>()); |
| 963 | } |
| 964 | #ifndef BOOST_NO_EXCEPTIONS |
| 965 | } |
| 966 | catch(const boost::math::rounding_error&){} // Just fall through, s is too large to round |
| 967 | catch(const std::overflow_error&){} |
| 968 | #endif |
| 969 | } |
| 970 | |
| 971 | if(fabs(s) < tools::root_epsilon<T>()) |
| 972 | { |
| 973 | result = -0.5f - constants::log_root_two_pi<T, Policy>() * s; |
| 974 | } |
| 975 | else if(s < 0) |
| 976 | { |
| 977 | std::swap(s, sc); |
| 978 | if(floor(sc/2) == sc/2) |
| 979 | result = 0; |
| 980 | else |
| 981 | { |
| 982 | if(s > max_factorial<T>::value) |
| 983 | { |
| 984 | T mult = boost::math::sin_pi(0.5f * sc, pol) * 2 * zeta_imp(s, sc, pol, tag); |
| 985 | result = boost::math::lgamma(s, pol); |
| 986 | result -= s * log(2 * constants::pi<T>()); |
| 987 | if(result > tools::log_max_value<T>()) |
| 988 | return sign(mult) * policies::raise_overflow_error<T>(function, 0, pol); |
| 989 | result = exp(result); |
| 990 | if(tools::max_value<T>() / fabs(mult) < result) |
| 991 | return boost::math::sign(mult) * policies::raise_overflow_error<T>(function, 0, pol); |
| 992 | result *= mult; |
| 993 | } |
| 994 | else |
| 995 | { |
| 996 | result = boost::math::sin_pi(0.5f * sc, pol) |
| 997 | * 2 * pow(2 * constants::pi<T>(), -s) |
| 998 | * boost::math::tgamma(s, pol) |
| 999 | * zeta_imp(s, sc, pol, tag); |
| 1000 | } |
| 1001 | } |
| 1002 | } |
| 1003 | else |
| 1004 | { |
| 1005 | result = zeta_imp_prec(s, sc, pol, tag); |
| 1006 | } |
| 1007 | return result; |
| 1008 | } |
| 1009 | |
| 1010 | template <class T, class Policy, class tag> |
| 1011 | struct zeta_initializer |
| 1012 | { |
| 1013 | struct init |
| 1014 | { |
| 1015 | init() |
| 1016 | { |
| 1017 | do_init(tag()); |
| 1018 | } |
| 1019 | static void do_init(const boost::integral_constant<int, 0>&){ boost::math::zeta(static_cast<T>(5), Policy()); } |
| 1020 | static void do_init(const boost::integral_constant<int, 53>&){ boost::math::zeta(static_cast<T>(5), Policy()); } |
| 1021 | static void do_init(const boost::integral_constant<int, 64>&) |
| 1022 | { |
| 1023 | boost::math::zeta(static_cast<T>(0.5), Policy()); |
| 1024 | boost::math::zeta(static_cast<T>(1.5), Policy()); |
| 1025 | boost::math::zeta(static_cast<T>(3.5), Policy()); |
| 1026 | boost::math::zeta(static_cast<T>(6.5), Policy()); |
| 1027 | boost::math::zeta(static_cast<T>(14.5), Policy()); |
| 1028 | boost::math::zeta(static_cast<T>(40.5), Policy()); |
| 1029 | |
| 1030 | boost::math::zeta(static_cast<T>(5), Policy()); |
| 1031 | } |
| 1032 | static void do_init(const boost::integral_constant<int, 113>&) |
| 1033 | { |
| 1034 | boost::math::zeta(static_cast<T>(0.5), Policy()); |
| 1035 | boost::math::zeta(static_cast<T>(1.5), Policy()); |
| 1036 | boost::math::zeta(static_cast<T>(3.5), Policy()); |
| 1037 | boost::math::zeta(static_cast<T>(5.5), Policy()); |
| 1038 | boost::math::zeta(static_cast<T>(9.5), Policy()); |
| 1039 | boost::math::zeta(static_cast<T>(16.5), Policy()); |
| 1040 | boost::math::zeta(static_cast<T>(25.5), Policy()); |
| 1041 | boost::math::zeta(static_cast<T>(70.5), Policy()); |
| 1042 | |
| 1043 | boost::math::zeta(static_cast<T>(5), Policy()); |
| 1044 | } |
| 1045 | void force_instantiate()const{} |
| 1046 | }; |
| 1047 | static const init initializer; |
| 1048 | static void force_instantiate() |
| 1049 | { |
| 1050 | initializer.force_instantiate(); |
| 1051 | } |
| 1052 | }; |
| 1053 | |
| 1054 | template <class T, class Policy, class tag> |
| 1055 | const typename zeta_initializer<T, Policy, tag>::init zeta_initializer<T, Policy, tag>::initializer; |
| 1056 | |
| 1057 | } // detail |
| 1058 | |
| 1059 | template <class T, class Policy> |
| 1060 | inline typename tools::promote_args<T>::type zeta(T s, const Policy&) |
| 1061 | { |
| 1062 | typedef typename tools::promote_args<T>::type result_type; |
| 1063 | typedef typename policies::evaluation<result_type, Policy>::type value_type; |
| 1064 | typedef typename policies::precision<result_type, Policy>::type precision_type; |
| 1065 | typedef typename policies::normalise< |
| 1066 | Policy, |
| 1067 | policies::promote_float<false>, |
| 1068 | policies::promote_double<false>, |
| 1069 | policies::discrete_quantile<>, |
| 1070 | policies::assert_undefined<> >::type forwarding_policy; |
| 1071 | typedef boost::integral_constant<int, |
| 1072 | precision_type::value <= 0 ? 0 : |
| 1073 | precision_type::value <= 53 ? 53 : |
| 1074 | precision_type::value <= 64 ? 64 : |
| 1075 | precision_type::value <= 113 ? 113 : 0 |
| 1076 | > tag_type; |
| 1077 | |
| 1078 | detail::zeta_initializer<value_type, forwarding_policy, tag_type>::force_instantiate(); |
| 1079 | |
| 1080 | return policies::checked_narrowing_cast<result_type, forwarding_policy>(detail::zeta_imp( |
| 1081 | static_cast<value_type>(s), |
| 1082 | static_cast<value_type>(1 - static_cast<value_type>(s)), |
| 1083 | forwarding_policy(), |
| 1084 | tag_type()), "boost::math::zeta<%1%>(%1%)" ); |
| 1085 | } |
| 1086 | |
| 1087 | template <class T> |
| 1088 | inline typename tools::promote_args<T>::type zeta(T s) |
| 1089 | { |
| 1090 | return zeta(s, policies::policy<>()); |
| 1091 | } |
| 1092 | |
| 1093 | }} // namespaces |
| 1094 | |
| 1095 | #endif // BOOST_MATH_ZETA_HPP |
| 1096 | |
| 1097 | |
| 1098 | |
| 1099 | |