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