owens_t.hpp 49 KB

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