owens_t.hpp 48 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877878879880881882883884885886887888889890891892893894895896897898899900901902903904905906907908909910911912913914915916917918919920921922923924925926927928929930931932933934935936937938939940941942943944945946947948949950951952953954955956957958959960961962963964965966967968969970971972973974975976977978979980981982983984985986987988989990991992993994995996997998999100010011002100310041005100610071008100910101011101210131014101510161017101810191020102110221023102410251026102710281029103010311032103310341035103610371038103910401041104210431044104510461047104810491050105110521053105410551056105710581059106010611062106310641065106610671068
  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/math/special_functions/erf.hpp>
  16. #include <boost/math/special_functions/expm1.hpp>
  17. #include <boost/math/tools/throw_exception.hpp>
  18. #include <boost/math/tools/assert.hpp>
  19. #include <boost/math/constants/constants.hpp>
  20. #include <boost/math/tools/big_constant.hpp>
  21. #include <stdexcept>
  22. #include <cmath>
  23. #ifdef _MSC_VER
  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. // LCOV_EXCL_START
  63. static const RealType hrange[] =
  64. { 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 };
  65. static const RealType arange[] = { 0.025f, 0.09f, 0.15f, 0.36f, 0.5f, 0.9f, 0.99999f };
  66. /*
  67. original select array from paper:
  68. 1, 1, 2,13,13,13,13,13,13,13,13,16,16,16, 9
  69. 1, 2, 2, 3, 3, 5, 5,14,14,15,15,16,16,16, 9
  70. 2, 2, 3, 3, 3, 5, 5,15,15,15,15,16,16,16,10
  71. 2, 2, 3, 5, 5, 5, 5, 7, 7,16,16,16,16,16,10
  72. 2, 3, 3, 5, 5, 6, 6, 8, 8,17,17,17,12,12,11
  73. 2, 3, 5, 5, 5, 6, 6, 8, 8,17,17,17,12,12,12
  74. 2, 3, 4, 4, 6, 6, 8, 8,17,17,17,17,17,12,12
  75. 2, 3, 4, 4, 6, 6,18,18,18,18,17,17,17,12,12
  76. */
  77. // subtract one because the array is written in FORTRAN in mind - in C arrays start @ zero
  78. static const unsigned short select[] =
  79. {
  80. 0, 0 , 1 , 12 ,12 , 12 , 12 , 12 , 12 , 12 , 12 , 15 , 15 , 15 , 8,
  81. 0 , 1 , 1 , 2 , 2 , 4 , 4 , 13 , 13 , 14 , 14 , 15 , 15 , 15 , 8,
  82. 1 , 1 , 2 , 2 , 2 , 4 , 4 , 14 , 14 , 14 , 14 , 15 , 15 , 15 , 9,
  83. 1 , 1 , 2 , 4 , 4 , 4 , 4 , 6 , 6 , 15 , 15 , 15 , 15 , 15 , 9,
  84. 1 , 2 , 2 , 4 , 4 , 5 , 5 , 7 , 7 , 16 ,16 , 16 , 11 , 11 , 10,
  85. 1 , 2 , 4 , 4 , 4 , 5 , 5 , 7 , 7 , 16 , 16 , 16 , 11 , 11 , 11,
  86. 1 , 2 , 3 , 3 , 5 , 5 , 7 , 7 , 16 , 16 , 16 , 16 , 16 , 11 , 11,
  87. 1 , 2 , 3 , 3 , 5 , 5 , 17 , 17 , 17 , 17 , 16 , 16 , 16 , 11 , 11
  88. };
  89. // LCOV_EXCL_STOP
  90. unsigned short ihint = 14, iaint = 7;
  91. for(unsigned short i = 0; i != 14; i++)
  92. {
  93. if( h <= hrange[i] )
  94. {
  95. ihint = i;
  96. break;
  97. }
  98. } // for(unsigned short i = 0; i != 14; i++)
  99. for(unsigned short i = 0; i != 7; i++)
  100. {
  101. if( a <= arange[i] )
  102. {
  103. iaint = i;
  104. break;
  105. }
  106. } // for(unsigned short i = 0; i != 7; i++)
  107. // interpret select array as 8x15 matrix
  108. BOOST_MATH_ASSERT(iaint * 15 + ihint < (int)(sizeof(select) / sizeof(select[0])));
  109. return select[iaint*15 + ihint];
  110. } // unsigned short owens_t_compute_code(const RealType h, const RealType a)
  111. template<typename RealType>
  112. inline unsigned short owens_t_get_order_imp(const unsigned short icode, RealType, const std::integral_constant<int, 53>&)
  113. {
  114. // LCOV_EXCL_START
  115. 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
  116. // LCOV_EXCL_STOP
  117. BOOST_MATH_ASSERT(icode<18);
  118. return ord[icode];
  119. } // unsigned short owens_t_get_order(const unsigned short icode, RealType, std::integral_constant<int, 53> const&)
  120. template<typename RealType>
  121. inline unsigned short owens_t_get_order_imp(const unsigned short icode, RealType, const std::integral_constant<int, 64>&)
  122. {
  123. // method ================>>> {1, 1, 1, 1, 1, 1, 1, 1, 2, 2, 2, 3, 4, 4, 4, 4, 5, 6}
  124. // LCOV_EXCL_START
  125. 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
  126. // LCOV_EXCL_STOP
  127. BOOST_MATH_ASSERT(icode<18);
  128. return ord[icode];
  129. } // unsigned short owens_t_get_order(const unsigned short icode, RealType, std::integral_constant<int, 64> const&)
  130. template<typename RealType, typename Policy>
  131. inline unsigned short owens_t_get_order(const unsigned short icode, RealType r, const Policy&)
  132. {
  133. typedef typename policies::precision<RealType, Policy>::type precision_type;
  134. typedef std::integral_constant<int,
  135. precision_type::value <= 0 ? 64 :
  136. precision_type::value <= 53 ? 53 : 64
  137. > tag_type;
  138. return owens_t_get_order_imp(icode, r, tag_type());
  139. }
  140. // compute the value of Owen's T function with method T1 from the reference paper
  141. template<typename RealType, typename Policy>
  142. inline RealType owens_t_T1(const RealType h, const RealType a, const unsigned short m, const Policy& pol)
  143. {
  144. BOOST_MATH_STD_USING
  145. using namespace boost::math::constants;
  146. const RealType hs = -h*h*half<RealType>();
  147. const RealType dhs = exp( hs );
  148. const RealType as = a*a;
  149. unsigned short j=1;
  150. RealType jj = 1;
  151. RealType aj = a * one_div_two_pi<RealType>();
  152. RealType dj = boost::math::expm1( hs, pol);
  153. RealType gj = hs*dhs;
  154. RealType val = atan( a ) * one_div_two_pi<RealType>();
  155. while( true )
  156. {
  157. val += dj*aj/jj;
  158. if( m <= j )
  159. break;
  160. j++;
  161. jj += static_cast<RealType>(2);
  162. aj *= as;
  163. dj = gj - dj;
  164. gj *= hs / static_cast<RealType>(j);
  165. } // while( true )
  166. return val;
  167. } // RealType owens_t_T1(const RealType h, const RealType a, const unsigned short m)
  168. // compute the value of Owen's T function with method T2 from the reference paper
  169. template<typename RealType, class Policy>
  170. 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&)
  171. {
  172. BOOST_MATH_STD_USING
  173. using namespace boost::math::constants;
  174. const unsigned short maxii = m+m+1;
  175. const RealType hs = h*h;
  176. const RealType as = -a*a;
  177. const RealType y = static_cast<RealType>(1) / hs;
  178. unsigned short ii = 1;
  179. RealType val = 0;
  180. RealType vi = a * exp( -ah*ah*half<RealType>() ) * one_div_root_two_pi<RealType>();
  181. RealType z = owens_t_znorm1(ah, pol)/h;
  182. while( true )
  183. {
  184. val += z;
  185. if( maxii <= ii )
  186. {
  187. val *= exp( -hs*half<RealType>() ) * one_div_root_two_pi<RealType>();
  188. break;
  189. } // if( maxii <= ii )
  190. z = y * ( vi - static_cast<RealType>(ii) * z );
  191. vi *= as;
  192. ii += 2;
  193. } // while( true )
  194. return val;
  195. } // RealType owens_t_T2(const RealType h, const RealType a, const unsigned short m, const RealType ah)
  196. // compute the value of Owen's T function with method T3 from the reference paper
  197. template<typename RealType, class Policy>
  198. inline RealType owens_t_T3_imp(const RealType h, const RealType a, const RealType ah, const std::integral_constant<int, 53>&, const Policy& pol)
  199. {
  200. BOOST_MATH_STD_USING
  201. using namespace boost::math::constants;
  202. const unsigned short m = 20;
  203. // LCOV_EXCL_START
  204. static const RealType c2[] =
  205. {
  206. static_cast<RealType>(0.99999999999999987510),
  207. static_cast<RealType>(-0.99999999999988796462), static_cast<RealType>(0.99999999998290743652),
  208. static_cast<RealType>(-0.99999999896282500134), static_cast<RealType>(0.99999996660459362918),
  209. static_cast<RealType>(-0.99999933986272476760), static_cast<RealType>(0.99999125611136965852),
  210. static_cast<RealType>(-0.99991777624463387686), static_cast<RealType>(0.99942835555870132569),
  211. static_cast<RealType>(-0.99697311720723000295), static_cast<RealType>(0.98751448037275303682),
  212. static_cast<RealType>(-0.95915857980572882813), static_cast<RealType>(0.89246305511006708555),
  213. static_cast<RealType>(-0.76893425990463999675), static_cast<RealType>(0.58893528468484693250),
  214. static_cast<RealType>(-0.38380345160440256652), static_cast<RealType>(0.20317601701045299653),
  215. static_cast<RealType>(-0.82813631607004984866E-01), static_cast<RealType>(0.24167984735759576523E-01),
  216. static_cast<RealType>(-0.44676566663971825242E-02), static_cast<RealType>(0.39141169402373836468E-03)
  217. };
  218. // LCOV_EXCL_STOP
  219. const RealType as = a*a;
  220. const RealType hs = h*h;
  221. const RealType y = static_cast<RealType>(1)/hs;
  222. RealType ii = 1;
  223. unsigned short i = 0;
  224. RealType vi = a * exp( -ah*ah*half<RealType>() ) * one_div_root_two_pi<RealType>();
  225. RealType zi = owens_t_znorm1(ah, pol)/h;
  226. RealType val = 0;
  227. while( true )
  228. {
  229. BOOST_MATH_ASSERT(i < 21);
  230. val += zi*c2[i];
  231. if( m <= i ) // if( m < i+1 )
  232. {
  233. val *= exp( -hs*half<RealType>() ) * one_div_root_two_pi<RealType>();
  234. break;
  235. } // if( m < i )
  236. zi = y * (ii*zi - vi);
  237. vi *= as;
  238. ii += 2;
  239. i++;
  240. } // while( true )
  241. return val;
  242. } // RealType owens_t_T3(const RealType h, const RealType a, const RealType ah)
  243. // compute the value of Owen's T function with method T3 from the reference paper
  244. template<class RealType, class Policy>
  245. inline RealType owens_t_T3_imp(const RealType h, const RealType a, const RealType ah, const std::integral_constant<int, 64>&, const Policy& pol)
  246. {
  247. BOOST_MATH_STD_USING
  248. using namespace boost::math::constants;
  249. const unsigned short m = 30;
  250. // LCOV_EXCL_START
  251. static const RealType c2[] =
  252. {
  253. BOOST_MATH_BIG_CONSTANT(RealType, 260, 0.99999999999999999999999729978162447266851932041876728736094298092917625009873),
  254. BOOST_MATH_BIG_CONSTANT(RealType, 260, -0.99999999999999999999467056379678391810626533251885323416799874878563998732905968),
  255. BOOST_MATH_BIG_CONSTANT(RealType, 260, 0.99999999999999999824849349313270659391127814689133077036298754586814091034842536),
  256. BOOST_MATH_BIG_CONSTANT(RealType, 260, -0.9999999999999997703859616213643405880166422891953033591551179153879839440241685),
  257. BOOST_MATH_BIG_CONSTANT(RealType, 260, 0.99999999999998394883415238173334565554173013941245103172035286759201504179038147),
  258. BOOST_MATH_BIG_CONSTANT(RealType, 260, -0.9999999999993063616095509371081203145247992197457263066869044528823599399470977),
  259. BOOST_MATH_BIG_CONSTANT(RealType, 260, 0.9999999999797336340409464429599229870590160411238245275855903767652432017766116267),
  260. BOOST_MATH_BIG_CONSTANT(RealType, 260, -0.999999999574958412069046680119051639753412378037565521359444170241346845522403274),
  261. BOOST_MATH_BIG_CONSTANT(RealType, 260, 0.9999999933226234193375324943920160947158239076786103108097456617750134812033362048),
  262. BOOST_MATH_BIG_CONSTANT(RealType, 260, -0.9999999188923242461073033481053037468263536806742737922476636768006622772762168467),
  263. BOOST_MATH_BIG_CONSTANT(RealType, 260, 0.9999992195143483674402853783549420883055129680082932629160081128947764415749728967),
  264. BOOST_MATH_BIG_CONSTANT(RealType, 260, -0.999993935137206712830997921913316971472227199741857386575097250553105958772041501),
  265. BOOST_MATH_BIG_CONSTANT(RealType, 260, 0.99996135597690552745362392866517133091672395614263398912807169603795088421057688716),
  266. BOOST_MATH_BIG_CONSTANT(RealType, 260, -0.99979556366513946026406788969630293820987757758641211293079784585126692672425362469),
  267. BOOST_MATH_BIG_CONSTANT(RealType, 260, 0.999092789629617100153486251423850590051366661947344315423226082520411961968929483),
  268. BOOST_MATH_BIG_CONSTANT(RealType, 260, -0.996593837411918202119308620432614600338157335862888580671450938858935084316004769854),
  269. BOOST_MATH_BIG_CONSTANT(RealType, 260, 0.98910017138386127038463510314625339359073956513420458166238478926511821146316469589567),
  270. BOOST_MATH_BIG_CONSTANT(RealType, 260, -0.970078558040693314521331982203762771512160168582494513347846407314584943870399016019),
  271. BOOST_MATH_BIG_CONSTANT(RealType, 260, 0.92911438683263187495758525500033707204091967947532160289872782771388170647150321633673),
  272. BOOST_MATH_BIG_CONSTANT(RealType, 260, -0.8542058695956156057286980736842905011429254735181323743367879525470479126968822863),
  273. BOOST_MATH_BIG_CONSTANT(RealType, 260, 0.73796526033030091233118357742803709382964420335559408722681794195743240930748630755),
  274. BOOST_MATH_BIG_CONSTANT(RealType, 260, -0.58523469882837394570128599003785154144164680587615878645171632791404210655891158),
  275. BOOST_MATH_BIG_CONSTANT(RealType, 260, 0.415997776145676306165661663581868460503874205343014196580122174949645271353372263),
  276. BOOST_MATH_BIG_CONSTANT(RealType, 260, -0.2588210875241943574388730510317252236407805082485246378222935376279663808416534365),
  277. BOOST_MATH_BIG_CONSTANT(RealType, 260, 0.1375535825163892648504646951500265585055789019410617565727090346559210218472356689),
  278. BOOST_MATH_BIG_CONSTANT(RealType, 260, -0.0607952766325955730493900985022020434830339794955745989150270485056436844239206648),
  279. BOOST_MATH_BIG_CONSTANT(RealType, 260, 0.0216337683299871528059836483840390514275488679530797294557060229266785853764115),
  280. BOOST_MATH_BIG_CONSTANT(RealType, 260, -0.00593405693455186729876995814181203900550014220428843483927218267309209471516256),
  281. BOOST_MATH_BIG_CONSTANT(RealType, 260, 0.0011743414818332946510474576182739210553333860106811865963485870668929503649964142),
  282. BOOST_MATH_BIG_CONSTANT(RealType, 260, -1.489155613350368934073453260689881330166342484405529981510694514036264969925132e-4),
  283. BOOST_MATH_BIG_CONSTANT(RealType, 260, 9.072354320794357587710929507988814669454281514268844884841547607134260303118208e-6)
  284. };
  285. // LCOV_EXCL_STOP
  286. const RealType as = a*a;
  287. const RealType hs = h*h;
  288. const RealType y = 1 / hs;
  289. RealType ii = 1;
  290. unsigned short i = 0;
  291. RealType vi = a * exp( -ah*ah*half<RealType>() ) * one_div_root_two_pi<RealType>();
  292. RealType zi = owens_t_znorm1(ah, pol)/h;
  293. RealType val = 0;
  294. while( true )
  295. {
  296. BOOST_MATH_ASSERT(i < 31);
  297. val += zi*c2[i];
  298. if( m <= i ) // if( m < i+1 )
  299. {
  300. val *= exp( -hs*half<RealType>() ) * one_div_root_two_pi<RealType>();
  301. break;
  302. } // if( m < i )
  303. zi = y * (ii*zi - vi);
  304. vi *= as;
  305. ii += 2;
  306. i++;
  307. } // while( true )
  308. return val;
  309. } // RealType owens_t_T3(const RealType h, const RealType a, const RealType ah)
  310. template<class RealType, class Policy>
  311. inline RealType owens_t_T3(const RealType h, const RealType a, const RealType ah, const Policy& pol)
  312. {
  313. typedef typename policies::precision<RealType, Policy>::type precision_type;
  314. typedef std::integral_constant<int,
  315. precision_type::value <= 0 ? 64 :
  316. precision_type::value <= 53 ? 53 : 64
  317. > tag_type;
  318. return owens_t_T3_imp(h, a, ah, tag_type(), pol);
  319. }
  320. // compute the value of Owen's T function with method T4 from the reference paper
  321. template<typename RealType>
  322. inline RealType owens_t_T4(const RealType h, const RealType a, const unsigned short m)
  323. {
  324. BOOST_MATH_STD_USING
  325. using namespace boost::math::constants;
  326. const unsigned short maxii = m+m+1;
  327. const RealType hs = h*h;
  328. const RealType as = -a*a;
  329. unsigned short ii = 1;
  330. RealType ai = a * exp( -hs*(static_cast<RealType>(1)-as)*half<RealType>() ) * one_div_two_pi<RealType>();
  331. RealType yi = 1;
  332. RealType val = 0;
  333. while( true )
  334. {
  335. val += ai*yi;
  336. if( maxii <= ii )
  337. break;
  338. ii += 2;
  339. yi = (static_cast<RealType>(1)-hs*yi) / static_cast<RealType>(ii);
  340. ai *= as;
  341. } // while( true )
  342. return val;
  343. } // RealType owens_t_T4(const RealType h, const RealType a, const unsigned short m)
  344. // compute the value of Owen's T function with method T5 from the reference paper
  345. template<typename RealType>
  346. inline RealType owens_t_T5_imp(const RealType h, const RealType a, const std::integral_constant<int, 53>&)
  347. {
  348. BOOST_MATH_STD_USING
  349. /*
  350. NOTICE:
  351. - The pts[] array contains the squares (!) of the abscissas, i.e. the roots of the Legendre
  352. polynomial P_n(x), instead of the plain roots as required in Gauss-Legendre
  353. quadrature, because T5(h,a,m) contains only x^2 terms.
  354. - The wts[] array contains the weights for Gauss-Legendre quadrature scaled with a factor
  355. of 1/(2*pi) according to T5(h,a,m).
  356. */
  357. const unsigned short m = 13;
  358. // LCOV_EXCL_START
  359. static const RealType pts[] = {
  360. static_cast<RealType>(0.35082039676451715489E-02),
  361. static_cast<RealType>(0.31279042338030753740E-01), static_cast<RealType>(0.85266826283219451090E-01),
  362. static_cast<RealType>(0.16245071730812277011), static_cast<RealType>(0.25851196049125434828),
  363. static_cast<RealType>(0.36807553840697533536), static_cast<RealType>(0.48501092905604697475),
  364. static_cast<RealType>(0.60277514152618576821), static_cast<RealType>(0.71477884217753226516),
  365. static_cast<RealType>(0.81475510988760098605), static_cast<RealType>(0.89711029755948965867),
  366. static_cast<RealType>(0.95723808085944261843), static_cast<RealType>(0.99178832974629703586) };
  367. static const RealType wts[] = {
  368. static_cast<RealType>(0.18831438115323502887E-01),
  369. static_cast<RealType>(0.18567086243977649478E-01), static_cast<RealType>(0.18042093461223385584E-01),
  370. static_cast<RealType>(0.17263829606398753364E-01), static_cast<RealType>(0.16243219975989856730E-01),
  371. static_cast<RealType>(0.14994592034116704829E-01), static_cast<RealType>(0.13535474469662088392E-01),
  372. static_cast<RealType>(0.11886351605820165233E-01), static_cast<RealType>(0.10070377242777431897E-01),
  373. static_cast<RealType>(0.81130545742299586629E-02), static_cast<RealType>(0.60419009528470238773E-02),
  374. static_cast<RealType>(0.38862217010742057883E-02), static_cast<RealType>(0.16793031084546090448E-02) };
  375. const RealType as = a*a;
  376. const RealType hs = -h*h*boost::math::constants::half<RealType>();
  377. // LCOV_EXCL_STOP
  378. RealType val = 0;
  379. for(unsigned short i = 0; i < m; ++i)
  380. {
  381. BOOST_MATH_ASSERT(i < 13);
  382. const RealType r = static_cast<RealType>(1) + as*pts[i];
  383. val += wts[i] * exp( hs*r ) / r;
  384. } // for(unsigned short i = 0; i < m; ++i)
  385. return val*a;
  386. } // RealType owens_t_T5(const RealType h, const RealType a)
  387. // compute the value of Owen's T function with method T5 from the reference paper
  388. template<typename RealType>
  389. inline RealType owens_t_T5_imp(const RealType h, const RealType a, const std::integral_constant<int, 64>&)
  390. {
  391. BOOST_MATH_STD_USING
  392. /*
  393. NOTICE:
  394. - The pts[] array contains the squares (!) of the abscissas, i.e. the roots of the Legendre
  395. polynomial P_n(x), instead of the plain roots as required in Gauss-Legendre
  396. quadrature, because T5(h,a,m) contains only x^2 terms.
  397. - The wts[] array contains the weights for Gauss-Legendre quadrature scaled with a factor
  398. of 1/(2*pi) according to T5(h,a,m).
  399. */
  400. const unsigned short m = 19;
  401. // LCOV_EXCL_START
  402. static const RealType pts[] = {
  403. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.0016634282895983227941),
  404. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.014904509242697054183),
  405. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.04103478879005817919),
  406. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.079359853513391511008),
  407. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.1288612130237615133),
  408. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.18822336642448518856),
  409. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.25586876186122962384),
  410. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.32999972011807857222),
  411. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.40864620815774761438),
  412. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.48971819306044782365),
  413. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.57106118513245543894),
  414. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.6505134942981533829),
  415. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.72596367859928091618),
  416. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.79540665919549865924),
  417. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.85699701386308739244),
  418. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.90909804422384697594),
  419. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.95032536436570154409),
  420. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.97958418733152273717),
  421. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.99610366384229088321)
  422. };
  423. static const RealType wts[] = {
  424. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.012975111395684900835),
  425. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.012888764187499150078),
  426. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.012716644398857307844),
  427. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.012459897461364705691),
  428. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.012120231988292330388),
  429. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.011699908404856841158),
  430. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.011201723906897224448),
  431. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.010628993848522759853),
  432. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.0099855296835573320047),
  433. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.0092756136096132857933),
  434. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.0085039700881139589055),
  435. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.0076757344408814561254),
  436. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.0067964187616556459109),
  437. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.005871875456524750363),
  438. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.0049082589542498110071),
  439. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.0039119870792519721409),
  440. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.0028897090921170700834),
  441. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.0018483371329504443947),
  442. BOOST_MATH_BIG_CONSTANT(RealType, 64, 0.00079623320100438873578)
  443. };
  444. // LCOV_EXCL_STOP
  445. const RealType as = a*a;
  446. const RealType hs = -h*h*boost::math::constants::half<RealType>();
  447. RealType val = 0;
  448. for(unsigned short i = 0; i < m; ++i)
  449. {
  450. BOOST_MATH_ASSERT(i < 19);
  451. const RealType r = 1 + as*pts[i];
  452. val += wts[i] * exp( hs*r ) / r;
  453. } // for(unsigned short i = 0; i < m; ++i)
  454. return val*a;
  455. } // RealType owens_t_T5(const RealType h, const RealType a)
  456. template<class RealType, class Policy>
  457. inline RealType owens_t_T5(const RealType h, const RealType a, const Policy&)
  458. {
  459. typedef typename policies::precision<RealType, Policy>::type precision_type;
  460. typedef std::integral_constant<int,
  461. precision_type::value <= 0 ? 64 :
  462. precision_type::value <= 53 ? 53 : 64
  463. > tag_type;
  464. return owens_t_T5_imp(h, a, tag_type());
  465. }
  466. // compute the value of Owen's T function with method T6 from the reference paper
  467. template<typename RealType, class Policy>
  468. inline RealType owens_t_T6(const RealType h, const RealType a, const Policy& pol)
  469. {
  470. BOOST_MATH_STD_USING
  471. using namespace boost::math::constants;
  472. const RealType normh = owens_t_znorm2(h, pol);
  473. const RealType y = static_cast<RealType>(1) - a;
  474. const RealType r = atan2(y, static_cast<RealType>(1 + a) );
  475. RealType val = normh * ( static_cast<RealType>(1) - normh ) * half<RealType>();
  476. if( r != 0 )
  477. val -= r * exp( -y*h*h*half<RealType>()/r ) * one_div_two_pi<RealType>();
  478. return val;
  479. } // RealType owens_t_T6(const RealType h, const RealType a, const unsigned short m)
  480. template <class T, class Policy>
  481. std::pair<T, T> owens_t_T1_accelerated(T h, T a, const Policy& pol)
  482. {
  483. //
  484. // This is the same series as T1, but:
  485. // * The Taylor series for atan has been combined with that for T1,
  486. // reducing but not eliminating cancellation error.
  487. // * The resulting alternating series is then accelerated using method 1
  488. // from H. Cohen, F. Rodriguez Villegas, D. Zagier,
  489. // "Convergence acceleration of alternating series", Bonn, (1991).
  490. //
  491. BOOST_MATH_STD_USING
  492. static const char* function = "boost::math::owens_t<%1%>(%1%, %1%)";
  493. T half_h_h = h * h / 2;
  494. T a_pow = a;
  495. T aa = a * a;
  496. T exp_term = exp(-h * h / 2);
  497. T one_minus_dj_sum = exp_term;
  498. T sum = a_pow * exp_term;
  499. T dj_pow = exp_term;
  500. T term = sum;
  501. T abs_err;
  502. int j = 1;
  503. //
  504. // Normally with this form of series acceleration we can calculate
  505. // up front how many terms will be required - based on the assumption
  506. // that each term decreases in size by a factor of 3. However,
  507. // that assumption does not apply here, as the underlying T1 series can
  508. // go quite strongly divergent in the early terms, before strongly
  509. // converging later. Various "guesstimates" have been tried to take account
  510. // of this, but they don't always work.... so instead set "n" to the
  511. // largest value that won't cause overflow later, and abort iteration
  512. // when the last accelerated term was small enough...
  513. //
  514. int n;
  515. #ifndef BOOST_MATH_NO_EXCEPTIONS
  516. try
  517. {
  518. #endif
  519. n = itrunc(T(tools::log_max_value<T>() / 6));
  520. #ifndef BOOST_MATH_NO_EXCEPTIONS
  521. }
  522. catch(...)
  523. {
  524. n = (std::numeric_limits<int>::max)();
  525. }
  526. #endif
  527. n = (std::min)(n, 1500);
  528. T d = pow(3 + sqrt(T(8)), T(n));
  529. d = (d + 1 / d) / 2;
  530. T b = -1;
  531. T c = -d;
  532. c = b - c;
  533. sum *= c;
  534. b = -n * n * b * 2;
  535. abs_err = ldexp(fabs(sum), -tools::digits<T>());
  536. while(j < n)
  537. {
  538. a_pow *= aa;
  539. dj_pow *= half_h_h / j;
  540. one_minus_dj_sum += dj_pow;
  541. term = one_minus_dj_sum * a_pow / (2 * j + 1);
  542. c = b - c;
  543. sum += c * term;
  544. abs_err += ldexp((std::max)(T(fabs(sum)), T(fabs(c*term))), -tools::digits<T>());
  545. b = (j + n) * (j - n) * b / ((j + T(0.5)) * (j + 1));
  546. ++j;
  547. //
  548. // Include an escape route to prevent calculating too many terms:
  549. //
  550. if((j > 10) && (fabs(sum * tools::epsilon<T>()) > fabs(c * term)))
  551. break;
  552. }
  553. abs_err += fabs(c * term);
  554. if(sum < 0) // sum must always be positive, if it's negative something really bad has happened:
  555. policies::raise_evaluation_error(function, 0, T(0), pol);
  556. return std::pair<T, T>((sum / d) / boost::math::constants::two_pi<T>(), abs_err / sum);
  557. }
  558. template<typename RealType, class Policy>
  559. 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&)
  560. {
  561. BOOST_MATH_STD_USING
  562. using namespace boost::math::constants;
  563. const unsigned short maxii = m+m+1;
  564. const RealType hs = h*h;
  565. const RealType as = -a*a;
  566. const RealType y = static_cast<RealType>(1) / hs;
  567. unsigned short ii = 1;
  568. RealType val = 0;
  569. RealType vi = a * exp( -ah*ah*half<RealType>() ) / root_two_pi<RealType>();
  570. RealType z = owens_t_znorm1(ah, pol)/h;
  571. RealType last_z = fabs(z);
  572. RealType lim = policies::get_epsilon<RealType, Policy>();
  573. while( true )
  574. {
  575. val += z;
  576. //
  577. // This series stops converging after a while, so put a limit
  578. // on how far we go before returning our best guess:
  579. //
  580. if((fabs(lim * val) > fabs(z)) || ((ii > maxii) && (fabs(z) > last_z)) || (z == 0))
  581. {
  582. val *= exp( -hs*half<RealType>() ) / root_two_pi<RealType>();
  583. break;
  584. } // if( maxii <= ii )
  585. last_z = fabs(z);
  586. z = y * ( vi - static_cast<RealType>(ii) * z );
  587. vi *= as;
  588. ii += 2;
  589. } // while( true )
  590. return val;
  591. } // RealType owens_t_T2(const RealType h, const RealType a, const unsigned short m, const RealType ah)
  592. template<typename RealType, class Policy>
  593. inline std::pair<RealType, RealType> owens_t_T2_accelerated(const RealType h, const RealType a, const RealType ah, const Policy& pol)
  594. {
  595. //
  596. // This is the same series as T2, but with acceleration applied.
  597. // Note that we have to be *very* careful to check that nothing bad
  598. // has happened during evaluation - this series will go divergent
  599. // and/or fail to alternate at a drop of a hat! :-(
  600. //
  601. BOOST_MATH_STD_USING
  602. using namespace boost::math::constants;
  603. const RealType hs = h*h;
  604. const RealType as = -a*a;
  605. const RealType y = static_cast<RealType>(1) / hs;
  606. unsigned short ii = 1;
  607. RealType val = 0;
  608. RealType vi = a * exp( -ah*ah*half<RealType>() ) / root_two_pi<RealType>();
  609. RealType z = boost::math::detail::owens_t_znorm1(ah, pol)/h;
  610. RealType last_z = fabs(z);
  611. //
  612. // Normally with this form of series acceleration we can calculate
  613. // up front how many terms will be required - based on the assumption
  614. // that each term decreases in size by a factor of 3. However,
  615. // that assumption does not apply here, as the underlying T1 series can
  616. // go quite strongly divergent in the early terms, before strongly
  617. // converging later. Various "guesstimates" have been tried to take account
  618. // of this, but they don't always work.... so instead set "n" to the
  619. // largest value that won't cause overflow later, and abort iteration
  620. // when the last accelerated term was small enough...
  621. //
  622. int n;
  623. #ifndef BOOST_MATH_NO_EXCEPTIONS
  624. try
  625. {
  626. #endif
  627. n = itrunc(RealType(tools::log_max_value<RealType>() / 6));
  628. #ifndef BOOST_MATH_NO_EXCEPTIONS
  629. }
  630. catch(...)
  631. {
  632. n = (std::numeric_limits<int>::max)();
  633. }
  634. #endif
  635. n = (std::min)(n, 1500);
  636. RealType d = pow(3 + sqrt(RealType(8)), RealType(n));
  637. d = (d + 1 / d) / 2;
  638. RealType b = -1;
  639. RealType c = -d;
  640. int s = 1;
  641. for(int k = 0; k < n; ++k)
  642. {
  643. //
  644. // Check for both convergence and whether the series has gone bad:
  645. //
  646. if(
  647. (fabs(z) > last_z) // Series has gone divergent, abort
  648. || (fabs(val) * tools::epsilon<RealType>() > fabs(c * s * z)) // Convergence!
  649. || (z * s < 0) // Series has stopped alternating - all bets are off - abort.
  650. )
  651. {
  652. break;
  653. }
  654. c = b - c;
  655. val += c * s * z;
  656. b = (k + n) * (k - n) * b / ((k + RealType(0.5)) * (k + 1));
  657. last_z = fabs(z);
  658. s = -s;
  659. z = y * ( vi - static_cast<RealType>(ii) * z );
  660. vi *= as;
  661. ii += 2;
  662. } // while( true )
  663. RealType err = fabs(c * z) / val;
  664. return std::pair<RealType, RealType>(val * exp( -hs*half<RealType>() ) / (d * root_two_pi<RealType>()), err);
  665. } // RealType owens_t_T2_accelerated(const RealType h, const RealType a, const RealType ah, const Policy&)
  666. template<typename RealType, typename Policy>
  667. inline RealType T4_mp(const RealType h, const RealType a, const Policy& pol)
  668. {
  669. BOOST_MATH_STD_USING
  670. const RealType hs = h*h;
  671. const RealType as = -a*a;
  672. unsigned short ii = 1;
  673. RealType ai = constants::one_div_two_pi<RealType>() * a * exp( -0.5*hs*(1.0-as) );
  674. RealType yi = 1.0;
  675. RealType val = 0.0;
  676. RealType lim = boost::math::policies::get_epsilon<RealType, Policy>();
  677. while( true )
  678. {
  679. RealType term = ai*yi;
  680. val += term;
  681. if((yi != 0) && (fabs(val * lim) > fabs(term)))
  682. break;
  683. ii += 2;
  684. yi = (1.0-hs*yi) / static_cast<RealType>(ii);
  685. ai *= as;
  686. if(ii > (std::min)(1500, (int)policies::get_max_series_iterations<Policy>()))
  687. policies::raise_evaluation_error("boost::math::owens_t<%1%>", 0, val, pol);
  688. } // while( true )
  689. return val;
  690. } // arg_type owens_t_T4(const arg_type h, const arg_type a, const unsigned short m)
  691. // This routine dispatches the call to one of six subroutines, depending on the values
  692. // of h and a.
  693. // preconditions: h >= 0, 0<=a<=1, ah=a*h
  694. //
  695. // Note there are different versions for different precisions....
  696. template<typename RealType, typename Policy>
  697. inline RealType owens_t_dispatch(const RealType h, const RealType a, const RealType ah, const Policy& pol, std::integral_constant<int, 64> const&)
  698. {
  699. // Simple main case for 64-bit precision or less, this is as per the Patefield-Tandy paper:
  700. BOOST_MATH_STD_USING
  701. //
  702. // Handle some special cases first, these are from
  703. // page 1077 of Owen's original paper:
  704. //
  705. if(h == 0)
  706. {
  707. return atan(a) * constants::one_div_two_pi<RealType>();
  708. }
  709. if(a == 0)
  710. {
  711. return 0;
  712. }
  713. if(a == 1)
  714. {
  715. return owens_t_znorm2(RealType(-h), pol) * owens_t_znorm2(h, pol) / 2;
  716. }
  717. // Rationale: when a>1 we call this routine with 1/a:
  718. BOOST_MATH_ASSERT(a <= 1);
  719. RealType val = 0; // avoid compiler warnings, 0 will be overwritten in any case
  720. const unsigned short icode = owens_t_compute_code(h, a);
  721. const unsigned short m = owens_t_get_order(icode, val /* just a dummy for the type */, pol);
  722. 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
  723. BOOST_MATH_ASSERT(icode < sizeof(meth) / sizeof(meth[0]));
  724. // determine the appropriate method, T1 ... T6
  725. switch( meth[icode] )
  726. {
  727. case 1: // T1
  728. val = owens_t_T1(h,a,m,pol);
  729. break;
  730. case 2: // T2
  731. typedef typename policies::precision<RealType, Policy>::type precision_type;
  732. typedef std::integral_constant<bool, (precision_type::value == 0) || (precision_type::value > 64)> tag_type;
  733. val = owens_t_T2(h, a, m, ah, pol, tag_type());
  734. break;
  735. case 3: // T3
  736. val = owens_t_T3(h,a,ah, pol);
  737. break;
  738. case 4: // T4
  739. val = owens_t_T4(h,a,m);
  740. break;
  741. case 5: // T5
  742. val = owens_t_T5(h,a, pol);
  743. break;
  744. case 6: // T6
  745. val = owens_t_T6(h,a, pol);
  746. break;
  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 std::integral_constant<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), pol) * owens_t_znorm2(h, pol) / 2;
  770. }
  771. if(a >= tools::max_value<RealType>())
  772. {
  773. return owens_t_znorm2(RealType(fabs(h)), pol);
  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_MATH_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_MATH_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_MATH_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_MATH_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_MATH_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_MATH_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_MATH_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_MATH_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_MATH_NO_EXCEPTIONS
  853. try
  854. {
  855. #endif
  856. return T4_mp(h, a, pol);
  857. #ifndef BOOST_MATH_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) < RealType(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, std::integral_constant<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 std::integral_constant<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, std::integral_constant<int, 64>());
  880. return owens_t_dispatch(h, a, ah, pol, std::integral_constant<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 std::integral_constant<int,
  888. precision_type::value <= 0 ? 0 :
  889. precision_type::value <= 64 ? 64 : 65
  890. > tag_type;
  891. return owens_t_dispatch(h, a, ah, pol, tag_type());
  892. }
  893. // compute Owen's T function, T(h,a), for arbitrary values of h and a
  894. template<typename RealType, class Policy>
  895. inline RealType owens_t(RealType h, RealType a, const Policy& pol)
  896. {
  897. BOOST_MATH_STD_USING
  898. // exploit that T(-h,a) == T(h,a)
  899. h = fabs(h);
  900. // Use equation (2) in the paper to remap the arguments
  901. // such that h>=0 and 0<=a<=1 for the call of the actual
  902. // computation routine.
  903. const RealType fabs_a = fabs(a);
  904. const RealType fabs_ah = fabs_a*h;
  905. RealType val = static_cast<RealType>(0.0f); // avoid compiler warnings, 0.0 will be overwritten in any case
  906. if(fabs_a <= 1)
  907. {
  908. val = owens_t_dispatch(h, fabs_a, fabs_ah, pol);
  909. } // if(fabs_a <= 1.0)
  910. else
  911. {
  912. if( h <= RealType(0.67) )
  913. {
  914. const RealType normh = owens_t_znorm1(h, pol);
  915. const RealType normah = owens_t_znorm1(fabs_ah, pol);
  916. val = static_cast<RealType>(1)/static_cast<RealType>(4) - normh*normah -
  917. owens_t_dispatch(fabs_ah, static_cast<RealType>(1 / fabs_a), h, pol);
  918. } // if( h <= 0.67 )
  919. else
  920. {
  921. const RealType normh = detail::owens_t_znorm2(h, pol);
  922. const RealType normah = detail::owens_t_znorm2(fabs_ah, pol);
  923. val = constants::half<RealType>()*(normh+normah) - normh*normah -
  924. owens_t_dispatch(fabs_ah, static_cast<RealType>(1 / fabs_a), h, pol);
  925. } // else [if( h <= 0.67 )]
  926. } // else [if(fabs_a <= 1)]
  927. // exploit that T(h,-a) == -T(h,a)
  928. if(a < 0)
  929. {
  930. return -val;
  931. } // if(a < 0)
  932. return val;
  933. } // RealType owens_t(RealType h, RealType a)
  934. } // namespace detail
  935. template <class T1, class T2, class Policy>
  936. inline typename tools::promote_args<T1, T2>::type owens_t(T1 h, T2 a, const Policy& pol)
  937. {
  938. typedef typename tools::promote_args<T1, T2>::type result_type;
  939. typedef typename policies::evaluation<result_type, Policy>::type value_type;
  940. 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%)");
  941. }
  942. template <class T1, class T2>
  943. inline typename tools::promote_args<T1, T2>::type owens_t(T1 h, T2 a)
  944. {
  945. return owens_t(h, a, policies::policy<>());
  946. }
  947. } // namespace math
  948. } // namespace boost
  949. #ifdef _MSC_VER
  950. #pragma warning(pop)
  951. #endif
  952. #endif
  953. // EOF