Mercurial > octave
view scripts/specfun/gammaincinv.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 | 662faf9de127 |
children | 6652d3823428 |
line wrap: on
line source
## Copyright (C) 2017 Michele Ginesi ## ## This file is part of Octave. ## ## Octave is free software: you can redistribute it and/or modify it ## under the terms of the GNU General Public License as published by ## the Free Software Foundation; either version 3 of the License, or ## (at your option) any later version. ## ## Octave is distributed in the hope that it will be useful, but ## WITHOUT ANY WARRANTY; without even the implied warranty of ## MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the ## GNU General Public License for more details. ## ## You should have received a copy of the GNU General Public License ## along with Octave; see the file COPYING. If not, see ## <https://www.gnu.org/licenses/>. ## Author: Michele Ginesi <michele.ginesi@gmail.com> ## -*- texinfo -*- ## @deftypefn {} {} gammaincinv (@var{y}, @var{a}) ## @deftypefnx {} {} gammaincinv (@var{y}, @var{a}, @var{tail}) ## Compute the inverse of the normalized incomplete gamma function. ## ## The normalized incomplete gamma function is defined as ## @tex ## $$ ## \gamma (x, a) = {1 \over {\Gamma (a)}}\displaystyle{\int_0^x t^{a-1} e^{-t} dt} ## $$ ## @end tex ## @ifnottex ## ## @example ## @group ## x ## 1 / ## gammainc (x, a) = --------- | exp (-t) t^(a-1) dt ## gamma (a) / ## t=0 ## @end group ## @end example ## ## @end ifnottex ## ## and @code{gammaincinv (gammainc (@var{x}, @var{a}), @var{a}) = @var{x}} ## for each non-negative value of @var{x}. If @var{a} is scalar then ## @code{gammaincinv (@var{y}, @var{a})} is returned for each element of ## @var{y} and vice versa. ## ## If neither @var{y} nor @var{a} is scalar then the sizes of @var{y} and ## @var{a} must agree, and @code{gammaincinv} is applied element-by-element. ## The variable @var{y} must be in the interval @math{[0,1]} while @var{a} must ## be real and positive. ## ## By default, @var{tail} is @qcode{"lower"} and the inverse of the incomplete ## gamma function integrated from 0 to @var{x} is computed. If @var{tail} is ## @qcode{"upper"}, then the complementary function integrated from @var{x} to ## infinity is inverted. ## ## The function is computed with Newton's method by solving ## @tex ## $$ ## y - \gamma (x, a) = 0 ## $$ ## @end tex ## @ifnottex ## ## @example ## @var{y} - gammainc (@var{x}, @var{a}) = 0 ## @end example ## ## @end ifnottex ## ## Reference: @nospell{A. Gil, J. Segura, and N. M. Temme}, @cite{Efficient and ## accurate algorithms for the computation and inversion of the incomplete ## gamma function ratios}, @nospell{SIAM J. Sci. Computing}, pp. A2965--A2981, ## Vol 34, 2012. ## ## @seealso{gammainc, gamma, gammaln} ## @end deftypefn function x = gammaincinv (y, a, tail = "lower") if (nargin < 2 || nargin > 3) print_usage (); endif [err, y, a] = common_size (y, a); if (err > 0) error ("gammaincinv: Y and A must be of common size or scalars"); endif if (iscomplex (y) || iscomplex (a)) error ("gammaincinv: all inputs must be real"); endif ## Remember original shape of data, but convert to column vector for calcs. orig_sz = size (y); y = y(:); a = a(:); if (any ((y < 0) | (y > 1))) error ("gammaincinv: Y must be in the range [0, 1]"); endif if (any (a <= 0)) error ("gammaincinv: A must be strictly positive"); endif ## If any of the arguments is single then the output should be as well. if (strcmp (class (y), "single") || strcmp (class (a), "single")) y = single (y); a = single (a); endif ## Convert to floating point if necessary if (isinteger (y)) y = double (y); endif if (isinteger (a)) a = double (a); endif ## Initialize output array x = zeros (size (y), class (y)); maxit = 20; tol = eps (class (y)); ## Special cases, a = 1 or y = 0, 1. if (strcmpi (tail, "lower")) x(a == 1) = - log (1 - y(a == 1)); x(y == 0) = 0; x(y == 1) = Inf; p = y; q = 1 - p; elseif (strcmpi (tail, "upper")) x(a == 1) = - log (y(a == 1)); x(y == 0) = Inf; x(y == 1) = 0; q = y; p = 1 - q; else error ("gammaincinv: invalid value for TAIL") endif todo = (a != 1) & (y != 0) & (y != 1); ## Case 1: p small. i_flag_1 = todo & (p < ((0.2 * (1 + a)) .^ a) ./ gamma (1 + a)); aa = a(i_flag_1); pp = p(i_flag_1); ## Initial guess. r = (pp .* gamma (1 + aa)) .^ (1 ./ aa); c2 = 1 ./ (aa + 1); c3 = (3 * aa + 5) ./ (2 * (aa + 1) .^2 .* (aa + 2)); c4 = (8 * aa .^ 2 + 33 * aa + 31) ./ (3 * (aa + 1) .^ 3 .* (aa + 2) .* ... (aa + 3)); c5 = (125 * aa .^ 4 + 1179 * aa .^ 3 + 3971 * aa.^2 + 5661 * aa + 2888) ... ./ (24 * (1 + aa) .^4 .* (aa + 2) .^ 2 .* (aa + 3) .* (aa + 4)); ## FIXME: Would polyval() be better here for more accuracy? x0 = r + c2 .* r .^ 2 + c3 .* r .^ 3 + c4 .* r .^4 + c5 .* r .^ 5; ## For this case we invert the lower version. F = @(p, a, x) p - gammainc (x, a, "lower"); JF = @(a, x) - exp (-gammaln (a) - x + (a - 1) .* log (x)); x(i_flag_1) = newton_method (F, JF, pp, aa, x0, tol, maxit); todo(i_flag_1) = false; ## Case 2: q small. i_flag_2 = (q < exp (- 0.5 * a) ./ gamma (1 + a)) & (a > 0) & (a < 10); i_flag_2 &= todo; aa = a(i_flag_2); qq = q(i_flag_2); ## Initial guess. x0 = (-log (qq) - gammaln (aa)); ## For this case, we invert the upper version. F = @(q, a, x) q - gammainc (x, a, "upper"); JF = @(a, x) exp (- gammaln (a) - x) .* x .^ (a - 1); x(i_flag_2) = newton_method (F, JF, qq, aa, x0, tol, maxit); todo(i_flag_2) = false; ## Case 3: a small. i_flag_3 = todo & ((a > 0) & (a < 1)); aa = a(i_flag_3); pp = p(i_flag_3); ## Initial guess xl = (pp .* gamma (aa + 1)) .^ (1 ./ aa); x0 = xl; ## For this case, we invert the lower version. F = @(p, a, x) p - gammainc (x, a, "lower"); JF = @(a, x) - exp (-gammaln (a) - x) .* x .^ (a - 1); x(i_flag_3) = newton_method (F, JF, pp, aa, x0, tol, maxit); todo(i_flag_3) = false; ## Case 4: a large. i_flag_4 = todo; aa = a(i_flag_4); qq = q(i_flag_4); ## Initial guess d = 1 ./ (9 * aa); t = 1 - d + sqrt (2) * erfcinv (2 * qq) .* sqrt (d); x0 = aa .* (t .^ 3); ## For this case, we invert the upper version. F = @(q, a, x) q - gammainc (x, a, "upper"); JF = @(a, x) exp (- gammaln (a) - x + (a - 1) .* log (x)); x(i_flag_4) = newton_method (F, JF, qq, aa, x0, tol, maxit); ## Restore original shape x = reshape (x, orig_sz); endfunction ## Subfunction: Newton's Method function x = newton_method (F, JF, y, a, x0, tol, maxit); l = numel (y); res = -F (y, a, x0) ./ JF (a, x0); todo = (abs (res) >= tol * abs (x0)); x = x0; it = 0; while (any (todo) && (it++ < maxit)) x(todo) += res(todo); res(todo) = -F (y(todo), a(todo), x(todo)) ./ JF (a(todo), x(todo)); todo = (abs (res) >= tol * abs (x)); endwhile x += res; endfunction %!test %! x = [1e-10, 1e-09, 1e-08, 1e-07]; %! a = [2, 3, 4]; %! [x, a] = ndgrid (x, a); %! xx = gammainc (gammaincinv (x, a), a); %! assert (xx, x, -3e-14); %!test %! x = [1e-10, 1e-09, 1e-08, 1e-07]; %! a = [2, 3, 4]; %! [x, a] = ndgrid (x, a); %! xx = gammainc (gammaincinv (x, a, "upper"), a, "upper"); %! assert (xx, x, -3e-14); %!test %! x = linspace (0, 1)'; %! a = [linspace(0.1, 1, 10), 2:5]; %! [x, a] = ndgrid (x, a); %! xx = gammainc (gammaincinv (x, a), a); %! assert (xx, x, -1e-13); %!test %! x = linspace (0, 1)'; %! a = [linspace(0.1, 1, 10), 2:5]; %! [x, a] = ndgrid (x, a); %! xx = gammainc (gammaincinv (x, a, "upper"), a, "upper"); %! assert (xx, x, -1e-13); ## Test the conservation of the input class %!assert (class (gammaincinv (0.5, 1)), "double") %!assert (class (gammaincinv (single (0.5), 1)), "single") %!assert (class (gammaincinv (0.5, single (1))), "single") %!assert (class (gammaincinv (int8 (0), 1)), "double") %!assert (class (gammaincinv (0.5, int8 (1))), "double") %!assert (class (gammaincinv (int8 (0), single (1))), "single") %!assert (class (gammaincinv (single (0.5), int8 (1))), "single") ## Test input validation %!error gammaincinv () %!error gammaincinv (1) %!error gammaincinv (1, 2, 3, 4) %!error <must be of common size or scalars> %! gammaincinv (ones (2,2), ones (1,2), 1); %!error <all inputs must be real> gammaincinv (0.5i, 1) %!error <all inputs must be real> gammaincinv (0, 1i) %!error <Y must be in the range \[0, 1\]> gammaincinv (-0.1,1) %!error <Y must be in the range \[0, 1\]> gammaincinv (1.1,1) %!error <Y must be in the range \[0, 1\]> %! y = ones (1, 1, 2); %! y(1,1,2) = -1; %! gammaincinv (y,1); %!error <A must be strictly positive> gammaincinv (0.5, 0) %!error <A must be strictly positive> %! a = ones (1, 1, 2); %! a(1,1,2) = 0; %! gammaincinv (1,a,1); %!error <invalid value for TAIL> gammaincinv (1,2, "foobar")