lambert_w.hpp 93 KB

12345678910111213141516171819202122232425262728293031323334353637383940414243444546474849505152535455565758596061626364656667686970717273747576777879808182838485868788899091929394959697989910010110210310410510610710810911011111211311411511611711811912012112212312412512612712812913013113213313413513613713813914014114214314414514614714814915015115215315415515615715815916016116216316416516616716816917017117217317417517617717817918018118218318418518618718818919019119219319419519619719819920020120220320420520620720820921021121221321421521621721821922022122222322422522622722822923023123223323423523623723823924024124224324424524624724824925025125225325425525625725825926026126226326426526626726826927027127227327427527627727827928028128228328428528628728828929029129229329429529629729829930030130230330430530630730830931031131231331431531631731831932032132232332432532632732832933033133233333433533633733833934034134234334434534634734834935035135235335435535635735835936036136236336436536636736836937037137237337437537637737837938038138238338438538638738838939039139239339439539639739839940040140240340440540640740840941041141241341441541641741841942042142242342442542642742842943043143243343443543643743843944044144244344444544644744844945045145245345445545645745845946046146246346446546646746846947047147247347447547647747847948048148248348448548648748848949049149249349449549649749849950050150250350450550650750850951051151251351451551651751851952052152252352452552652752852953053153253353453553653753853954054154254354454554654754854955055155255355455555655755855956056156256356456556656756856957057157257357457557657757857958058158258358458558658758858959059159259359459559659759859960060160260360460560660760860961061161261361461561661761861962062162262362462562662762862963063163263363463563663763863964064164264364464564664764864965065165265365465565665765865966066166266366466566666766866967067167267367467567667767867968068168268368468568668768868969069169269369469569669769869970070170270370470570670770870971071171271371471571671771871972072172272372472572672772872973073173273373473573673773873974074174274374474574674774874975075175275375475575675775875976076176276376476576676776876977077177277377477577677777877978078178278378478578678778878979079179279379479579679779879980080180280380480580680780880981081181281381481581681781881982082182282382482582682782882983083183283383483583683783883984084184284384484584684784884985085185285385485585685785885986086186286386486586686786886987087187287387487587687787887988088188288388488588688788888989089189289389489589689789889990090190290390490590690790890991091191291391491591691791891992092192292392492592692792892993093193293393493593693793893994094194294394494594694794894995095195295395495595695795895996096196296396496596696796896997097197297397497597697797897998098198298398498598698798898999099199299399499599699799899910001001100210031004100510061007100810091010101110121013101410151016101710181019102010211022102310241025102610271028102910301031103210331034103510361037103810391040104110421043104410451046104710481049105010511052105310541055105610571058105910601061106210631064106510661067106810691070107110721073107410751076107710781079108010811082108310841085108610871088108910901091109210931094109510961097109810991100110111021103110411051106110711081109111011111112111311141115111611171118111911201121112211231124112511261127112811291130113111321133113411351136113711381139114011411142114311441145114611471148114911501151115211531154115511561157115811591160116111621163116411651166116711681169117011711172117311741175117611771178117911801181118211831184118511861187118811891190119111921193119411951196119711981199120012011202120312041205120612071208120912101211121212131214121512161217121812191220122112221223122412251226122712281229123012311232123312341235123612371238123912401241124212431244124512461247124812491250125112521253125412551256125712581259126012611262126312641265126612671268126912701271127212731274127512761277127812791280128112821283128412851286128712881289129012911292129312941295129612971298129913001301130213031304130513061307130813091310131113121313131413151316131713181319132013211322132313241325132613271328132913301331133213331334133513361337133813391340134113421343134413451346134713481349135013511352135313541355135613571358135913601361136213631364136513661367136813691370137113721373137413751376137713781379138013811382138313841385138613871388138913901391139213931394139513961397139813991400140114021403140414051406140714081409141014111412141314141415141614171418141914201421142214231424142514261427142814291430143114321433143414351436143714381439144014411442144314441445144614471448144914501451145214531454145514561457145814591460146114621463146414651466146714681469147014711472147314741475147614771478147914801481148214831484148514861487148814891490149114921493149414951496149714981499150015011502150315041505150615071508150915101511151215131514151515161517151815191520152115221523152415251526152715281529153015311532153315341535153615371538153915401541154215431544154515461547154815491550155115521553155415551556155715581559156015611562156315641565156615671568156915701571157215731574157515761577157815791580158115821583158415851586158715881589159015911592159315941595159615971598159916001601160216031604160516061607160816091610161116121613161416151616161716181619162016211622162316241625162616271628162916301631163216331634163516361637163816391640164116421643164416451646164716481649165016511652165316541655165616571658165916601661166216631664166516661667166816691670167116721673167416751676167716781679168016811682168316841685168616871688168916901691169216931694169516961697169816991700170117021703170417051706170717081709171017111712171317141715171617171718171917201721172217231724172517261727172817291730173117321733173417351736173717381739174017411742174317441745174617471748174917501751175217531754175517561757175817591760176117621763176417651766176717681769177017711772177317741775177617771778177917801781178217831784178517861787178817891790179117921793179417951796179717981799180018011802180318041805180618071808180918101811181218131814181518161817181818191820182118221823182418251826182718281829183018311832183318341835183618371838183918401841184218431844184518461847184818491850185118521853185418551856185718581859186018611862186318641865186618671868186918701871187218731874187518761877187818791880188118821883188418851886188718881889189018911892189318941895189618971898189919001901190219031904190519061907190819091910191119121913191419151916191719181919192019211922192319241925192619271928192919301931193219331934193519361937193819391940194119421943194419451946194719481949195019511952195319541955195619571958195919601961196219631964196519661967196819691970197119721973197419751976197719781979198019811982198319841985198619871988198919901991199219931994199519961997199819992000200120022003200420052006200720082009201020112012201320142015201620172018201920202021202220232024202520262027202820292030203120322033203420352036203720382039204020412042204320442045204620472048204920502051205220532054205520562057205820592060206120622063206420652066206720682069207020712072207320742075207620772078207920802081208220832084208520862087208820892090209120922093209420952096209720982099210021012102210321042105210621072108210921102111211221132114211521162117211821192120212121222123212421252126212721282129213021312132213321342135213621372138213921402141214221432144214521462147214821492150215121522153215421552156215721582159216021612162216321642165216621672168216921702171217221732174217521762177217821792180218121822183218421852186
  1. // Copyright John Maddock 2017.
  2. // Copyright Paul A. Bristow 2016, 2017, 2018.
  3. // Copyright Nicholas Thompson 2018
  4. // Distributed under the Boost Software License, Version 1.0.
  5. // (See accompanying file LICENSE_1_0.txt or
  6. // copy at http ://www.boost.org/LICENSE_1_0.txt).
  7. #ifndef BOOST_MATH_SF_LAMBERT_W_HPP
  8. #define BOOST_MATH_SF_LAMBERT_W_HPP
  9. #ifdef _MSC_VER
  10. #pragma warning(disable : 4127)
  11. #endif
  12. /*
  13. Implementation of an algorithm for the Lambert W0 and W-1 real-only functions.
  14. This code is based in part on the algorithm by
  15. Toshio Fukushima,
  16. "Precise and fast computation of Lambert W-functions without transcendental function evaluations",
  17. J.Comp.Appl.Math. 244 (2013) 77-89,
  18. and on a C/C++ version by Darko Veberic, darko.veberic@ijs.si
  19. based on the Fukushima algorithm and Toshio Fukushima's FORTRAN version of his algorithm.
  20. First derivative of Lambert_w is derived from
  21. Princeton Companion to Applied Mathematics, 'The Lambert-W function', Section 1.3: Series and Generating Functions.
  22. */
  23. /*
  24. TODO revise this list of macros.
  25. Some macros that will show some (or much) diagnostic values if #defined.
  26. //[boost_math_instrument_lambert_w_macros
  27. // #define-able macros
  28. BOOST_MATH_INSTRUMENT_LAMBERT_W_HALLEY // Halley refinement diagnostics.
  29. BOOST_MATH_INSTRUMENT_LAMBERT_W_PRECISION // Precision.
  30. BOOST_MATH_INSTRUMENT_LAMBERT_WM1 // W1 branch diagnostics.
  31. BOOST_MATH_INSTRUMENT_LAMBERT_WM1_HALLEY // Halley refinement diagnostics only for W-1 branch.
  32. BOOST_MATH_INSTRUMENT_LAMBERT_WM1_TINY // K > 64, z > -1.0264389699511303e-26
  33. BOOST_MATH_INSTRUMENT_LAMBERT_WM1_LOOKUP // Show results from W-1 lookup table.
  34. BOOST_MATH_INSTRUMENT_LAMBERT_W_SCHROEDER // Schroeder refinement diagnostics.
  35. BOOST_MATH_INSTRUMENT_LAMBERT_W_TERMS // Number of terms used for near-singularity series.
  36. BOOST_MATH_INSTRUMENT_LAMBERT_W_SINGULARITY_SERIES // Show evaluation of series near branch singularity.
  37. BOOST_MATH_INSTRUMENT_LAMBERT_W_SMALL_Z_SERIES
  38. BOOST_MATH_INSTRUMENT_LAMBERT_W_SMALL_Z_SERIES_ITERATIONS // Show evaluation of series for small z.
  39. //] [/boost_math_instrument_lambert_w_macros]
  40. */
  41. #include <boost/math/policies/error_handling.hpp>
  42. #include <boost/math/policies/policy.hpp>
  43. #include <boost/math/tools/promotion.hpp>
  44. #include <boost/math/special_functions/fpclassify.hpp>
  45. #include <boost/math/special_functions/log1p.hpp> // for log (1 + x)
  46. #include <boost/math/constants/constants.hpp> // For exp_minus_one == 3.67879441171442321595523770161460867e-01.
  47. #include <boost/math/special_functions/pow.hpp> // powers with compile time exponent, used in arbitrary precision code.
  48. #include <boost/math/tools/series.hpp> // series functor.
  49. //#include <boost/math/tools/polynomial.hpp> // polynomial.
  50. #include <boost/math/tools/rational.hpp> // evaluate_polynomial.
  51. #include <boost/mpl/int.hpp>
  52. #include <boost/type_traits/is_integral.hpp>
  53. #include <boost/math/tools/precision.hpp> // boost::math::tools::max_value().
  54. #include <boost/math/tools/big_constant.hpp>
  55. #include <limits>
  56. #include <cmath>
  57. #include <limits>
  58. #include <exception>
  59. // Needed for testing and diagnostics only.
  60. #include <iostream>
  61. #include <typeinfo>
  62. #include <boost/math/special_functions/next.hpp> // For float_distance.
  63. typedef double lookup_t; // Type for lookup table (double or float, or even long double?)
  64. //#include "J:\Cpp\Misc\lambert_w_lookup_table_generator\lambert_w_lookup_table.ipp"
  65. // #include "lambert_w_lookup_table.ipp" // Boost.Math version.
  66. #include <boost/math/special_functions/detail/lambert_w_lookup_table.ipp>
  67. #if defined(__GNUC__) && defined(BOOST_MATH_USE_FLOAT128)
  68. //
  69. // This is the only way we can avoid
  70. // warning: non-standard suffix on floating constant [-Wpedantic]
  71. // when building with -Wall -pedantic. Neither __extension__
  72. // nor #pragma dianostic ignored work :(
  73. //
  74. #pragma GCC system_header
  75. #endif
  76. namespace boost {
  77. namespace math {
  78. namespace lambert_w_detail {
  79. //! \brief Applies a single Halley step to make a better estimate of Lambert W.
  80. //! \details Used the simplified formulae obtained from
  81. //! http://www.wolframalpha.com/input/?i=%5B2(z+exp(z)-w)+d%2Fdx+(z+exp(z)-w)%5D+%2F+%5B2+(d%2Fdx+(z+exp(z)-w))%5E2+-+(z+exp(z)-w)+d%5E2%2Fdx%5E2+(z+exp(z)-w)%5D
  82. //! [2(z exp(z)-w) d/dx (z exp(z)-w)] / [2 (d/dx (z exp(z)-w))^2 - (z exp(z)-w) d^2/dx^2 (z exp(z)-w)]
  83. //! \tparam T floating-point (or fixed-point) type.
  84. //! \param w_est Lambert W estimate.
  85. //! \param z Argument z for Lambert_w function.
  86. //! \returns New estimate of Lambert W, hopefully improved.
  87. //!
  88. template <class T>
  89. inline T lambert_w_halley_step(T w_est, const T z)
  90. {
  91. BOOST_MATH_STD_USING
  92. T e = exp(w_est);
  93. w_est -= 2 * (w_est + 1) * (e * w_est - z) / (z * (w_est + 2) + e * (w_est * (w_est + 2) + 2));
  94. return w_est;
  95. } // template <class T> lambert_w_halley_step(T w_est, T z)
  96. //! \brief Halley iterate to refine Lambert_w estimate,
  97. //! taking at least one Halley_step.
  98. //! Repeat Halley steps until the *last step* had fewer than half the digits wrong,
  99. //! the step we've just taken should have been sufficient to have completed the iteration.
  100. //! \tparam T floating-point (or fixed-point) type.
  101. //! \param z Argument z for Lambert_w function.
  102. //! \param w_est Lambert w estimate.
  103. template <class T>
  104. inline
  105. T lambert_w_halley_iterate(T w_est, const T z)
  106. {
  107. BOOST_MATH_STD_USING
  108. static const T max_diff = boost::math::tools::root_epsilon<T>() * fabs(w_est);
  109. T w_new = lambert_w_halley_step(w_est, z);
  110. T diff = fabs(w_est - w_new);
  111. while (diff > max_diff)
  112. {
  113. w_est = w_new;
  114. w_new = lambert_w_halley_step(w_est, z);
  115. diff = fabs(w_est - w_new);
  116. }
  117. return w_new;
  118. } // template <class T> lambert_w_halley_iterate(T w_est, T z)
  119. // Two Halley function versions that either
  120. // single step (if mpl::false_) or iterate (if mpl::true_).
  121. // Selected at compile-time using parameter 3.
  122. template <class T>
  123. inline
  124. T lambert_w_maybe_halley_iterate(T z, T w, mpl::false_ const&)
  125. {
  126. return lambert_w_halley_step(z, w); // Single step.
  127. }
  128. template <class T>
  129. inline
  130. T lambert_w_maybe_halley_iterate(T z, T w, mpl::true_ const&)
  131. {
  132. return lambert_w_halley_iterate(z, w); // Iterate steps.
  133. }
  134. //! maybe_reduce_to_double function,
  135. //! Two versions that have a compile-time option to
  136. //! reduce argument z to double precision (if mpl::true_).
  137. //! Version is selected at compile-time using parameter 2.
  138. template <class T>
  139. inline
  140. double maybe_reduce_to_double(const T& z, const mpl::true_&)
  141. {
  142. return static_cast<double>(z); // Reduce to double precision.
  143. }
  144. template <class T>
  145. inline
  146. T maybe_reduce_to_double(const T& z, const mpl::false_&)
  147. { // Don't reduce to double.
  148. return z;
  149. }
  150. template <class T>
  151. inline
  152. double must_reduce_to_double(const T& z, const mpl::true_&)
  153. {
  154. return static_cast<double>(z); // Reduce to double precision.
  155. }
  156. template <class T>
  157. inline
  158. double must_reduce_to_double(const T& z, const mpl::false_&)
  159. { // try a lexical_cast and hope for the best:
  160. return boost::lexical_cast<double>(z);
  161. }
  162. //! \brief Schroeder method, fifth-order update formula,
  163. //! \details See T. Fukushima page 80-81, and
  164. //! A. Householder, The Numerical Treatment of a Single Nonlinear Equation,
  165. //! McGraw-Hill, New York, 1970, section 4.4.
  166. //! Fukushima algorithm switches to @c schroeder_update after pre-computed bisections,
  167. //! chosen to ensure that the result will be achieve the +/- 10 epsilon target.
  168. //! \param w Lambert w estimate from bisection or series.
  169. //! \param y bracketing value from bisection.
  170. //! \returns Refined estimate of Lambert w.
  171. // Schroeder refinement, called unless NOT required by precision policy.
  172. template<typename T>
  173. inline
  174. T schroeder_update(const T w, const T y)
  175. {
  176. // Compute derivatives using 5th order Schroeder refinement.
  177. // Since this is the final step, it will always use the highest precision type T.
  178. // Example of Call:
  179. // result = schroeder_update(w, y);
  180. //where
  181. // w is estimate of Lambert W (from bisection or series).
  182. // y is z * e^-w.
  183. BOOST_MATH_STD_USING // Aid argument dependent lookup of abs.
  184. #ifdef BOOST_MATH_INSTRUMENT_LAMBERT_W_SCHROEDER
  185. std::streamsize saved_precision = std::cout.precision(std::numeric_limits<T>::max_digits10);
  186. using boost::math::float_distance;
  187. T fd = float_distance<T>(w, y);
  188. std::cout << "Schroder ";
  189. if (abs(fd) < 214748000.)
  190. {
  191. std::cout << " Distance = "<< static_cast<int>(fd);
  192. }
  193. else
  194. {
  195. std::cout << "Difference w - y = " << (w - y) << ".";
  196. }
  197. std::cout << std::endl;
  198. #endif // BOOST_MATH_INSTRUMENT_LAMBERT_W_SCHROEDER
  199. // Fukushima equation 18, page 6.
  200. const T f0 = w - y; // f0 = w - y.
  201. const T f1 = 1 + y; // f1 = df/dW
  202. const T f00 = f0 * f0;
  203. const T f11 = f1 * f1;
  204. const T f0y = f0 * y;
  205. const T result =
  206. w - 4 * f0 * (6 * f1 * (f11 + f0y) + f00 * y) /
  207. (f11 * (24 * f11 + 36 * f0y) +
  208. f00 * (6 * y * y + 8 * f1 * y + f0y)); // Fukushima Page 81, equation 21 from equation 20.
  209. #ifdef BOOST_MATH_INSTRUMENT_LAMBERT_W_SCHROEDER
  210. std::cout << "Schroeder refined " << w << " " << y << ", difference " << w-y << ", change " << w - result << ", to result " << result << std::endl;
  211. std::cout.precision(saved_precision); // Restore.
  212. #endif // BOOST_MATH_INSTRUMENT_LAMBERT_W_SCHROEDER
  213. return result;
  214. } // template<typename T = double> T schroeder_update(const T w, const T y)
  215. //! \brief Series expansion used near the singularity/branch point z = -exp(-1) = -3.6787944.
  216. //! Wolfram InverseSeries[Series[sqrt[2(p Exp[1 + p] + 1)], {p,-1, 20}]]
  217. //! Wolfram command used to obtain 40 series terms at 50 decimal digit precision was
  218. //! N[InverseSeries[Series[Sqrt[2(p Exp[1 + p] + 1)], { p,-1,40 }]], 50]
  219. //! -1+p-p^2/3+(11 p^3)/72-(43 p^4)/540+(769 p^5)/17280-(221 p^6)/8505+(680863 p^7)/43545600 ...
  220. //! Decimal values of specifications for built-in floating-point types below
  221. //! are at least 21 digits precision == max_digits10 for long double.
  222. //! Longer decimal digits strings are rationals evaluated using Wolfram.
  223. template<typename T>
  224. T lambert_w_singularity_series(const T p)
  225. {
  226. #ifdef BOOST_MATH_INSTRUMENT_LAMBERT_W_SINGULARITY_SERIES
  227. std::size_t saved_precision = std::cout.precision(3);
  228. std::cout << "Singularity_series Lambert_w p argument = " << p << std::endl;
  229. std::cout
  230. //<< "Argument Type = " << typeid(T).name()
  231. //<< ", max_digits10 = " << std::numeric_limits<T>::max_digits10
  232. //<< ", epsilon = " << std::numeric_limits<T>::epsilon()
  233. << std::endl;
  234. std::cout.precision(saved_precision);
  235. #endif // BOOST_MATH_INSTRUMENT_LAMBERT_W_SINGULARITY_SERIES
  236. static const T q[] =
  237. {
  238. -static_cast<T>(1), // j0
  239. +T(1), // j1
  240. -T(1) / 3, // 1/3 j2
  241. +T(11) / 72, // 0.152777777777777778, // 11/72 j3
  242. -T(43) / 540, // 0.0796296296296296296, // 43/540 j4
  243. +T(769) / 17280, // 0.0445023148148148148, j5
  244. -T(221) / 8505, // 0.0259847148736037625, j6
  245. //+T(0.0156356325323339212L), // j7
  246. //+T(0.015635632532333921222810111699000587889476778365667L), // j7 from Wolfram N[680863/43545600, 50]
  247. +T(680863uLL) / 43545600uLL, // +0.0156356325323339212, j7
  248. //-T(0.00961689202429943171L), // j8
  249. -T(1963uLL) / 204120uLL, // 0.00961689202429943171, j8
  250. //-T(0.0096168920242994317068391142465216539290613364687439L), // j8 from Wolfram N[1963/204120, 50]
  251. +T(226287557uLL) / 37623398400uLL, // 0.00601454325295611786, j9
  252. -T(5776369uLL) / 1515591000uLL, // 0.00381129803489199923, j10
  253. //+T(0.00244087799114398267L), j11 0.0024408779911439826658968585286437530215699919795550
  254. +T(169709463197uLL) / 69528040243200uLL, // j11
  255. // -T(0.00157693034468678425L), // j12 -0.0015769303446867842539234095399314115973161850314723
  256. -T(1118511313uLL) / 709296588000uLL, // j12
  257. +T(667874164916771uLL) / 650782456676352000uLL, // j13
  258. //+T(0.00102626332050760715L), // j13 0.0010262633205076071544375481533906861056468041465973
  259. -T(500525573uLL) / 744761417400uLL, // j14
  260. // -T(0.000672061631156136204L), j14
  261. //+T(1003663334225097487uLL) / 234281684403486720000uLL, // j15 0.00044247306181462090993020760858473726479232802068800 error C2177: constant too big
  262. //+T(0.000442473061814620910L, // j15
  263. BOOST_MATH_BIG_CONSTANT(T, 64, +0.000442473061814620910), // j15
  264. // -T(0.000292677224729627445L), // j16
  265. BOOST_MATH_BIG_CONSTANT(T, 64, -0.000292677224729627445), // j16
  266. //+T(0.000194387276054539318L), // j17
  267. BOOST_MATH_BIG_CONSTANT(T, 64, 0.000194387276054539318), // j17
  268. //-T(0.000129574266852748819L), // j18
  269. BOOST_MATH_BIG_CONSTANT(T, 64, -0.000129574266852748819), // j18
  270. //+T(0.0000866503580520812717L), // j19 N[+1150497127780071399782389/13277465363600276402995200000, 50] 0.000086650358052081271660451590462390293190597827783288
  271. BOOST_MATH_BIG_CONSTANT(T, 64, +0.0000866503580520812717), // j19
  272. //-T(0.0000581136075044138168L) // j20 N[2853534237182741069/49102686267859224000000, 50] 0.000058113607504413816772205464778828177256611844221913
  273. // -T(2853534237182741069uLL) / 49102686267859224000000uLL // j20 // error C2177: constant too big,
  274. // so must use BOOST_MATH_BIG_CONSTANT(T, ) format in hope of using suffix Q for quad or decimal digits string for others.
  275. //-T(0.000058113607504413816772205464778828177256611844221913L), // j20 N[2853534237182741069/49102686267859224000000, 50] 0.000058113607504413816772205464778828177256611844221913
  276. BOOST_MATH_BIG_CONSTANT(T, 113, -0.000058113607504413816772205464778828177256611844221913) // j20 - last used by Fukushima
  277. // More terms don't seem to give any improvement (worse in fact) and are not use for many z values.
  278. //BOOST_MATH_BIG_CONSTANT(T, +0.000039076684867439051635395583044527492132109160553593), // j21
  279. //BOOST_MATH_BIG_CONSTANT(T, -0.000026338064747231098738584082718649443078703982217219), // j22
  280. //BOOST_MATH_BIG_CONSTANT(T, +0.000017790345805079585400736282075184540383274460464169), // j23
  281. //BOOST_MATH_BIG_CONSTANT(T, -0.000012040352739559976942274116578992585158113153190354), // j24
  282. //BOOST_MATH_BIG_CONSTANT(T, +8.1635319824966121713827512573558687050675701559448E-6), // j25
  283. //BOOST_MATH_BIG_CONSTANT(T, -5.5442032085673591366657251660804575198155559225316E-6) // j26
  284. // -T(5.5442032085673591366657251660804575198155559225316E-6L) // j26
  285. // 21 to 26 Added for long double.
  286. }; // static const T q[]
  287. /*
  288. // Temporary copy of original double values for comparison; these are reproduced well.
  289. static const T q[] =
  290. {
  291. -1L, // j0
  292. +1L, // j1
  293. -0.333333333333333333L, // 1/3 j2
  294. +0.152777777777777778L, // 11/72 j3
  295. -0.0796296296296296296L, // 43/540
  296. +0.0445023148148148148L,
  297. -0.0259847148736037625L,
  298. +0.0156356325323339212L,
  299. -0.00961689202429943171L,
  300. +0.00601454325295611786L,
  301. -0.00381129803489199923L,
  302. +0.00244087799114398267L,
  303. -0.00157693034468678425L,
  304. +0.00102626332050760715L,
  305. -0.000672061631156136204L,
  306. +0.000442473061814620910L,
  307. -0.000292677224729627445L,
  308. +0.000194387276054539318L,
  309. -0.000129574266852748819L,
  310. +0.0000866503580520812717L,
  311. -0.0000581136075044138168L // j20
  312. };
  313. */
  314. // Decide how many series terms to use, increasing as z approaches the singularity,
  315. // balancing run-time versus computational noise from round-off.
  316. // In practice, we truncate the series expansion at a certain order.
  317. // If the order is too large, not only does the amount of computation increase,
  318. // but also the round-off errors accumulate.
  319. // See Fukushima equation 35, page 85 for logic of choice of number of series terms.
  320. BOOST_MATH_STD_USING // Aid argument dependent lookup (ADL) of abs.
  321. const T absp = abs(p);
  322. #ifdef BOOST_MATH_INSTRUMENT_LAMBERT_W_TERMS
  323. {
  324. int terms = 20; // Default to using all terms.
  325. if (absp < 0.001150)
  326. { // Very near singularity.
  327. terms = 6;
  328. }
  329. else if (absp < 0.0766)
  330. { // Near singularity.
  331. terms = 10;
  332. }
  333. std::streamsize saved_precision = std::cout.precision(3);
  334. std::cout << "abs(p) = " << absp << ", terms = " << terms << std::endl;
  335. std::cout.precision(saved_precision);
  336. }
  337. #endif // BOOST_MATH_INSTRUMENT_LAMBERT_W_TERMS
  338. if (absp < 0.01159)
  339. { // Only 6 near-singularity series terms are useful.
  340. return
  341. -1 +
  342. p * (1 +
  343. p * (q[2] +
  344. p * (q[3] +
  345. p * (q[4] +
  346. p * (q[5] +
  347. p * q[6]
  348. )))));
  349. }
  350. else if (absp < 0.0766) // Use 10 near-singularity series terms.
  351. { // Use 10 near-singularity series terms.
  352. return
  353. -1 +
  354. p * (1 +
  355. p * (q[2] +
  356. p * (q[3] +
  357. p * (q[4] +
  358. p * (q[5] +
  359. p * (q[6] +
  360. p * (q[7] +
  361. p * (q[8] +
  362. p * (q[9] +
  363. p * q[10]
  364. )))))))));
  365. }
  366. else
  367. { // Use all 20 near-singularity series terms.
  368. return
  369. -1 +
  370. p * (1 +
  371. p * (q[2] +
  372. p * (q[3] +
  373. p * (q[4] +
  374. p * (q[5] +
  375. p * (q[6] +
  376. p * (q[7] +
  377. p * (q[8] +
  378. p * (q[9] +
  379. p * (q[10] +
  380. p * (q[11] +
  381. p * (q[12] +
  382. p * (q[13] +
  383. p * (q[14] +
  384. p * (q[15] +
  385. p * (q[16] +
  386. p * (q[17] +
  387. p * (q[18] +
  388. p * (q[19] +
  389. p * q[20] // Last Fukushima term.
  390. )))))))))))))))))));
  391. // + // more terms for more precise T: long double ...
  392. //// but makes almost no difference, so don't use more terms?
  393. // p*q[21] +
  394. // p*q[22] +
  395. // p*q[23] +
  396. // p*q[24] +
  397. // p*q[25]
  398. // )))))))))))))))))));
  399. }
  400. } // template<typename T = double> T lambert_w_singularity_series(const T p)
  401. /////////////////////////////////////////////////////////////////////////////////////////////
  402. //! \brief Series expansion used near zero (abs(z) < 0.05).
  403. //! \details
  404. //! Coefficients of the inverted series expansion of the Lambert W function around z = 0.
  405. //! Tosio Fukushima always uses all 17 terms of a Taylor series computed using Wolfram with
  406. //! InverseSeries[Series[z Exp[z],{z,0,17}]]
  407. //! Tosio Fukushima / Journal of Computational and Applied Mathematics 244 (2013) page 86.
  408. //! Decimal values of specifications for built-in floating-point types below
  409. //! are 21 digits precision == max_digits10 for long double.
  410. //! Care! Some coefficients might overflow some fixed_point types.
  411. //! This version is intended to allow use by user-defined types
  412. //! like Boost.Multiprecision quad and cpp_dec_float types.
  413. //! The three specializations below for built-in float, double
  414. //! (and perhaps long double) will be chosen in preference for these types.
  415. //! This version uses rationals computed by Wolfram as far as possible,
  416. //! limited by maximum size of uLL integers.
  417. //! For higher term, uses decimal digit strings computed by Wolfram up to the maximum possible using uLL rationals,
  418. //! and then higher coefficients are computed as necessary using function lambert_w0_small_z_series_term
  419. //! until the precision required by the policy is achieved.
  420. //! InverseSeries[Series[z Exp[z],{z,0,34}]] also computed.
  421. // Series evaluation for LambertW(z) as z -> 0.
  422. // See http://functions.wolfram.com/ElementaryFunctions/ProductLog/06/01/01/0003/
  423. // http://functions.wolfram.com/ElementaryFunctions/ProductLog/06/01/01/0003/MainEq1.L.gif
  424. //! \brief lambert_w0_small_z uses a tag_type to select a variant depending on the size of the type.
  425. //! The Lambert W is computed by lambert_w0_small_z for small z.
  426. //! The cutoff for z smallness determined by Tosio Fukushima by trial and error is (abs(z) < 0.05),
  427. //! but the optimum might be a function of the size of the type of z.
  428. //! \details
  429. //! The tag_type selection is based on the value @c std::numeric_limits<T>::max_digits10.
  430. //! This allows distinguishing between long double types that commonly vary between 64 and 80-bits,
  431. //! and also compilers that have a float type using 64 bits and/or long double using 128-bits.
  432. //! It assumes that max_digits10 is defined correctly or this might fail to make the correct selection.
  433. //! causing very small differences in computing lambert_w that would be very difficult to detect and diagnose.
  434. //! Cannot switch on @c std::numeric_limits<>::max() because comparison values may overflow the compiler limit.
  435. //! Cannot switch on @c std::numeric_limits<long double>::max_exponent10()
  436. //! because both 80 and 128 bit floating-point implementations use 11 bits for the exponent.
  437. //! So must rely on @c std::numeric_limits<long double>::max_digits10.
  438. //! Specialization of float zero series expansion used for small z (abs(z) < 0.05).
  439. //! Specializations of lambert_w0_small_z for built-in types.
  440. //! These specializations should be chosen in preference to T version.
  441. //! For example: lambert_w0_small_z(0.001F) should use the float version.
  442. //! (Parameter Policy is not used by built-in types when all terms are used during an inline computation,
  443. //! but for the tag_type selection to work, they all must include Policy in their signature.
  444. // Forward declaration of variants of lambert_w0_small_z.
  445. template <class T, class Policy>
  446. T lambert_w0_small_z(T x, const Policy&, boost::mpl::int_<0> const&); // for float (32-bit) type.
  447. template <class T, class Policy>
  448. T lambert_w0_small_z(T x, const Policy&, boost::mpl::int_<1> const&); // for double (64-bit) type.
  449. template <class T, class Policy>
  450. T lambert_w0_small_z(T x, const Policy&, boost::mpl::int_<2> const&); // for long double (double extended 80-bit) type.
  451. template <class T, class Policy>
  452. T lambert_w0_small_z(T x, const Policy&, boost::mpl::int_<3> const&); // for long double (128-bit) type.
  453. template <class T, class Policy>
  454. T lambert_w0_small_z(T x, const Policy&, boost::mpl::int_<4> const&); // for float128 quadmath Q type.
  455. template <class T, class Policy>
  456. T lambert_w0_small_z(T x, const Policy&, boost::mpl::int_<5> const&); // Generic multiprecision T.
  457. // Set tag_type depending on max_digits10.
  458. template <class T, class Policy>
  459. T lambert_w0_small_z(T x, const Policy& pol)
  460. { //std::numeric_limits<T>::max_digits10 == 36 ? 3 : // 128-bit long double.
  461. typedef boost::mpl::int_
  462. <
  463. std::numeric_limits<T>::is_specialized == 0 ? 5 :
  464. #ifndef BOOST_NO_CXX11_NUMERIC_LIMITS
  465. std::numeric_limits<T>::max_digits10 <= 9 ? 0 : // for float 32-bit.
  466. std::numeric_limits<T>::max_digits10 <= 17 ? 1 : // for double 64-bit.
  467. std::numeric_limits<T>::max_digits10 <= 22 ? 2 : // for 80-bit double extended.
  468. std::numeric_limits<T>::max_digits10 < 37 ? 4 // for both 128-bit long double (3) and 128-bit quad suffix Q type (4).
  469. #else
  470. std::numeric_limits<T>::radix != 2 ? 5 :
  471. std::numeric_limits<T>::digits <= 24 ? 0 : // for float 32-bit.
  472. std::numeric_limits<T>::digits <= 53 ? 1 : // for double 64-bit.
  473. std::numeric_limits<T>::digits <= 64 ? 2 : // for 80-bit double extended.
  474. std::numeric_limits<T>::digits <= 113 ? 4 // for both 128-bit long double (3) and 128-bit quad suffix Q type (4).
  475. #endif
  476. : 5 // All Generic multiprecision types.
  477. > tag_type;
  478. // std::cout << "\ntag type = " << tag_type << std::endl; // error C2275: 'tag_type': illegal use of this type as an expression.
  479. return lambert_w0_small_z(x, pol, tag_type());
  480. } // template <class T> T lambert_w0_small_z(T x)
  481. //! Specialization of float (32-bit) series expansion used for small z (abs(z) < 0.05).
  482. // Only 9 Coefficients are computed to 21 decimal digits precision, ample for 32-bit float used by most platforms.
  483. // Taylor series coefficients used are computed by Wolfram to 50 decimal digits using instruction
  484. // N[InverseSeries[Series[z Exp[z],{z,0,34}]],50],
  485. // as proposed by Tosio Fukushima and implemented by Darko Veberic.
  486. template <class T, class Policy>
  487. T lambert_w0_small_z(T z, const Policy&, boost::mpl::int_<0> const&)
  488. {
  489. #ifdef BOOST_MATH_INSTRUMENT_LAMBERT_W_SMALL_Z_SERIES
  490. std::streamsize prec = std::cout.precision(std::numeric_limits<T>::max_digits10); // Save.
  491. std::cout << "\ntag_type 0 float lambert_w0_small_z called with z = " << z << " using " << 9 << " terms of precision "
  492. << std::numeric_limits<float>::max_digits10 << " decimal digits. " << std::endl;
  493. #endif // BOOST_MATH_INSTRUMENT_LAMBERT_W_SMALL_Z_SERIES
  494. T result =
  495. z * (1 - // j1 z^1 term = 1
  496. z * (1 - // j2 z^2 term = -1
  497. z * (static_cast<float>(3uLL) / 2uLL - // 3/2 // j3 z^3 term = 1.5.
  498. z * (2.6666666666666666667F - // 8/3 // j4
  499. z * (5.2083333333333333333F - // -125/24 // j5
  500. z * (10.8F - // j6
  501. z * (23.343055555555555556F - // j7
  502. z * (52.012698412698412698F - // j8
  503. z * 118.62522321428571429F)))))))); // j9
  504. #ifdef BOOST_MATH_INSTRUMENT_LAMBERT_W_SMALL_Z_SERIES
  505. std::cout << "return w = " << result << std::endl;
  506. std::cout.precision(prec); // Restore.
  507. #endif // BOOST_MATH_INSTRUMENT_LAMBERT_W_SMALL_Z_SERIES
  508. return result;
  509. } // template <class T> T lambert_w0_small_z(T x, mpl::int_<0> const&)
  510. //! Specialization of double (64-bit double) series expansion used for small z (abs(z) < 0.05).
  511. // 17 Coefficients are computed to 21 decimal digits precision suitable for 64-bit double used by most platforms.
  512. // Taylor series coefficients used are computed by Wolfram to 50 decimal digits using instruction
  513. // N[InverseSeries[Series[z Exp[z],{z,0,34}]],50], as proposed by Tosio Fukushima and implemented by Veberic.
  514. template <class T, class Policy>
  515. T lambert_w0_small_z(const T z, const Policy&, boost::mpl::int_<1> const&)
  516. {
  517. #ifdef BOOST_MATH_INSTRUMENT_LAMBERT_W_SMALL_Z_SERIES
  518. std::streamsize prec = std::cout.precision(std::numeric_limits<T>::max_digits10); // Save.
  519. std::cout << "\ntag_type 1 double lambert_w0_small_z called with z = " << z << " using " << 17 << " terms of precision, "
  520. << std::numeric_limits<double>::max_digits10 << " decimal digits. " << std::endl;
  521. #endif // BOOST_MATH_INSTRUMENT_LAMBERT_W_SMALL_Z_SERIES
  522. T result =
  523. z * (1. - // j1 z^1
  524. z * (1. - // j2 z^2
  525. z * (1.5 - // 3/2 // j3 z^3
  526. z * (2.6666666666666666667 - // 8/3 // j4
  527. z * (5.2083333333333333333 - // -125/24 // j5
  528. z * (10.8 - // j6
  529. z * (23.343055555555555556 - // j7
  530. z * (52.012698412698412698 - // j8
  531. z * (118.62522321428571429 - // j9
  532. z * (275.57319223985890653 - // j10
  533. z * (649.78717234347442681 - // j11
  534. z * (1551.1605194805194805 - // j12
  535. z * (3741.4497029592385495 - // j13
  536. z * (9104.5002411580189358 - // j14
  537. z * (22324.308512706601434 - // j15
  538. z * (55103.621972903835338 - // j16
  539. z * 136808.86090394293563)))))))))))))))); // j17 z^17
  540. #ifdef BOOST_MATH_INSTRUMENT_LAMBERT_W_SMALL_Z_SERIES
  541. std::cout << "return w = " << result << std::endl;
  542. std::cout.precision(prec); // Restore.
  543. #endif // BOOST_MATH_INSTRUMENT_LAMBERT_W_SMALL_Z_SERIES
  544. return result;
  545. } // T lambert_w0_small_z(const T z, boost::mpl::int_<1> const&)
  546. //! Specialization of long double (80-bit double extended) series expansion used for small z (abs(z) < 0.05).
  547. // 21 Coefficients are computed to 21 decimal digits precision suitable for 80-bit long double used by some
  548. // platforms including GCC and Clang when generating for Intel X86 floating-point processors with 80-bit operations enabled (the default).
  549. // (This is NOT used by Microsoft Visual Studio where double and long always both use only 64-bit type.
  550. // Nor used for 128-bit float128.)
  551. template <class T, class Policy>
  552. T lambert_w0_small_z(const T z, const Policy&, boost::mpl::int_<2> const&)
  553. {
  554. #ifdef BOOST_MATH_INSTRUMENT_LAMBERT_W_SMALL_Z_SERIES
  555. std::streamsize precision = std::cout.precision(std::numeric_limits<T>::max_digits10); // Save.
  556. std::cout << "\ntag_type 2 long double (80-bit double extended) lambert_w0_small_z called with z = " << z << " using " << 21 << " terms of precision, "
  557. << std::numeric_limits<long double>::max_digits10 << " decimal digits. " << std::endl;
  558. #endif // BOOST_MATH_INSTRUMENT_LAMBERT_W_SMALL_Z_SERIES
  559. // T result =
  560. // z * (1.L - // j1 z^1
  561. // z * (1.L - // j2 z^2
  562. // z * (1.5L - // 3/2 // j3
  563. // z * (2.6666666666666666667L - // 8/3 // j4
  564. // z * (5.2083333333333333333L - // -125/24 // j5
  565. // z * (10.800000000000000000L - // j6
  566. // z * (23.343055555555555556L - // j7
  567. // z * (52.012698412698412698L - // j8
  568. // z * (118.62522321428571429L - // j9
  569. // z * (275.57319223985890653L - // j10
  570. // z * (649.78717234347442681L - // j11
  571. // z * (1551.1605194805194805L - // j12
  572. // z * (3741.4497029592385495L - // j13
  573. // z * (9104.5002411580189358L - // j14
  574. // z * (22324.308512706601434L - // j15
  575. // z * (55103.621972903835338L - // j16
  576. // z * (136808.86090394293563L - // j17 z^17 last term used by Fukushima double.
  577. // z * (341422.050665838363317L - // z^18
  578. // z * (855992.9659966075514633L - // z^19
  579. // z * (2.154990206091088289321e6L - // z^20
  580. // z * 5.4455529223144624316423e6L // z^21
  581. // ))))))))))))))))))));
  582. //
  583. T result =
  584. z * (1.L - // z j1
  585. z * (1.L - // z^2
  586. z * (1.500000000000000000000000000000000L - // z^3
  587. z * (2.666666666666666666666666666666666L - // z ^ 4
  588. z * (5.208333333333333333333333333333333L - // z ^ 5
  589. z * (10.80000000000000000000000000000000L - // z ^ 6
  590. z * (23.34305555555555555555555555555555L - // z ^ 7
  591. z * (52.01269841269841269841269841269841L - // z ^ 8
  592. z * (118.6252232142857142857142857142857L - // z ^ 9
  593. z * (275.5731922398589065255731922398589L - // z ^ 10
  594. z * (649.7871723434744268077601410934744L - // z ^ 11
  595. z * (1551.160519480519480519480519480519L - // z ^ 12
  596. z * (3741.449702959238549516327294105071L - //z ^ 13
  597. z * (9104.500241158018935796713574491352L - // z ^ 14
  598. z * (22324.308512706601434280005708577137L - // z ^ 15
  599. z * (55103.621972903835337697771560205422L - // z ^ 16
  600. z * (136808.86090394293563342215789305736L - // z ^ 17
  601. z * (341422.05066583836331735491399356945L - // z^18
  602. z * (855992.9659966075514633630250633224L - // z^19
  603. z * (2.154990206091088289321708745358647e6L // z^20 distance -5 without term 20
  604. ))))))))))))))))))));
  605. #ifdef BOOST_MATH_INSTRUMENT_LAMBERT_W_SMALL_Z_SERIES
  606. std::cout << "return w = " << result << std::endl;
  607. std::cout.precision(precision); // Restore.
  608. #endif // BOOST_MATH_INSTRUMENT_LAMBERT_W_SMALL_Z_SERIES
  609. return result;
  610. } // long double lambert_w0_small_z(const T z, boost::mpl::int_<1> const&)
  611. //! Specialization of 128-bit long double series expansion used for small z (abs(z) < 0.05).
  612. // 34 Taylor series coefficients used are computed by Wolfram to 50 decimal digits using instruction
  613. // N[InverseSeries[Series[z Exp[z],{z,0,34}]],50],
  614. // and are suffixed by L as they are assumed of type long double.
  615. // (This is NOT used for 128-bit quad boost::multiprecision::float128 type which required a suffix Q
  616. // nor multiprecision type cpp_bin_float_quad that can only be initialised at full precision of the type
  617. // constructed with a decimal digit string like "2.6666666666666666666666666666666666666666666666667".)
  618. template <class T, class Policy>
  619. T lambert_w0_small_z(const T z, const Policy&, boost::mpl::int_<3> const&)
  620. {
  621. #ifdef BOOST_MATH_INSTRUMENT_LAMBERT_W_SMALL_Z_SERIES
  622. std::streamsize precision = std::cout.precision(std::numeric_limits<T>::max_digits10); // Save.
  623. std::cout << "\ntag_type 3 long double (128-bit) lambert_w0_small_z called with z = " << z << " using " << 17 << " terms of precision, "
  624. << std::numeric_limits<double>::max_digits10 << " decimal digits. " << std::endl;
  625. #endif // BOOST_MATH_INSTRUMENT_LAMBERT_W_SMALL_Z_SERIES
  626. T result =
  627. z * (1.L - // j1
  628. z * (1.L - // j2
  629. z * (1.5L - // 3/2 // j3
  630. z * (2.6666666666666666666666666666666666L - // 8/3 // j4
  631. z * (5.2052083333333333333333333333333333L - // -125/24 // j5
  632. z * (10.800000000000000000000000000000000L - // j6
  633. z * (23.343055555555555555555555555555555L - // j7
  634. z * (52.0126984126984126984126984126984126L - // j8
  635. z * (118.625223214285714285714285714285714L - // j9
  636. z * (275.57319223985890652557319223985890L - // * z ^ 10 - // j10
  637. z * (649.78717234347442680776014109347442680776014109347L - // j11
  638. z * (1551.1605194805194805194805194805194805194805194805L - // j12
  639. z * (3741.4497029592385495163272941050718828496606274384L - // j13
  640. z * (9104.5002411580189357967135744913522691300469078247L - // j14
  641. z * (22324.308512706601434280005708577137148565719994291L - // j15
  642. z * (55103.621972903835337697771560205422639285073147507L - // j16
  643. z * 136808.86090394293563342215789305736395683485630576L // j17
  644. ))))))))))))))));
  645. #ifdef BOOST_MATH_INSTRUMENT_LAMBERT_W_SMALL_Z_SERIES
  646. std::cout << "return w = " << result << std::endl;
  647. std::cout.precision(precision); // Restore.
  648. #endif // BOOST_MATH_INSTRUMENT_LAMBERT_W_SMALL_Z_SERIES
  649. return result;
  650. } // T lambert_w0_small_z(const T z, boost::mpl::int_<3> const&)
  651. //! Specialization of 128-bit quad series expansion used for small z (abs(z) < 0.05).
  652. // 34 Taylor series coefficients used were computed by Wolfram to 50 decimal digits using instruction
  653. // N[InverseSeries[Series[z Exp[z],{z,0,34}]],50],
  654. // and are suffixed by Q as they are assumed of type quad.
  655. // This could be used for 128-bit quad (which requires a suffix Q for full precision).
  656. // But experiments with GCC 7.2.0 show that while this gives full 128-bit precision
  657. // when the -f-ext-numeric-literals option is in force and the libquadmath library available,
  658. // over the range -0.049 to +0.049,
  659. // it is slightly slower than getting a double approximation followed by a single Halley step.
  660. #ifdef BOOST_HAS_FLOAT128
  661. template <class T, class Policy>
  662. T lambert_w0_small_z(const T z, const Policy&, boost::mpl::int_<4> const&)
  663. {
  664. #ifdef BOOST_MATH_INSTRUMENT_LAMBERT_W_SMALL_Z_SERIES
  665. std::streamsize precision = std::cout.precision(std::numeric_limits<T>::max_digits10); // Save.
  666. std::cout << "\ntag_type 4 128-bit quad float128 lambert_w0_small_z called with z = " << z << " using " << 34 << " terms of precision, "
  667. << std::numeric_limits<float128>::max_digits10 << " max decimal digits." << std::endl;
  668. #endif // BOOST_MATH_INSTRUMENT_LAMBERT_W_SMALL_Z_SERIES
  669. T result =
  670. z * (1.Q - // z j1
  671. z * (1.Q - // z^2
  672. z * (1.500000000000000000000000000000000Q - // z^3
  673. z * (2.666666666666666666666666666666666Q - // z ^ 4
  674. z * (5.208333333333333333333333333333333Q - // z ^ 5
  675. z * (10.80000000000000000000000000000000Q - // z ^ 6
  676. z * (23.34305555555555555555555555555555Q - // z ^ 7
  677. z * (52.01269841269841269841269841269841Q - // z ^ 8
  678. z * (118.6252232142857142857142857142857Q - // z ^ 9
  679. z * (275.5731922398589065255731922398589Q - // z ^ 10
  680. z * (649.7871723434744268077601410934744Q - // z ^ 11
  681. z * (1551.160519480519480519480519480519Q - // z ^ 12
  682. z * (3741.449702959238549516327294105071Q - //z ^ 13
  683. z * (9104.500241158018935796713574491352Q - // z ^ 14
  684. z * (22324.308512706601434280005708577137Q - // z ^ 15
  685. z * (55103.621972903835337697771560205422Q - // z ^ 16
  686. z * (136808.86090394293563342215789305736Q - // z ^ 17
  687. z * (341422.05066583836331735491399356945Q - // z^18
  688. z * (855992.9659966075514633630250633224Q - // z^19
  689. z * (2.154990206091088289321708745358647e6Q - // 20
  690. z * (5.445552922314462431642316420035073e6Q - // 21
  691. z * (1.380733000216662949061923813184508e7Q - // 22
  692. z * (3.511704498513923292853869855945334e7Q - // 23
  693. z * (8.956800256102797693072819557780090e7Q - // 24
  694. z * (2.290416846187949813964782641734774e8Q - // 25
  695. z * (5.871035041171798492020292225245235e8Q - // 26
  696. z * (1.508256053857792919641317138812957e9Q - // 27
  697. z * (3.882630161293188940385873468413841e9Q - // 28
  698. z * (1.001394313665482968013913601565723e10Q - // 29
  699. z * (2.587356736265760638992878359024929e10Q - // 30
  700. z * (6.696209709358073856946120522333454e10Q - // 31
  701. z * (1.735711659599198077777078238043644e11Q - // 32
  702. z * (4.505680465642353886756098108484670e11Q - // 33
  703. z * (1.171223178256487391904047636564823e12Q //z^34
  704. ))))))))))))))))))))))))))))))))));
  705. #ifdef BOOST_MATH_INSTRUMENT_LAMBERT_W_SMALL_Z_SERIES
  706. std::cout << "return w = " << result << std::endl;
  707. std::cout.precision(precision); // Restore.
  708. #endif // BOOST_MATH_INSTRUMENT_LAMBERT_W_SMALL_Z_SERIES
  709. return result;
  710. } // T lambert_w0_small_z(const T z, boost::mpl::int_<4> const&) float128
  711. #else
  712. template <class T, class Policy>
  713. inline T lambert_w0_small_z(const T z, const Policy& pol, boost::mpl::int_<4> const&)
  714. {
  715. return lambert_w0_small_z(z, pol, boost::mpl::int_<5>());
  716. }
  717. #endif // BOOST_HAS_FLOAT128
  718. //! Series functor to compute series term using pow and factorial.
  719. //! \details Functor is called after evaluating polynomial with the coefficients as rationals below.
  720. template <class T>
  721. struct lambert_w0_small_z_series_term
  722. {
  723. typedef T result_type;
  724. //! \param _z Lambert W argument z.
  725. //! \param -term -pow<18>(z) / 6402373705728000uLL
  726. //! \param _k number of terms == initially 18
  727. // Note *after* evaluating N terms, its internal state has k = N and term = (-1)^N z^N.
  728. lambert_w0_small_z_series_term(T _z, T _term, int _k)
  729. : k(_k), z(_z), term(_term) { }
  730. T operator()()
  731. { // Called by sum_series until needs precision set by factor (policy::get_epsilon).
  732. using std::pow;
  733. ++k;
  734. term *= -z / k;
  735. //T t = pow(z, k) * pow(T(k), -1 + k) / factorial<T>(k); // (z^k * k(k-1)^k) / k!
  736. T result = term * pow(T(k), -1 + k); // term * k^(k-1)
  737. // std::cout << " k = " << k << ", term = " << term << ", result = " << result << std::endl;
  738. return result; //
  739. }
  740. private:
  741. int k;
  742. T z;
  743. T term;
  744. }; // template <class T> struct lambert_w0_small_z_series_term
  745. //! Generic variant for T a User-defined types like Boost.Multiprecision.
  746. template <class T, class Policy>
  747. inline T lambert_w0_small_z(T z, const Policy& pol, boost::mpl::int_<5> const&)
  748. {
  749. #ifdef BOOST_MATH_INSTRUMENT_LAMBERT_W_SMALL_Z_SERIES
  750. std::streamsize precision = std::cout.precision(std::numeric_limits<T>::max_digits10); // Save.
  751. std::cout << "Generic lambert_w0_small_z called with z = " << z << " using as many terms needed for precision." << std::endl;
  752. std::cout << "Argument z is of type " << typeid(T).name() << std::endl;
  753. #endif // BOOST_MATH_INSTRUMENT_LAMBERT_W_SMALL_Z_SERIES
  754. // First several terms of the series are tabulated and evaluated as a polynomial:
  755. // this will save us a bunch of expensive calls to pow.
  756. // Then our series functor is initialized "as if" it had already reached term 18,
  757. // enough evaluation of built-in 64-bit double and float (and 80-bit long double?) types.
  758. // Coefficients should be stored such that the coefficients for the x^i terms are in poly[i].
  759. static const T coeff[] =
  760. {
  761. 0, // z^0 Care: zeroth term needed by tools::evaluate_polynomial, but not in the Wolfram equation, so indexes are one different!
  762. 1, // z^1 term.
  763. -1, // z^2 term
  764. static_cast<T>(3uLL) / 2uLL, // z^3 term.
  765. -static_cast<T>(8uLL) / 3uLL, // z^4
  766. static_cast<T>(125uLL) / 24uLL, // z^5
  767. -static_cast<T>(54uLL) / 5uLL, // z^6
  768. static_cast<T>(16807uLL) / 720uLL, // z^7
  769. -static_cast<T>(16384uLL) / 315uLL, // z^8
  770. static_cast<T>(531441uLL) / 4480uLL, // z^9
  771. -static_cast<T>(156250uLL) / 567uLL, // z^10
  772. static_cast<T>(2357947691uLL) / 3628800uLL, // z^11
  773. -static_cast<T>(2985984uLL) / 1925uLL, // z^12
  774. static_cast<T>(1792160394037uLL) / 479001600uLL, // z^13
  775. -static_cast<T>(7909306972uLL) / 868725uLL, // z^14
  776. static_cast<T>(320361328125uLL) / 14350336uLL, // z^15
  777. -static_cast<T>(35184372088832uLL) / 638512875uLL, // z^16
  778. static_cast<T>(2862423051509815793uLL) / 20922789888000uLL, // z^17 term
  779. -static_cast<T>(5083731656658uLL) / 14889875uLL,
  780. // z^18 term. = 136808.86090394293563342215789305735851647769682393
  781. // z^18 is biggest that can be computed as rational using the largest possible uLL integers,
  782. // so higher terms cannot be potentially compiler-computed as uLL rationals.
  783. // Wolfram (5083731656658 z ^ 18) / 14889875 or
  784. // -341422.05066583836331735491399356945575432970390954 z^18
  785. // See note below calling the functor to compute another term,
  786. // sufficient for 80-bit long double precision.
  787. // Wolfram -341422.05066583836331735491399356945575432970390954 z^19 term.
  788. // (5480386857784802185939 z^19)/6402373705728000
  789. // But now this variant is not used to compute long double
  790. // as specializations are provided above.
  791. }; // static const T coeff[]
  792. /*
  793. Table of 19 computed coefficients:
  794. #0 0
  795. #1 1
  796. #2 -1
  797. #3 1.5
  798. #4 -2.6666666666666666666666666666666665382713370408509
  799. #5 5.2083333333333333333333333333333330765426740817019
  800. #6 -10.800000000000000000000000000000000616297582203915
  801. #7 23.343055555555555555555555555555555076212991619177
  802. #8 -52.012698412698412698412698412698412659282693193402
  803. #9 118.62522321428571428571428571428571146835390992496
  804. #10 -275.57319223985890652557319223985891400375196748314
  805. #11 649.7871723434744268077601410934743969785223845882
  806. #12 -1551.1605194805194805194805194805194947599566007429
  807. #13 3741.4497029592385495163272941050719510009019331763
  808. #14 -9104.5002411580189357967135744913524243896052869184
  809. #15 22324.308512706601434280005708577137322392070452582
  810. #16 -55103.621972903835337697771560205423203318720697224
  811. #17 136808.86090394293563342215789305735851647769682393
  812. 136808.86090394293563342215789305735851647769682393 == Exactly same as Wolfram computed value.
  813. #18 -341422.05066583836331735491399356947486381600607416
  814. 341422.05066583836331735491399356945575432970390954 z^19 Wolfram value differs at 36 decimal digit, as expected.
  815. */
  816. using boost::math::policies::get_epsilon; // for type T.
  817. using boost::math::tools::sum_series;
  818. using boost::math::tools::evaluate_polynomial;
  819. // http://www.boost.org/doc/libs/release/libs/math/doc/html/math_toolkit/roots/rational.html
  820. // std::streamsize prec = std::cout.precision(std::numeric_limits <T>::max_digits10);
  821. T result = evaluate_polynomial(coeff, z);
  822. // template <std::size_t N, class T, class V>
  823. // V evaluate_polynomial(const T(&poly)[N], const V& val);
  824. // Size of coeff found from N
  825. //std::cout << "evaluate_polynomial(coeff, z); == " << result << std::endl;
  826. //std::cout << "result = " << result << std::endl;
  827. // It's an artefact of the way I wrote the functor: *after* evaluating N
  828. // terms, its internal state has k = N and term = (-1)^N z^N. So after
  829. // evaluating 18 terms, we initialize the functor to the term we've just
  830. // evaluated, and then when it's called, it increments itself to the next term.
  831. // So 18!is 6402373705728000, which is where that comes from.
  832. // The 19th coefficient of the polynomial is actually, 19 ^ 18 / 19!=
  833. // 104127350297911241532841 / 121645100408832000 which after removing GCDs
  834. // reduces down to Wolfram rational 5480386857784802185939 / 6402373705728000.
  835. // Wolfram z^19 term +(5480386857784802185939 z^19) /6402373705728000
  836. // +855992.96599660755146336302506332246623424823099755 z^19
  837. //! Evaluate Functor.
  838. lambert_w0_small_z_series_term<T> s(z, -pow<18>(z) / 6402373705728000uLL, 18);
  839. // Temporary to list the coefficients.
  840. //std::cout << " Table of coefficients" << std::endl;
  841. //std::streamsize saved_precision = std::cout.precision(50);
  842. //for (size_t i = 0; i != 19; i++)
  843. //{
  844. // std::cout << "#" << i << " " << coeff[i] << std::endl;
  845. //}
  846. //std::cout.precision(saved_precision);
  847. boost::uintmax_t max_iter = policies::get_max_series_iterations<Policy>(); // Max iterations from policy.
  848. #ifdef BOOST_MATH_INSTRUMENT_LAMBERT_W_SMALL_Z_SERIES
  849. std::cout << "max iter from policy = " << max_iter << std::endl;
  850. // // max iter from policy = 1000000 is default.
  851. #endif // BOOST_MATH_INSTRUMENT_LAMBERT_W_SMALL_Z_SERIES
  852. result = sum_series(s, get_epsilon<T, Policy>(), max_iter, result);
  853. // result == evaluate_polynomial.
  854. //sum_series(Functor& func, int bits, boost::uintmax_t& max_terms, const U& init_value)
  855. // std::cout << "sum_series(s, get_epsilon<T, Policy>(), max_iter, result); = " << result << std::endl;
  856. //T epsilon = get_epsilon<T, Policy>();
  857. //std::cout << "epilson from policy = " << epsilon << std::endl;
  858. // epilson from policy = 1.93e-34 for T == quad
  859. // 5.35e-51 for t = cpp_bin_float_50
  860. // std::cout << " get eps = " << get_epsilon<T, Policy>() << std::endl; // quad eps = 1.93e-34, bin_float_50 eps = 5.35e-51
  861. policies::check_series_iterations<T>("boost::math::lambert_w0_small_z<%1%>(%1%)", max_iter, pol);
  862. #ifdef BOOST_MATH_INSTRUMENT_LAMBERT_W_SMALL_Z_SERIES_ITERATIONS
  863. std::cout << "z = " << z << " needed " << max_iter << " iterations." << std::endl;
  864. std::cout.precision(prec); // Restore.
  865. #endif // BOOST_MATH_INSTRUMENT_LAMBERT_W_SMALL_Z_SERIES_ITERATIONS
  866. return result;
  867. } // template <class T, class Policy> inline T lambert_w0_small_z_series(T z, const Policy& pol)
  868. // Approximate lambert_w0 (used for z values that are outside range of lookup table or rational functions)
  869. // Corless equation 4.19, page 349, and Chapeau-Blondeau equation 20, page 2162.
  870. template <typename T>
  871. inline
  872. T lambert_w0_approx(T z)
  873. {
  874. BOOST_MATH_STD_USING
  875. T lz = log(z);
  876. T llz = log(lz);
  877. T w = lz - llz + (llz / lz); // Corless equation 4.19, page 349, and Chapeau-Blondeau equation 20, page 2162.
  878. return w;
  879. // std::cout << "w max " << max_w << std::endl; // double 703.227
  880. }
  881. //////////////////////////////////////////////////////////////////////////////////////////
  882. //! \brief Lambert_w0 implementations for float, double and higher precisions.
  883. //! 3rd parameter used to select which version is used.
  884. //! /details Rational polynomials are provided for several range of argument z.
  885. //! For very small values of z, and for z very near the branch singularity at -e^-1 (~= -0.367879),
  886. //! two other series functions are used.
  887. //! float precision polynomials are used for 32-bit (usually float) precision (for speed)
  888. //! double precision polynomials are used for 64-bit (usually double) precision.
  889. //! For higher precisions, a 64-bit double approximation is computed first,
  890. //! and then refined using Halley interations.
  891. template <class T>
  892. inline T get_near_singularity_param(T z)
  893. {
  894. BOOST_MATH_STD_USING
  895. const T p2 = 2 * (boost::math::constants::e<T>() * z + 1);
  896. const T p = sqrt(p2);
  897. return p;
  898. }
  899. inline float get_near_singularity_param(float z)
  900. {
  901. return static_cast<float>(get_near_singularity_param((double)z));
  902. }
  903. inline double get_near_singularity_param(double z)
  904. {
  905. return static_cast<double>(get_near_singularity_param((long double)z));
  906. }
  907. // Forward declarations:
  908. //template <class T, class Policy> T lambert_w0_small_z(T z, const Policy& pol);
  909. //template <class T, class Policy>
  910. //T lambert_w0_imp(T w, const Policy& pol, const mpl::int_<0>&); // 32 bit usually float.
  911. //template <class T, class Policy>
  912. //T lambert_w0_imp(T w, const Policy& pol, const mpl::int_<1>&); // 64 bit usually double.
  913. //template <class T, class Policy>
  914. //T lambert_w0_imp(T w, const Policy& pol, const mpl::int_<2>&); // 80-bit long double.
  915. template <class T>
  916. T lambert_w_positive_rational_float(T z)
  917. {
  918. BOOST_MATH_STD_USING
  919. if (z < 2)
  920. {
  921. if (z < 0.5)
  922. { // 0.05 < z < 0.5
  923. // Maximum Deviation Found: 2.993e-08
  924. // Expected Error Term : 2.993e-08
  925. // Maximum Relative Change in Control Points : 7.555e-04 Y offset : -8.196592331e-01
  926. static const T Y = 8.196592331e-01f;
  927. static const T P[] = {
  928. 1.803388345e-01f,
  929. -4.820256838e-01f,
  930. -1.068349741e+00f,
  931. -3.506624319e-02f,
  932. };
  933. static const T Q[] = {
  934. 1.000000000e+00f,
  935. 2.871703469e+00f,
  936. 1.690949264e+00f,
  937. };
  938. return z * (Y + boost::math::tools::evaluate_polynomial(P, z) / boost::math::tools::evaluate_polynomial(Q, z));
  939. }
  940. else
  941. { // 0.5 < z < 2
  942. // Max error in interpolated form: 1.018e-08
  943. static const T Y = 5.503368378e-01f;
  944. static const T P[] = {
  945. 4.493332766e-01f,
  946. 2.543432707e-01f,
  947. -4.808788799e-01f,
  948. -1.244425316e-01f,
  949. };
  950. static const T Q[] = {
  951. 1.000000000e+00f,
  952. 2.780661241e+00f,
  953. 1.830840318e+00f,
  954. 2.407221031e-01f,
  955. };
  956. return z * (Y + boost::math::tools::evaluate_rational(P, Q, z));
  957. }
  958. }
  959. else if (z < 6)
  960. {
  961. // 2 < z < 6
  962. // Max error in interpolated form: 2.944e-08
  963. static const T Y = 1.162393570e+00f;
  964. static const T P[] = {
  965. -1.144183394e+00f,
  966. -4.712732855e-01f,
  967. 1.563162512e-01f,
  968. 1.434010911e-02f,
  969. };
  970. static const T Q[] = {
  971. 1.000000000e+00f,
  972. 1.192626340e+00f,
  973. 2.295580708e-01f,
  974. 5.477869455e-03f,
  975. };
  976. return Y + boost::math::tools::evaluate_rational(P, Q, z);
  977. }
  978. else if (z < 18)
  979. {
  980. // 6 < z < 18
  981. // Max error in interpolated form: 5.893e-08
  982. static const T Y = 1.809371948e+00f;
  983. static const T P[] = {
  984. -1.689291769e+00f,
  985. -3.337812742e-01f,
  986. 3.151434873e-02f,
  987. 1.134178734e-03f,
  988. };
  989. static const T Q[] = {
  990. 1.000000000e+00f,
  991. 5.716915685e-01f,
  992. 4.489521292e-02f,
  993. 4.076716763e-04f,
  994. };
  995. return Y + boost::math::tools::evaluate_rational(P, Q, z);
  996. }
  997. else if (z < 9897.12905874) // 2.8 < log(z) < 9.2
  998. {
  999. // Max error in interpolated form: 1.771e-08
  1000. static const T Y = -1.402973175e+00f;
  1001. static const T P[] = {
  1002. 1.966174312e+00f,
  1003. 2.350864728e-01f,
  1004. -5.098074353e-02f,
  1005. -1.054818339e-02f,
  1006. };
  1007. static const T Q[] = {
  1008. 1.000000000e+00f,
  1009. 4.388208264e-01f,
  1010. 8.316639634e-02f,
  1011. 3.397187918e-03f,
  1012. -1.321489743e-05f,
  1013. };
  1014. T log_w = log(z);
  1015. return log_w + Y + boost::math::tools::evaluate_polynomial(P, log_w) / boost::math::tools::evaluate_polynomial(Q, log_w);
  1016. }
  1017. else if (z < 7.896296e+13) // 9.2 < log(z) <= 32
  1018. {
  1019. // Max error in interpolated form: 5.821e-08
  1020. static const T Y = -2.735729218e+00f;
  1021. static const T P[] = {
  1022. 3.424903470e+00f,
  1023. 7.525631787e-02f,
  1024. -1.427309584e-02f,
  1025. -1.435974178e-05f,
  1026. };
  1027. static const T Q[] = {
  1028. 1.000000000e+00f,
  1029. 2.514005579e-01f,
  1030. 6.118994652e-03f,
  1031. -1.357889535e-05f,
  1032. 7.312865624e-08f,
  1033. };
  1034. T log_w = log(z);
  1035. return log_w + Y + boost::math::tools::evaluate_polynomial(P, log_w) / boost::math::tools::evaluate_polynomial(Q, log_w);
  1036. }
  1037. else // 32 < log(z) < 100
  1038. {
  1039. // Max error in interpolated form: 1.491e-08
  1040. static const T Y = -4.012863159e+00f;
  1041. static const T P[] = {
  1042. 4.431629226e+00f,
  1043. 2.756690487e-01f,
  1044. -2.992956930e-03f,
  1045. -4.912259384e-05f,
  1046. };
  1047. static const T Q[] = {
  1048. 1.000000000e+00f,
  1049. 2.015434591e-01f,
  1050. 4.949426142e-03f,
  1051. 1.609659944e-05f,
  1052. -5.111523436e-09f,
  1053. };
  1054. T log_w = log(z);
  1055. return log_w + Y + boost::math::tools::evaluate_polynomial(P, log_w) / boost::math::tools::evaluate_polynomial(Q, log_w);
  1056. }
  1057. }
  1058. template <class T, class Policy>
  1059. T lambert_w_negative_rational_float(T z, const Policy& pol)
  1060. {
  1061. BOOST_MATH_STD_USING
  1062. if (z > -0.27)
  1063. {
  1064. if (z < -0.051)
  1065. {
  1066. // -0.27 < z < -0.051
  1067. // Max error in interpolated form: 5.080e-08
  1068. static const T Y = 1.255809784e+00f;
  1069. static const T P[] = {
  1070. -2.558083412e-01f,
  1071. -2.306524098e+00f,
  1072. -5.630887033e+00f,
  1073. -3.803974556e+00f,
  1074. };
  1075. static const T Q[] = {
  1076. 1.000000000e+00f,
  1077. 5.107680783e+00f,
  1078. 7.914062868e+00f,
  1079. 3.501498501e+00f,
  1080. };
  1081. return z * (Y + boost::math::tools::evaluate_rational(P, Q, z));
  1082. }
  1083. else
  1084. {
  1085. // Very small z so use a series function.
  1086. return lambert_w0_small_z(z, pol);
  1087. }
  1088. }
  1089. else if (z > -0.3578794411714423215955237701)
  1090. { // Very close to branch singularity.
  1091. // Max error in interpolated form: 5.269e-08
  1092. static const T Y = 1.220928431e-01f;
  1093. static const T P[] = {
  1094. -1.221787446e-01f,
  1095. -6.816155875e+00f,
  1096. 7.144582035e+01f,
  1097. 1.128444390e+03f,
  1098. };
  1099. static const T Q[] = {
  1100. 1.000000000e+00f,
  1101. 6.480326790e+01f,
  1102. 1.869145243e+02f,
  1103. -1.361804274e+03f,
  1104. 1.117826726e+03f,
  1105. };
  1106. T d = z + 0.367879441171442321595523770161460867445811f;
  1107. return -d / (Y + boost::math::tools::evaluate_polynomial(P, d) / boost::math::tools::evaluate_polynomial(Q, d));
  1108. }
  1109. else
  1110. {
  1111. // z is very close (within 0.01) of the singularity at e^-1.
  1112. return lambert_w_singularity_series(get_near_singularity_param(z));
  1113. }
  1114. }
  1115. //! Lambert_w0 @b 'float' implementation, selected when T is 32-bit precision.
  1116. template <class T, class Policy>
  1117. inline T lambert_w0_imp(T z, const Policy& pol, const mpl::int_<1>&)
  1118. {
  1119. static const char* function = "boost::math::lambert_w0<%1%>"; // For error messages.
  1120. BOOST_MATH_STD_USING // Aid ADL of std functions.
  1121. if ((boost::math::isnan)(z))
  1122. {
  1123. return boost::math::policies::raise_domain_error<T>(function, "Expected a value > -e^-1 (-0.367879...) but got %1%.", z, pol);
  1124. }
  1125. if ((boost::math::isinf)(z))
  1126. {
  1127. return boost::math::policies::raise_overflow_error<T>(function, "Expected a finite value but got %1%.", z, pol);
  1128. }
  1129. if (z >= 0.05) // Fukushima switch point.
  1130. // if (z >= 0.045) // 34 terms makes 128-bit 'exact' below 0.045.
  1131. { // Normal ranges using several rational polynomials.
  1132. return lambert_w_positive_rational_float(z);
  1133. }
  1134. else if (z <= -0.3678794411714423215955237701614608674458111310f)
  1135. {
  1136. if (z < -0.3678794411714423215955237701614608674458111310f)
  1137. return boost::math::policies::raise_domain_error<T>(function, "Expected z >= -e^-1 (-0.367879...) but got %1%.", z, pol);
  1138. return -1;
  1139. }
  1140. else // z < 0.05
  1141. {
  1142. return lambert_w_negative_rational_float(z, pol);
  1143. }
  1144. } // T lambert_w0_imp(T z, const Policy& pol, const mpl::int_<1>&) for 32-bit usually float.
  1145. template <class T>
  1146. T lambert_w_positive_rational_double(T z)
  1147. {
  1148. BOOST_MATH_STD_USING
  1149. if (z < 2)
  1150. {
  1151. if (z < 0.5)
  1152. {
  1153. // Max error in interpolated form: 2.255e-17
  1154. static const T offset = 8.19659233093261719e-01;
  1155. static const T P[] = {
  1156. 1.80340766906685177e-01,
  1157. 3.28178241493119307e-01,
  1158. -2.19153620687139706e+00,
  1159. -7.24750929074563990e+00,
  1160. -7.28395876262524204e+00,
  1161. -2.57417169492512916e+00,
  1162. -2.31606948888704503e-01
  1163. };
  1164. static const T Q[] = {
  1165. 1.00000000000000000e+00,
  1166. 7.36482529307436604e+00,
  1167. 2.03686007856430677e+01,
  1168. 2.62864592096657307e+01,
  1169. 1.59742041380858333e+01,
  1170. 4.03760534788374589e+00,
  1171. 2.91327346750475362e-01
  1172. };
  1173. return z * (offset + boost::math::tools::evaluate_polynomial(P, z) / boost::math::tools::evaluate_polynomial(Q, z));
  1174. }
  1175. else
  1176. {
  1177. // Max error in interpolated form: 3.806e-18
  1178. static const T offset = 5.50335884094238281e-01;
  1179. static const T P[] = {
  1180. 4.49664083944098322e-01,
  1181. 1.90417666196776909e+00,
  1182. 1.99951368798255994e+00,
  1183. -6.91217310299270265e-01,
  1184. -1.88533935998617058e+00,
  1185. -7.96743968047750836e-01,
  1186. -1.02891726031055254e-01,
  1187. -3.09156013592636568e-03
  1188. };
  1189. static const T Q[] = {
  1190. 1.00000000000000000e+00,
  1191. 6.45854489419584014e+00,
  1192. 1.54739232422116048e+01,
  1193. 1.72606164253337843e+01,
  1194. 9.29427055609544096e+00,
  1195. 2.29040824649748117e+00,
  1196. 2.21610620995418981e-01,
  1197. 5.70597669908194213e-03
  1198. };
  1199. return z * (offset + boost::math::tools::evaluate_rational(P, Q, z));
  1200. }
  1201. }
  1202. else if (z < 6)
  1203. {
  1204. // 2 < z < 6
  1205. // Max error in interpolated form: 1.216e-17
  1206. static const T Y = 1.16239356994628906e+00;
  1207. static const T P[] = {
  1208. -1.16230494982099475e+00,
  1209. -3.38528144432561136e+00,
  1210. -2.55653717293161565e+00,
  1211. -3.06755172989214189e-01,
  1212. 1.73149743765268289e-01,
  1213. 3.76906042860014206e-02,
  1214. 1.84552217624706666e-03,
  1215. 1.69434126904822116e-05,
  1216. };
  1217. static const T Q[] = {
  1218. 1.00000000000000000e+00,
  1219. 3.77187616711220819e+00,
  1220. 4.58799960260143701e+00,
  1221. 2.24101228462292447e+00,
  1222. 4.54794195426212385e-01,
  1223. 3.60761772095963982e-02,
  1224. 9.25176499518388571e-04,
  1225. 4.43611344705509378e-06,
  1226. };
  1227. return Y + boost::math::tools::evaluate_rational(P, Q, z);
  1228. }
  1229. else if (z < 18)
  1230. {
  1231. // 6 < z < 18
  1232. // Max error in interpolated form: 1.985e-19
  1233. static const T offset = 1.80937194824218750e+00;
  1234. static const T P[] =
  1235. {
  1236. -1.80690935424793635e+00,
  1237. -3.66995929380314602e+00,
  1238. -1.93842957940149781e+00,
  1239. -2.94269984375794040e-01,
  1240. 1.81224710627677778e-03,
  1241. 2.48166798603547447e-03,
  1242. 1.15806592415397245e-04,
  1243. 1.43105573216815533e-06,
  1244. 3.47281483428369604e-09
  1245. };
  1246. static const T Q[] = {
  1247. 1.00000000000000000e+00,
  1248. 2.57319080723908597e+00,
  1249. 1.96724528442680658e+00,
  1250. 5.84501352882650722e-01,
  1251. 7.37152837939206240e-02,
  1252. 3.97368430940416778e-03,
  1253. 8.54941838187085088e-05,
  1254. 6.05713225608426678e-07,
  1255. 8.17517283816615732e-10
  1256. };
  1257. return offset + boost::math::tools::evaluate_rational(P, Q, z);
  1258. }
  1259. else if (z < 9897.12905874) // 2.8 < log(z) < 9.2
  1260. {
  1261. // Max error in interpolated form: 1.195e-18
  1262. static const T Y = -1.40297317504882812e+00;
  1263. static const T P[] = {
  1264. 1.97011826279311924e+00,
  1265. 1.05639945701546704e+00,
  1266. 3.33434529073196304e-01,
  1267. 3.34619153200386816e-02,
  1268. -5.36238353781326675e-03,
  1269. -2.43901294871308604e-03,
  1270. -2.13762095619085404e-04,
  1271. -4.85531936495542274e-06,
  1272. -2.02473518491905386e-08,
  1273. };
  1274. static const T Q[] = {
  1275. 1.00000000000000000e+00,
  1276. 8.60107275833921618e-01,
  1277. 4.10420467985504373e-01,
  1278. 1.18444884081994841e-01,
  1279. 2.16966505556021046e-02,
  1280. 2.24529766630769097e-03,
  1281. 9.82045090226437614e-05,
  1282. 1.36363515125489502e-06,
  1283. 3.44200749053237945e-09,
  1284. };
  1285. T log_w = log(z);
  1286. return log_w + Y + boost::math::tools::evaluate_rational(P, Q, log_w);
  1287. }
  1288. else if (z < 7.896296e+13) // 9.2 < log(z) <= 32
  1289. {
  1290. // Max error in interpolated form: 6.529e-18
  1291. static const T Y = -2.73572921752929688e+00;
  1292. static const T P[] = {
  1293. 3.30547638424076217e+00,
  1294. 1.64050071277550167e+00,
  1295. 4.57149576470736039e-01,
  1296. 4.03821227745424840e-02,
  1297. -4.99664976882514362e-04,
  1298. -1.28527893803052956e-04,
  1299. -2.95470325373338738e-06,
  1300. -1.76662025550202762e-08,
  1301. -1.98721972463709290e-11,
  1302. };
  1303. static const T Q[] = {
  1304. 1.00000000000000000e+00,
  1305. 6.91472559412458759e-01,
  1306. 2.48154578891676774e-01,
  1307. 4.60893578284335263e-02,
  1308. 3.60207838982301946e-03,
  1309. 1.13001153242430471e-04,
  1310. 1.33690948263488455e-06,
  1311. 4.97253225968548872e-09,
  1312. 3.39460723731970550e-12,
  1313. };
  1314. T log_w = log(z);
  1315. return log_w + Y + boost::math::tools::evaluate_rational(P, Q, log_w);
  1316. }
  1317. else if (z < 2.6881171e+43) // 32 < log(z) < 100
  1318. {
  1319. // Max error in interpolated form: 2.015e-18
  1320. static const T Y = -4.01286315917968750e+00;
  1321. static const T P[] = {
  1322. 5.07714858354309672e+00,
  1323. -3.32994414518701458e+00,
  1324. -8.61170416909864451e-01,
  1325. -4.01139705309486142e-02,
  1326. -1.85374201771834585e-04,
  1327. 1.08824145844270666e-05,
  1328. 1.17216905810452396e-07,
  1329. 2.97998248101385990e-10,
  1330. 1.42294856434176682e-13,
  1331. };
  1332. static const T Q[] = {
  1333. 1.00000000000000000e+00,
  1334. -4.85840770639861485e-01,
  1335. -3.18714850604827580e-01,
  1336. -3.20966129264610534e-02,
  1337. -1.06276178044267895e-03,
  1338. -1.33597828642644955e-05,
  1339. -6.27900905346219472e-08,
  1340. -9.35271498075378319e-11,
  1341. -2.60648331090076845e-14,
  1342. };
  1343. T log_w = log(z);
  1344. return log_w + Y + boost::math::tools::evaluate_rational(P, Q, log_w);
  1345. }
  1346. else // 100 < log(z) < 710
  1347. {
  1348. // Max error in interpolated form: 5.277e-18
  1349. static const T Y = -5.70115661621093750e+00;
  1350. static const T P[] = {
  1351. 6.42275660145116698e+00,
  1352. 1.33047964073367945e+00,
  1353. 6.72008923401652816e-02,
  1354. 1.16444069958125895e-03,
  1355. 7.06966760237470501e-06,
  1356. 5.48974896149039165e-09,
  1357. -7.00379652018853621e-11,
  1358. -1.89247635913659556e-13,
  1359. -1.55898770790170598e-16,
  1360. -4.06109208815303157e-20,
  1361. -2.21552699006496737e-24,
  1362. };
  1363. static const T Q[] = {
  1364. 1.00000000000000000e+00,
  1365. 3.34498588416632854e-01,
  1366. 2.51519862456384983e-02,
  1367. 6.81223810622416254e-04,
  1368. 7.94450897106903537e-06,
  1369. 4.30675039872881342e-08,
  1370. 1.10667669458467617e-10,
  1371. 1.31012240694192289e-13,
  1372. 6.53282047177727125e-17,
  1373. 1.11775518708172009e-20,
  1374. 3.78250395617836059e-25,
  1375. };
  1376. T log_w = log(z);
  1377. return log_w + Y + boost::math::tools::evaluate_rational(P, Q, log_w);
  1378. }
  1379. }
  1380. template <class T, class Policy>
  1381. T lambert_w_negative_rational_double(T z, const Policy& pol)
  1382. {
  1383. BOOST_MATH_STD_USING
  1384. if (z > -0.1)
  1385. {
  1386. if (z < -0.051)
  1387. {
  1388. // -0.1 < z < -0.051
  1389. // Maximum Deviation Found: 4.402e-22
  1390. // Expected Error Term : 4.240e-22
  1391. // Maximum Relative Change in Control Points : 4.115e-03
  1392. static const T Y = 1.08633995056152344e+00;
  1393. static const T P[] = {
  1394. -8.63399505615014331e-02,
  1395. -1.64303871814816464e+00,
  1396. -7.71247913918273738e+00,
  1397. -1.41014495545382454e+01,
  1398. -1.02269079949257616e+01,
  1399. -2.17236002836306691e+00,
  1400. };
  1401. static const T Q[] = {
  1402. 1.00000000000000000e+00,
  1403. 7.44775406945739243e+00,
  1404. 2.04392643087266541e+01,
  1405. 2.51001961077774193e+01,
  1406. 1.31256080849023319e+01,
  1407. 2.11640324843601588e+00,
  1408. };
  1409. return z * (Y + boost::math::tools::evaluate_rational(P, Q, z));
  1410. }
  1411. else
  1412. {
  1413. // Very small z > 0.051:
  1414. return lambert_w0_small_z(z, pol);
  1415. }
  1416. }
  1417. else if (z > -0.2)
  1418. {
  1419. // -0.2 < z < -0.1
  1420. // Maximum Deviation Found: 2.898e-20
  1421. // Expected Error Term : 2.873e-20
  1422. // Maximum Relative Change in Control Points : 3.779e-04
  1423. static const T Y = 1.20359611511230469e+00;
  1424. static const T P[] = {
  1425. -2.03596115108465635e-01,
  1426. -2.95029082937201859e+00,
  1427. -1.54287922188671648e+01,
  1428. -3.81185809571116965e+01,
  1429. -4.66384358235575985e+01,
  1430. -2.59282069989642468e+01,
  1431. -4.70140451266553279e+00,
  1432. };
  1433. static const T Q[] = {
  1434. 1.00000000000000000e+00,
  1435. 9.57921436074599929e+00,
  1436. 3.60988119290234377e+01,
  1437. 6.73977699505546007e+01,
  1438. 6.41104992068148823e+01,
  1439. 2.82060127225153607e+01,
  1440. 4.10677610657724330e+00,
  1441. };
  1442. return z * (Y + boost::math::tools::evaluate_rational(P, Q, z));
  1443. }
  1444. else if (z > -0.3178794411714423215955237)
  1445. {
  1446. // Max error in interpolated form: 6.996e-18
  1447. static const T Y = 3.49680423736572266e-01;
  1448. static const T P[] = {
  1449. -3.49729841718749014e-01,
  1450. -6.28207407760709028e+01,
  1451. -2.57226178029669171e+03,
  1452. -2.50271008623093747e+04,
  1453. 1.11949239154711388e+05,
  1454. 1.85684566607844318e+06,
  1455. 4.80802490427638643e+06,
  1456. 2.76624752134636406e+06,
  1457. };
  1458. static const T Q[] = {
  1459. 1.00000000000000000e+00,
  1460. 1.82717661215113000e+02,
  1461. 8.00121119810280100e+03,
  1462. 1.06073266717010129e+05,
  1463. 3.22848993926057721e+05,
  1464. -8.05684814514171256e+05,
  1465. -2.59223192927265737e+06,
  1466. -5.61719645211570871e+05,
  1467. 6.27765369292636844e+04,
  1468. };
  1469. T d = z + 0.367879441171442321595523770161460867445811;
  1470. return -d / (Y + boost::math::tools::evaluate_polynomial(P, d) / boost::math::tools::evaluate_polynomial(Q, d));
  1471. }
  1472. else if (z > -0.3578794411714423215955237701)
  1473. {
  1474. // Max error in interpolated form: 1.404e-17
  1475. static const T Y = 5.00126481056213379e-02;
  1476. static const T P[] = {
  1477. -5.00173570682372162e-02,
  1478. -4.44242461870072044e+01,
  1479. -9.51185533619946042e+03,
  1480. -5.88605699015429386e+05,
  1481. -1.90760843597427751e+06,
  1482. 5.79797663818311404e+08,
  1483. 1.11383352508459134e+10,
  1484. 5.67791253678716467e+10,
  1485. 6.32694500716584572e+10,
  1486. };
  1487. static const T Q[] = {
  1488. 1.00000000000000000e+00,
  1489. 9.08910517489981551e+02,
  1490. 2.10170163753340133e+05,
  1491. 1.67858612416470327e+07,
  1492. 4.90435561733227953e+08,
  1493. 4.54978142622939917e+09,
  1494. 2.87716585708739168e+09,
  1495. -4.59414247951143131e+10,
  1496. -1.72845216404874299e+10,
  1497. };
  1498. T d = z + 0.36787944117144232159552377016146086744581113103176804;
  1499. return -d / (Y + boost::math::tools::evaluate_polynomial(P, d) / boost::math::tools::evaluate_polynomial(Q, d));
  1500. }
  1501. else
  1502. { // z is very close (within 0.01) of the singularity at -e^-1,
  1503. // so use a series expansion from R. M. Corless et al.
  1504. const T p2 = 2 * (boost::math::constants::e<T>() * z + 1);
  1505. const T p = sqrt(p2);
  1506. return lambert_w_detail::lambert_w_singularity_series(p);
  1507. }
  1508. }
  1509. //! Lambert_w0 @b 'double' implementation, selected when T is 64-bit precision.
  1510. template <class T, class Policy>
  1511. inline T lambert_w0_imp(T z, const Policy& pol, const mpl::int_<2>&)
  1512. {
  1513. static const char* function = "boost::math::lambert_w0<%1%>";
  1514. BOOST_MATH_STD_USING // Aid ADL of std functions.
  1515. // Detect unusual case of 32-bit double with a wider/64-bit long double
  1516. BOOST_STATIC_ASSERT_MSG(std::numeric_limits<double>::digits >= 53,
  1517. "Our double precision coefficients will be truncated, "
  1518. "please file a bug report with details of your platform's floating point types "
  1519. "- or possibly edit the coefficients to have "
  1520. "an appropriate size-suffix for 64-bit floats on your platform - L?");
  1521. if ((boost::math::isnan)(z))
  1522. {
  1523. return boost::math::policies::raise_domain_error<T>(function, "Expected a value > -e^-1 (-0.367879...) but got %1%.", z, pol);
  1524. }
  1525. if ((boost::math::isinf)(z))
  1526. {
  1527. return boost::math::policies::raise_overflow_error<T>(function, "Expected a finite value but got %1%.", z, pol);
  1528. }
  1529. if (z >= 0.05)
  1530. {
  1531. return lambert_w_positive_rational_double(z);
  1532. }
  1533. else if (z <= -0.36787944117144232159552377016146086744581113103176804) // Precision is max_digits10(cpp_bin_float_50).
  1534. {
  1535. if (z < -0.36787944117144232159552377016146086744581113103176804)
  1536. {
  1537. return boost::math::policies::raise_domain_error<T>(function, "Expected z >= -e^-1 (-0.367879...) but got %1%.", z, pol);
  1538. }
  1539. return -1;
  1540. }
  1541. else
  1542. {
  1543. return lambert_w_negative_rational_double(z, pol);
  1544. }
  1545. } // T lambert_w0_imp(T z, const Policy& pol, const mpl::int_<2>&) 64-bit precision, usually double.
  1546. //! lambert_W0 implementation for extended precision types including
  1547. //! long double (80-bit and 128-bit), ???
  1548. //! quad float128, Boost.Multiprecision types like cpp_bin_float_quad, cpp_bin_float_50...
  1549. template <class T, class Policy>
  1550. inline T lambert_w0_imp(T z, const Policy& pol, const mpl::int_<0>&)
  1551. {
  1552. static const char* function = "boost::math::lambert_w0<%1%>";
  1553. BOOST_MATH_STD_USING // Aid ADL of std functions.
  1554. // Filter out special cases first:
  1555. if ((boost::math::isnan)(z))
  1556. {
  1557. return boost::math::policies::raise_domain_error<T>(function, "Expected z >= -e^-1 (-0.367879...) but got %1%.", z, pol);
  1558. }
  1559. if (fabs(z) <= 0.05f)
  1560. {
  1561. // Very small z:
  1562. return lambert_w0_small_z(z, pol);
  1563. }
  1564. if (z > (std::numeric_limits<double>::max)())
  1565. {
  1566. if ((boost::math::isinf)(z))
  1567. {
  1568. return policies::raise_overflow_error<T>(function, 0, pol);
  1569. // Or might return infinity if available else max_value,
  1570. // but other Boost.Math special functions raise overflow.
  1571. }
  1572. // z is larger than the largest double, so cannot use the polynomial to get an approximation,
  1573. // so use the asymptotic approximation and Halley iterate:
  1574. T w = lambert_w0_approx(z); // Make an inline function as also used elsewhere.
  1575. //T lz = log(z);
  1576. //T llz = log(lz);
  1577. //T w = lz - llz + (llz / lz); // Corless equation 4.19, page 349, and Chapeau-Blondeau equation 20, page 2162.
  1578. return lambert_w_halley_iterate(w, z);
  1579. }
  1580. if (z < -0.3578794411714423215955237701)
  1581. { // Very close to branch point so rational polynomials are not usable.
  1582. if (z <= -boost::math::constants::exp_minus_one<T>())
  1583. {
  1584. if (z == -boost::math::constants::exp_minus_one<T>())
  1585. { // Exactly at the branch point singularity.
  1586. return -1;
  1587. }
  1588. return boost::math::policies::raise_domain_error<T>(function, "Expected z >= -e^-1 (-0.367879...) but got %1%.", z, pol);
  1589. }
  1590. // z is very close (within 0.01) of the branch singularity at -e^-1
  1591. // so use a series approximation proposed by Corless et al.
  1592. const T p2 = 2 * (boost::math::constants::e<T>() * z + 1);
  1593. const T p = sqrt(p2);
  1594. T w = lambert_w_detail::lambert_w_singularity_series(p);
  1595. return lambert_w_halley_iterate(w, z);
  1596. }
  1597. // Phew! If we get here we are in the normal range of the function,
  1598. // so get a double precision approximation first, then iterate to full precision of T.
  1599. // We define a tag_type that is:
  1600. // mpl::true_ if there are so many digits precision wanted that iteration is necessary.
  1601. // mpl::false_ if a single Halley step is sufficient.
  1602. typedef typename policies::precision<T, Policy>::type precision_type;
  1603. typedef mpl::bool_<
  1604. (precision_type::value == 0) || (precision_type::value > 113) ?
  1605. true // Unknown at compile-time, variable/arbitrary, or more than float128 or cpp_bin_quad 128-bit precision.
  1606. : false // float, double, float128, cpp_bin_quad 128-bit, so single Halley step.
  1607. > tag_type;
  1608. // For speed, we also cast z to type double when that is possible
  1609. // if (boost::is_constructible<double, T>() == true).
  1610. T w = lambert_w0_imp(maybe_reduce_to_double(z, boost::is_constructible<double, T>()), pol, mpl::int_<2>());
  1611. return lambert_w_maybe_halley_iterate(w, z, tag_type());
  1612. } // T lambert_w0_imp(T z, const Policy& pol, const mpl::int_<0>&) all extended precision types.
  1613. // Lambert w-1 implementation
  1614. // ==============================================================================================
  1615. //! Lambert W for W-1 branch, -max(z) < z <= -1/e.
  1616. // TODO is -max(z) allowed?
  1617. template<typename T, class Policy>
  1618. T lambert_wm1_imp(const T z, const Policy& pol)
  1619. {
  1620. // Catch providing an integer value as parameter x to lambert_w, for example, lambert_w(1).
  1621. // Need to ensure it is a floating-point type (of the desired type, float 1.F, double 1., or long double 1.L),
  1622. // or static_casted integer, for example: static_cast<float>(1) or static_cast<cpp_dec_float_50>(1).
  1623. // Want to allow fixed_point types too, so do not just test for floating-point.
  1624. // Integral types should be promoted to double by user Lambert w functions.
  1625. // If integral type provided to user function lambert_w0 or lambert_wm1,
  1626. // then should already have been promoted to double.
  1627. BOOST_STATIC_ASSERT_MSG(!boost::is_integral<T>::value,
  1628. "Must be floating-point or fixed type (not integer type), for example: lambert_wm1(1.), not lambert_wm1(1)!");
  1629. BOOST_MATH_STD_USING // Aid argument dependent lookup (ADL) of abs.
  1630. const char* function = "boost::math::lambert_wm1<RealType>(<RealType>)"; // Used for error messages.
  1631. // Check for edge and corner cases first:
  1632. if ((boost::math::isnan)(z))
  1633. {
  1634. return policies::raise_domain_error(function,
  1635. "Argument z is NaN!",
  1636. z, pol);
  1637. } // isnan
  1638. if ((boost::math::isinf)(z))
  1639. {
  1640. return policies::raise_domain_error(function,
  1641. "Argument z is infinite!",
  1642. z, pol);
  1643. } // isinf
  1644. if (z == static_cast<T>(0))
  1645. { // z is exactly zero so return -std::numeric_limits<T>::infinity();
  1646. if (std::numeric_limits<T>::has_infinity)
  1647. {
  1648. return -std::numeric_limits<T>::infinity();
  1649. }
  1650. else
  1651. {
  1652. return -tools::max_value<T>();
  1653. }
  1654. }
  1655. if (std::numeric_limits<T>::has_denorm)
  1656. { // All real types except arbitrary precision.
  1657. if (!(boost::math::isnormal)(z))
  1658. { // Almost zero - might also just return infinity like z == 0 or max_value?
  1659. return policies::raise_overflow_error(function,
  1660. "Argument z = %1% is denormalized! (must be z > (std::numeric_limits<RealType>::min)() or z == 0)",
  1661. z, pol);
  1662. }
  1663. }
  1664. if (z > static_cast<T>(0))
  1665. { //
  1666. return policies::raise_domain_error(function,
  1667. "Argument z = %1% is out of range (z <= 0) for Lambert W-1 branch! (Try Lambert W0 branch?)",
  1668. z, pol);
  1669. }
  1670. if (z > -boost::math::tools::min_value<T>())
  1671. { // z is denormalized, so cannot be computed.
  1672. // -std::numeric_limits<T>::min() is smallest for type T,
  1673. // for example, for double: lambert_wm1(-2.2250738585072014e-308) = -714.96865723796634
  1674. return policies::raise_overflow_error(function,
  1675. "Argument z = %1% is too small (z < -std::numeric_limits<T>::min so denormalized) for Lambert W-1 branch!",
  1676. z, pol);
  1677. }
  1678. if (z == -boost::math::constants::exp_minus_one<T>()) // == singularity/branch point z = -exp(-1) = -3.6787944.
  1679. { // At singularity, so return exactly -1.
  1680. return -static_cast<T>(1);
  1681. }
  1682. // z is too negative for the W-1 (or W0) branch.
  1683. if (z < -boost::math::constants::exp_minus_one<T>()) // > singularity/branch point z = -exp(-1) = -3.6787944.
  1684. {
  1685. return policies::raise_domain_error(function,
  1686. "Argument z = %1% is out of range (z < -exp(-1) = -3.6787944... <= 0) for Lambert W-1 (or W0) branch!",
  1687. z, pol);
  1688. }
  1689. if (z < static_cast<T>(-0.35))
  1690. { // Close to singularity/branch point z = -0.3678794411714423215955237701614608727 but on W-1 branch.
  1691. const T p2 = 2 * (boost::math::constants::e<T>() * z + 1);
  1692. if (p2 == 0)
  1693. { // At the singularity at branch point.
  1694. return -1;
  1695. }
  1696. if (p2 > 0)
  1697. {
  1698. T w_series = lambert_w_singularity_series(T(-sqrt(p2)));
  1699. if (boost::math::tools::digits<T>() > 53)
  1700. { // Multiprecision, so try a Halley refinement.
  1701. w_series = lambert_w_detail::lambert_w_halley_iterate(w_series, z);
  1702. #ifdef BOOST_MATH_INSTRUMENT_LAMBERT_WM1_NOT_BUILTIN
  1703. std::streamsize saved_precision = std::cout.precision(std::numeric_limits<T>::max_digits10);
  1704. std::cout << "Lambert W-1 Halley updated to " << w_series << std::endl;
  1705. std::cout.precision(saved_precision);
  1706. #endif // BOOST_MATH_INSTRUMENT_LAMBERT_WM1_NOT_BUILTIN
  1707. }
  1708. return w_series;
  1709. }
  1710. // Should not get here.
  1711. return policies::raise_domain_error(function,
  1712. "Argument z = %1% is out of range for Lambert W-1 branch. (Should not get here - please report!)",
  1713. z, pol);
  1714. } // if (z < -0.35)
  1715. using lambert_w_lookup::wm1es;
  1716. using lambert_w_lookup::wm1zs;
  1717. using lambert_w_lookup::noof_wm1zs; // size == 64
  1718. // std::cout <<" Wm1zs[63] (== G[64]) = " << " " << wm1zs[63] << std::endl; // Wm1zs[63] (== G[64]) = -1.0264389699511283e-26
  1719. // Check that z argument value is not smaller than lookup_table G[64]
  1720. // std::cout << "(z > wm1zs[63]) = " << std::boolalpha << (z > wm1zs[63]) << std::endl;
  1721. if (z >= wm1zs[63]) // wm1zs[63] = -1.0264389699511282259046957018510946438e-26L W = 64.00000000000000000
  1722. { // z >= -1.0264389699511303e-26 (but z != 0 and z >= std::numeric_limits<T>::min() and so NOT denormalized).
  1723. // Some info on Lambert W-1 values for extreme values of z.
  1724. // std::streamsize saved_precision = std::cout.precision(std::numeric_limits<T>::max_digits10);
  1725. // std::cout << "-std::numeric_limits<float>::min() = " << -(std::numeric_limits<float>::min)() << std::endl;
  1726. // std::cout << "-std::numeric_limits<double>::min() = " << -(std::numeric_limits<double>::min)() << std::endl;
  1727. // -std::numeric_limits<float>::min() = -1.1754943508222875e-38
  1728. // -std::numeric_limits<double>::min() = -2.2250738585072014e-308
  1729. // N[productlog(-1, -1.1754943508222875 * 10^-38 ), 50] = -91.856775324595479509567756730093823993834155027858
  1730. // N[productlog(-1, -2.2250738585072014e-308 * 10^-308 ), 50] = -1424.8544521230553853558132180518404363617968042942
  1731. // N[productlog(-1, -1.4325445274604020119111357113179868158* 10^-27), 37] = -65.99999999999999999999999999999999955
  1732. // R.M.Corless, G.H.Gonnet, D.E.G.Hare, D.J.Jeffrey, and D.E.Knuth,
  1733. // On the Lambert W function, Adv.Comput.Math., vol. 5, pp. 329, 1996.
  1734. // Francois Chapeau-Blondeau and Abdelilah Monir
  1735. // Numerical Evaluation of the Lambert W Function
  1736. // IEEE Transactions On Signal Processing, VOL. 50, NO. 9, Sep 2002
  1737. // https://pdfs.semanticscholar.org/7a5a/76a9369586dd0dd34dda156d8f2779d1fd59.pdf
  1738. // Estimate Lambert W using ln(-z) ...
  1739. // This is roughly the power of ten * ln(10) ~= 2.3. n ~= 10^n
  1740. // and improve by adding a second term -ln(ln(-z))
  1741. T guess; // bisect lowest possible Gk[=64] (for lookup_t type)
  1742. T lz = log(-z);
  1743. T llz = log(-lz);
  1744. guess = lz - llz + (llz / lz); // Chapeau-Blondeau equation 20, page 2162.
  1745. #ifdef BOOST_MATH_INSTRUMENT_LAMBERT_WM1_TINY
  1746. std::streamsize saved_precision = std::cout.precision(std::numeric_limits<T>::max_digits10);
  1747. std::cout << "z = " << z << ", guess = " << guess << ", ln(-z) = " << lz << ", ln(-ln(-z) = " << llz << ", llz/lz = " << (llz / lz) << std::endl;
  1748. // z = -1.0000000000000001e-30, guess = -73.312782616731482, ln(-z) = -69.077552789821368, ln(-ln(-z) = 4.2352298269101114, llz/lz = -0.061311231447304194
  1749. // z = -9.9999999999999999e-91, guess = -212.56650048504233, ln(-z) = -207.23265836946410, ln(-ln(-z) = 5.3338421155782205, llz/lz = -0.025738424423764311
  1750. // >z = -2.2250738585072014e-308, guess = -714.95942238244606, ln(-z) = -708.39641853226408, ln(-ln(-z) = 6.5630038501819854, llz/lz = -0.0092645920821846622
  1751. int d10 = policies::digits_base10<T, Policy>(); // policy template parameter digits10
  1752. int d2 = policies::digits<T, Policy>(); // digits base 2 from policy.
  1753. std::cout << "digits10 = " << d10 << ", digits2 = " << d2 // For example: digits10 = 1, digits2 = 5
  1754. << std::endl;
  1755. std::cout.precision(saved_precision);
  1756. #endif // BOOST_MATH_INSTRUMENT_LAMBERT_WM1_TINY
  1757. if (policies::digits<T, Policy>() < 12)
  1758. { // For the worst case near w = 64, the error in the 'guess' is ~0.008, ratio ~ 0.0001 or 1 in 10,000 digits 10 ~= 4, or digits2 ~= 12.
  1759. return guess;
  1760. }
  1761. T result = lambert_w_detail::lambert_w_halley_iterate(guess, z);
  1762. return result;
  1763. // Was Fukushima
  1764. // G[k=64] == g[63] == -1.02643897e-26
  1765. //return policies::raise_domain_error(function,
  1766. // "Argument z = %1% is too small (< -1.02643897e-26) ! (Should not occur, please report.",
  1767. // z, pol);
  1768. } // Z too small so use approximation and Halley.
  1769. // Else Use a lookup table to find the nearest integer part of Lambert W-1 as starting point for Bisection.
  1770. if (boost::math::tools::digits<T>() > 53)
  1771. { // T is more precise than 64-bit double (or long double, or ?),
  1772. // so compute an approximate value using only one Schroeder refinement,
  1773. // (avoiding any double-precision Halley refinement from policy double_digits2<50> 53 - 3 = 50
  1774. // because are next going to use Halley refinement at full/high precision using this as an approximation).
  1775. using boost::math::policies::precision;
  1776. using boost::math::policies::digits10;
  1777. using boost::math::policies::digits2;
  1778. using boost::math::policies::policy;
  1779. // Compute a 50-bit precision approximate W0 in a double (no Halley refinement).
  1780. T double_approx(static_cast<T>(lambert_wm1_imp(must_reduce_to_double(z, boost::is_constructible<double, T>()), policy<digits2<50> >())));
  1781. #ifdef BOOST_MATH_INSTRUMENT_LAMBERT_WM1_NOT_BUILTIN
  1782. std::streamsize saved_precision = std::cout.precision(std::numeric_limits<T>::max_digits10);
  1783. std::cout << "Lambert_wm1 Argument Type " << typeid(T).name() << " approximation double = " << double_approx << std::endl;
  1784. std::cout.precision(saved_precision);
  1785. #endif // BOOST_MATH_INSTRUMENT_LAMBERT_WM1
  1786. // Perform additional Halley refinement(s) to ensure that
  1787. // get a near as possible to correct result (usually +/- one epsilon).
  1788. T result = lambert_w_halley_iterate(double_approx, z);
  1789. #ifdef BOOST_MATH_INSTRUMENT_LAMBERT_WM1
  1790. std::streamsize saved_precision = std::cout.precision(std::numeric_limits<T>::max_digits10);
  1791. std::cout << "Result " << typeid(T).name() << " precision Halley refinement = " << result << std::endl;
  1792. std::cout.precision(saved_precision);
  1793. #endif // BOOST_MATH_INSTRUMENT_LAMBERT_WM1
  1794. return result;
  1795. } // digits > 53 - higher precision than double.
  1796. else // T is double or less precision.
  1797. { // Use a lookup table to find the nearest integer part of Lambert W as starting point for Bisection.
  1798. using namespace boost::math::lambert_w_detail::lambert_w_lookup;
  1799. // Bracketing sequence n = (2, 4, 8, 16, 32, 64) for W-1 branch. (0 is -infinity)
  1800. // Since z is probably quite small, start with lowest n (=2).
  1801. int n = 2;
  1802. if (wm1zs[n - 1] > z)
  1803. {
  1804. goto bisect;
  1805. }
  1806. for (int j = 1; j <= 5; ++j)
  1807. {
  1808. n *= 2;
  1809. if (wm1zs[n - 1] > z)
  1810. {
  1811. goto overshot;
  1812. }
  1813. }
  1814. // else z < g[63] == -1.0264389699511303e-26, so Lambert W-1 integer part > 64.
  1815. // This should not now occur (should be caught by test and code above) so should be a logic_error?
  1816. return policies::raise_domain_error(function,
  1817. "Argument z = %1% is too small (< -1.026439e-26) (logic error - please report!)",
  1818. z, pol);
  1819. overshot:
  1820. {
  1821. int nh = n / 2;
  1822. for (int j = 1; j <= 5; ++j)
  1823. {
  1824. nh /= 2; // halve step size.
  1825. if (nh <= 0)
  1826. {
  1827. break; // goto bisect;
  1828. }
  1829. if (wm1zs[n - nh - 1] > z)
  1830. {
  1831. n -= nh;
  1832. }
  1833. }
  1834. }
  1835. bisect:
  1836. --n;
  1837. // g[n] now holds lambert W of floor integer n and g[n+1] the ceil part;
  1838. // these are used as initial values for bisection.
  1839. #ifdef BOOST_MATH_INSTRUMENT_LAMBERT_WM1_LOOKUP
  1840. std::streamsize saved_precision = std::cout.precision(std::numeric_limits<T>::max_digits10);
  1841. std::cout << "Result lookup W-1(" << z << ") bisection between wm1zs[" << n - 1 << "] = " << wm1zs[n - 1] << " and wm1zs[" << n << "] = " << wm1zs[n]
  1842. << ", bisect mean = " << (wm1zs[n - 1] + wm1zs[n]) / 2 << std::endl;
  1843. std::cout.precision(saved_precision);
  1844. #endif // BOOST_MATH_INSTRUMENT_LAMBERT_WM1_LOOKUP
  1845. // Compute bisections is the number of bisections computed from n,
  1846. // such that a single application of the fifth-order Schroeder update formula
  1847. // after the bisections is enough to evaluate Lambert W-1 with (near?) 53-bit accuracy.
  1848. // Fukushima established these by trial and error?
  1849. int bisections = 11; // Assume maximum number of bisections will be needed (most common case).
  1850. if (n >= 8)
  1851. {
  1852. bisections = 8;
  1853. }
  1854. else if (n >= 3)
  1855. {
  1856. bisections = 9;
  1857. }
  1858. else if (n >= 2)
  1859. {
  1860. bisections = 10;
  1861. }
  1862. // Bracketing, Fukushima section 2.3, page 82:
  1863. // (Avoiding using exponential function for speed).
  1864. // Only use @c lookup_t precision, default double, for bisection (again for speed),
  1865. // and use later Halley refinement for higher precisions.
  1866. using lambert_w_lookup::halves;
  1867. using lambert_w_lookup::sqrtwm1s;
  1868. typedef typename mpl::if_c<boost::is_constructible<lookup_t, T>::value, lookup_t, T>::type calc_type;
  1869. calc_type w = -static_cast<calc_type>(n); // Equation 25,
  1870. calc_type y = static_cast<calc_type>(z * wm1es[n - 1]); // Equation 26,
  1871. // Perform the bisections fractional bisections for necessary precision.
  1872. for (int j = 0; j < bisections; ++j)
  1873. { // Equation 27.
  1874. calc_type wj = w - halves[j]; // Subtract 1/2, 1/4, 1/8 ...
  1875. calc_type yj = y * sqrtwm1s[j]; // Multiply by sqrt(1/e), ...
  1876. if (wj < yj)
  1877. {
  1878. w = wj;
  1879. y = yj;
  1880. }
  1881. } // for j
  1882. return static_cast<T>(schroeder_update(w, y)); // Schroeder 5th order method refinement.
  1883. // else // Perform additional Halley refinement(s) to ensure that
  1884. // // get a near as possible to correct result (usually +/- epsilon).
  1885. // {
  1886. // // result = lambert_w_halley_iterate(result, z);
  1887. // result = lambert_w_halley_step(result, z); // Just one Halley step should be enough.
  1888. //#ifdef BOOST_MATH_INSTRUMENT_LAMBERT_WM1_HALLEY
  1889. // std::streamsize saved_precision = std::cout.precision(std::numeric_limits<T>::max_digits10);
  1890. // std::cout << "Halley refinement estimate = " << result << std::endl;
  1891. // std::cout.precision(saved_precision);
  1892. //#endif // BOOST_MATH_INSTRUMENT_LAMBERT_W1_HALLEY
  1893. // return result; // Halley
  1894. // } // Schroeder or Schroeder and Halley.
  1895. }
  1896. } // template<typename T = double> T lambert_wm1_imp(const T z)
  1897. } // namespace lambert_w_detail
  1898. ///////////////////////////// User Lambert w functions. //////////////////////////////
  1899. //! Lambert W0 using User-defined policy.
  1900. template <class T, class Policy>
  1901. inline
  1902. typename boost::math::tools::promote_args<T>::type
  1903. lambert_w0(T z, const Policy& pol)
  1904. {
  1905. // Promote integer or expression template arguments to double,
  1906. // without doing any other internal promotion like float to double.
  1907. typedef typename tools::promote_args<T>::type result_type;
  1908. // Work out what precision has been selected,
  1909. // based on the Policy and the number type.
  1910. typedef typename policies::precision<result_type, Policy>::type precision_type;
  1911. // and then select the correct implementation based on that precision (not the type T):
  1912. typedef mpl::int_<
  1913. (precision_type::value == 0) || (precision_type::value > 53) ?
  1914. 0 // either variable precision (0), or greater than 64-bit precision.
  1915. : (precision_type::value <= 24) ? 1 // 32-bit (probably float) precision.
  1916. : 2 // 64-bit (probably double) precision.
  1917. > tag_type;
  1918. return lambert_w_detail::lambert_w0_imp(result_type(z), pol, tag_type()); //
  1919. } // lambert_w0(T z, const Policy& pol)
  1920. //! Lambert W0 using default policy.
  1921. template <class T>
  1922. inline
  1923. typename tools::promote_args<T>::type
  1924. lambert_w0(T z)
  1925. {
  1926. // Promote integer or expression template arguments to double,
  1927. // without doing any other internal promotion like float to double.
  1928. typedef typename tools::promote_args<T>::type result_type;
  1929. // Work out what precision has been selected, based on the Policy and the number type.
  1930. // For the default policy version, we want the *default policy* precision for T.
  1931. typedef typename policies::precision<result_type, policies::policy<> >::type precision_type;
  1932. // and then select the correct implementation based on that (not the type T):
  1933. typedef mpl::int_<
  1934. (precision_type::value == 0) || (precision_type::value > 53) ?
  1935. 0 // either variable precision (0), or greater than 64-bit precision.
  1936. : (precision_type::value <= 24) ? 1 // 32-bit (probably float) precision.
  1937. : 2 // 64-bit (probably double) precision.
  1938. > tag_type;
  1939. return lambert_w_detail::lambert_w0_imp(result_type(z), policies::policy<>(), tag_type());
  1940. } // lambert_w0(T z) using default policy.
  1941. //! W-1 branch (-max(z) < z <= -1/e).
  1942. //! Lambert W-1 using User-defined policy.
  1943. template <class T, class Policy>
  1944. inline
  1945. typename tools::promote_args<T>::type
  1946. lambert_wm1(T z, const Policy& pol)
  1947. {
  1948. // Promote integer or expression template arguments to double,
  1949. // without doing any other internal promotion like float to double.
  1950. typedef typename tools::promote_args<T>::type result_type;
  1951. return lambert_w_detail::lambert_wm1_imp(result_type(z), pol); //
  1952. }
  1953. //! Lambert W-1 using default policy.
  1954. template <class T>
  1955. inline
  1956. typename tools::promote_args<T>::type
  1957. lambert_wm1(T z)
  1958. {
  1959. typedef typename tools::promote_args<T>::type result_type;
  1960. return lambert_w_detail::lambert_wm1_imp(result_type(z), policies::policy<>());
  1961. } // lambert_wm1(T z)
  1962. // First derivative of Lambert W0 and W-1.
  1963. template <class T, class Policy>
  1964. inline typename tools::promote_args<T>::type
  1965. lambert_w0_prime(T z, const Policy& pol)
  1966. {
  1967. typedef typename tools::promote_args<T>::type result_type;
  1968. using std::numeric_limits;
  1969. if (z == 0)
  1970. {
  1971. return static_cast<result_type>(1);
  1972. }
  1973. // This is the sensible choice if we regard the Lambert-W function as complex analytic.
  1974. // Of course on the real line, it's just undefined.
  1975. if (z == - boost::math::constants::exp_minus_one<result_type>())
  1976. {
  1977. return numeric_limits<result_type>::has_infinity ? numeric_limits<result_type>::infinity() : boost::math::tools::max_value<result_type>();
  1978. }
  1979. // if z < -1/e, we'll let lambert_w0 do the error handling:
  1980. result_type w = lambert_w0(result_type(z), pol);
  1981. // If w ~ -1, then presumably this can get inaccurate.
  1982. // Is there an accurate way to evaluate 1 + W(-1/e + eps)?
  1983. // Yes: This is discussed in the Princeton Companion to Applied Mathematics,
  1984. // 'The Lambert-W function', Section 1.3: Series and Generating Functions.
  1985. // 1 + W(-1/e + x) ~ sqrt(2ex).
  1986. // Nick is not convinced this formula is more accurate than the naive one.
  1987. // However, for z != -1/e, we never get rounded to w = -1 in any precision I've tested (up to cpp_bin_float_100).
  1988. return w / (z * (1 + w));
  1989. } // lambert_w0_prime(T z)
  1990. template <class T>
  1991. inline typename tools::promote_args<T>::type
  1992. lambert_w0_prime(T z)
  1993. {
  1994. return lambert_w0_prime(z, policies::policy<>());
  1995. }
  1996. template <class T, class Policy>
  1997. inline typename tools::promote_args<T>::type
  1998. lambert_wm1_prime(T z, const Policy& pol)
  1999. {
  2000. using std::numeric_limits;
  2001. typedef typename tools::promote_args<T>::type result_type;
  2002. //if (z == 0)
  2003. //{
  2004. // return static_cast<result_type>(1);
  2005. //}
  2006. //if (z == - boost::math::constants::exp_minus_one<result_type>())
  2007. if (z == 0 || z == - boost::math::constants::exp_minus_one<result_type>())
  2008. {
  2009. return numeric_limits<result_type>::has_infinity ? -numeric_limits<result_type>::infinity() : -boost::math::tools::max_value<result_type>();
  2010. }
  2011. result_type w = lambert_wm1(z, pol);
  2012. return w/(z*(1+w));
  2013. } // lambert_wm1_prime(T z)
  2014. template <class T>
  2015. inline typename tools::promote_args<T>::type
  2016. lambert_wm1_prime(T z)
  2017. {
  2018. return lambert_wm1_prime(z, policies::policy<>());
  2019. }
  2020. }} //boost::math namespaces
  2021. #endif // #ifdef BOOST_MATH_SF_LAMBERT_W_HPP