Mercurial > octave
annotate scripts/specfun/betaincinv.m @ 29665:5f37ef6e7114
maint: merge stable to default.
author | Rik <rik@octave.org> |
---|---|
date | Fri, 14 May 2021 08:12:42 -0700 |
parents | ce4436d2b206 2a1f57067fbf |
children | 26bb2cbf6da2 |
rev | line source |
---|---|
27923
bd51beb6205e
update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents:
27919
diff
changeset
|
1 ######################################################################## |
bd51beb6205e
update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents:
27919
diff
changeset
|
2 ## |
29358
0a5b15007766
update Octave Project Developers copyright for the new year
John W. Eaton <jwe@octave.org>
parents:
27984
diff
changeset
|
3 ## Copyright (C) 2017-2021 The Octave Project Developers |
27918
b442ec6dda5c
use centralized file for copyright info for individual contributors
John W. Eaton <jwe@octave.org>
parents:
26376
diff
changeset
|
4 ## |
27923
bd51beb6205e
update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents:
27919
diff
changeset
|
5 ## See the file COPYRIGHT.md in the top-level directory of this |
bd51beb6205e
update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents:
27919
diff
changeset
|
6 ## distribution or <https://octave.org/copyright/>. |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
7 ## |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
8 ## This file is part of Octave. |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
9 ## |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
10 ## Octave is free software: you can redistribute it and/or modify it |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
11 ## under the terms of the GNU General Public License as published by |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
12 ## the Free Software Foundation, either version 3 of the License, or |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
13 ## (at your option) any later version. |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
14 ## |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
15 ## Octave is distributed in the hope that it will be useful, but |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
16 ## WITHOUT ANY WARRANTY; without even the implied warranty of |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
17 ## MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
18 ## GNU General Public License for more details. |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
19 ## |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
20 ## You should have received a copy of the GNU General Public License |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
21 ## along with Octave; see the file COPYING. If not, see |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
22 ## <https://www.gnu.org/licenses/>. |
27923
bd51beb6205e
update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents:
27919
diff
changeset
|
23 ## |
bd51beb6205e
update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents:
27919
diff
changeset
|
24 ######################################################################## |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
25 |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
26 ## -*- texinfo -*- |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
27 ## @deftypefn {} {} betaincinv (@var{y}, @var{a}, @var{b}) |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
28 ## @deftypefnx {} {} betaincinv (@var{y}, @var{a}, @var{b}, "lower") |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
29 ## @deftypefnx {} {} betaincinv (@var{y}, @var{a}, @var{b}, "upper") |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
30 ## Compute the inverse of the normalized incomplete beta function. |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
31 ## |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
32 ## The normalized incomplete beta function is defined as |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
33 ## @tex |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
34 ## $$ |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
35 ## I_x (a, b) = {1 \over {B(a,b)}} \displaystyle{\int_0^x t^{a-1} (1-t)^{b-1} dt} |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
36 ## $$ |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
37 ## @end tex |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
38 ## @ifnottex |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
39 ## |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
40 ## @example |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
41 ## @group |
29664
2a1f57067fbf
betainc.m, betaincinv.m: Correct non-Tex definition of beta incomplete integral.
Rik <rik@octave.org>
parents:
29656
diff
changeset
|
42 ## x |
2a1f57067fbf
betainc.m, betaincinv.m: Correct non-Tex definition of beta incomplete integral.
Rik <rik@octave.org>
parents:
29656
diff
changeset
|
43 ## / |
2a1f57067fbf
betainc.m, betaincinv.m: Correct non-Tex definition of beta incomplete integral.
Rik <rik@octave.org>
parents:
29656
diff
changeset
|
44 ## 1 | |
2a1f57067fbf
betainc.m, betaincinv.m: Correct non-Tex definition of beta incomplete integral.
Rik <rik@octave.org>
parents:
29656
diff
changeset
|
45 ## I_x (a, b) = ---------- | t^(a-1) (1-t)^(b-1) dt |
2a1f57067fbf
betainc.m, betaincinv.m: Correct non-Tex definition of beta incomplete integral.
Rik <rik@octave.org>
parents:
29656
diff
changeset
|
46 ## beta (a,b) | |
2a1f57067fbf
betainc.m, betaincinv.m: Correct non-Tex definition of beta incomplete integral.
Rik <rik@octave.org>
parents:
29656
diff
changeset
|
47 ## / |
2a1f57067fbf
betainc.m, betaincinv.m: Correct non-Tex definition of beta incomplete integral.
Rik <rik@octave.org>
parents:
29656
diff
changeset
|
48 ## 0 |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
49 ## @end group |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
50 ## @end example |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
51 ## |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
52 ## @end ifnottex |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
53 ## |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
54 ## If two inputs are scalar, then @code{betaincinv (@var{y}, @var{a}, @var{b})} |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
55 ## is returned for each of the other inputs. |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
56 ## |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
57 ## If two or more inputs are not scalar, the sizes of them must agree, and |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
58 ## @code{betaincinv} is applied element-by-element. |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
59 ## |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
60 ## The variable @var{y} must be in the interval [0,1], while @var{a} and |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
61 ## @var{b} must be real and strictly positive. |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
62 ## |
25579
07c2c42f457e
doc: Miscellaneous documentation fixes all over the manual (bug #54288).
Rik <rik@octave.org>
parents:
25054
diff
changeset
|
63 ## By default, @var{tail} is @qcode{"lower"} and the inverse of the incomplete |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
64 ## beta function integrated from 0 to @var{x} is computed. If @var{tail} is |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
65 ## @qcode{"upper"} then the complementary function integrated from @var{x} to 1 |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
66 ## is inverted. |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
67 ## |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
68 ## The function is computed by standard Newton's method, by solving |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
69 ## @tex |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
70 ## $$ |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
71 ## y - I_x (a, b) = 0 |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
72 ## $$ |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
73 ## @end tex |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
74 ## @ifnottex |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
75 ## |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
76 ## @example |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
77 ## @var{y} - betainc (@var{x}, @var{a}, @var{b}) = 0 |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
78 ## @end example |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
79 ## |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
80 ## @end ifnottex |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
81 ## |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
82 ## @seealso{betainc, beta, betaln} |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
83 ## @end deftypefn |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
84 |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
85 function x = betaincinv (y, a, b, tail = "lower") |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
86 |
29656
7f5bd197fea6
maint: use std::size_t in more instances (bug #60471)
Rik <rik@octave.org>
parents:
29616
diff
changeset
|
87 if (nargin < 3) |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
88 print_usage (); |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
89 endif |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
90 |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
91 [err, y, a, b] = common_size (y, a, b); |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
92 if (err > 0) |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
93 error ("betaincinv: Y, A, and B must be of common size or scalars"); |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
94 endif |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
95 |
29656
7f5bd197fea6
maint: use std::size_t in more instances (bug #60471)
Rik <rik@octave.org>
parents:
29616
diff
changeset
|
96 if (! (isfloat (y) && isfloat (a) && isfloat (b) |
7f5bd197fea6
maint: use std::size_t in more instances (bug #60471)
Rik <rik@octave.org>
parents:
29616
diff
changeset
|
97 && isreal (y) && isreal (a) && isreal (b))) |
7f5bd197fea6
maint: use std::size_t in more instances (bug #60471)
Rik <rik@octave.org>
parents:
29616
diff
changeset
|
98 error ("betaincinv: Y, A, and B must be real, floating point values"); |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
99 endif |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
100 |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
101 ## Remember original shape of data, but convert to column vector for calcs. |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
102 orig_sz = size (y); |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
103 y = y(:); |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
104 a = a(:); |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
105 b = b(:); |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
106 |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
107 if (any ((y < 0) | (y > 1))) |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
108 error ("betaincinv: Y must be in the range [0, 1]"); |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
109 endif |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
110 |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
111 if (any (a <= 0)) |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
112 error ("betaincinv: A must be strictly positive"); |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
113 endif |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
114 |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
115 if (any (b <= 0)) |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
116 error ("betaincinv: B must be strictly positive"); |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
117 endif |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
118 |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
119 ## If any of the arguments is single then the output should be as well. |
29656
7f5bd197fea6
maint: use std::size_t in more instances (bug #60471)
Rik <rik@octave.org>
parents:
29616
diff
changeset
|
120 if (isa (y, "single") || isa (a, "single") || isa (b, "single")) |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
121 y = single (y); |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
122 a = single (a); |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
123 b = single (b); |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
124 endif |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
125 |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
126 ## Initialize output array |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
127 x = zeros (size (y), class (y)); |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
128 |
29656
7f5bd197fea6
maint: use std::size_t in more instances (bug #60471)
Rik <rik@octave.org>
parents:
29616
diff
changeset
|
129 ## Parameters for the Newton's method search |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
130 maxit = 20; |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
131 tol = eps (class (y)); |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
132 |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
133 if (strcmpi (tail, "lower")) |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
134 p = y; |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
135 q = 1 - y; |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
136 x(y == 0) = 0; |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
137 x(y == 1) = 1; |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
138 elseif (strcmpi (tail, "upper")) |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
139 p = 1 - y; |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
140 q = y; |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
141 x(y == 0) = 1; |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
142 x(y == 1) = 0; |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
143 else |
28928
ae7ce8358953
maint: Add semicolon to end of all warning() and error() invocations.
Rik <rik@octave.org>
parents:
28912
diff
changeset
|
144 error ("betaincinv: invalid value for TAIL"); |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
145 endif |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
146 |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
147 ## Special values have been already computed. |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
148 todo = (y != 0) & (y != 1); |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
149 |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
150 ## We will invert the lower version for p < 0.5 and the upper otherwise. |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
151 i_low = (p < 0.5); |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
152 i_upp = (! i_low); |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
153 |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
154 idx = todo & i_low; |
26268
6dd232798997
maint: Remove useless ';' from end of for, if, while, etc. statements.
Rik <rik@octave.org>
parents:
25579
diff
changeset
|
155 if (any (idx)) |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
156 n = nnz (idx); |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
157 ## Function and derivative of the lower version. |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
158 F = @(x, a, b, y) y - betainc (x, a, b); |
29656
7f5bd197fea6
maint: use std::size_t in more instances (bug #60471)
Rik <rik@octave.org>
parents:
29616
diff
changeset
|
159 JF = @(x, a, b) -real (exp ((a-1) .* log (x) + (b-1) .* log1p (-x) + ... |
7f5bd197fea6
maint: use std::size_t in more instances (bug #60471)
Rik <rik@octave.org>
parents:
29616
diff
changeset
|
160 gammaln (a+b) - gammaln (a) - gammaln (b))); |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
161 |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
162 ## Compute the initial guess with a bisection method of 10 steps. |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
163 x0 = bisection_method (F, zeros (n,1), ones (n,1), ... |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
164 a(i_low), b(i_low), p(i_low), 10); |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
165 |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
166 ## Use Newton's method to iteratively find solution. |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
167 x(i_low) = newton_method (F, JF, x0, a(i_low), b(i_low), p(i_low), ... |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
168 tol, maxit); |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
169 endif |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
170 |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
171 idx = todo & i_upp; |
26268
6dd232798997
maint: Remove useless ';' from end of for, if, while, etc. statements.
Rik <rik@octave.org>
parents:
25579
diff
changeset
|
172 if (any (idx)) |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
173 n = nnz (idx); |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
174 ## Function and derivative of the upper version. |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
175 F = @(x, a, b, y) y - betainc (x, a, b, "upper"); |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
176 JF = @(x, a, b) real (exp ((a-1) .* log (x) + (b-1) .* log1p (-x) + ... |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
177 gammaln (a+b) - gammaln (a) - gammaln (b))); |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
178 |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
179 ## Compute the initial guess with a bisection method of 10 steps. |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
180 x0 = bisection_method (F, zeros (n,1), ones (n,1), ... |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
181 a(i_upp), b(i_upp), q(i_upp), 10); |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
182 |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
183 ## Use Newton's method to iteratively find solution. |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
184 x(i_upp) = newton_method (F, JF, x0, a(i_upp), b(i_upp), q(i_upp), ... |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
185 tol, maxit); |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
186 endif |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
187 |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
188 ## Restore original shape |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
189 x = reshape (x, orig_sz); |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
190 |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
191 endfunction |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
192 |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
193 |
27984
b09432b20a84
maint: Remove special cases of old version control keywords in code base.
Rik <rik@octave.org>
parents:
27978
diff
changeset
|
194 ## subfunctions: Bisection and Newton Methods |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
195 function xc = bisection_method (F, xl, xr, a, b, y, maxit) |
28945
6e460773bdda
maint: Use newlines after "function" and before "endfunction" for clarity.
Rik <rik@octave.org>
parents:
28928
diff
changeset
|
196 |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
197 F_l = F (xl, a, b, y); |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
198 F_r = F (xr, a, b, y); |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
199 for it = 1:maxit |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
200 xc = (xl + xr) / 2; |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
201 F_c = F (xc, a, b, y); |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
202 flag_l = ((F_c .* F_r) < 0); |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
203 flag_r = ((F_c .* F_l) < 0); |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
204 flag_c = (F_c == 0); |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
205 xl(flag_l) = xc(flag_l); |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
206 xr(flag_r) = xc(flag_r); |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
207 xl(flag_c) = xr(flag_c) = xc(flag_c); |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
208 F_l(flag_l) = F_c(flag_l); |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
209 F_r(flag_r) = F_c(flag_r); |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
210 F_l(flag_c) = F_r(flag_c) = 0; |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
211 endfor |
28945
6e460773bdda
maint: Use newlines after "function" and before "endfunction" for clarity.
Rik <rik@octave.org>
parents:
28928
diff
changeset
|
212 |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
213 endfunction |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
214 |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
215 function x = newton_method (F, JF, x0, a, b, y, tol, maxit); |
28945
6e460773bdda
maint: Use newlines after "function" and before "endfunction" for clarity.
Rik <rik@octave.org>
parents:
28928
diff
changeset
|
216 |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
217 res = -F (x0, a, b, y) ./ JF (x0, a, b); |
29656
7f5bd197fea6
maint: use std::size_t in more instances (bug #60471)
Rik <rik@octave.org>
parents:
29616
diff
changeset
|
218 todo = (abs (res) >= tol * abs (x0)); |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
219 x = x0; |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
220 it = 0; |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
221 while (any (todo) && (it < maxit)) |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
222 it++; |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
223 x(todo) += res(todo); |
29616
69b6b783a8ab
betaincinv.m: Correctly handle inputs very close to 1.0 (bug #60528)
Rik <rik@octave.org>
parents:
29614
diff
changeset
|
224 ## Avoid intermediate values outside betainc() range of [0, 1], bug #60528 |
69b6b783a8ab
betaincinv.m: Correctly handle inputs very close to 1.0 (bug #60528)
Rik <rik@octave.org>
parents:
29614
diff
changeset
|
225 x(x(todo) < 0) = eps; |
69b6b783a8ab
betaincinv.m: Correctly handle inputs very close to 1.0 (bug #60528)
Rik <rik@octave.org>
parents:
29614
diff
changeset
|
226 x(x(todo) > 1) = 1-eps; |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
227 res(todo) = -F(x(todo), a(todo), b(todo), y(todo)) ... |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
228 ./ JF (x(todo), a(todo), b(todo)); |
29616
69b6b783a8ab
betaincinv.m: Correctly handle inputs very close to 1.0 (bug #60528)
Rik <rik@octave.org>
parents:
29614
diff
changeset
|
229 todo = (abs (res) >= tol * abs (x)); |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
230 endwhile |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
231 x += res; |
28945
6e460773bdda
maint: Use newlines after "function" and before "endfunction" for clarity.
Rik <rik@octave.org>
parents:
28928
diff
changeset
|
232 |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
233 endfunction |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
234 |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
235 |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
236 %!test |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
237 %! x = linspace (0.1, 0.9, 11); |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
238 %! a = [2, 3, 4]; |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
239 %! [x,a,b] = ndgrid (x,a,a); |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
240 %! xx = betaincinv (betainc (x, a, b), a, b); |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
241 %! assert (xx, x, 3e-15); |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
242 |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
243 %!test |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
244 %! x = linspace (0.1, 0.9, 11); |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
245 %! a = [2, 3, 4]; |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
246 %! [x,a,b] = ndgrid (x,a,a); |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
247 %! xx = betaincinv (betainc (x, a, b, "upper"), a, b, "upper"); |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
248 %! assert (xx, x, 3e-15); |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
249 |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
250 %!test |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
251 %! x = linspace (0.1, 0.9, 11); |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
252 %! a = [0.1:0.1:1]; |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
253 %! [x,a,b] = ndgrid (x,a,a); |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
254 %! xx = betaincinv (betainc (x, a, b), a, b); |
24996
f80e68529bbf
test: relax tolerances on specfun tests to pass on i386 (bug #53437)
Mike Miller <mtmiller@octave.org>
parents:
24927
diff
changeset
|
255 %! assert (xx, x, 5e-15); |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
256 |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
257 %!test |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
258 %! x = linspace (0.1, 0.9, 11); |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
259 %! a = [0.1:0.1:1]; |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
260 %! [x,a,b] = ndgrid (x,a,a); |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
261 %! xx = betaincinv (betainc (x, a, b, "upper"), a, b, "upper"); |
24996
f80e68529bbf
test: relax tolerances on specfun tests to pass on i386 (bug #53437)
Mike Miller <mtmiller@octave.org>
parents:
24927
diff
changeset
|
262 %! assert (xx, x, 5e-15); |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
263 |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
264 ## Test the conservation of the input class |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
265 %!assert (class (betaincinv (0.5, 1, 1)), "double") |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
266 %!assert (class (betaincinv (single (0.5), 1, 1)), "single") |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
267 %!assert (class (betaincinv (0.5, single (1), 1)), "single") |
29656
7f5bd197fea6
maint: use std::size_t in more instances (bug #60471)
Rik <rik@octave.org>
parents:
29616
diff
changeset
|
268 %!assert (class (betaincinv (0.5, 1, single (1))), "single") |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
269 |
29616
69b6b783a8ab
betaincinv.m: Correctly handle inputs very close to 1.0 (bug #60528)
Rik <rik@octave.org>
parents:
29614
diff
changeset
|
270 %!assert <*60528> (betaincinv (1e-6, 1, 3), 3.3333344444450657e-07, 5*eps) |
69b6b783a8ab
betaincinv.m: Correctly handle inputs very close to 1.0 (bug #60528)
Rik <rik@octave.org>
parents:
29614
diff
changeset
|
271 %!assert <*60528> (betaincinv (1-1e-6, 3, 1), 0.999999666666556, 5*eps) |
29614
26ba91f0eea7
betaincinv.m: Correctly handle small inputs (bug #60528)
Rik <rik@octave.org>
parents:
29358
diff
changeset
|
272 |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
273 ## Test input validation |
28896
90fea9cc9caa
test: Add expected error message <Invalid call> to BIST tests for nargin.
Rik <rik@octave.org>
parents:
28789
diff
changeset
|
274 %!error <Invalid call> betaincinv () |
90fea9cc9caa
test: Add expected error message <Invalid call> to BIST tests for nargin.
Rik <rik@octave.org>
parents:
28789
diff
changeset
|
275 %!error <Invalid call> betaincinv (1) |
90fea9cc9caa
test: Add expected error message <Invalid call> to BIST tests for nargin.
Rik <rik@octave.org>
parents:
28789
diff
changeset
|
276 %!error <Invalid call> betaincinv (1,2) |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
277 %!error <must be of common size or scalars> |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
278 %! betaincinv (ones (2,2), ones (1,2), 1); |
29656
7f5bd197fea6
maint: use std::size_t in more instances (bug #60471)
Rik <rik@octave.org>
parents:
29616
diff
changeset
|
279 %!error <must be .* floating point> betaincinv ('a', 1, 2) |
7f5bd197fea6
maint: use std::size_t in more instances (bug #60471)
Rik <rik@octave.org>
parents:
29616
diff
changeset
|
280 %!error <must be .* floating point> betaincinv (0, int8 (1), 1) |
7f5bd197fea6
maint: use std::size_t in more instances (bug #60471)
Rik <rik@octave.org>
parents:
29616
diff
changeset
|
281 %!error <must be .* floating point> betaincinv (0, 1, true) |
7f5bd197fea6
maint: use std::size_t in more instances (bug #60471)
Rik <rik@octave.org>
parents:
29616
diff
changeset
|
282 %!error <must be real> betaincinv (0.5i, 1, 2) |
7f5bd197fea6
maint: use std::size_t in more instances (bug #60471)
Rik <rik@octave.org>
parents:
29616
diff
changeset
|
283 %!error <must be real> betaincinv (0, 1i, 1) |
7f5bd197fea6
maint: use std::size_t in more instances (bug #60471)
Rik <rik@octave.org>
parents:
29616
diff
changeset
|
284 %!error <must be real> betaincinv (0, 1, 1i) |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
285 %!error <Y must be in the range \[0, 1\]> betaincinv (-0.1,1,1) |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
286 %!error <Y must be in the range \[0, 1\]> betaincinv (1.1,1,1) |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
287 %!error <Y must be in the range \[0, 1\]> |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
288 %! y = ones (1, 1, 2); |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
289 %! y(1,1,2) = -1; |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
290 %! betaincinv (y,1,1); |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
291 %!error <A must be strictly positive> betaincinv (0.5,0,1) |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
292 %!error <A must be strictly positive> |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
293 %! a = ones (1, 1, 2); |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
294 %! a(1,1,2) = 0; |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
295 %! betaincinv (1,a,1); |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
296 %!error <B must be strictly positive> betaincinv (0.5,1,0) |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
297 %!error <B must be strictly positive> |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
298 %! b = ones (1, 1, 2); |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
299 %! b(1,1,2) = 0; |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
300 %! betaincinv (1,1,b); |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24907
diff
changeset
|
301 %!error <invalid value for TAIL> betaincinv (1,2,3, "foobar") |