Mercurial > octave
annotate scripts/specfun/betainc.m @ 24927:c280560d9c96 stable
Overhaul special functions modified by GSOC2018 project.
* NEWS: Add note about new functions added. Add note explaining the changes
done to existing functions.
* __betainc__.cc: Renamed from __betainc_lentz__.cc. Use standard GPL v3
copyright block. Add missing #include "dNDArray.h". Add one-line texinfo
documentation for internal function. Remove fourth input argument to function.
Use is_single_type() to decide whether to operate with FloatNDArray or NDArray.
Delete temporary variables x_arg_s, a_arg_s, b_arg_s used in input validation.
Rename len_x, len_a, len_b to numel_[xab] for clarity. Remove input validation
that numel match (internal function, no need). Rename function output to
retval according to Octave coding conventions. Use more spacing and newlines
for readability of code.
* __expint__.cc: Renamed from __expint_lentz__.cc. Use standard GPL v3
copyright block. Cull list of #includes to just the necessary ones. Use
Complex and FloatComplex typedefs defined by Octave. Add one-line texinfo
documentation for internal function. Remove second input argument to function.
Use is_single_type() to decide whether to operate with FloatComplexNDArray or
ComplexNDArray. Delete temporary variables x_arg_s used in input validation.
Rename len_x to numel_x for clarity. Use constructor with dim_vector and
scalar value rather than fill() after creating array. Rename function output
to retval according to Octave coding conventions. Use more spacing and
newlines for readability of code.
* __gammainc__.cc: Renamed from __gammainc_lentz__.cc. Use standard GPL v3
copyright block. Add one-line texinfo documentation for internal function.
Remove third input argument to function. Remove input validation that numel
match (internal function, no need). Use is_single_type() to decide whether to
operate with FloatNDArray or NDArray. Delete temporary variables x_arg, a_arg
used in input validation. Use constructor with dim_vector and scalar value
rather than fill() after creating array. Rename function output to retval
according to Octave coding conventions. Use more spacing and newlines for
readability of code.
* __betainc_lentz__.cc, __expint_lentz__.cc, __gammainc_lentz__.cc: Removed.
* libinterp/corefcn/module.mk: Add renamed functions to build system.
* betainc.m: Use Octave standard GPL block. Rewrite parts of docstring.
Don't use array brackets around single output of function. Remove isscalar
checks on inputs because common_size() function will already handle it.
Use capital variable names in error messages to match documentation as displayed
in terminal. Reshape all inputs in to column vectors quickly so that input
validation tests that depend on all/any will pass with N-D arrays. Add
comments to code. Check for specific error messages in input validation BIST
tests.
* betaincinv.m: Use Octave standard GPL block. Rewrite parts of docstring.
Don't use array brackets around single output of function. Remove isscalar
checks on inputs because common_size() function will already handle it.
Use capital variable names in error messages to match documentation as displayed
in terminal. Reshape all inputs in to column vectors quickly so that input
validation tests that depend on all/any will pass with N-D arrays. Put most
common case of tail ("lower") first in if/elseif trees. Call functions directly
with function handle rather than using unnecessary feval() call. Use numel
in preference to length. Rename variable i_miss to todo for clarity. Add
comments to code. Check for specific error messages in input validation BIST
tests.
* cosint.m: Use Octave standard GPL block. Rewrite parts of docstring.
Don't use array brackets around single output of function. Remove isscalar
checks on inputs because common_size() function will already handle it.
Add input validation check for isnumeric value. Convert integer classes to
double before proceeding. Rename i_miss to todo for clarity. Use isinf to
detect both -Inf and +Inf rather than separate tests. Use ++it in while
loop conditional to shorten loop blocks. Add Input validation BIST tests.
* expint.m: Remove Sylvain who was not actually an author on this file.
Rewrite parts of docstring. Rename variable sparse_x to orig_sparse.
Eliminate temporary variables res_tmp, x_s_tmp, ssum_tmp. Rename i_miss to
todo. Use Octave coding conventions throughout. Add comments to code.
* gammainc.m: Use Octave standard GPL block. Rewrite parts of docstring.
Remove isscalar checks on inputs because common_size() function will already
handle it. Use capital variable names in error messages to match documentation
as displayed in terminal. Reshape all inputs in to column vectors quickly so
that input validation tests that depend on all/any will pass with N-D arrays.
Put most common case of tail ("lower") first in if/elseif trees. Add input
validation of tail. Rename variable ii to idx for clarity. Rename variable
i_done to todo and switch polarity so that the '!' operator is not required
every time the variable is updated. Use indexing and direct assignment to
update todo rather than logical operator '&' which is slower. Use tolower on
tail variable and then switch strcmpi calls to strcmp. Reformat %!test blocks
in to %!assert blocks to be more compact. Check for specific error messages in
input validation BIST tests.
* gammaincinv.m: Use Octave standard GPL block. Rewrite parts of docstring.
Don't use array brackets around single output of function. Remove isscalar
checks on inputs because common_size() function will already handle it.
Add input validation check for iscomplex value. Use capital variable names in
error messages to match documentation as displayed in terminal. Reshape all
inputs in to column vectors quickly so that input validation tests that depend
on all/any will pass with N-D arrays. Rename i_miss to todo. Use numel
in preference to length. Call functions directly with function handle rather
than using unnecessary feval() call.Rename variable i_miss to todo for clarity.
Use it++ in while loop conditional to shorten loop blocks. Add comments to
code. Check for specific error messages in input validation BIST tests.
Add input validation BIST tests for all error messages.
* sinint.m: Use Octave standard GPL block. Rewrite parts of docstring.
Add input validation for isnumeric. Convert integers to double for
calculation. Reshape input to column vector. Rename variable sz to orig_sz
for clarity. rename i_miss to todo. Reformat BIST tests to mak them more
compact. Add input validation BIST tests.
author | Rik <rik@octave.org> |
---|---|
date | Mon, 19 Mar 2018 10:01:48 -0700 |
parents | 40ab8be7d7ec |
children | f80e68529bbf |
rev | line source |
---|---|
24921
6b33ee8aad0f
Fixed authors in betainc
Michele Ginesi <michele.ginesi@gmail.com>
parents:
24919
diff
changeset
|
1 ## Copyright (C) 2018 Stefan Schlögl |
6b33ee8aad0f
Fixed authors in betainc
Michele Ginesi <michele.ginesi@gmail.com>
parents:
24919
diff
changeset
|
2 ## Copyright (C) 2018 Michele Ginesi |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
3 ## |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
4 ## 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
|
5 ## |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
6 ## 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
|
7 ## under the terms of the GNU General Public License as published by |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
8 ## the Free Software Foundation; either version 3 of the License, or |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
9 ## (at your option) any later version. |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
10 ## |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
11 ## 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
|
12 ## WITHOUT ANY WARRANTY; without even the implied warranty of |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
13 ## MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
14 ## GNU General Public License for more details. |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
15 ## |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
16 ## 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
|
17 ## along with Octave; see the file COPYING. If not, see |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
18 ## <https://www.gnu.org/licenses/>. |
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 ## -*- texinfo -*- |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
21 ## @deftypefn {} {} betainc (@var{x}, @var{a}, @var{b}) |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
22 ## @deftypefnx {} {} betainc (@var{x}, @var{a}, @var{b}, @var{tail}) |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
23 ## Compute the incomplete beta function. |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
24 ## |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
25 ## This is defined as |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
26 ## @tex |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
27 ## $$ |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
28 ## 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
|
29 ## $$ |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
30 ## @end tex |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
31 ## @ifnottex |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
32 ## |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
33 ## @example |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
34 ## @group |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
35 ## x |
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 ## | |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
38 ## I_x (a, b) = | 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
|
39 ## | |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
40 ## / |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
41 ## 0 |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
42 ## @end group |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
43 ## @end example |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
44 ## |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
45 ## @end ifnottex |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
46 ## |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
47 ## with real @var{x} in the range [0,1]. The inputs @var{a} and @var{b} must |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
48 ## be real and strictly positive (> 0). If one of the inputs is not a scalar |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
49 ## then the other inputs must be scalar or of compatible dimensions. |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
50 ## |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
51 ## By default, @var{tail} is @qcode{"lower"} and the incomplete beta function |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
52 ## integrated from 0 to @var{x} is computed. If @var{tail} is @qcode{"upper"} |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
53 ## then the complementary function integrated from @var{x} to 1 is calculated. |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
54 ## The two choices are related by |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
55 ## |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
56 ## betainc (@var{x}, @var{a}, @var{b}, @qcode{"upper"}) = |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
57 ## 1 - betainc (@var{x}, @var{a}, @var{b}, @qcode{"lower"}). |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
58 ## |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
59 ## @code{betainc} uses a more sophisticated algorithm than subtraction to |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
60 ## get numerically accurate results when the @qcode{"lower"} value is small. |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
61 ## |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
62 ## Reference: @nospell{A. Cuyt, V. Brevik Petersen, B. Verdonk, H. Waadeland, |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
63 ## W.B. Jones}, @cite{Handbook of Continued Fractions for Special Functions}, |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
64 ## ch.@: 18. |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
65 ## |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
66 ## @seealso{beta, betaincinv, betaln} |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
67 ## @end deftypefn |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
68 |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
69 function y = betainc (x, 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
|
70 |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
71 if (nargin < 3 || nargin > 4) |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
72 print_usage (); |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
73 endif |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
74 |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
75 [err, x, a, b] = common_size (x, a, b); |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
76 if (err > 0) |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
77 error ("betainc: X, A, and B must be of common size or scalars"); |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
78 endif |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
79 |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
80 if (iscomplex (x) || iscomplex (a) || iscomplex (b)) |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
81 error ("betainc: all inputs must be real"); |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
82 endif |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
83 |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
84 ## 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:
24923
diff
changeset
|
85 orig_sz = size (x); |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
86 x = x(:); |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
87 a = a(:); |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
88 b = b(:); |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
89 |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
90 if (any ((x < 0) | (x > 1))) |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
91 error ("betainc: X 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
|
92 endif |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
93 |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
94 if (any (a <= 0)) |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
95 error ("betainc: A must be strictly positive"); |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
96 endif |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
97 |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
98 if (any (b <= 0)) |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
99 error ("betainc: 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
|
100 endif |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
101 |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
102 ## If any of the arguments is single then the output should be as well. |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
103 if (strcmp (class (x), "single") || strcmp (class (a), "single") |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
104 || strcmp (class (b), "single")) |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
105 a = single (a); |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
106 b = single (b); |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
107 x = single (x); |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
108 endif |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
109 |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
110 ## Convert to floating point if necessary |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
111 if (isinteger (x)) |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
112 y = double (x); |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
113 endif |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
114 if (isinteger (a)) |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
115 a = double (a); |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
116 endif |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
117 if (isinteger (b)) |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
118 b = double (b); |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
119 endif |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
120 |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
121 ## Initialize output array |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
122 y = zeros (size (x), class (x)); |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
123 |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
124 ## In the following, we use the fact that the continued fraction Octave uses |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
125 ## is more efficient when x <= a / (a + b). Moreover, to compute the upper |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
126 ## version, which is defined as I_x(a,b,"upper") = 1 - I_x(a,b) we use the |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
127 ## property I_x(a,b) + I_(1-x) (b,a) = 1. |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
128 |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
129 if (strcmpi (tail, "lower")) |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
130 fflag = (x > a./(a+b)); |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
131 x(fflag) = 1 - x(fflag); |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
132 [a(fflag), b(fflag)] = deal (b(fflag), a(fflag)); |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
133 elseif (strcmpi (tail, "upper")) |
24923
40ab8be7d7ec
Fixed style in specfun scripts
Michele Ginesi <michele.ginesi@gmail.com>
parents:
24921
diff
changeset
|
134 fflag = (x < (a ./ (a + b))); |
40ab8be7d7ec
Fixed style in specfun scripts
Michele Ginesi <michele.ginesi@gmail.com>
parents:
24921
diff
changeset
|
135 x(! fflag) = 1 - x(! fflag); |
40ab8be7d7ec
Fixed style in specfun scripts
Michele Ginesi <michele.ginesi@gmail.com>
parents:
24921
diff
changeset
|
136 [a(! fflag), b(! fflag)] = deal (b(! fflag), a(! fflag)); |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
137 else |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
138 error ("betainc: 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
|
139 endif |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
140 |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
141 f = zeros (size (x), class (x)); |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
142 |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
143 ## Continued fractions: CPVWJ, formula 18.5.20, modified Lentz algorithm |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
144 ## implemented in a separate .cc file. This particular continued fraction |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
145 ## gives (B(a,b) * I_x(a,b)) / (x^a * (1-x)^b). |
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:
24923
diff
changeset
|
147 f = __betainc__ (x, a, b); |
24919
ed6f6bbed604
betainc: vectorized the Lentz's algorithm
Michele Ginesi <michele.ginesi@gmail.com>
parents:
24907
diff
changeset
|
148 |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
149 ## Divide continued fraction by B(a,b) / (x^a * (1-x)^b) to obtain I_x(a,b). |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
150 y = a .* log (x) + b .* log1p (-x) + gammaln (a + b) - ... |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
151 gammaln (a) - gammaln (b) + log (f); |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
152 y = real (exp (y)); |
24923
40ab8be7d7ec
Fixed style in specfun scripts
Michele Ginesi <michele.ginesi@gmail.com>
parents:
24921
diff
changeset
|
153 y(fflag) = 1 - y(fflag); |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
154 |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
155 ## Restore original shape |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
156 y = reshape (y, orig_sz); |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
157 |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
158 endfunction |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
159 |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
160 |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
161 ## Double precision |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
162 %!test |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
163 %! a = [1, 1.5, 2, 3]; |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
164 %! b = [4, 3, 2, 1]; |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
165 %! v1 = betainc (1, a, b); |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
166 %! v2 = [1,1,1,1]; |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
167 %! x = [.2, .4, .6, .8]; |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
168 %! v3 = betainc (x, a, b); |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
169 %! v4 = 1 - betainc (1-x, b, a); |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
170 %! assert (v1, v2, sqrt (eps)); |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
171 %! assert (v3, v4, sqrt (eps)); |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
172 |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
173 ## Single precision |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
174 %!test |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
175 %! a = single ([1, 1.5, 2, 3]); |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
176 %! b = single ([4, 3, 2, 1]); |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
177 %! v1 = betainc (1, a, b); |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
178 %! v2 = single ([1,1,1,1]); |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
179 %! x = single ([.2, .4, .6, .8]); |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
180 %! v3 = betainc (x, a, b); |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
181 %! v4 = 1 - betainc (1-x, b, a); |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
182 %! assert (v1, v2, sqrt (eps ("single"))); |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
183 %! assert (v3, v4, sqrt (eps ("single"))); |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
184 |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
185 ## Mixed double/single precision |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
186 %!test |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
187 %! a = single ([1, 1.5, 2, 3]); |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
188 %! b = [4, 3, 2, 1]; |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
189 %! v1 = betainc (1,a,b); |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
190 %! v2 = single ([1,1,1,1]); |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
191 %! x = [.2, .4, .6, .8]; |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
192 %! v3 = betainc (x, a, b); |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
193 %! v4 = 1-betainc (1.-x, b, a); |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
194 %! assert (v1, v2, sqrt (eps ("single"))); |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
195 %! assert (v3, v4, sqrt (eps ("single"))); |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
196 |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
197 %!test <51157> |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
198 %! y = betainc ([0.00780;0.00782;0.00784],250.005,49750.995); |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
199 %! y_ex = [0.999999999999989; 0.999999999999992; 0.999999999999995]; |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
200 %! assert (y, y_ex, -1e-14); |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
201 |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
202 %!assert (betainc (0.001, 20, 30), 2.750687665855991e-47, -3e-14); |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
203 %!assert (betainc (0.0001, 20, 30), 2.819953178893307e-67, -3e-14); |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
204 %!assert (betainc (0.99, 20, 30, "upper"), 1.5671643161872703e-47, -3e-14); |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
205 %!assert (betainc (0.999, 20, 30, "upper"), 1.850806276141535e-77, -3e-14); |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
206 %!assert (betainc (0.5, 200, 300), 0.9999964565197356, -1e-15); |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
207 %!assert (betainc (0.5, 200, 300, "upper"), 3.54348026439253e-06, -1e-13); |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
208 |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
209 ## Test trivial values |
24907
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
210 %!test |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
211 %! [a,b] = ndgrid (linspace (1e-4, 100, 20), linspace (1e-4, 100, 20)); |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
212 %! assert (betainc (0, a, b), zeros (20)); |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
213 %! assert (betainc (1, a, b), ones (20)); |
24907
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 ## Test input validation |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
216 %!error betainc () |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
217 %!error betainc (1) |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
218 %!error betainc (1,2) |
bd89440407aa
Incomplete beta function moved to a .m file, fixing accuracy and
Michele Ginesi <michele.ginesi@gmail.com>
parents:
diff
changeset
|
219 %!error betainc (1,2,3,4,5) |
24927
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
220 %!error <must be of common size or scalars> betainc (ones (2,2), ones (1,2), 1) |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
221 %!error <all inputs must be real> betainc (0.5i, 1, 2) |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
222 %!error <all inputs must be real> betainc (0, 1i, 1) |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
223 %!error <all inputs must be real> betainc (0, 1, 1i) |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
224 %!error <X must be in the range \[0, 1\]> betainc (-0.1,1,1) |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
225 %!error <X must be in the range \[0, 1\]> betainc (1.1,1,1) |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
226 %!error <X must be in the range \[0, 1\]> |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
227 %! x = ones (1, 1, 2); |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
228 %! x(1,1,2) = -1; |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
229 %! betainc (x,1,1); |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
230 %!error <A must be strictly positive> betainc (0.5,0,1) |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
231 %!error <A must be strictly positive> |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
232 %! a = ones (1, 1, 2); |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
233 %! a(1,1,2) = 0; |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
234 %! betainc (1,a,1); |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
235 %!error <B must be strictly positive> betainc (0.5,1,0) |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
236 %!error <B must be strictly positive> |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
237 %! b = ones (1, 1, 2); |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
238 %! b(1,1,2) = 0; |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
239 %! betainc (1,1,b); |
c280560d9c96
Overhaul special functions modified by GSOC2018 project.
Rik <rik@octave.org>
parents:
24923
diff
changeset
|
240 %!error <invalid value for TAIL> betainc (1,2,3, "foobar") |