Chris@16
|
1 // Copyright John Maddock 2006.
|
Chris@16
|
2 // Copyright Paul A. Bristow 2007
|
Chris@16
|
3 // Use, modification and distribution are subject to the
|
Chris@16
|
4 // Boost Software License, Version 1.0. (See accompanying file
|
Chris@16
|
5 // LICENSE_1_0.txt or copy at http://www.boost.org/LICENSE_1_0.txt)
|
Chris@16
|
6
|
Chris@16
|
7 #ifndef BOOST_MATH_SPECIAL_FUNCTIONS_IBETA_INVERSE_HPP
|
Chris@16
|
8 #define BOOST_MATH_SPECIAL_FUNCTIONS_IBETA_INVERSE_HPP
|
Chris@16
|
9
|
Chris@16
|
10 #ifdef _MSC_VER
|
Chris@16
|
11 #pragma once
|
Chris@16
|
12 #endif
|
Chris@16
|
13
|
Chris@16
|
14 #include <boost/math/special_functions/beta.hpp>
|
Chris@16
|
15 #include <boost/math/special_functions/erf.hpp>
|
Chris@16
|
16 #include <boost/math/tools/roots.hpp>
|
Chris@16
|
17 #include <boost/math/special_functions/detail/t_distribution_inv.hpp>
|
Chris@16
|
18
|
Chris@16
|
19 namespace boost{ namespace math{ namespace detail{
|
Chris@16
|
20
|
Chris@16
|
21 //
|
Chris@16
|
22 // Helper object used by root finding
|
Chris@16
|
23 // code to convert eta to x.
|
Chris@16
|
24 //
|
Chris@16
|
25 template <class T>
|
Chris@16
|
26 struct temme_root_finder
|
Chris@16
|
27 {
|
Chris@16
|
28 temme_root_finder(const T t_, const T a_) : t(t_), a(a_) {}
|
Chris@16
|
29
|
Chris@16
|
30 boost::math::tuple<T, T> operator()(T x)
|
Chris@16
|
31 {
|
Chris@16
|
32 BOOST_MATH_STD_USING // ADL of std names
|
Chris@16
|
33
|
Chris@16
|
34 T y = 1 - x;
|
Chris@16
|
35 if(y == 0)
|
Chris@16
|
36 {
|
Chris@16
|
37 T big = tools::max_value<T>() / 4;
|
Chris@16
|
38 return boost::math::make_tuple(static_cast<T>(-big), static_cast<T>(-big));
|
Chris@16
|
39 }
|
Chris@16
|
40 if(x == 0)
|
Chris@16
|
41 {
|
Chris@16
|
42 T big = tools::max_value<T>() / 4;
|
Chris@16
|
43 return boost::math::make_tuple(static_cast<T>(-big), big);
|
Chris@16
|
44 }
|
Chris@16
|
45 T f = log(x) + a * log(y) + t;
|
Chris@16
|
46 T f1 = (1 / x) - (a / (y));
|
Chris@16
|
47 return boost::math::make_tuple(f, f1);
|
Chris@16
|
48 }
|
Chris@16
|
49 private:
|
Chris@16
|
50 T t, a;
|
Chris@16
|
51 };
|
Chris@16
|
52 //
|
Chris@16
|
53 // See:
|
Chris@16
|
54 // "Asymptotic Inversion of the Incomplete Beta Function"
|
Chris@16
|
55 // N.M. Temme
|
Chris@16
|
56 // Journal of Computation and Applied Mathematics 41 (1992) 145-157.
|
Chris@16
|
57 // Section 2.
|
Chris@16
|
58 //
|
Chris@16
|
59 template <class T, class Policy>
|
Chris@16
|
60 T temme_method_1_ibeta_inverse(T a, T b, T z, const Policy& pol)
|
Chris@16
|
61 {
|
Chris@16
|
62 BOOST_MATH_STD_USING // ADL of std names
|
Chris@16
|
63
|
Chris@16
|
64 const T r2 = sqrt(T(2));
|
Chris@16
|
65 //
|
Chris@16
|
66 // get the first approximation for eta from the inverse
|
Chris@16
|
67 // error function (Eq: 2.9 and 2.10).
|
Chris@16
|
68 //
|
Chris@16
|
69 T eta0 = boost::math::erfc_inv(2 * z, pol);
|
Chris@16
|
70 eta0 /= -sqrt(a / 2);
|
Chris@16
|
71
|
Chris@16
|
72 T terms[4] = { eta0 };
|
Chris@16
|
73 T workspace[7];
|
Chris@16
|
74 //
|
Chris@16
|
75 // calculate powers:
|
Chris@16
|
76 //
|
Chris@16
|
77 T B = b - a;
|
Chris@16
|
78 T B_2 = B * B;
|
Chris@16
|
79 T B_3 = B_2 * B;
|
Chris@16
|
80 //
|
Chris@16
|
81 // Calculate correction terms:
|
Chris@16
|
82 //
|
Chris@16
|
83
|
Chris@16
|
84 // See eq following 2.15:
|
Chris@16
|
85 workspace[0] = -B * r2 / 2;
|
Chris@16
|
86 workspace[1] = (1 - 2 * B) / 8;
|
Chris@16
|
87 workspace[2] = -(B * r2 / 48);
|
Chris@16
|
88 workspace[3] = T(-1) / 192;
|
Chris@16
|
89 workspace[4] = -B * r2 / 3840;
|
Chris@16
|
90 terms[1] = tools::evaluate_polynomial(workspace, eta0, 5);
|
Chris@16
|
91 // Eq Following 2.17:
|
Chris@16
|
92 workspace[0] = B * r2 * (3 * B - 2) / 12;
|
Chris@16
|
93 workspace[1] = (20 * B_2 - 12 * B + 1) / 128;
|
Chris@16
|
94 workspace[2] = B * r2 * (20 * B - 1) / 960;
|
Chris@16
|
95 workspace[3] = (16 * B_2 + 30 * B - 15) / 4608;
|
Chris@16
|
96 workspace[4] = B * r2 * (21 * B + 32) / 53760;
|
Chris@16
|
97 workspace[5] = (-32 * B_2 + 63) / 368640;
|
Chris@16
|
98 workspace[6] = -B * r2 * (120 * B + 17) / 25804480;
|
Chris@16
|
99 terms[2] = tools::evaluate_polynomial(workspace, eta0, 7);
|
Chris@16
|
100 // Eq Following 2.17:
|
Chris@16
|
101 workspace[0] = B * r2 * (-75 * B_2 + 80 * B - 16) / 480;
|
Chris@16
|
102 workspace[1] = (-1080 * B_3 + 868 * B_2 - 90 * B - 45) / 9216;
|
Chris@16
|
103 workspace[2] = B * r2 * (-1190 * B_2 + 84 * B + 373) / 53760;
|
Chris@16
|
104 workspace[3] = (-2240 * B_3 - 2508 * B_2 + 2100 * B - 165) / 368640;
|
Chris@16
|
105 terms[3] = tools::evaluate_polynomial(workspace, eta0, 4);
|
Chris@16
|
106 //
|
Chris@16
|
107 // Bring them together to get a final estimate for eta:
|
Chris@16
|
108 //
|
Chris@16
|
109 T eta = tools::evaluate_polynomial(terms, T(1/a), 4);
|
Chris@16
|
110 //
|
Chris@16
|
111 // now we need to convert eta to x, by solving the appropriate
|
Chris@16
|
112 // quadratic equation:
|
Chris@16
|
113 //
|
Chris@16
|
114 T eta_2 = eta * eta;
|
Chris@16
|
115 T c = -exp(-eta_2 / 2);
|
Chris@16
|
116 T x;
|
Chris@16
|
117 if(eta_2 == 0)
|
Chris@16
|
118 x = 0.5;
|
Chris@16
|
119 else
|
Chris@16
|
120 x = (1 + eta * sqrt((1 + c) / eta_2)) / 2;
|
Chris@16
|
121
|
Chris@16
|
122 BOOST_ASSERT(x >= 0);
|
Chris@16
|
123 BOOST_ASSERT(x <= 1);
|
Chris@16
|
124 BOOST_ASSERT(eta * (x - 0.5) >= 0);
|
Chris@16
|
125 #ifdef BOOST_INSTRUMENT
|
Chris@16
|
126 std::cout << "Estimating x with Temme method 1: " << x << std::endl;
|
Chris@16
|
127 #endif
|
Chris@16
|
128 return x;
|
Chris@16
|
129 }
|
Chris@16
|
130 //
|
Chris@16
|
131 // See:
|
Chris@16
|
132 // "Asymptotic Inversion of the Incomplete Beta Function"
|
Chris@16
|
133 // N.M. Temme
|
Chris@16
|
134 // Journal of Computation and Applied Mathematics 41 (1992) 145-157.
|
Chris@16
|
135 // Section 3.
|
Chris@16
|
136 //
|
Chris@16
|
137 template <class T, class Policy>
|
Chris@16
|
138 T temme_method_2_ibeta_inverse(T /*a*/, T /*b*/, T z, T r, T theta, const Policy& pol)
|
Chris@16
|
139 {
|
Chris@16
|
140 BOOST_MATH_STD_USING // ADL of std names
|
Chris@16
|
141
|
Chris@16
|
142 //
|
Chris@16
|
143 // Get first estimate for eta, see Eq 3.9 and 3.10,
|
Chris@16
|
144 // but note there is a typo in Eq 3.10:
|
Chris@16
|
145 //
|
Chris@16
|
146 T eta0 = boost::math::erfc_inv(2 * z, pol);
|
Chris@16
|
147 eta0 /= -sqrt(r / 2);
|
Chris@16
|
148
|
Chris@16
|
149 T s = sin(theta);
|
Chris@16
|
150 T c = cos(theta);
|
Chris@16
|
151 //
|
Chris@16
|
152 // Now we need to purturb eta0 to get eta, which we do by
|
Chris@16
|
153 // evaluating the polynomial in 1/r at the bottom of page 151,
|
Chris@16
|
154 // to do this we first need the error terms e1, e2 e3
|
Chris@16
|
155 // which we'll fill into the array "terms". Since these
|
Chris@16
|
156 // terms are themselves polynomials, we'll need another
|
Chris@16
|
157 // array "workspace" to calculate those...
|
Chris@16
|
158 //
|
Chris@16
|
159 T terms[4] = { eta0 };
|
Chris@16
|
160 T workspace[6];
|
Chris@16
|
161 //
|
Chris@16
|
162 // some powers of sin(theta)cos(theta) that we'll need later:
|
Chris@16
|
163 //
|
Chris@16
|
164 T sc = s * c;
|
Chris@16
|
165 T sc_2 = sc * sc;
|
Chris@16
|
166 T sc_3 = sc_2 * sc;
|
Chris@16
|
167 T sc_4 = sc_2 * sc_2;
|
Chris@16
|
168 T sc_5 = sc_2 * sc_3;
|
Chris@16
|
169 T sc_6 = sc_3 * sc_3;
|
Chris@16
|
170 T sc_7 = sc_4 * sc_3;
|
Chris@16
|
171 //
|
Chris@16
|
172 // Calculate e1 and put it in terms[1], see the middle of page 151:
|
Chris@16
|
173 //
|
Chris@16
|
174 workspace[0] = (2 * s * s - 1) / (3 * s * c);
|
Chris@16
|
175 static const BOOST_MATH_INT_TABLE_TYPE(T, int) co1[] = { -1, -5, 5 };
|
Chris@16
|
176 workspace[1] = -tools::evaluate_even_polynomial(co1, s, 3) / (36 * sc_2);
|
Chris@16
|
177 static const BOOST_MATH_INT_TABLE_TYPE(T, int) co2[] = { 1, 21, -69, 46 };
|
Chris@16
|
178 workspace[2] = tools::evaluate_even_polynomial(co2, s, 4) / (1620 * sc_3);
|
Chris@16
|
179 static const BOOST_MATH_INT_TABLE_TYPE(T, int) co3[] = { 7, -2, 33, -62, 31 };
|
Chris@16
|
180 workspace[3] = -tools::evaluate_even_polynomial(co3, s, 5) / (6480 * sc_4);
|
Chris@16
|
181 static const BOOST_MATH_INT_TABLE_TYPE(T, int) co4[] = { 25, -52, -17, 88, -115, 46 };
|
Chris@16
|
182 workspace[4] = tools::evaluate_even_polynomial(co4, s, 6) / (90720 * sc_5);
|
Chris@16
|
183 terms[1] = tools::evaluate_polynomial(workspace, eta0, 5);
|
Chris@16
|
184 //
|
Chris@16
|
185 // Now evaluate e2 and put it in terms[2]:
|
Chris@16
|
186 //
|
Chris@16
|
187 static const BOOST_MATH_INT_TABLE_TYPE(T, int) co5[] = { 7, 12, -78, 52 };
|
Chris@16
|
188 workspace[0] = -tools::evaluate_even_polynomial(co5, s, 4) / (405 * sc_3);
|
Chris@16
|
189 static const BOOST_MATH_INT_TABLE_TYPE(T, int) co6[] = { -7, 2, 183, -370, 185 };
|
Chris@16
|
190 workspace[1] = tools::evaluate_even_polynomial(co6, s, 5) / (2592 * sc_4);
|
Chris@16
|
191 static const BOOST_MATH_INT_TABLE_TYPE(T, int) co7[] = { -533, 776, -1835, 10240, -13525, 5410 };
|
Chris@16
|
192 workspace[2] = -tools::evaluate_even_polynomial(co7, s, 6) / (204120 * sc_5);
|
Chris@16
|
193 static const BOOST_MATH_INT_TABLE_TYPE(T, int) co8[] = { -1579, 3747, -3372, -15821, 45588, -45213, 15071 };
|
Chris@16
|
194 workspace[3] = -tools::evaluate_even_polynomial(co8, s, 7) / (2099520 * sc_6);
|
Chris@16
|
195 terms[2] = tools::evaluate_polynomial(workspace, eta0, 4);
|
Chris@16
|
196 //
|
Chris@16
|
197 // And e3, and put it in terms[3]:
|
Chris@16
|
198 //
|
Chris@16
|
199 static const BOOST_MATH_INT_TABLE_TYPE(T, int) co9[] = {449, -1259, -769, 6686, -9260, 3704 };
|
Chris@16
|
200 workspace[0] = tools::evaluate_even_polynomial(co9, s, 6) / (102060 * sc_5);
|
Chris@16
|
201 static const BOOST_MATH_INT_TABLE_TYPE(T, int) co10[] = { 63149, -151557, 140052, -727469, 2239932, -2251437, 750479 };
|
Chris@16
|
202 workspace[1] = -tools::evaluate_even_polynomial(co10, s, 7) / (20995200 * sc_6);
|
Chris@16
|
203 static const BOOST_MATH_INT_TABLE_TYPE(T, int) co11[] = { 29233, -78755, 105222, 146879, -1602610, 3195183, -2554139, 729754 };
|
Chris@16
|
204 workspace[2] = tools::evaluate_even_polynomial(co11, s, 8) / (36741600 * sc_7);
|
Chris@16
|
205 terms[3] = tools::evaluate_polynomial(workspace, eta0, 3);
|
Chris@16
|
206 //
|
Chris@16
|
207 // Bring the correction terms together to evaluate eta,
|
Chris@16
|
208 // this is the last equation on page 151:
|
Chris@16
|
209 //
|
Chris@16
|
210 T eta = tools::evaluate_polynomial(terms, T(1/r), 4);
|
Chris@16
|
211 //
|
Chris@16
|
212 // Now that we have eta we need to back solve for x,
|
Chris@16
|
213 // we seek the value of x that gives eta in Eq 3.2.
|
Chris@16
|
214 // The two methods used are described in section 5.
|
Chris@16
|
215 //
|
Chris@16
|
216 // Begin by defining a few variables we'll need later:
|
Chris@16
|
217 //
|
Chris@16
|
218 T x;
|
Chris@16
|
219 T s_2 = s * s;
|
Chris@16
|
220 T c_2 = c * c;
|
Chris@16
|
221 T alpha = c / s;
|
Chris@16
|
222 alpha *= alpha;
|
Chris@16
|
223 T lu = (-(eta * eta) / (2 * s_2) + log(s_2) + c_2 * log(c_2) / s_2);
|
Chris@16
|
224 //
|
Chris@16
|
225 // Temme doesn't specify what value to switch on here,
|
Chris@16
|
226 // but this seems to work pretty well:
|
Chris@16
|
227 //
|
Chris@16
|
228 if(fabs(eta) < 0.7)
|
Chris@16
|
229 {
|
Chris@16
|
230 //
|
Chris@16
|
231 // Small eta use the expansion Temme gives in the second equation
|
Chris@16
|
232 // of section 5, it's a polynomial in eta:
|
Chris@16
|
233 //
|
Chris@16
|
234 workspace[0] = s * s;
|
Chris@16
|
235 workspace[1] = s * c;
|
Chris@16
|
236 workspace[2] = (1 - 2 * workspace[0]) / 3;
|
Chris@16
|
237 static const BOOST_MATH_INT_TABLE_TYPE(T, int) co12[] = { 1, -13, 13 };
|
Chris@16
|
238 workspace[3] = tools::evaluate_polynomial(co12, workspace[0], 3) / (36 * s * c);
|
Chris@16
|
239 static const BOOST_MATH_INT_TABLE_TYPE(T, int) co13[] = { 1, 21, -69, 46 };
|
Chris@16
|
240 workspace[4] = tools::evaluate_polynomial(co13, workspace[0], 4) / (270 * workspace[0] * c * c);
|
Chris@16
|
241 x = tools::evaluate_polynomial(workspace, eta, 5);
|
Chris@16
|
242 #ifdef BOOST_INSTRUMENT
|
Chris@16
|
243 std::cout << "Estimating x with Temme method 2 (small eta): " << x << std::endl;
|
Chris@16
|
244 #endif
|
Chris@16
|
245 }
|
Chris@16
|
246 else
|
Chris@16
|
247 {
|
Chris@16
|
248 //
|
Chris@16
|
249 // If eta is large we need to solve Eq 3.2 more directly,
|
Chris@16
|
250 // begin by getting an initial approximation for x from
|
Chris@16
|
251 // the last equation on page 155, this is a polynomial in u:
|
Chris@16
|
252 //
|
Chris@16
|
253 T u = exp(lu);
|
Chris@16
|
254 workspace[0] = u;
|
Chris@16
|
255 workspace[1] = alpha;
|
Chris@16
|
256 workspace[2] = 0;
|
Chris@16
|
257 workspace[3] = 3 * alpha * (3 * alpha + 1) / 6;
|
Chris@16
|
258 workspace[4] = 4 * alpha * (4 * alpha + 1) * (4 * alpha + 2) / 24;
|
Chris@16
|
259 workspace[5] = 5 * alpha * (5 * alpha + 1) * (5 * alpha + 2) * (5 * alpha + 3) / 120;
|
Chris@16
|
260 x = tools::evaluate_polynomial(workspace, u, 6);
|
Chris@16
|
261 //
|
Chris@16
|
262 // At this point we may or may not have the right answer, Eq-3.2 has
|
Chris@16
|
263 // two solutions for x for any given eta, however the mapping in 3.2
|
Chris@16
|
264 // is 1:1 with the sign of eta and x-sin^2(theta) being the same.
|
Chris@16
|
265 // So we can check if we have the right root of 3.2, and if not
|
Chris@16
|
266 // switch x for 1-x. This transformation is motivated by the fact
|
Chris@16
|
267 // that the distribution is *almost* symetric so 1-x will be in the right
|
Chris@16
|
268 // ball park for the solution:
|
Chris@16
|
269 //
|
Chris@16
|
270 if((x - s_2) * eta < 0)
|
Chris@16
|
271 x = 1 - x;
|
Chris@16
|
272 #ifdef BOOST_INSTRUMENT
|
Chris@16
|
273 std::cout << "Estimating x with Temme method 2 (large eta): " << x << std::endl;
|
Chris@16
|
274 #endif
|
Chris@16
|
275 }
|
Chris@16
|
276 //
|
Chris@16
|
277 // The final step is a few Newton-Raphson iterations to
|
Chris@16
|
278 // clean up our approximation for x, this is pretty cheap
|
Chris@16
|
279 // in general, and very cheap compared to an incomplete beta
|
Chris@16
|
280 // evaluation. The limits set on x come from the observation
|
Chris@16
|
281 // that the sign of eta and x-sin^2(theta) are the same.
|
Chris@16
|
282 //
|
Chris@16
|
283 T lower, upper;
|
Chris@16
|
284 if(eta < 0)
|
Chris@16
|
285 {
|
Chris@16
|
286 lower = 0;
|
Chris@16
|
287 upper = s_2;
|
Chris@16
|
288 }
|
Chris@16
|
289 else
|
Chris@16
|
290 {
|
Chris@16
|
291 lower = s_2;
|
Chris@16
|
292 upper = 1;
|
Chris@16
|
293 }
|
Chris@16
|
294 //
|
Chris@16
|
295 // If our initial approximation is out of bounds then bisect:
|
Chris@16
|
296 //
|
Chris@16
|
297 if((x < lower) || (x > upper))
|
Chris@16
|
298 x = (lower+upper) / 2;
|
Chris@16
|
299 //
|
Chris@16
|
300 // And iterate:
|
Chris@16
|
301 //
|
Chris@16
|
302 x = tools::newton_raphson_iterate(
|
Chris@16
|
303 temme_root_finder<T>(-lu, alpha), x, lower, upper, policies::digits<T, Policy>() / 2);
|
Chris@16
|
304
|
Chris@16
|
305 return x;
|
Chris@16
|
306 }
|
Chris@16
|
307 //
|
Chris@16
|
308 // See:
|
Chris@16
|
309 // "Asymptotic Inversion of the Incomplete Beta Function"
|
Chris@16
|
310 // N.M. Temme
|
Chris@16
|
311 // Journal of Computation and Applied Mathematics 41 (1992) 145-157.
|
Chris@16
|
312 // Section 4.
|
Chris@16
|
313 //
|
Chris@16
|
314 template <class T, class Policy>
|
Chris@16
|
315 T temme_method_3_ibeta_inverse(T a, T b, T p, T q, const Policy& pol)
|
Chris@16
|
316 {
|
Chris@16
|
317 BOOST_MATH_STD_USING // ADL of std names
|
Chris@16
|
318
|
Chris@16
|
319 //
|
Chris@16
|
320 // Begin by getting an initial approximation for the quantity
|
Chris@16
|
321 // eta from the dominant part of the incomplete beta:
|
Chris@16
|
322 //
|
Chris@16
|
323 T eta0;
|
Chris@16
|
324 if(p < q)
|
Chris@16
|
325 eta0 = boost::math::gamma_q_inv(b, p, pol);
|
Chris@16
|
326 else
|
Chris@16
|
327 eta0 = boost::math::gamma_p_inv(b, q, pol);
|
Chris@16
|
328 eta0 /= a;
|
Chris@16
|
329 //
|
Chris@16
|
330 // Define the variables and powers we'll need later on:
|
Chris@16
|
331 //
|
Chris@16
|
332 T mu = b / a;
|
Chris@16
|
333 T w = sqrt(1 + mu);
|
Chris@16
|
334 T w_2 = w * w;
|
Chris@16
|
335 T w_3 = w_2 * w;
|
Chris@16
|
336 T w_4 = w_2 * w_2;
|
Chris@16
|
337 T w_5 = w_3 * w_2;
|
Chris@16
|
338 T w_6 = w_3 * w_3;
|
Chris@16
|
339 T w_7 = w_4 * w_3;
|
Chris@16
|
340 T w_8 = w_4 * w_4;
|
Chris@16
|
341 T w_9 = w_5 * w_4;
|
Chris@16
|
342 T w_10 = w_5 * w_5;
|
Chris@16
|
343 T d = eta0 - mu;
|
Chris@16
|
344 T d_2 = d * d;
|
Chris@16
|
345 T d_3 = d_2 * d;
|
Chris@16
|
346 T d_4 = d_2 * d_2;
|
Chris@16
|
347 T w1 = w + 1;
|
Chris@16
|
348 T w1_2 = w1 * w1;
|
Chris@16
|
349 T w1_3 = w1 * w1_2;
|
Chris@16
|
350 T w1_4 = w1_2 * w1_2;
|
Chris@16
|
351 //
|
Chris@16
|
352 // Now we need to compute the purturbation error terms that
|
Chris@16
|
353 // convert eta0 to eta, these are all polynomials of polynomials.
|
Chris@16
|
354 // Probably these should be re-written to use tabulated data
|
Chris@16
|
355 // (see examples above), but it's less of a win in this case as we
|
Chris@16
|
356 // need to calculate the individual powers for the denominator terms
|
Chris@16
|
357 // anyway, so we might as well use them for the numerator-polynomials
|
Chris@16
|
358 // as well....
|
Chris@16
|
359 //
|
Chris@16
|
360 // Refer to p154-p155 for the details of these expansions:
|
Chris@16
|
361 //
|
Chris@16
|
362 T e1 = (w + 2) * (w - 1) / (3 * w);
|
Chris@16
|
363 e1 += (w_3 + 9 * w_2 + 21 * w + 5) * d / (36 * w_2 * w1);
|
Chris@16
|
364 e1 -= (w_4 - 13 * w_3 + 69 * w_2 + 167 * w + 46) * d_2 / (1620 * w1_2 * w_3);
|
Chris@16
|
365 e1 -= (7 * w_5 + 21 * w_4 + 70 * w_3 + 26 * w_2 - 93 * w - 31) * d_3 / (6480 * w1_3 * w_4);
|
Chris@16
|
366 e1 -= (75 * w_6 + 202 * w_5 + 188 * w_4 - 888 * w_3 - 1345 * w_2 + 118 * w + 138) * d_4 / (272160 * w1_4 * w_5);
|
Chris@16
|
367
|
Chris@16
|
368 T e2 = (28 * w_4 + 131 * w_3 + 402 * w_2 + 581 * w + 208) * (w - 1) / (1620 * w1 * w_3);
|
Chris@16
|
369 e2 -= (35 * w_6 - 154 * w_5 - 623 * w_4 - 1636 * w_3 - 3983 * w_2 - 3514 * w - 925) * d / (12960 * w1_2 * w_4);
|
Chris@16
|
370 e2 -= (2132 * w_7 + 7915 * w_6 + 16821 * w_5 + 35066 * w_4 + 87490 * w_3 + 141183 * w_2 + 95993 * w + 21640) * d_2 / (816480 * w_5 * w1_3);
|
Chris@16
|
371 e2 -= (11053 * w_8 + 53308 * w_7 + 117010 * w_6 + 163924 * w_5 + 116188 * w_4 - 258428 * w_3 - 677042 * w_2 - 481940 * w - 105497) * d_3 / (14696640 * w1_4 * w_6);
|
Chris@16
|
372
|
Chris@16
|
373 T e3 = -((3592 * w_7 + 8375 * w_6 - 1323 * w_5 - 29198 * w_4 - 89578 * w_3 - 154413 * w_2 - 116063 * w - 29632) * (w - 1)) / (816480 * w_5 * w1_2);
|
Chris@16
|
374 e3 -= (442043 * w_9 + 2054169 * w_8 + 3803094 * w_7 + 3470754 * w_6 + 2141568 * w_5 - 2393568 * w_4 - 19904934 * w_3 - 34714674 * w_2 - 23128299 * w - 5253353) * d / (146966400 * w_6 * w1_3);
|
Chris@16
|
375 e3 -= (116932 * w_10 + 819281 * w_9 + 2378172 * w_8 + 4341330 * w_7 + 6806004 * w_6 + 10622748 * w_5 + 18739500 * w_4 + 30651894 * w_3 + 30869976 * w_2 + 15431867 * w + 2919016) * d_2 / (146966400 * w1_4 * w_7);
|
Chris@16
|
376 //
|
Chris@16
|
377 // Combine eta0 and the error terms to compute eta (Second eqaution p155):
|
Chris@16
|
378 //
|
Chris@16
|
379 T eta = eta0 + e1 / a + e2 / (a * a) + e3 / (a * a * a);
|
Chris@16
|
380 //
|
Chris@16
|
381 // Now we need to solve Eq 4.2 to obtain x. For any given value of
|
Chris@16
|
382 // eta there are two solutions to this equation, and since the distribtion
|
Chris@16
|
383 // may be very skewed, these are not related by x ~ 1-x we used when
|
Chris@16
|
384 // implementing section 3 above. However we know that:
|
Chris@16
|
385 //
|
Chris@16
|
386 // cross < x <= 1 ; iff eta < mu
|
Chris@16
|
387 // x == cross ; iff eta == mu
|
Chris@16
|
388 // 0 <= x < cross ; iff eta > mu
|
Chris@16
|
389 //
|
Chris@16
|
390 // Where cross == 1 / (1 + mu)
|
Chris@16
|
391 // Many thanks to Prof Temme for clarifying this point.
|
Chris@16
|
392 //
|
Chris@16
|
393 // Therefore we'll just jump straight into Newton iterations
|
Chris@16
|
394 // to solve Eq 4.2 using these bounds, and simple bisection
|
Chris@16
|
395 // as the first guess, in practice this converges pretty quickly
|
Chris@16
|
396 // and we only need a few digits correct anyway:
|
Chris@16
|
397 //
|
Chris@16
|
398 if(eta <= 0)
|
Chris@16
|
399 eta = tools::min_value<T>();
|
Chris@16
|
400 T u = eta - mu * log(eta) + (1 + mu) * log(1 + mu) - mu;
|
Chris@16
|
401 T cross = 1 / (1 + mu);
|
Chris@16
|
402 T lower = eta < mu ? cross : 0;
|
Chris@16
|
403 T upper = eta < mu ? 1 : cross;
|
Chris@16
|
404 T x = (lower + upper) / 2;
|
Chris@16
|
405 x = tools::newton_raphson_iterate(
|
Chris@16
|
406 temme_root_finder<T>(u, mu), x, lower, upper, policies::digits<T, Policy>() / 2);
|
Chris@16
|
407 #ifdef BOOST_INSTRUMENT
|
Chris@16
|
408 std::cout << "Estimating x with Temme method 3: " << x << std::endl;
|
Chris@16
|
409 #endif
|
Chris@16
|
410 return x;
|
Chris@16
|
411 }
|
Chris@16
|
412
|
Chris@16
|
413 template <class T, class Policy>
|
Chris@16
|
414 struct ibeta_roots
|
Chris@16
|
415 {
|
Chris@16
|
416 ibeta_roots(T _a, T _b, T t, bool inv = false)
|
Chris@16
|
417 : a(_a), b(_b), target(t), invert(inv) {}
|
Chris@16
|
418
|
Chris@16
|
419 boost::math::tuple<T, T, T> operator()(T x)
|
Chris@16
|
420 {
|
Chris@16
|
421 BOOST_MATH_STD_USING // ADL of std names
|
Chris@16
|
422
|
Chris@16
|
423 BOOST_FPU_EXCEPTION_GUARD
|
Chris@16
|
424
|
Chris@16
|
425 T f1;
|
Chris@16
|
426 T y = 1 - x;
|
Chris@16
|
427 T f = ibeta_imp(a, b, x, Policy(), invert, true, &f1) - target;
|
Chris@16
|
428 if(invert)
|
Chris@16
|
429 f1 = -f1;
|
Chris@16
|
430 if(y == 0)
|
Chris@16
|
431 y = tools::min_value<T>() * 64;
|
Chris@16
|
432 if(x == 0)
|
Chris@16
|
433 x = tools::min_value<T>() * 64;
|
Chris@16
|
434
|
Chris@16
|
435 T f2 = f1 * (-y * a + (b - 2) * x + 1);
|
Chris@16
|
436 if(fabs(f2) < y * x * tools::max_value<T>())
|
Chris@16
|
437 f2 /= (y * x);
|
Chris@16
|
438 if(invert)
|
Chris@16
|
439 f2 = -f2;
|
Chris@16
|
440
|
Chris@16
|
441 // make sure we don't have a zero derivative:
|
Chris@16
|
442 if(f1 == 0)
|
Chris@16
|
443 f1 = (invert ? -1 : 1) * tools::min_value<T>() * 64;
|
Chris@16
|
444
|
Chris@16
|
445 return boost::math::make_tuple(f, f1, f2);
|
Chris@16
|
446 }
|
Chris@16
|
447 private:
|
Chris@16
|
448 T a, b, target;
|
Chris@16
|
449 bool invert;
|
Chris@16
|
450 };
|
Chris@16
|
451
|
Chris@16
|
452 template <class T, class Policy>
|
Chris@16
|
453 T ibeta_inv_imp(T a, T b, T p, T q, const Policy& pol, T* py)
|
Chris@16
|
454 {
|
Chris@16
|
455 BOOST_MATH_STD_USING // For ADL of math functions.
|
Chris@16
|
456
|
Chris@16
|
457 //
|
Chris@16
|
458 // The flag invert is set to true if we swap a for b and p for q,
|
Chris@16
|
459 // in which case the result has to be subtracted from 1:
|
Chris@16
|
460 //
|
Chris@16
|
461 bool invert = false;
|
Chris@16
|
462 //
|
Chris@16
|
463 // Handle trivial cases first:
|
Chris@16
|
464 //
|
Chris@16
|
465 if(q == 0)
|
Chris@16
|
466 {
|
Chris@16
|
467 if(py) *py = 0;
|
Chris@16
|
468 return 1;
|
Chris@16
|
469 }
|
Chris@16
|
470 else if(p == 0)
|
Chris@16
|
471 {
|
Chris@16
|
472 if(py) *py = 1;
|
Chris@16
|
473 return 0;
|
Chris@16
|
474 }
|
Chris@16
|
475 else if(a == 1)
|
Chris@16
|
476 {
|
Chris@16
|
477 if(b == 1)
|
Chris@16
|
478 {
|
Chris@16
|
479 if(py) *py = 1 - p;
|
Chris@16
|
480 return p;
|
Chris@16
|
481 }
|
Chris@16
|
482 // Change things around so we can handle as b == 1 special case below:
|
Chris@16
|
483 std::swap(a, b);
|
Chris@16
|
484 std::swap(p, q);
|
Chris@16
|
485 invert = true;
|
Chris@16
|
486 }
|
Chris@16
|
487 //
|
Chris@16
|
488 // Depending upon which approximation method we use, we may end up
|
Chris@16
|
489 // calculating either x or y initially (where y = 1-x):
|
Chris@16
|
490 //
|
Chris@16
|
491 T x = 0; // Set to a safe zero to avoid a
|
Chris@16
|
492 // MSVC 2005 warning C4701: potentially uninitialized local variable 'x' used
|
Chris@16
|
493 // But code inspection appears to ensure that x IS assigned whatever the code path.
|
Chris@16
|
494 T y;
|
Chris@16
|
495
|
Chris@16
|
496 // For some of the methods we can put tighter bounds
|
Chris@16
|
497 // on the result than simply [0,1]:
|
Chris@16
|
498 //
|
Chris@16
|
499 T lower = 0;
|
Chris@16
|
500 T upper = 1;
|
Chris@16
|
501 //
|
Chris@16
|
502 // Student's T with b = 0.5 gets handled as a special case, swap
|
Chris@16
|
503 // around if the arguments are in the "wrong" order:
|
Chris@16
|
504 //
|
Chris@16
|
505 if(a == 0.5f)
|
Chris@16
|
506 {
|
Chris@16
|
507 if(b == 0.5f)
|
Chris@16
|
508 {
|
Chris@16
|
509 x = sin(p * constants::half_pi<T>());
|
Chris@16
|
510 x *= x;
|
Chris@16
|
511 if(py)
|
Chris@16
|
512 {
|
Chris@16
|
513 *py = sin(q * constants::half_pi<T>());
|
Chris@16
|
514 *py *= *py;
|
Chris@16
|
515 }
|
Chris@16
|
516 return x;
|
Chris@16
|
517 }
|
Chris@16
|
518 else if(b > 0.5f)
|
Chris@16
|
519 {
|
Chris@16
|
520 std::swap(a, b);
|
Chris@16
|
521 std::swap(p, q);
|
Chris@16
|
522 invert = !invert;
|
Chris@16
|
523 }
|
Chris@16
|
524 }
|
Chris@16
|
525 //
|
Chris@16
|
526 // Select calculation method for the initial estimate:
|
Chris@16
|
527 //
|
Chris@16
|
528 if((b == 0.5f) && (a >= 0.5f) && (p != 1))
|
Chris@16
|
529 {
|
Chris@16
|
530 //
|
Chris@16
|
531 // We have a Student's T distribution:
|
Chris@16
|
532 x = find_ibeta_inv_from_t_dist(a, p, q, &y, pol);
|
Chris@16
|
533 }
|
Chris@16
|
534 else if(b == 1)
|
Chris@16
|
535 {
|
Chris@16
|
536 if(p < q)
|
Chris@16
|
537 {
|
Chris@16
|
538 if(a > 1)
|
Chris@16
|
539 {
|
Chris@16
|
540 x = pow(p, 1 / a);
|
Chris@101
|
541 y = -boost::math::expm1(log(p) / a, pol);
|
Chris@16
|
542 }
|
Chris@16
|
543 else
|
Chris@16
|
544 {
|
Chris@16
|
545 x = pow(p, 1 / a);
|
Chris@16
|
546 y = 1 - x;
|
Chris@16
|
547 }
|
Chris@16
|
548 }
|
Chris@16
|
549 else
|
Chris@16
|
550 {
|
Chris@101
|
551 x = exp(boost::math::log1p(-q, pol) / a);
|
Chris@101
|
552 y = -boost::math::expm1(boost::math::log1p(-q, pol) / a, pol);
|
Chris@16
|
553 }
|
Chris@16
|
554 if(invert)
|
Chris@16
|
555 std::swap(x, y);
|
Chris@16
|
556 if(py)
|
Chris@16
|
557 *py = y;
|
Chris@16
|
558 return x;
|
Chris@16
|
559 }
|
Chris@16
|
560 else if(a + b > 5)
|
Chris@16
|
561 {
|
Chris@16
|
562 //
|
Chris@16
|
563 // When a+b is large then we can use one of Prof Temme's
|
Chris@16
|
564 // asymptotic expansions, begin by swapping things around
|
Chris@16
|
565 // so that p < 0.5, we do this to avoid cancellations errors
|
Chris@16
|
566 // when p is large.
|
Chris@16
|
567 //
|
Chris@16
|
568 if(p > 0.5)
|
Chris@16
|
569 {
|
Chris@16
|
570 std::swap(a, b);
|
Chris@16
|
571 std::swap(p, q);
|
Chris@16
|
572 invert = !invert;
|
Chris@16
|
573 }
|
Chris@16
|
574 T minv = (std::min)(a, b);
|
Chris@16
|
575 T maxv = (std::max)(a, b);
|
Chris@16
|
576 if((sqrt(minv) > (maxv - minv)) && (minv > 5))
|
Chris@16
|
577 {
|
Chris@16
|
578 //
|
Chris@16
|
579 // When a and b differ by a small amount
|
Chris@16
|
580 // the curve is quite symmetrical and we can use an error
|
Chris@16
|
581 // function to approximate the inverse. This is the cheapest
|
Chris@16
|
582 // of the three Temme expantions, and the calculated value
|
Chris@16
|
583 // for x will never be much larger than p, so we don't have
|
Chris@16
|
584 // to worry about cancellation as long as p is small.
|
Chris@16
|
585 //
|
Chris@16
|
586 x = temme_method_1_ibeta_inverse(a, b, p, pol);
|
Chris@16
|
587 y = 1 - x;
|
Chris@16
|
588 }
|
Chris@16
|
589 else
|
Chris@16
|
590 {
|
Chris@16
|
591 T r = a + b;
|
Chris@16
|
592 T theta = asin(sqrt(a / r));
|
Chris@16
|
593 T lambda = minv / r;
|
Chris@16
|
594 if((lambda >= 0.2) && (lambda <= 0.8) && (r >= 10))
|
Chris@16
|
595 {
|
Chris@16
|
596 //
|
Chris@16
|
597 // The second error function case is the next cheapest
|
Chris@16
|
598 // to use, it brakes down when the result is likely to be
|
Chris@16
|
599 // very small, if a+b is also small, but we can use a
|
Chris@16
|
600 // cheaper expansion there in any case. As before x won't
|
Chris@16
|
601 // be much larger than p, so as long as p is small we should
|
Chris@16
|
602 // be free of cancellation error.
|
Chris@16
|
603 //
|
Chris@16
|
604 T ppa = pow(p, 1/a);
|
Chris@16
|
605 if((ppa < 0.0025) && (a + b < 200))
|
Chris@16
|
606 {
|
Chris@16
|
607 x = ppa * pow(a * boost::math::beta(a, b, pol), 1/a);
|
Chris@16
|
608 }
|
Chris@16
|
609 else
|
Chris@16
|
610 x = temme_method_2_ibeta_inverse(a, b, p, r, theta, pol);
|
Chris@16
|
611 y = 1 - x;
|
Chris@16
|
612 }
|
Chris@16
|
613 else
|
Chris@16
|
614 {
|
Chris@16
|
615 //
|
Chris@16
|
616 // If we get here then a and b are very different in magnitude
|
Chris@16
|
617 // and we need to use the third of Temme's methods which
|
Chris@16
|
618 // involves inverting the incomplete gamma. This is much more
|
Chris@16
|
619 // expensive than the other methods. We also can only use this
|
Chris@16
|
620 // method when a > b, which can lead to cancellation errors
|
Chris@16
|
621 // if we really want y (as we will when x is close to 1), so
|
Chris@16
|
622 // a different expansion is used in that case.
|
Chris@16
|
623 //
|
Chris@16
|
624 if(a < b)
|
Chris@16
|
625 {
|
Chris@16
|
626 std::swap(a, b);
|
Chris@16
|
627 std::swap(p, q);
|
Chris@16
|
628 invert = !invert;
|
Chris@16
|
629 }
|
Chris@16
|
630 //
|
Chris@16
|
631 // Try and compute the easy way first:
|
Chris@16
|
632 //
|
Chris@16
|
633 T bet = 0;
|
Chris@16
|
634 if(b < 2)
|
Chris@16
|
635 bet = boost::math::beta(a, b, pol);
|
Chris@16
|
636 if(bet != 0)
|
Chris@16
|
637 {
|
Chris@16
|
638 y = pow(b * q * bet, 1/b);
|
Chris@16
|
639 x = 1 - y;
|
Chris@16
|
640 }
|
Chris@16
|
641 else
|
Chris@16
|
642 y = 1;
|
Chris@16
|
643 if(y > 1e-5)
|
Chris@16
|
644 {
|
Chris@16
|
645 x = temme_method_3_ibeta_inverse(a, b, p, q, pol);
|
Chris@16
|
646 y = 1 - x;
|
Chris@16
|
647 }
|
Chris@16
|
648 }
|
Chris@16
|
649 }
|
Chris@16
|
650 }
|
Chris@16
|
651 else if((a < 1) && (b < 1))
|
Chris@16
|
652 {
|
Chris@16
|
653 //
|
Chris@16
|
654 // Both a and b less than 1,
|
Chris@16
|
655 // there is a point of inflection at xs:
|
Chris@16
|
656 //
|
Chris@16
|
657 T xs = (1 - a) / (2 - a - b);
|
Chris@16
|
658 //
|
Chris@16
|
659 // Now we need to ensure that we start our iteration from the
|
Chris@16
|
660 // right side of the inflection point:
|
Chris@16
|
661 //
|
Chris@16
|
662 T fs = boost::math::ibeta(a, b, xs, pol) - p;
|
Chris@16
|
663 if(fabs(fs) / p < tools::epsilon<T>() * 3)
|
Chris@16
|
664 {
|
Chris@16
|
665 // The result is at the point of inflection, best just return it:
|
Chris@16
|
666 *py = invert ? xs : 1 - xs;
|
Chris@16
|
667 return invert ? 1-xs : xs;
|
Chris@16
|
668 }
|
Chris@16
|
669 if(fs < 0)
|
Chris@16
|
670 {
|
Chris@16
|
671 std::swap(a, b);
|
Chris@16
|
672 std::swap(p, q);
|
Chris@16
|
673 invert = !invert;
|
Chris@16
|
674 xs = 1 - xs;
|
Chris@16
|
675 }
|
Chris@16
|
676 T xg = pow(a * p * boost::math::beta(a, b, pol), 1/a);
|
Chris@16
|
677 x = xg / (1 + xg);
|
Chris@16
|
678 y = 1 / (1 + xg);
|
Chris@16
|
679 //
|
Chris@16
|
680 // And finally we know that our result is below the inflection
|
Chris@16
|
681 // point, so set an upper limit on our search:
|
Chris@16
|
682 //
|
Chris@16
|
683 if(x > xs)
|
Chris@16
|
684 x = xs;
|
Chris@16
|
685 upper = xs;
|
Chris@16
|
686 }
|
Chris@16
|
687 else if((a > 1) && (b > 1))
|
Chris@16
|
688 {
|
Chris@16
|
689 //
|
Chris@16
|
690 // Small a and b, both greater than 1,
|
Chris@16
|
691 // there is a point of inflection at xs,
|
Chris@16
|
692 // and it's complement is xs2, we must always
|
Chris@16
|
693 // start our iteration from the right side of the
|
Chris@16
|
694 // point of inflection.
|
Chris@16
|
695 //
|
Chris@16
|
696 T xs = (a - 1) / (a + b - 2);
|
Chris@16
|
697 T xs2 = (b - 1) / (a + b - 2);
|
Chris@16
|
698 T ps = boost::math::ibeta(a, b, xs, pol) - p;
|
Chris@16
|
699
|
Chris@16
|
700 if(ps < 0)
|
Chris@16
|
701 {
|
Chris@16
|
702 std::swap(a, b);
|
Chris@16
|
703 std::swap(p, q);
|
Chris@16
|
704 std::swap(xs, xs2);
|
Chris@16
|
705 invert = !invert;
|
Chris@16
|
706 }
|
Chris@16
|
707 //
|
Chris@16
|
708 // Estimate x and y, using expm1 to get a good estimate
|
Chris@16
|
709 // for y when it's very small:
|
Chris@16
|
710 //
|
Chris@16
|
711 T lx = log(p * a * boost::math::beta(a, b, pol)) / a;
|
Chris@16
|
712 x = exp(lx);
|
Chris@16
|
713 y = x < 0.9 ? T(1 - x) : (T)(-boost::math::expm1(lx, pol));
|
Chris@16
|
714
|
Chris@16
|
715 if((b < a) && (x < 0.2))
|
Chris@16
|
716 {
|
Chris@16
|
717 //
|
Chris@16
|
718 // Under a limited range of circumstances we can improve
|
Chris@16
|
719 // our estimate for x, frankly it's clear if this has much effect!
|
Chris@16
|
720 //
|
Chris@16
|
721 T ap1 = a - 1;
|
Chris@16
|
722 T bm1 = b - 1;
|
Chris@16
|
723 T a_2 = a * a;
|
Chris@16
|
724 T a_3 = a * a_2;
|
Chris@16
|
725 T b_2 = b * b;
|
Chris@16
|
726 T terms[5] = { 0, 1 };
|
Chris@16
|
727 terms[2] = bm1 / ap1;
|
Chris@16
|
728 ap1 *= ap1;
|
Chris@16
|
729 terms[3] = bm1 * (3 * a * b + 5 * b + a_2 - a - 4) / (2 * (a + 2) * ap1);
|
Chris@16
|
730 ap1 *= (a + 1);
|
Chris@16
|
731 terms[4] = bm1 * (33 * a * b_2 + 31 * b_2 + 8 * a_2 * b_2 - 30 * a * b - 47 * b + 11 * a_2 * b + 6 * a_3 * b + 18 + 4 * a - a_3 + a_2 * a_2 - 10 * a_2)
|
Chris@16
|
732 / (3 * (a + 3) * (a + 2) * ap1);
|
Chris@16
|
733 x = tools::evaluate_polynomial(terms, x, 5);
|
Chris@16
|
734 }
|
Chris@16
|
735 //
|
Chris@16
|
736 // And finally we know that our result is below the inflection
|
Chris@16
|
737 // point, so set an upper limit on our search:
|
Chris@16
|
738 //
|
Chris@16
|
739 if(x > xs)
|
Chris@16
|
740 x = xs;
|
Chris@16
|
741 upper = xs;
|
Chris@16
|
742 }
|
Chris@16
|
743 else /*if((a <= 1) != (b <= 1))*/
|
Chris@16
|
744 {
|
Chris@16
|
745 //
|
Chris@16
|
746 // If all else fails we get here, only one of a and b
|
Chris@16
|
747 // is above 1, and a+b is small. Start by swapping
|
Chris@16
|
748 // things around so that we have a concave curve with b > a
|
Chris@16
|
749 // and no points of inflection in [0,1]. As long as we expect
|
Chris@16
|
750 // x to be small then we can use the simple (and cheap) power
|
Chris@16
|
751 // term to estimate x, but when we expect x to be large then
|
Chris@16
|
752 // this greatly underestimates x and leaves us trying to
|
Chris@16
|
753 // iterate "round the corner" which may take almost forever...
|
Chris@16
|
754 //
|
Chris@16
|
755 // We could use Temme's inverse gamma function case in that case,
|
Chris@16
|
756 // this works really rather well (albeit expensively) even though
|
Chris@16
|
757 // strictly speaking we're outside it's defined range.
|
Chris@16
|
758 //
|
Chris@16
|
759 // However it's expensive to compute, and an alternative approach
|
Chris@16
|
760 // which models the curve as a distorted quarter circle is much
|
Chris@16
|
761 // cheaper to compute, and still keeps the number of iterations
|
Chris@16
|
762 // required down to a reasonable level. With thanks to Prof Temme
|
Chris@16
|
763 // for this suggestion.
|
Chris@16
|
764 //
|
Chris@16
|
765 if(b < a)
|
Chris@16
|
766 {
|
Chris@16
|
767 std::swap(a, b);
|
Chris@16
|
768 std::swap(p, q);
|
Chris@16
|
769 invert = !invert;
|
Chris@16
|
770 }
|
Chris@16
|
771 if(pow(p, 1/a) < 0.5)
|
Chris@16
|
772 {
|
Chris@16
|
773 x = pow(p * a * boost::math::beta(a, b, pol), 1 / a);
|
Chris@16
|
774 if(x == 0)
|
Chris@16
|
775 x = boost::math::tools::min_value<T>();
|
Chris@16
|
776 y = 1 - x;
|
Chris@16
|
777 }
|
Chris@16
|
778 else /*if(pow(q, 1/b) < 0.1)*/
|
Chris@16
|
779 {
|
Chris@16
|
780 // model a distorted quarter circle:
|
Chris@16
|
781 y = pow(1 - pow(p, b * boost::math::beta(a, b, pol)), 1/b);
|
Chris@16
|
782 if(y == 0)
|
Chris@16
|
783 y = boost::math::tools::min_value<T>();
|
Chris@16
|
784 x = 1 - y;
|
Chris@16
|
785 }
|
Chris@16
|
786 }
|
Chris@16
|
787
|
Chris@16
|
788 //
|
Chris@16
|
789 // Now we have a guess for x (and for y) we can set things up for
|
Chris@16
|
790 // iteration. If x > 0.5 it pays to swap things round:
|
Chris@16
|
791 //
|
Chris@16
|
792 if(x > 0.5)
|
Chris@16
|
793 {
|
Chris@16
|
794 std::swap(a, b);
|
Chris@16
|
795 std::swap(p, q);
|
Chris@16
|
796 std::swap(x, y);
|
Chris@16
|
797 invert = !invert;
|
Chris@16
|
798 T l = 1 - upper;
|
Chris@16
|
799 T u = 1 - lower;
|
Chris@16
|
800 lower = l;
|
Chris@16
|
801 upper = u;
|
Chris@16
|
802 }
|
Chris@16
|
803 //
|
Chris@16
|
804 // lower bound for our search:
|
Chris@16
|
805 //
|
Chris@16
|
806 // We're not interested in denormalised answers as these tend to
|
Chris@16
|
807 // these tend to take up lots of iterations, given that we can't get
|
Chris@16
|
808 // accurate derivatives in this area (they tend to be infinite).
|
Chris@16
|
809 //
|
Chris@16
|
810 if(lower == 0)
|
Chris@16
|
811 {
|
Chris@16
|
812 if(invert && (py == 0))
|
Chris@16
|
813 {
|
Chris@16
|
814 //
|
Chris@16
|
815 // We're not interested in answers smaller than machine epsilon:
|
Chris@16
|
816 //
|
Chris@16
|
817 lower = boost::math::tools::epsilon<T>();
|
Chris@16
|
818 if(x < lower)
|
Chris@16
|
819 x = lower;
|
Chris@16
|
820 }
|
Chris@16
|
821 else
|
Chris@16
|
822 lower = boost::math::tools::min_value<T>();
|
Chris@16
|
823 if(x < lower)
|
Chris@16
|
824 x = lower;
|
Chris@16
|
825 }
|
Chris@16
|
826 //
|
Chris@16
|
827 // Figure out how many digits to iterate towards:
|
Chris@16
|
828 //
|
Chris@16
|
829 int digits = boost::math::policies::digits<T, Policy>() / 2;
|
Chris@16
|
830 if((x < 1e-50) && ((a < 1) || (b < 1)))
|
Chris@16
|
831 {
|
Chris@16
|
832 //
|
Chris@16
|
833 // If we're in a region where the first derivative is very
|
Chris@16
|
834 // large, then we have to take care that the root-finder
|
Chris@16
|
835 // doesn't terminate prematurely. We'll bump the precision
|
Chris@16
|
836 // up to avoid this, but we have to take care not to set the
|
Chris@16
|
837 // precision too high or the last few iterations will just
|
Chris@16
|
838 // thrash around and convergence may be slow in this case.
|
Chris@16
|
839 // Try 3/4 of machine epsilon:
|
Chris@16
|
840 //
|
Chris@16
|
841 digits *= 3;
|
Chris@16
|
842 digits /= 2;
|
Chris@16
|
843 }
|
Chris@16
|
844 //
|
Chris@16
|
845 // Now iterate, we can use either p or q as the target here
|
Chris@16
|
846 // depending on which is smaller:
|
Chris@16
|
847 //
|
Chris@16
|
848 boost::uintmax_t max_iter = policies::get_max_root_iterations<Policy>();
|
Chris@16
|
849 x = boost::math::tools::halley_iterate(
|
Chris@16
|
850 boost::math::detail::ibeta_roots<T, Policy>(a, b, (p < q ? p : q), (p < q ? false : true)), x, lower, upper, digits, max_iter);
|
Chris@16
|
851 policies::check_root_iterations<T>("boost::math::ibeta<%1%>(%1%, %1%, %1%)", max_iter, pol);
|
Chris@16
|
852 //
|
Chris@16
|
853 // We don't really want these asserts here, but they are useful for sanity
|
Chris@16
|
854 // checking that we have the limits right, uncomment if you suspect bugs *only*.
|
Chris@16
|
855 //
|
Chris@16
|
856 //BOOST_ASSERT(x != upper);
|
Chris@16
|
857 //BOOST_ASSERT((x != lower) || (x == boost::math::tools::min_value<T>()) || (x == boost::math::tools::epsilon<T>()));
|
Chris@16
|
858 //
|
Chris@16
|
859 // Tidy up, if we "lower" was too high then zero is the best answer we have:
|
Chris@16
|
860 //
|
Chris@16
|
861 if(x == lower)
|
Chris@16
|
862 x = 0;
|
Chris@16
|
863 if(py)
|
Chris@16
|
864 *py = invert ? x : 1 - x;
|
Chris@16
|
865 return invert ? 1-x : x;
|
Chris@16
|
866 }
|
Chris@16
|
867
|
Chris@16
|
868 } // namespace detail
|
Chris@16
|
869
|
Chris@16
|
870 template <class T1, class T2, class T3, class T4, class Policy>
|
Chris@16
|
871 inline typename tools::promote_args<T1, T2, T3, T4>::type
|
Chris@16
|
872 ibeta_inv(T1 a, T2 b, T3 p, T4* py, const Policy& pol)
|
Chris@16
|
873 {
|
Chris@16
|
874 static const char* function = "boost::math::ibeta_inv<%1%>(%1%,%1%,%1%)";
|
Chris@16
|
875 BOOST_FPU_EXCEPTION_GUARD
|
Chris@16
|
876 typedef typename tools::promote_args<T1, T2, T3, T4>::type result_type;
|
Chris@16
|
877 typedef typename policies::evaluation<result_type, Policy>::type value_type;
|
Chris@16
|
878 typedef typename policies::normalise<
|
Chris@16
|
879 Policy,
|
Chris@16
|
880 policies::promote_float<false>,
|
Chris@16
|
881 policies::promote_double<false>,
|
Chris@16
|
882 policies::discrete_quantile<>,
|
Chris@16
|
883 policies::assert_undefined<> >::type forwarding_policy;
|
Chris@16
|
884
|
Chris@16
|
885 if(a <= 0)
|
Chris@16
|
886 return policies::raise_domain_error<result_type>(function, "The argument a to the incomplete beta function inverse must be greater than zero (got a=%1%).", a, pol);
|
Chris@16
|
887 if(b <= 0)
|
Chris@16
|
888 return policies::raise_domain_error<result_type>(function, "The argument b to the incomplete beta function inverse must be greater than zero (got b=%1%).", b, pol);
|
Chris@16
|
889 if((p < 0) || (p > 1))
|
Chris@16
|
890 return policies::raise_domain_error<result_type>(function, "Argument p outside the range [0,1] in the incomplete beta function inverse (got p=%1%).", p, pol);
|
Chris@16
|
891
|
Chris@16
|
892 value_type rx, ry;
|
Chris@16
|
893
|
Chris@16
|
894 rx = detail::ibeta_inv_imp(
|
Chris@16
|
895 static_cast<value_type>(a),
|
Chris@16
|
896 static_cast<value_type>(b),
|
Chris@16
|
897 static_cast<value_type>(p),
|
Chris@16
|
898 static_cast<value_type>(1 - p),
|
Chris@16
|
899 forwarding_policy(), &ry);
|
Chris@16
|
900
|
Chris@16
|
901 if(py) *py = policies::checked_narrowing_cast<T4, forwarding_policy>(ry, function);
|
Chris@16
|
902 return policies::checked_narrowing_cast<result_type, forwarding_policy>(rx, function);
|
Chris@16
|
903 }
|
Chris@16
|
904
|
Chris@16
|
905 template <class T1, class T2, class T3, class T4>
|
Chris@16
|
906 inline typename tools::promote_args<T1, T2, T3, T4>::type
|
Chris@16
|
907 ibeta_inv(T1 a, T2 b, T3 p, T4* py)
|
Chris@16
|
908 {
|
Chris@16
|
909 return ibeta_inv(a, b, p, py, policies::policy<>());
|
Chris@16
|
910 }
|
Chris@16
|
911
|
Chris@16
|
912 template <class T1, class T2, class T3>
|
Chris@16
|
913 inline typename tools::promote_args<T1, T2, T3>::type
|
Chris@16
|
914 ibeta_inv(T1 a, T2 b, T3 p)
|
Chris@16
|
915 {
|
Chris@16
|
916 typedef typename tools::promote_args<T1, T2, T3>::type result_type;
|
Chris@16
|
917 return ibeta_inv(a, b, p, static_cast<result_type*>(0), policies::policy<>());
|
Chris@16
|
918 }
|
Chris@16
|
919
|
Chris@16
|
920 template <class T1, class T2, class T3, class Policy>
|
Chris@16
|
921 inline typename tools::promote_args<T1, T2, T3>::type
|
Chris@16
|
922 ibeta_inv(T1 a, T2 b, T3 p, const Policy& pol)
|
Chris@16
|
923 {
|
Chris@16
|
924 typedef typename tools::promote_args<T1, T2, T3>::type result_type;
|
Chris@16
|
925 return ibeta_inv(a, b, p, static_cast<result_type*>(0), pol);
|
Chris@16
|
926 }
|
Chris@16
|
927
|
Chris@16
|
928 template <class T1, class T2, class T3, class T4, class Policy>
|
Chris@16
|
929 inline typename tools::promote_args<T1, T2, T3, T4>::type
|
Chris@16
|
930 ibetac_inv(T1 a, T2 b, T3 q, T4* py, const Policy& pol)
|
Chris@16
|
931 {
|
Chris@16
|
932 static const char* function = "boost::math::ibetac_inv<%1%>(%1%,%1%,%1%)";
|
Chris@16
|
933 BOOST_FPU_EXCEPTION_GUARD
|
Chris@16
|
934 typedef typename tools::promote_args<T1, T2, T3, T4>::type result_type;
|
Chris@16
|
935 typedef typename policies::evaluation<result_type, Policy>::type value_type;
|
Chris@16
|
936 typedef typename policies::normalise<
|
Chris@16
|
937 Policy,
|
Chris@16
|
938 policies::promote_float<false>,
|
Chris@16
|
939 policies::promote_double<false>,
|
Chris@16
|
940 policies::discrete_quantile<>,
|
Chris@16
|
941 policies::assert_undefined<> >::type forwarding_policy;
|
Chris@16
|
942
|
Chris@16
|
943 if(a <= 0)
|
Chris@101
|
944 return policies::raise_domain_error<result_type>(function, "The argument a to the incomplete beta function inverse must be greater than zero (got a=%1%).", a, pol);
|
Chris@16
|
945 if(b <= 0)
|
Chris@101
|
946 return policies::raise_domain_error<result_type>(function, "The argument b to the incomplete beta function inverse must be greater than zero (got b=%1%).", b, pol);
|
Chris@16
|
947 if((q < 0) || (q > 1))
|
Chris@101
|
948 return policies::raise_domain_error<result_type>(function, "Argument q outside the range [0,1] in the incomplete beta function inverse (got q=%1%).", q, pol);
|
Chris@16
|
949
|
Chris@16
|
950 value_type rx, ry;
|
Chris@16
|
951
|
Chris@16
|
952 rx = detail::ibeta_inv_imp(
|
Chris@16
|
953 static_cast<value_type>(a),
|
Chris@16
|
954 static_cast<value_type>(b),
|
Chris@16
|
955 static_cast<value_type>(1 - q),
|
Chris@16
|
956 static_cast<value_type>(q),
|
Chris@16
|
957 forwarding_policy(), &ry);
|
Chris@16
|
958
|
Chris@16
|
959 if(py) *py = policies::checked_narrowing_cast<T4, forwarding_policy>(ry, function);
|
Chris@16
|
960 return policies::checked_narrowing_cast<result_type, forwarding_policy>(rx, function);
|
Chris@16
|
961 }
|
Chris@16
|
962
|
Chris@16
|
963 template <class T1, class T2, class T3, class T4>
|
Chris@16
|
964 inline typename tools::promote_args<T1, T2, T3, T4>::type
|
Chris@16
|
965 ibetac_inv(T1 a, T2 b, T3 q, T4* py)
|
Chris@16
|
966 {
|
Chris@16
|
967 return ibetac_inv(a, b, q, py, policies::policy<>());
|
Chris@16
|
968 }
|
Chris@16
|
969
|
Chris@16
|
970 template <class RT1, class RT2, class RT3>
|
Chris@16
|
971 inline typename tools::promote_args<RT1, RT2, RT3>::type
|
Chris@16
|
972 ibetac_inv(RT1 a, RT2 b, RT3 q)
|
Chris@16
|
973 {
|
Chris@16
|
974 typedef typename tools::promote_args<RT1, RT2, RT3>::type result_type;
|
Chris@16
|
975 return ibetac_inv(a, b, q, static_cast<result_type*>(0), policies::policy<>());
|
Chris@16
|
976 }
|
Chris@16
|
977
|
Chris@16
|
978 template <class RT1, class RT2, class RT3, class Policy>
|
Chris@16
|
979 inline typename tools::promote_args<RT1, RT2, RT3>::type
|
Chris@16
|
980 ibetac_inv(RT1 a, RT2 b, RT3 q, const Policy& pol)
|
Chris@16
|
981 {
|
Chris@16
|
982 typedef typename tools::promote_args<RT1, RT2, RT3>::type result_type;
|
Chris@16
|
983 return ibetac_inv(a, b, q, static_cast<result_type*>(0), pol);
|
Chris@16
|
984 }
|
Chris@16
|
985
|
Chris@16
|
986 } // namespace math
|
Chris@16
|
987 } // namespace boost
|
Chris@16
|
988
|
Chris@16
|
989 #endif // BOOST_MATH_SPECIAL_FUNCTIONS_IGAMMA_INVERSE_HPP
|
Chris@16
|
990
|
Chris@16
|
991
|
Chris@16
|
992
|
Chris@16
|
993
|