owens_t.hpp 47 KB

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