owens_t.hpp 49 KB

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