owens_t.hpp 49 KB

1234567891011121314151617181920212223242526272829303132333435363738394041424344454647484950515253545556575859606162636465666768697071727374757677787980818283848586878889909192939495969798991001011021031041051061071081091101111121131141151161171181191201211221231241251261271281291301311321331341351361371381391401411421431441451461471481491501511521531541551561571581591601611621631641651661671681691701711721731741751761771781791801811821831841851861871881891901911921931941951961971981992002012022032042052062072082092102112122132142152162172182192202212222232242252262272282292302312322332342352362372382392402412422432442452462472482492502512522532542552562572582592602612622632642652662672682692702712722732742752762772782792802812822832842852862872882892902912922932942952962972982993003013023033043053063073083093103113123133143153163173183193203213223233243253263273283293303313323333343353363373383393403413423433443453463473483493503513523533543553563573583593603613623633643653663673683693703713723733743753763773783793803813823833843853863873883893903913923933943953963973983994004014024034044054064074084094104114124134144154164174184194204214224234244254264274284294304314324334344354364374384394404414424434444454464474484494504514524534544554564574584594604614624634644654664674684694704714724734744754764774784794804814824834844854864874884894904914924934944954964974984995005015025035045055065075085095105115125135145155165175185195205215225235245255265275285295305315325335345355365375385395405415425435445455465475485495505515525535545555565575585595605615625635645655665675685695705715725735745755765775785795805815825835845855865875885895905915925935945955965975985996006016026036046056066076086096106116126136146156166176186196206216226236246256266276286296306316326336346356366376386396406416426436446456466476486496506516526536546556566576586596606616626636646656666676686696706716726736746756766776786796806816826836846856866876886896906916926936946956966976986997007017027037047057067077087097107117127137147157167177187197207217227237247257267277287297307317327337347357367377387397407417427437447457467477487497507517527537547557567577587597607617627637647657667677687697707717727737747757767777787797807817827837847857867877887897907917927937947957967977987998008018028038048058068078088098108118128138148158168178188198208218228238248258268278288298308318328338348358368378388398408418428438448458468478488498508518528538548558568578588598608618628638648658668678688698708718728738748758768778788798808818828838848858868878888898908918928938948958968978988999009019029039049059069079089099109119129139149159169179189199209219229239249259269279289299309319329339349359369379389399409419429439449459469479489499509519529539549559569579589599609619629639649659669679689699709719729739749759769779789799809819829839849859869879889899909919929939949959969979989991000100110021003100410051006100710081009101010111012101310141015101610171018101910201021102210231024102510261027102810291030103110321033103410351036103710381039104010411042104310441045104610471048104910501051105210531054105510561057105810591060106110621063106410651066106710681069107010711072107310741075107610771078107910801081108210831084108510861087108810891090109110921093109410951096109710981099110011011102110311041105110611071108110911101111
  1. // Copyright Benjamin Sobotta 2012
  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. #ifndef BOOST_OWENS_T_HPP
  6. #define BOOST_OWENS_T_HPP
  7. // Reference:
  8. // Mike Patefield, David Tandy
  9. // FAST AND ACCURATE CALCULATION OF OWEN'S T-FUNCTION
  10. // Journal of Statistical Software, 5 (5), 1-25
  11. #ifdef _MSC_VER
  12. # pragma once
  13. #endif
  14. #include <boost/math/special_functions/math_fwd.hpp>
  15. #include <boost/config/no_tr1/cmath.hpp>
  16. #include <boost/math/special_functions/erf.hpp>
  17. #include <boost/math/special_functions/expm1.hpp>
  18. #include <boost/throw_exception.hpp>
  19. #include <boost/assert.hpp>
  20. #include <boost/math/constants/constants.hpp>
  21. #include <boost/math/tools/big_constant.hpp>
  22. #include <stdexcept>
  23. #ifdef BOOST_MSVC
  24. #pragma warning(push)
  25. #pragma warning(disable:4127)
  26. #endif
  27. #if defined(__GNUC__) && defined(BOOST_MATH_USE_FLOAT128)
  28. //
  29. // This is the only way we can avoid
  30. // warning: non-standard suffix on floating constant [-Wpedantic]
  31. // when building with -Wall -pedantic. Neither __extension__
  32. // nor #pragma dianostic ignored work :(
  33. //
  34. #pragma GCC system_header
  35. #endif
  36. namespace boost
  37. {
  38. namespace math
  39. {
  40. namespace detail
  41. {
  42. // owens_t_znorm1(x) = P(-oo<Z<=x)-0.5 with Z being normally distributed.
  43. template<typename RealType>
  44. inline RealType owens_t_znorm1(const RealType x)
  45. {
  46. using namespace boost::math::constants;
  47. return erf(x*one_div_root_two<RealType>())*half<RealType>();
  48. } // RealType owens_t_znorm1(const RealType x)
  49. // owens_t_znorm2(x) = P(x<=Z<oo) with Z being normally distributed.
  50. template<typename RealType>
  51. inline RealType owens_t_znorm2(const RealType x)
  52. {
  53. using namespace boost::math::constants;
  54. return erfc(x*one_div_root_two<RealType>())*half<RealType>();
  55. } // RealType owens_t_znorm2(const RealType x)
  56. // Auxiliary function, it computes an array key that is used to determine
  57. // the specific computation method for Owen's T and the order thereof
  58. // used in owens_t_dispatch.
  59. template<typename RealType>
  60. inline unsigned short owens_t_compute_code(const RealType h, const RealType a)
  61. {
  62. static const RealType hrange[] =
  63. { 0.02f, 0.06f, 0.09f, 0.125f, 0.26f, 0.4f, 0.6f, 1.6f, 1.7f, 2.33f, 2.4f, 3.36f, 3.4f, 4.8f };
  64. static const RealType arange[] = { 0.025f, 0.09f, 0.15f, 0.36f, 0.5f, 0.9f, 0.99999f };
  65. /*
  66. original select array from paper:
  67. 1, 1, 2,13,13,13,13,13,13,13,13,16,16,16, 9
  68. 1, 2, 2, 3, 3, 5, 5,14,14,15,15,16,16,16, 9
  69. 2, 2, 3, 3, 3, 5, 5,15,15,15,15,16,16,16,10
  70. 2, 2, 3, 5, 5, 5, 5, 7, 7,16,16,16,16,16,10
  71. 2, 3, 3, 5, 5, 6, 6, 8, 8,17,17,17,12,12,11
  72. 2, 3, 5, 5, 5, 6, 6, 8, 8,17,17,17,12,12,12
  73. 2, 3, 4, 4, 6, 6, 8, 8,17,17,17,17,17,12,12
  74. 2, 3, 4, 4, 6, 6,18,18,18,18,17,17,17,12,12
  75. */
  76. // subtract one because the array is written in FORTRAN in mind - in C arrays start @ zero
  77. static const unsigned short select[] =
  78. {
  79. 0, 0 , 1 , 12 ,12 , 12 , 12 , 12 , 12 , 12 , 12 , 15 , 15 , 15 , 8,
  80. 0 , 1 , 1 , 2 , 2 , 4 , 4 , 13 , 13 , 14 , 14 , 15 , 15 , 15 , 8,
  81. 1 , 1 , 2 , 2 , 2 , 4 , 4 , 14 , 14 , 14 , 14 , 15 , 15 , 15 , 9,
  82. 1 , 1 , 2 , 4 , 4 , 4 , 4 , 6 , 6 , 15 , 15 , 15 , 15 , 15 , 9,
  83. 1 , 2 , 2 , 4 , 4 , 5 , 5 , 7 , 7 , 16 ,16 , 16 , 11 , 11 , 10,
  84. 1 , 2 , 4 , 4 , 4 , 5 , 5 , 7 , 7 , 16 , 16 , 16 , 11 , 11 , 11,
  85. 1 , 2 , 3 , 3 , 5 , 5 , 7 , 7 , 16 , 16 , 16 , 16 , 16 , 11 , 11,
  86. 1 , 2 , 3 , 3 , 5 , 5 , 17 , 17 , 17 , 17 , 16 , 16 , 16 , 11 , 11
  87. };
  88. unsigned short ihint = 14, iaint = 7;
  89. for(unsigned short i = 0; i != 14; i++)
  90. {
  91. if( h <= hrange[i] )
  92. {
  93. ihint = i;
  94. break;
  95. }
  96. } // for(unsigned short i = 0; i != 14; i++)
  97. for(unsigned short i = 0; i != 7; i++)
  98. {
  99. if( a <= arange[i] )
  100. {
  101. iaint = i;
  102. break;
  103. }
  104. } // for(unsigned short i = 0; i != 7; i++)
  105. // interprete select array as 8x15 matrix
  106. return select[iaint*15 + ihint];
  107. } // unsigned short owens_t_compute_code(const RealType h, const RealType a)
  108. template<typename RealType>
  109. inline unsigned short owens_t_get_order_imp(const unsigned short icode, RealType, const mpl::int_<53>&)
  110. {
  111. static const unsigned short ord[] = {2, 3, 4, 5, 7, 10, 12, 18, 10, 20, 30, 0, 4, 7, 8, 20, 0, 0}; // 18 entries
  112. BOOST_ASSERT(icode<18);
  113. return ord[icode];
  114. } // unsigned short owens_t_get_order(const unsigned short icode, RealType, mpl::int<53> const&)
  115. template<typename RealType>
  116. inline unsigned short owens_t_get_order_imp(const unsigned short icode, RealType, const mpl::int_<64>&)
  117. {
  118. // method ================>>> {1, 1, 1, 1, 1, 1, 1, 1, 2, 2, 2, 3, 4, 4, 4, 4, 5, 6}
  119. static const unsigned short ord[] = {3, 4, 5, 6, 8, 11, 13, 19, 10, 20, 30, 0, 7, 10, 11, 23, 0, 0}; // 18 entries
  120. BOOST_ASSERT(icode<18);
  121. return ord[icode];
  122. } // unsigned short owens_t_get_order(const unsigned short icode, RealType, mpl::int<64> const&)
  123. template<typename RealType, typename Policy>
  124. inline unsigned short owens_t_get_order(const unsigned short icode, RealType r, const Policy&)
  125. {
  126. typedef typename policies::precision<RealType, Policy>::type precision_type;
  127. typedef typename mpl::if_<
  128. mpl::or_<
  129. mpl::less_equal<precision_type, mpl::int_<0> >,
  130. mpl::greater<precision_type, mpl::int_<53> >
  131. >,
  132. mpl::int_<64>,
  133. mpl::int_<53>
  134. >::type tag_type;
  135. return owens_t_get_order_imp(icode, r, tag_type());
  136. }
  137. // compute the value of Owen's T function with method T1 from the reference paper
  138. template<typename RealType, typename Policy>
  139. inline RealType owens_t_T1(const RealType h, const RealType a, const unsigned short m, const Policy& pol)
  140. {
  141. BOOST_MATH_STD_USING
  142. using namespace boost::math::constants;
  143. const RealType hs = -h*h*half<RealType>();
  144. const RealType dhs = exp( hs );
  145. const RealType as = a*a;
  146. unsigned short j=1;
  147. RealType jj = 1;
  148. RealType aj = a * one_div_two_pi<RealType>();
  149. RealType dj = boost::math::expm1( hs, pol);
  150. RealType gj = hs*dhs;
  151. RealType val = atan( a ) * one_div_two_pi<RealType>();
  152. while( true )
  153. {
  154. val += dj*aj/jj;
  155. if( m <= j )
  156. break;
  157. j++;
  158. jj += static_cast<RealType>(2);
  159. aj *= as;
  160. dj = gj - dj;
  161. gj *= hs / static_cast<RealType>(j);
  162. } // while( true )
  163. return val;
  164. } // RealType owens_t_T1(const RealType h, const RealType a, const unsigned short m)
  165. // compute the value of Owen's T function with method T2 from the reference paper
  166. template<typename RealType, class Policy>
  167. inline RealType owens_t_T2(const RealType h, const RealType a, const unsigned short m, const RealType ah, const Policy&, const mpl::false_&)
  168. {
  169. BOOST_MATH_STD_USING
  170. using namespace boost::math::constants;
  171. const unsigned short maxii = m+m+1;
  172. const RealType hs = h*h;
  173. const RealType as = -a*a;
  174. const RealType y = static_cast<RealType>(1) / hs;
  175. unsigned short ii = 1;
  176. RealType val = 0;
  177. RealType vi = a * exp( -ah*ah*half<RealType>() ) * one_div_root_two_pi<RealType>();
  178. RealType z = owens_t_znorm1(ah)/h;
  179. while( true )
  180. {
  181. val += z;
  182. if( maxii <= ii )
  183. {
  184. val *= exp( -hs*half<RealType>() ) * one_div_root_two_pi<RealType>();
  185. break;
  186. } // if( maxii <= ii )
  187. z = y * ( vi - static_cast<RealType>(ii) * z );
  188. vi *= as;
  189. ii += 2;
  190. } // while( true )
  191. return val;
  192. } // RealType owens_t_T2(const RealType h, const RealType a, const unsigned short m, const RealType ah)
  193. // compute the value of Owen's T function with method T3 from the reference paper
  194. template<typename RealType>
  195. inline RealType owens_t_T3_imp(const RealType h, const RealType a, const RealType ah, const mpl::int_<53>&)
  196. {
  197. BOOST_MATH_STD_USING
  198. using namespace boost::math::constants;
  199. const unsigned short m = 20;
  200. static const RealType c2[] =
  201. {
  202. static_cast<RealType>(0.99999999999999987510),
  203. static_cast<RealType>(-0.99999999999988796462), static_cast<RealType>(0.99999999998290743652),
  204. static_cast<RealType>(-0.99999999896282500134), static_cast<RealType>(0.99999996660459362918),
  205. static_cast<RealType>(-0.99999933986272476760), static_cast<RealType>(0.99999125611136965852),
  206. static_cast<RealType>(-0.99991777624463387686), static_cast<RealType>(0.99942835555870132569),
  207. static_cast<RealType>(-0.99697311720723000295), static_cast<RealType>(0.98751448037275303682),
  208. static_cast<RealType>(-0.95915857980572882813), static_cast<RealType>(0.89246305511006708555),
  209. static_cast<RealType>(-0.76893425990463999675), static_cast<RealType>(0.58893528468484693250),
  210. static_cast<RealType>(-0.38380345160440256652), static_cast<RealType>(0.20317601701045299653),
  211. static_cast<RealType>(-0.82813631607004984866E-01), static_cast<RealType>(0.24167984735759576523E-01),
  212. static_cast<RealType>(-0.44676566663971825242E-02), static_cast<RealType>(0.39141169402373836468E-03)
  213. };
  214. const RealType as = a*a;
  215. const RealType hs = h*h;
  216. const RealType y = static_cast<RealType>(1)/hs;
  217. RealType ii = 1;
  218. unsigned short i = 0;
  219. RealType vi = a * exp( -ah*ah*half<RealType>() ) * one_div_root_two_pi<RealType>();
  220. RealType zi = owens_t_znorm1(ah)/h;
  221. RealType val = 0;
  222. while( true )
  223. {
  224. BOOST_ASSERT(i < 21);
  225. val += zi*c2[i];
  226. if( m <= i ) // if( m < i+1 )
  227. {
  228. val *= exp( -hs*half<RealType>() ) * one_div_root_two_pi<RealType>();
  229. break;
  230. } // if( m < i )
  231. zi = y * (ii*zi - vi);
  232. vi *= as;
  233. ii += 2;
  234. i++;
  235. } // while( true )
  236. return val;
  237. } // RealType owens_t_T3(const RealType h, const RealType a, const RealType ah)
  238. // compute the value of Owen's T function with method T3 from the reference paper
  239. template<class RealType>
  240. inline RealType owens_t_T3_imp(const RealType h, const RealType a, const RealType ah, const mpl::int_<64>&)
  241. {
  242. BOOST_MATH_STD_USING
  243. using namespace boost::math::constants;
  244. const unsigned short m = 30;
  245. static const RealType c2[] =
  246. {
  247. BOOST_MATH_BIG_CONSTANT(RealType, 260, 0.99999999999999999999999729978162447266851932041876728736094298092917625009873),
  248. BOOST_MATH_BIG_CONSTANT(RealType, 260, -0.99999999999999999999467056379678391810626533251885323416799874878563998732905968),
  249. BOOST_MATH_BIG_CONSTANT(RealType, 260, 0.99999999999999999824849349313270659391127814689133077036298754586814091034842536),
  250. BOOST_MATH_BIG_CONSTANT(RealType, 260, -0.9999999999999997703859616213643405880166422891953033591551179153879839440241685),
  251. BOOST_MATH_BIG_CONSTANT(RealType, 260, 0.99999999999998394883415238173334565554173013941245103172035286759201504179038147),
  252. BOOST_MATH_BIG_CONSTANT(RealType, 260, -0.9999999999993063616095509371081203145247992197457263066869044528823599399470977),
  253. BOOST_MATH_BIG_CONSTANT(RealType, 260, 0.9999999999797336340409464429599229870590160411238245275855903767652432017766116267),
  254. BOOST_MATH_BIG_CONSTANT(RealType, 260, -0.999999999574958412069046680119051639753412378037565521359444170241346845522403274),
  255. BOOST_MATH_BIG_CONSTANT(RealType, 260, 0.9999999933226234193375324943920160947158239076786103108097456617750134812033362048),
  256. BOOST_MATH_BIG_CONSTANT(RealType, 260, -0.9999999188923242461073033481053037468263536806742737922476636768006622772762168467),
  257. BOOST_MATH_BIG_CONSTANT(RealType, 260, 0.9999992195143483674402853783549420883055129680082932629160081128947764415749728967),
  258. BOOST_MATH_BIG_CONSTANT(RealType, 260, -0.999993935137206712830997921913316971472227199741857386575097250553105958772041501),
  259. BOOST_MATH_BIG_CONSTANT(RealType, 260, 0.99996135597690552745362392866517133091672395614263398912807169603795088421057688716),
  260. BOOST_MATH_BIG_CONSTANT(RealType, 260, -0.99979556366513946026406788969630293820987757758641211293079784585126692672425362469),
  261. BOOST_MATH_BIG_CONSTANT(RealType, 260, 0.999092789629617100153486251423850590051366661947344315423226082520411961968929483),
  262. BOOST_MATH_BIG_CONSTANT(RealType, 260, -0.996593837411918202119308620432614600338157335862888580671450938858935084316004769854),
  263. BOOST_MATH_BIG_CONSTANT(RealType, 260, 0.98910017138386127038463510314625339359073956513420458166238478926511821146316469589567),
  264. BOOST_MATH_BIG_CONSTANT(RealType, 260, -0.970078558040693314521331982203762771512160168582494513347846407314584943870399016019),
  265. BOOST_MATH_BIG_CONSTANT(RealType, 260, 0.92911438683263187495758525500033707204091967947532160289872782771388170647150321633673),
  266. BOOST_MATH_BIG_CONSTANT(RealType, 260, -0.8542058695956156057286980736842905011429254735181323743367879525470479126968822863),
  267. BOOST_MATH_BIG_CONSTANT(RealType, 260, 0.73796526033030091233118357742803709382964420335559408722681794195743240930748630755),
  268. BOOST_MATH_BIG_CONSTANT(RealType, 260, -0.58523469882837394570128599003785154144164680587615878645171632791404210655891158),
  269. BOOST_MATH_BIG_CONSTANT(RealType, 260, 0.415997776145676306165661663581868460503874205343014196580122174949645271353372263),
  270. BOOST_MATH_BIG_CONSTANT(RealType, 260, -0.2588210875241943574388730510317252236407805082485246378222935376279663808416534365),
  271. BOOST_MATH_BIG_CONSTANT(RealType, 260, 0.1375535825163892648504646951500265585055789019410617565727090346559210218472356689),
  272. BOOST_MATH_BIG_CONSTANT(RealType, 260, -0.0607952766325955730493900985022020434830339794955745989150270485056436844239206648),
  273. BOOST_MATH_BIG_CONSTANT(RealType, 260, 0.0216337683299871528059836483840390514275488679530797294557060229266785853764115),
  274. BOOST_MATH_BIG_CONSTANT(RealType, 260, -0.00593405693455186729876995814181203900550014220428843483927218267309209471516256),
  275. BOOST_MATH_BIG_CONSTANT(RealType, 260, 0.0011743414818332946510474576182739210553333860106811865963485870668929503649964142),
  276. BOOST_MATH_BIG_CONSTANT(RealType, 260, -1.489155613350368934073453260689881330166342484405529981510694514036264969925132e-4),
  277. BOOST_MATH_BIG_CONSTANT(RealType, 260, 9.072354320794357587710929507988814669454281514268844884841547607134260303118208e-6)
  278. };
  279. const RealType as = a*a;
  280. const RealType hs = h*h;
  281. const RealType y = 1 / hs;
  282. RealType ii = 1;
  283. unsigned short i = 0;
  284. RealType vi = a * exp( -ah*ah*half<RealType>() ) * one_div_root_two_pi<RealType>();
  285. RealType zi = owens_t_znorm1(ah)/h;
  286. RealType val = 0;
  287. while( true )
  288. {
  289. BOOST_ASSERT(i < 31);
  290. val += zi*c2[i];
  291. if( m <= i ) // if( m < i+1 )
  292. {
  293. val *= exp( -hs*half<RealType>() ) * one_div_root_two_pi<RealType>();
  294. break;
  295. } // if( m < i )
  296. zi = y * (ii*zi - vi);
  297. vi *= as;
  298. ii += 2;
  299. i++;
  300. } // while( true )
  301. return val;
  302. } // RealType owens_t_T3(const RealType h, const RealType a, const RealType ah)
  303. template<class RealType, class Policy>
  304. inline RealType owens_t_T3(const RealType h, const RealType a, const RealType ah, const Policy&)
  305. {
  306. typedef typename policies::precision<RealType, Policy>::type precision_type;
  307. typedef typename mpl::if_<
  308. mpl::or_<
  309. mpl::less_equal<precision_type, mpl::int_<0> >,
  310. mpl::greater<precision_type, mpl::int_<53> >
  311. >,
  312. mpl::int_<64>,
  313. mpl::int_<53>
  314. >::type tag_type;
  315. return owens_t_T3_imp(h, a, ah, tag_type());
  316. }
  317. // compute the value of Owen's T function with method T4 from the reference paper
  318. template<typename RealType>
  319. inline RealType owens_t_T4(const RealType h, const RealType a, const unsigned short m)
  320. {
  321. BOOST_MATH_STD_USING
  322. using namespace boost::math::constants;
  323. const unsigned short maxii = m+m+1;
  324. const RealType hs = h*h;
  325. const RealType as = -a*a;
  326. unsigned short ii = 1;
  327. RealType ai = a * exp( -hs*(static_cast<RealType>(1)-as)*half<RealType>() ) * one_div_two_pi<RealType>();
  328. RealType yi = 1;
  329. RealType val = 0;
  330. while( true )
  331. {
  332. val += ai*yi;
  333. if( maxii <= ii )
  334. break;
  335. ii += 2;
  336. yi = (static_cast<RealType>(1)-hs*yi) / static_cast<RealType>(ii);
  337. ai *= as;
  338. } // while( true )
  339. return val;
  340. } // RealType owens_t_T4(const RealType h, const RealType a, const unsigned short m)
  341. // compute the value of Owen's T function with method T5 from the reference paper
  342. template<typename RealType>
  343. inline RealType owens_t_T5_imp(const RealType h, const RealType a, const mpl::int_<53>&)
  344. {
  345. BOOST_MATH_STD_USING
  346. /*
  347. NOTICE:
  348. - The pts[] array contains the squares (!) of the abscissas, i.e. the roots of the Legendre
  349. polynomial P_n(x), instead of the plain roots as required in Gauss-Legendre
  350. quadrature, because T5(h,a,m) contains only x^2 terms.
  351. - The wts[] array contains the weights for Gauss-Legendre quadrature scaled with a factor
  352. of 1/(2*pi) according to T5(h,a,m).
  353. */
  354. const unsigned short m = 13;
  355. static const RealType pts[] = {
  356. static_cast<RealType>(0.35082039676451715489E-02),
  357. static_cast<RealType>(0.31279042338030753740E-01), static_cast<RealType>(0.85266826283219451090E-01),
  358. static_cast<RealType>(0.16245071730812277011), static_cast<RealType>(0.25851196049125434828),
  359. static_cast<RealType>(0.36807553840697533536), static_cast<RealType>(0.48501092905604697475),
  360. static_cast<RealType>(0.60277514152618576821), static_cast<RealType>(0.71477884217753226516),
  361. static_cast<RealType>(0.81475510988760098605), static_cast<RealType>(0.89711029755948965867),
  362. static_cast<RealType>(0.95723808085944261843), static_cast<RealType>(0.99178832974629703586) };
  363. static const RealType wts[] = {
  364. static_cast<RealType>(0.18831438115323502887E-01),
  365. static_cast<RealType>(0.18567086243977649478E-01), static_cast<RealType>(0.18042093461223385584E-01),
  366. static_cast<RealType>(0.17263829606398753364E-01), static_cast<RealType>(0.16243219975989856730E-01),
  367. static_cast<RealType>(0.14994592034116704829E-01), static_cast<RealType>(0.13535474469662088392E-01),
  368. static_cast<RealType>(0.11886351605820165233E-01), static_cast<RealType>(0.10070377242777431897E-01),
  369. static_cast<RealType>(0.81130545742299586629E-02), static_cast<RealType>(0.60419009528470238773E-02),
  370. static_cast<RealType>(0.38862217010742057883E-02), static_cast<RealType>(0.16793031084546090448E-02) };
  371. const RealType as = a*a;
  372. const RealType hs = -h*h*boost::math::constants::half<RealType>();
  373. RealType val = 0;
  374. for(unsigned short i = 0; i < m; ++i)
  375. {
  376. BOOST_ASSERT(i < 13);
  377. const RealType r = static_cast<RealType>(1) + as*pts[i];
  378. val += wts[i] * exp( hs*r ) / r;
  379. } // for(unsigned short i = 0; i < m; ++i)
  380. return val*a;
  381. } // RealType owens_t_T5(const RealType h, const RealType a)
  382. // compute the value of Owen's T function with method T5 from the reference paper
  383. template<typename RealType>
  384. inline RealType owens_t_T5_imp(const RealType h, const RealType a, const mpl::int_<64>&)
  385. {
  386. BOOST_MATH_STD_USING
  387. /*
  388. NOTICE:
  389. - The pts[] array contains the squares (!) of the abscissas, i.e. the roots of the Legendre
  390. polynomial P_n(x), instead of the plain roots as required in Gauss-Legendre
  391. quadrature, because T5(h,a,m) contains only x^2 terms.
  392. - The wts[] array contains the weights for Gauss-Legendre quadrature scaled with a factor
  393. of 1/(2*pi) according to T5(h,a,m).
  394. */
  395. const unsigned short m = 19;
  396. static const RealType pts[] = {
  397. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.0016634282895983227941),
  398. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.014904509242697054183),
  399. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.04103478879005817919),
  400. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.079359853513391511008),
  401. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.1288612130237615133),
  402. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.18822336642448518856),
  403. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.25586876186122962384),
  404. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.32999972011807857222),
  405. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.40864620815774761438),
  406. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.48971819306044782365),
  407. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.57106118513245543894),
  408. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.6505134942981533829),
  409. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.72596367859928091618),
  410. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.79540665919549865924),
  411. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.85699701386308739244),
  412. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.90909804422384697594),
  413. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.95032536436570154409),
  414. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.97958418733152273717),
  415. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.99610366384229088321)
  416. };
  417. static const RealType wts[] = {
  418. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.012975111395684900835),
  419. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.012888764187499150078),
  420. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.012716644398857307844),
  421. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.012459897461364705691),
  422. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.012120231988292330388),
  423. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.011699908404856841158),
  424. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.011201723906897224448),
  425. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.010628993848522759853),
  426. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.0099855296835573320047),
  427. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.0092756136096132857933),
  428. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.0085039700881139589055),
  429. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.0076757344408814561254),
  430. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.0067964187616556459109),
  431. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.005871875456524750363),
  432. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.0049082589542498110071),
  433. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.0039119870792519721409),
  434. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.0028897090921170700834),
  435. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.0018483371329504443947),
  436. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.00079623320100438873578)
  437. };
  438. const RealType as = a*a;
  439. const RealType hs = -h*h*boost::math::constants::half<RealType>();
  440. RealType val = 0;
  441. for(unsigned short i = 0; i < m; ++i)
  442. {
  443. BOOST_ASSERT(i < 19);
  444. const RealType r = 1 + as*pts[i];
  445. val += wts[i] * exp( hs*r ) / r;
  446. } // for(unsigned short i = 0; i < m; ++i)
  447. return val*a;
  448. } // RealType owens_t_T5(const RealType h, const RealType a)
  449. template<class RealType, class Policy>
  450. inline RealType owens_t_T5(const RealType h, const RealType a, const Policy&)
  451. {
  452. typedef typename policies::precision<RealType, Policy>::type precision_type;
  453. typedef typename mpl::if_<
  454. mpl::or_<
  455. mpl::less_equal<precision_type, mpl::int_<0> >,
  456. mpl::greater<precision_type, mpl::int_<53> >
  457. >,
  458. mpl::int_<64>,
  459. mpl::int_<53>
  460. >::type tag_type;
  461. return owens_t_T5_imp(h, a, tag_type());
  462. }
  463. // compute the value of Owen's T function with method T6 from the reference paper
  464. template<typename RealType>
  465. inline RealType owens_t_T6(const RealType h, const RealType a)
  466. {
  467. BOOST_MATH_STD_USING
  468. using namespace boost::math::constants;
  469. const RealType normh = owens_t_znorm2( h );
  470. const RealType y = static_cast<RealType>(1) - a;
  471. const RealType r = atan2(y, static_cast<RealType>(1 + a) );
  472. RealType val = normh * ( static_cast<RealType>(1) - normh ) * half<RealType>();
  473. if( r != 0 )
  474. val -= r * exp( -y*h*h*half<RealType>()/r ) * one_div_two_pi<RealType>();
  475. return val;
  476. } // RealType owens_t_T6(const RealType h, const RealType a, const unsigned short m)
  477. template <class T, class Policy>
  478. std::pair<T, T> owens_t_T1_accelerated(T h, T a, const Policy& pol)
  479. {
  480. //
  481. // This is the same series as T1, but:
  482. // * The Taylor series for atan has been combined with that for T1,
  483. // reducing but not eliminating cancellation error.
  484. // * The resulting alternating series is then accelerated using method 1
  485. // from H. Cohen, F. Rodriguez Villegas, D. Zagier,
  486. // "Convergence acceleration of alternating series", Bonn, (1991).
  487. //
  488. BOOST_MATH_STD_USING
  489. static const char* function = "boost::math::owens_t<%1%>(%1%, %1%)";
  490. T half_h_h = h * h / 2;
  491. T a_pow = a;
  492. T aa = a * a;
  493. T exp_term = exp(-h * h / 2);
  494. T one_minus_dj_sum = exp_term;
  495. T sum = a_pow * exp_term;
  496. T dj_pow = exp_term;
  497. T term = sum;
  498. T abs_err;
  499. int j = 1;
  500. //
  501. // Normally with this form of series acceleration we can calculate
  502. // up front how many terms will be required - based on the assumption
  503. // that each term decreases in size by a factor of 3. However,
  504. // that assumption does not apply here, as the underlying T1 series can
  505. // go quite strongly divergent in the early terms, before strongly
  506. // converging later. Various "guestimates" have been tried to take account
  507. // of this, but they don't always work.... so instead set "n" to the
  508. // largest value that won't cause overflow later, and abort iteration
  509. // when the last accelerated term was small enough...
  510. //
  511. int n;
  512. #ifndef BOOST_NO_EXCEPTIONS
  513. try
  514. {
  515. #endif
  516. n = itrunc(T(tools::log_max_value<T>() / 6));
  517. #ifndef BOOST_NO_EXCEPTIONS
  518. }
  519. catch(...)
  520. {
  521. n = (std::numeric_limits<int>::max)();
  522. }
  523. #endif
  524. n = (std::min)(n, 1500);
  525. T d = pow(3 + sqrt(T(8)), n);
  526. d = (d + 1 / d) / 2;
  527. T b = -1;
  528. T c = -d;
  529. c = b - c;
  530. sum *= c;
  531. b = -n * n * b * 2;
  532. abs_err = ldexp(fabs(sum), -tools::digits<T>());
  533. while(j < n)
  534. {
  535. a_pow *= aa;
  536. dj_pow *= half_h_h / j;
  537. one_minus_dj_sum += dj_pow;
  538. term = one_minus_dj_sum * a_pow / (2 * j + 1);
  539. c = b - c;
  540. sum += c * term;
  541. abs_err += ldexp((std::max)(T(fabs(sum)), T(fabs(c*term))), -tools::digits<T>());
  542. b = (j + n) * (j - n) * b / ((j + T(0.5)) * (j + 1));
  543. ++j;
  544. //
  545. // Include an escape route to prevent calculating too many terms:
  546. //
  547. if((j > 10) && (fabs(sum * tools::epsilon<T>()) > fabs(c * term)))
  548. break;
  549. }
  550. abs_err += fabs(c * term);
  551. if(sum < 0) // sum must always be positive, if it's negative something really bad has happend:
  552. policies::raise_evaluation_error(function, 0, T(0), pol);
  553. return std::pair<T, T>((sum / d) / boost::math::constants::two_pi<T>(), abs_err / sum);
  554. }
  555. template<typename RealType, class Policy>
  556. inline RealType owens_t_T2(const RealType h, const RealType a, const unsigned short m, const RealType ah, const Policy&, const mpl::true_&)
  557. {
  558. BOOST_MATH_STD_USING
  559. using namespace boost::math::constants;
  560. const unsigned short maxii = m+m+1;
  561. const RealType hs = h*h;
  562. const RealType as = -a*a;
  563. const RealType y = static_cast<RealType>(1) / hs;
  564. unsigned short ii = 1;
  565. RealType val = 0;
  566. RealType vi = a * exp( -ah*ah*half<RealType>() ) / root_two_pi<RealType>();
  567. RealType z = owens_t_znorm1(ah)/h;
  568. RealType last_z = fabs(z);
  569. RealType lim = policies::get_epsilon<RealType, Policy>();
  570. while( true )
  571. {
  572. val += z;
  573. //
  574. // This series stops converging after a while, so put a limit
  575. // on how far we go before returning our best guess:
  576. //
  577. if((fabs(lim * val) > fabs(z)) || ((ii > maxii) && (fabs(z) > last_z)) || (z == 0))
  578. {
  579. val *= exp( -hs*half<RealType>() ) / root_two_pi<RealType>();
  580. break;
  581. } // if( maxii <= ii )
  582. last_z = fabs(z);
  583. z = y * ( vi - static_cast<RealType>(ii) * z );
  584. vi *= as;
  585. ii += 2;
  586. } // while( true )
  587. return val;
  588. } // RealType owens_t_T2(const RealType h, const RealType a, const unsigned short m, const RealType ah)
  589. template<typename RealType, class Policy>
  590. inline std::pair<RealType, RealType> owens_t_T2_accelerated(const RealType h, const RealType a, const RealType ah, const Policy&)
  591. {
  592. //
  593. // This is the same series as T2, but with acceleration applied.
  594. // Note that we have to be *very* careful to check that nothing bad
  595. // has happened during evaluation - this series will go divergent
  596. // and/or fail to alternate at a drop of a hat! :-(
  597. //
  598. BOOST_MATH_STD_USING
  599. using namespace boost::math::constants;
  600. const RealType hs = h*h;
  601. const RealType as = -a*a;
  602. const RealType y = static_cast<RealType>(1) / hs;
  603. unsigned short ii = 1;
  604. RealType val = 0;
  605. RealType vi = a * exp( -ah*ah*half<RealType>() ) / root_two_pi<RealType>();
  606. RealType z = boost::math::detail::owens_t_znorm1(ah)/h;
  607. RealType last_z = fabs(z);
  608. //
  609. // Normally with this form of series acceleration we can calculate
  610. // up front how many terms will be required - based on the assumption
  611. // that each term decreases in size by a factor of 3. However,
  612. // that assumption does not apply here, as the underlying T1 series can
  613. // go quite strongly divergent in the early terms, before strongly
  614. // converging later. Various "guestimates" have been tried to take account
  615. // of this, but they don't always work.... so instead set "n" to the
  616. // largest value that won't cause overflow later, and abort iteration
  617. // when the last accelerated term was small enough...
  618. //
  619. int n;
  620. #ifndef BOOST_NO_EXCEPTIONS
  621. try
  622. {
  623. #endif
  624. n = itrunc(RealType(tools::log_max_value<RealType>() / 6));
  625. #ifndef BOOST_NO_EXCEPTIONS
  626. }
  627. catch(...)
  628. {
  629. n = (std::numeric_limits<int>::max)();
  630. }
  631. #endif
  632. n = (std::min)(n, 1500);
  633. RealType d = pow(3 + sqrt(RealType(8)), n);
  634. d = (d + 1 / d) / 2;
  635. RealType b = -1;
  636. RealType c = -d;
  637. int s = 1;
  638. for(int k = 0; k < n; ++k)
  639. {
  640. //
  641. // Check for both convergence and whether the series has gone bad:
  642. //
  643. if(
  644. (fabs(z) > last_z) // Series has gone divergent, abort
  645. || (fabs(val) * tools::epsilon<RealType>() > fabs(c * s * z)) // Convergence!
  646. || (z * s < 0) // Series has stopped alternating - all bets are off - abort.
  647. )
  648. {
  649. break;
  650. }
  651. c = b - c;
  652. val += c * s * z;
  653. b = (k + n) * (k - n) * b / ((k + RealType(0.5)) * (k + 1));
  654. last_z = fabs(z);
  655. s = -s;
  656. z = y * ( vi - static_cast<RealType>(ii) * z );
  657. vi *= as;
  658. ii += 2;
  659. } // while( true )
  660. RealType err = fabs(c * z) / val;
  661. return std::pair<RealType, RealType>(val * exp( -hs*half<RealType>() ) / (d * root_two_pi<RealType>()), err);
  662. } // RealType owens_t_T2_accelerated(const RealType h, const RealType a, const RealType ah, const Policy&)
  663. template<typename RealType, typename Policy>
  664. inline RealType T4_mp(const RealType h, const RealType a, const Policy& pol)
  665. {
  666. BOOST_MATH_STD_USING
  667. const RealType hs = h*h;
  668. const RealType as = -a*a;
  669. unsigned short ii = 1;
  670. RealType ai = constants::one_div_two_pi<RealType>() * a * exp( -0.5*hs*(1.0-as) );
  671. RealType yi = 1.0;
  672. RealType val = 0.0;
  673. RealType lim = boost::math::policies::get_epsilon<RealType, Policy>();
  674. while( true )
  675. {
  676. RealType term = ai*yi;
  677. val += term;
  678. if((yi != 0) && (fabs(val * lim) > fabs(term)))
  679. break;
  680. ii += 2;
  681. yi = (1.0-hs*yi) / static_cast<RealType>(ii);
  682. ai *= as;
  683. if(ii > (std::min)(1500, (int)policies::get_max_series_iterations<Policy>()))
  684. policies::raise_evaluation_error("boost::math::owens_t<%1%>", 0, val, pol);
  685. } // while( true )
  686. return val;
  687. } // arg_type owens_t_T4(const arg_type h, const arg_type a, const unsigned short m)
  688. // This routine dispatches the call to one of six subroutines, depending on the values
  689. // of h and a.
  690. // preconditions: h >= 0, 0<=a<=1, ah=a*h
  691. //
  692. // Note there are different versions for different precisions....
  693. template<typename RealType, typename Policy>
  694. inline RealType owens_t_dispatch(const RealType h, const RealType a, const RealType ah, const Policy& pol, mpl::int_<64> const&)
  695. {
  696. // Simple main case for 64-bit precision or less, this is as per the Patefield-Tandy paper:
  697. BOOST_MATH_STD_USING
  698. //
  699. // Handle some special cases first, these are from
  700. // page 1077 of Owen's original paper:
  701. //
  702. if(h == 0)
  703. {
  704. return atan(a) * constants::one_div_two_pi<RealType>();
  705. }
  706. if(a == 0)
  707. {
  708. return 0;
  709. }
  710. if(a == 1)
  711. {
  712. return owens_t_znorm2(RealType(-h)) * owens_t_znorm2(h) / 2;
  713. }
  714. if(a >= tools::max_value<RealType>())
  715. {
  716. return owens_t_znorm2(RealType(fabs(h)));
  717. }
  718. RealType val = 0; // avoid compiler warnings, 0 will be overwritten in any case
  719. const unsigned short icode = owens_t_compute_code(h, a);
  720. const unsigned short m = owens_t_get_order(icode, val /* just a dummy for the type */, pol);
  721. static const unsigned short meth[] = {1, 1, 1, 1, 1, 1, 1, 1, 2, 2, 2, 3, 4, 4, 4, 4, 5, 6}; // 18 entries
  722. // determine the appropriate method, T1 ... T6
  723. switch( meth[icode] )
  724. {
  725. case 1: // T1
  726. val = owens_t_T1(h,a,m,pol);
  727. break;
  728. case 2: // T2
  729. typedef typename policies::precision<RealType, Policy>::type precision_type;
  730. typedef mpl::bool_<(precision_type::value == 0) || (precision_type::value > 64)> tag_type;
  731. val = owens_t_T2(h, a, m, ah, pol, tag_type());
  732. break;
  733. case 3: // T3
  734. val = owens_t_T3(h,a,ah, pol);
  735. break;
  736. case 4: // T4
  737. val = owens_t_T4(h,a,m);
  738. break;
  739. case 5: // T5
  740. val = owens_t_T5(h,a, pol);
  741. break;
  742. case 6: // T6
  743. val = owens_t_T6(h,a);
  744. break;
  745. default:
  746. BOOST_THROW_EXCEPTION(std::logic_error("selection routine in Owen's T function failed"));
  747. }
  748. return val;
  749. }
  750. template<typename RealType, typename Policy>
  751. inline RealType owens_t_dispatch(const RealType h, const RealType a, const RealType ah, const Policy& pol, const mpl::int_<65>&)
  752. {
  753. // Arbitrary precision version:
  754. BOOST_MATH_STD_USING
  755. //
  756. // Handle some special cases first, these are from
  757. // page 1077 of Owen's original paper:
  758. //
  759. if(h == 0)
  760. {
  761. return atan(a) * constants::one_div_two_pi<RealType>();
  762. }
  763. if(a == 0)
  764. {
  765. return 0;
  766. }
  767. if(a == 1)
  768. {
  769. return owens_t_znorm2(RealType(-h)) * owens_t_znorm2(h) / 2;
  770. }
  771. if(a >= tools::max_value<RealType>())
  772. {
  773. return owens_t_znorm2(RealType(fabs(h)));
  774. }
  775. // Attempt arbitrary precision code, this will throw if it goes wrong:
  776. typedef typename boost::math::policies::normalise<Policy, boost::math::policies::evaluation_error<> >::type forwarding_policy;
  777. std::pair<RealType, RealType> p1(0, tools::max_value<RealType>()), p2(0, tools::max_value<RealType>());
  778. RealType target_precision = policies::get_epsilon<RealType, Policy>() * 1000;
  779. bool have_t1(false), have_t2(false);
  780. if(ah < 3)
  781. {
  782. #ifndef BOOST_NO_EXCEPTIONS
  783. try
  784. {
  785. #endif
  786. have_t1 = true;
  787. p1 = owens_t_T1_accelerated(h, a, forwarding_policy());
  788. if(p1.second < target_precision)
  789. return p1.first;
  790. #ifndef BOOST_NO_EXCEPTIONS
  791. }
  792. catch(const boost::math::evaluation_error&){} // T1 may fail and throw, that's OK
  793. #endif
  794. }
  795. if(ah > 1)
  796. {
  797. #ifndef BOOST_NO_EXCEPTIONS
  798. try
  799. {
  800. #endif
  801. have_t2 = true;
  802. p2 = owens_t_T2_accelerated(h, a, ah, forwarding_policy());
  803. if(p2.second < target_precision)
  804. return p2.first;
  805. #ifndef BOOST_NO_EXCEPTIONS
  806. }
  807. catch(const boost::math::evaluation_error&){} // T2 may fail and throw, that's OK
  808. #endif
  809. }
  810. //
  811. // If we haven't tried T1 yet, do it now - sometimes it succeeds and the number of iterations
  812. // is fairly low compared to T4.
  813. //
  814. if(!have_t1)
  815. {
  816. #ifndef BOOST_NO_EXCEPTIONS
  817. try
  818. {
  819. #endif
  820. have_t1 = true;
  821. p1 = owens_t_T1_accelerated(h, a, forwarding_policy());
  822. if(p1.second < target_precision)
  823. return p1.first;
  824. #ifndef BOOST_NO_EXCEPTIONS
  825. }
  826. catch(const boost::math::evaluation_error&){} // T1 may fail and throw, that's OK
  827. #endif
  828. }
  829. //
  830. // If we haven't tried T2 yet, do it now - sometimes it succeeds and the number of iterations
  831. // is fairly low compared to T4.
  832. //
  833. if(!have_t2)
  834. {
  835. #ifndef BOOST_NO_EXCEPTIONS
  836. try
  837. {
  838. #endif
  839. have_t2 = true;
  840. p2 = owens_t_T2_accelerated(h, a, ah, forwarding_policy());
  841. if(p2.second < target_precision)
  842. return p2.first;
  843. #ifndef BOOST_NO_EXCEPTIONS
  844. }
  845. catch(const boost::math::evaluation_error&){} // T2 may fail and throw, that's OK
  846. #endif
  847. }
  848. //
  849. // OK, nothing left to do but try the most expensive option which is T4,
  850. // this is often slow to converge, but when it does converge it tends to
  851. // be accurate:
  852. #ifndef BOOST_NO_EXCEPTIONS
  853. try
  854. {
  855. #endif
  856. return T4_mp(h, a, pol);
  857. #ifndef BOOST_NO_EXCEPTIONS
  858. }
  859. catch(const boost::math::evaluation_error&){} // T4 may fail and throw, that's OK
  860. #endif
  861. //
  862. // Now look back at the results from T1 and T2 and see if either gave better
  863. // results than we could get from the 64-bit precision versions.
  864. //
  865. if((std::min)(p1.second, p2.second) < 1e-20)
  866. {
  867. return p1.second < p2.second ? p1.first : p2.first;
  868. }
  869. //
  870. // We give up - no arbitrary precision versions succeeded!
  871. //
  872. return owens_t_dispatch(h, a, ah, pol, mpl::int_<64>());
  873. } // RealType owens_t_dispatch(RealType h, RealType a, RealType ah)
  874. template<typename RealType, typename Policy>
  875. inline RealType owens_t_dispatch(const RealType h, const RealType a, const RealType ah, const Policy& pol, const mpl::int_<0>&)
  876. {
  877. // We don't know what the precision is until runtime:
  878. if(tools::digits<RealType>() <= 64)
  879. return owens_t_dispatch(h, a, ah, pol, mpl::int_<64>());
  880. return owens_t_dispatch(h, a, ah, pol, mpl::int_<65>());
  881. }
  882. template<typename RealType, typename Policy>
  883. inline RealType owens_t_dispatch(const RealType h, const RealType a, const RealType ah, const Policy& pol)
  884. {
  885. // Figure out the precision and forward to the correct version:
  886. typedef typename policies::precision<RealType, Policy>::type precision_type;
  887. typedef typename mpl::if_c<
  888. precision_type::value == 0,
  889. mpl::int_<0>,
  890. typename mpl::if_c<
  891. precision_type::value <= 64,
  892. mpl::int_<64>,
  893. mpl::int_<65>
  894. >::type
  895. >::type tag_type;
  896. return owens_t_dispatch(h, a, ah, pol, tag_type());
  897. }
  898. // compute Owen's T function, T(h,a), for arbitrary values of h and a
  899. template<typename RealType, class Policy>
  900. inline RealType owens_t(RealType h, RealType a, const Policy& pol)
  901. {
  902. BOOST_MATH_STD_USING
  903. // exploit that T(-h,a) == T(h,a)
  904. h = fabs(h);
  905. // Use equation (2) in the paper to remap the arguments
  906. // such that h>=0 and 0<=a<=1 for the call of the actual
  907. // computation routine.
  908. const RealType fabs_a = fabs(a);
  909. const RealType fabs_ah = fabs_a*h;
  910. RealType val = 0.0; // avoid compiler warnings, 0.0 will be overwritten in any case
  911. if(fabs_a <= 1)
  912. {
  913. val = owens_t_dispatch(h, fabs_a, fabs_ah, pol);
  914. } // if(fabs_a <= 1.0)
  915. else
  916. {
  917. if( h <= 0.67 )
  918. {
  919. const RealType normh = owens_t_znorm1(h);
  920. const RealType normah = owens_t_znorm1(fabs_ah);
  921. val = static_cast<RealType>(1)/static_cast<RealType>(4) - normh*normah -
  922. owens_t_dispatch(fabs_ah, static_cast<RealType>(1 / fabs_a), h, pol);
  923. } // if( h <= 0.67 )
  924. else
  925. {
  926. const RealType normh = detail::owens_t_znorm2(h);
  927. const RealType normah = detail::owens_t_znorm2(fabs_ah);
  928. val = constants::half<RealType>()*(normh+normah) - normh*normah -
  929. owens_t_dispatch(fabs_ah, static_cast<RealType>(1 / fabs_a), h, pol);
  930. } // else [if( h <= 0.67 )]
  931. } // else [if(fabs_a <= 1)]
  932. // exploit that T(h,-a) == -T(h,a)
  933. if(a < 0)
  934. {
  935. return -val;
  936. } // if(a < 0)
  937. return val;
  938. } // RealType owens_t(RealType h, RealType a)
  939. template <class T, class Policy, class tag>
  940. struct owens_t_initializer
  941. {
  942. struct init
  943. {
  944. init()
  945. {
  946. do_init(tag());
  947. }
  948. template <int N>
  949. static void do_init(const mpl::int_<N>&){}
  950. static void do_init(const mpl::int_<64>&)
  951. {
  952. boost::math::owens_t(static_cast<T>(7), static_cast<T>(0.96875), Policy());
  953. boost::math::owens_t(static_cast<T>(2), static_cast<T>(0.5), Policy());
  954. }
  955. void force_instantiate()const{}
  956. };
  957. static const init initializer;
  958. static void force_instantiate()
  959. {
  960. initializer.force_instantiate();
  961. }
  962. };
  963. template <class T, class Policy, class tag>
  964. const typename owens_t_initializer<T, Policy, tag>::init owens_t_initializer<T, Policy, tag>::initializer;
  965. } // namespace detail
  966. template <class T1, class T2, class Policy>
  967. inline typename tools::promote_args<T1, T2>::type owens_t(T1 h, T2 a, const Policy& pol)
  968. {
  969. typedef typename tools::promote_args<T1, T2>::type result_type;
  970. typedef typename policies::evaluation<result_type, Policy>::type value_type;
  971. typedef typename policies::precision<value_type, Policy>::type precision_type;
  972. typedef typename mpl::if_c<
  973. precision_type::value == 0,
  974. mpl::int_<0>,
  975. typename mpl::if_c<
  976. precision_type::value <= 64,
  977. mpl::int_<64>,
  978. mpl::int_<65>
  979. >::type
  980. >::type tag_type;
  981. detail::owens_t_initializer<result_type, Policy, tag_type>::force_instantiate();
  982. return policies::checked_narrowing_cast<result_type, Policy>(detail::owens_t(static_cast<value_type>(h), static_cast<value_type>(a), pol), "boost::math::owens_t<%1%>(%1%,%1%)");
  983. }
  984. template <class T1, class T2>
  985. inline typename tools::promote_args<T1, T2>::type owens_t(T1 h, T2 a)
  986. {
  987. return owens_t(h, a, policies::policy<>());
  988. }
  989. } // namespace math
  990. } // namespace boost
  991. #ifdef BOOST_MSVC
  992. #pragma warning(pop)
  993. #endif
  994. #endif
  995. // EOF