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