view scripts/specfun/betainc.m @ 30875:5d3faba0342e

doc: Ensure documentation lists output argument when it exists for all m-files. For new users of Octave it is best to show explicit calling forms in the documentation and to show a return argument when it exists. * bp-table.cc, shift.m, accumarray.m, accumdim.m, bincoeff.m, bitcmp.m, bitget.m, bitset.m, blkdiag.m, celldisp.m, cplxpair.m, dblquad.m, flip.m, fliplr.m, flipud.m, idivide.m, int2str.m, interpft.m, logspace.m, num2str.m, polyarea.m, postpad.m, prepad.m, randi.m, repmat.m, rng.m, rot90.m, rotdim.m, structfun.m, triplequad.m, uibuttongroup.m, uicontrol.m, uipanel.m, uipushtool.m, uitoggletool.m, uitoolbar.m, waitforbuttonpress.m, help.m, __additional_help_message__.m, hsv.m, im2double.m, im2frame.m, javachk.m, usejava.m, argnames.m, char.m, formula.m, inline.m, __vectorize__.m, findstr.m, flipdim.m, strmatch.m, vectorize.m, commutation_matrix.m, cond.m, cross.m, duplication_matrix.m, expm.m, orth.m, rank.m, rref.m, trace.m, vech.m, cast.m, compare_versions.m, delete.m, dir.m, fileattrib.m, grabcode.m, gunzip.m, inputname.m, license.m, list_primes.m, ls.m, mexext.m, movefile.m, namelengthmax.m, nargoutchk.m, nthargout.m, substruct.m, swapbytes.m, ver.m, verLessThan.m, what.m, fminunc.m, fsolve.m, fzero.m, optimget.m, __fdjac__.m, matlabroot.m, savepath.m, campos.m, camroll.m, camtarget.m, camup.m, camva.m, camzoom.m, clabel.m, diffuse.m, legend.m, orient.m, rticks.m, specular.m, thetaticks.m, xlim.m, xtickangle.m, xticklabels.m, xticks.m, ylim.m, ytickangle.m, yticklabels.m, yticks.m, zlim.m, ztickangle.m, zticklabels.m, zticks.m, ellipsoid.m, isocolors.m, isonormals.m, stairs.m, surfnorm.m, __actual_axis_position__.m, __pltopt__.m, close.m, graphics_toolkit.m, pan.m, print.m, printd.m, __ghostscript__.m, __gnuplot_print__.m, __opengl_print__.m, rotate3d.m, subplot.m, zoom.m, compan.m, conv.m, poly.m, polyaffine.m, polyder.m, polyint.m, polyout.m, polyreduce.m, polyvalm.m, roots.m, prefdir.m, prefsfile.m, profexplore.m, profexport.m, profshow.m, powerset.m, unique.m, arch_rnd.m, arma_rnd.m, autoreg_matrix.m, bartlett.m, blackman.m, detrend.m, durbinlevinson.m, fftconv.m, fftfilt.m, fftshift.m, fractdiff.m, hamming.m, hanning.m, hurst.m, ifftshift.m, rectangle_lw.m, rectangle_sw.m, triangle_lw.m, sinc.m, sinetone.m, sinewave.m, spectral_adf.m, spectral_xdf.m, spencer.m, ilu.m, __sprand__.m, sprand.m, sprandn.m, sprandsym.m, treelayout.m, beta.m, betainc.m, betaincinv.m, betaln.m, cosint.m, expint.m, factorial.m, gammainc.m, gammaincinv.m, lcm.m, nthroot.m, perms.m, reallog.m, realpow.m, realsqrt.m, sinint.m, hadamard.m, hankel.m, hilb.m, invhilb.m, magic.m, pascal.m, rosser.m, toeplitz.m, vander.m, wilkinson.m, center.m, corr.m, cov.m, discrete_cdf.m, discrete_inv.m, discrete_pdf.m, discrete_rnd.m, empirical_cdf.m, empirical_inv.m, empirical_pdf.m, empirical_rnd.m, kendall.m, kurtosis.m, mad.m, mean.m, meansq.m, median.m, mode.m, moment.m, range.m, ranks.m, run_count.m, skewness.m, spearman.m, statistics.m, std.m, base2dec.m, bin2dec.m, blanks.m, cstrcat.m, deblank.m, dec2base.m, dec2bin.m, dec2hex.m, hex2dec.m, index.m, regexptranslate.m, rindex.m, strcat.m, strjust.m, strtrim.m, strtrunc.m, substr.m, untabify.m, __have_feature__.m, __prog_output_assert__.m, __run_test_suite__.m, example.m, fail.m, asctime.m, calendar.m, ctime.m, date.m, etime.m: Add return arguments to @deftypefn macros where they were missing. Rename variables in functions (particularly generic "retval") to match documentation. Rename some return variables for (hopefully) better clarity (e.g., 'ax' to 'hax' to indicate it is a graphics handle to an axes object).
author Rik <rik@octave.org>
date Wed, 30 Mar 2022 20:40:27 -0700
parents 796f54d4ddbf
children 74f8854958d2
line wrap: on
line source

########################################################################
##
## Copyright (C) 2018-2022 The Octave Project Developers
##
## See the file COPYRIGHT.md in the top-level directory of this
## distribution or <https://octave.org/copyright/>.
##
## 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/>.
##
########################################################################

## -*- texinfo -*-
## @deftypefn  {} {@var{I} =} betainc (@var{x}, @var{a}, @var{b})
## @deftypefnx {} {@var{I} =} betainc (@var{x}, @var{a}, @var{b}, @var{tail})
## Compute the incomplete beta function.
##
## This is defined as
## @tex
## $$
## I_x (a, b) = {1 \over {B(a,b)}} \displaystyle{\int_0^x t^{a-1} (1-t)^{b-1} dt}
## $$
## @end tex
## @ifnottex
##
## @example
## @group
##                           x
##                          /
##                  1       |
## I_x (a, b) = ----------  | t^(a-1) (1-t)^(b-1) dt
##              beta (a,b)  |
##                          /
##                         0
## @end group
## @end example
##
## @end ifnottex
##
## with real @var{x} in the range [0,1].  The inputs @var{a} and @var{b} must
## be real and strictly positive (> 0).  If one of the inputs is not a scalar
## then the other inputs must be scalar or of compatible dimensions.
##
## By default, @var{tail} is @qcode{"lower"} and the incomplete beta function
## integrated from 0 to @var{x} is computed.  If @var{tail} is @qcode{"upper"}
## then the complementary function integrated from @var{x} to 1 is calculated.
## The two choices are related by
##
## betainc (@var{x}, @var{a}, @var{b}, @qcode{"upper"}) =
## 1 - betainc (@var{x}, @var{a}, @var{b}, @qcode{"lower"}).
##
## @code{betainc} uses a more sophisticated algorithm than subtraction to
## get numerically accurate results when the @qcode{"lower"} value is small.
##
## Reference: @nospell{A. Cuyt, V. Brevik Petersen, B. Verdonk, H. Waadeland,
## W.B. Jones}, @cite{Handbook of Continued Fractions for Special Functions},
## ch.@: 18.
##
## @seealso{beta, betaincinv, betaln}
## @end deftypefn

function I = betainc (x, a, b, tail = "lower")

  if (nargin < 3)
    print_usage ();
  endif

  [err, x, a, b] = common_size (x, a, b);
  if (err > 0)
    error ("betainc: X, A, and B must be of common size or scalars");
  endif

  if (iscomplex (x) || iscomplex (a) || iscomplex (b))
    error ("betainc: all inputs must be real");
  endif

  ## Remember original shape of data, but convert to column vector for calcs.
  orig_sz = size (x);
  x = x(:);
  a = a(:);
  b = b(:);

  if (any ((x < 0) | (x > 1)))
    error ("betainc: X must be in the range [0, 1]");
  endif

  if (any (a <= 0))
    error ("betainc: A must be strictly positive");
  endif

  if (any (b <= 0))
    error ("betainc: B must be strictly positive");
  endif

  ## If any of the arguments is single then the output should be as well.
  if (strcmp (class (x), "single") || strcmp (class (a), "single")
      || strcmp (class (b), "single"))
    a = single (a);
    b = single (b);
    x = single (x);
  endif

  ## Convert to floating point if necessary
  if (isinteger (x))
    I = double (x);
  endif
  if (isinteger (a))
    a = double (a);
  endif
  if (isinteger (b))
    b = double (b);
  endif

  ## Initialize output array
  I = zeros (size (x), class (x));

  ## Trivial cases (long code here trades memory for speed)
  a_one = (a == 1);
  b_one = (b == 1);
  a_b_one = a_one & b_one;
  a_not_one = ! a_one;
  b_not_one = ! b_one;
  non_trivial = a_not_one & b_not_one;
  a_one &= b_not_one;
  b_one &= a_not_one;

  if (strcmpi (tail, "lower"))
    I(a_b_one) = x(a_b_one);
    I(a_one) = 1 - (1 - x(a_one)) .^ b(a_one);
    I(b_one) = x(b_one) .^ a(b_one);
  elseif (strcmpi (tail, "upper"))
    I(a_b_one) = 1 - x(a_b_one);
    I(a_one) = (1 - x(a_one)) .^ b(a_one);
    I(b_one) = 1 - x(b_one) .^ a(b_one);
  endif

  ## Non-Trivial cases
  ## In the following, we use the fact that the continued fraction Octave uses
  ## is more efficient when x <= a / (a + b).  Moreover, to compute the upper
  ## version, which is defined as I_x(a,b,"upper") = 1 - I_x(a,b) we use the
  ## property I_x(a,b) + I_(1-x) (b,a) = 1.

  x = x(non_trivial);
  a = a(non_trivial);
  b = b(non_trivial);

  if (strcmpi (tail, "lower"))
    fflag = (x > a./(a + b));
    x(fflag) = 1 - x(fflag);
    [a(fflag), b(fflag)] = deal (b(fflag), a(fflag));
  elseif (strcmpi (tail, "upper"))
    fflag = (x < (a ./ (a + b)));
    x(! fflag) = 1 - x(! fflag);
    [a(! fflag), b(! fflag)] = deal (b(! fflag), a(! fflag));
  else
    error ("betainc: invalid value for TAIL");
  endif

  f = zeros (size (x), class (x));

  ## Continued fractions: CPVWJ, formula 18.5.20, modified Lentz algorithm
  ## implemented in a separate .cc file.  This particular continued fraction
  ## gives (B(a,b) * I_x(a,b)) / (x^a * (1-x)^b).

  f = __betainc__ (x, a, b);

  ## Divide continued fraction by B(a,b) / (x^a * (1-x)^b) to obtain I_x(a,b).
  y_nt = a .* log (x) + b .* log1p (-x) ...
         + (gammaln (a + b) - gammaln (a) - gammaln (b)) + log (f);
  y_nt = real (exp (y_nt));
  y_nt(fflag) = 1 - y_nt(fflag);

  I(non_trivial) = y_nt;

  ## Restore original shape
  I = reshape (I, orig_sz);

endfunction


## Double precision
%!test
%! a = [1, 1.5, 2, 3];
%! b = [4, 3, 2, 1];
%! v1 = betainc (1, a, b);
%! v2 = [1,1,1,1];
%! x = [.2, .4, .6, .8];
%! v3 = betainc (x, a, b);
%! v4 = 1 - betainc (1-x, b, a);
%! assert (v1, v2, sqrt (eps));
%! assert (v3, v4, sqrt (eps));

## Single precision
%!test
%! a = single ([1, 1.5, 2, 3]);
%! b = single ([4, 3, 2, 1]);
%! v1 = betainc (1, a, b);
%! v2 = single ([1,1,1,1]);
%! x = single ([.2, .4, .6, .8]);
%! v3 = betainc (x, a, b);
%! v4 = 1 - betainc (1-x, b, a);
%! assert (v1, v2, sqrt (eps ("single")));
%! assert (v3, v4, sqrt (eps ("single")));

## Mixed double/single precision
%!test
%! a = single ([1, 1.5, 2, 3]);
%! b = [4, 3, 2, 1];
%! v1 = betainc (1,a,b);
%! v2 = single ([1,1,1,1]);
%! x = [.2, .4, .6, .8];
%! v3 = betainc (x, a, b);
%! v4 = 1 - betainc (1. - x, b, a);
%! assert (v1, v2, sqrt (eps ("single")));
%! assert (v3, v4, sqrt (eps ("single")));

%!test <*51157>
%! y = betainc ([0.00780;0.00782;0.00784],250.005,49750.995);
%! y_ex = [0.999999999999989; 0.999999999999992; 0.999999999999995];
%! assert (y, y_ex, -1e-14);

%!assert (betainc (0.001, 20, 30), 2.750687665855991e-47, -3e-14)
%!assert (betainc (0.0001, 20, 30), 2.819953178893307e-67, -7e-14)
%!assert <*54383> (betainc (0.99, 20, 30, "upper"),
%!                 1.5671643161872703e-47, -7e-14)
%!assert (betainc (0.999, 20, 30, "upper"), 1.850806276141535e-77, -7e-14)
%!assert (betainc (0.5, 200, 300), 0.9999964565197356, -1e-15)
%!assert (betainc (0.5, 200, 300, "upper"), 3.54348026439253e-06, -3e-13)

## Test trivial values
%!test
%! [a,b] = ndgrid (linspace (1e-4, 100, 20), linspace (1e-4, 100, 20));
%! assert (betainc (0, a, b), zeros (20));
%! assert (betainc (1, a, b), ones (20));

%!test <*34405>
%! assert (betainc (NaN, 1, 2), NaN);
%! assert (betainc (0.5, 1, Inf), 1);

## Test input validation
%!error <Invalid call> betainc ()
%!error <Invalid call> betainc (1)
%!error <Invalid call> betainc (1,2)
%!error <must be of common size or scalars> betainc (ones (2,2), ones (1,2), 1)
%!error <all inputs must be real> betainc (0.5i, 1, 2)
%!error <all inputs must be real> betainc (0, 1i, 1)
%!error <all inputs must be real> betainc (0, 1, 1i)
%!error <X must be in the range \[0, 1\]> betainc (-0.1,1,1)
%!error <X must be in the range \[0, 1\]> betainc (1.1,1,1)
%!error <X must be in the range \[0, 1\]>
%! x = ones (1, 1, 2);
%! x(1,1,2) = -1;
%! betainc (x,1,1);
%!error <A must be strictly positive> betainc (0.5,0,1)
%!error <A must be strictly positive>
%! a = ones (1, 1, 2);
%! a(1,1,2) = 0;
%! betainc (1,a,1);
%!error <B must be strictly positive> betainc (0.5,1,0)
%!error <B must be strictly positive>
%! b = ones (1, 1, 2);
%! b(1,1,2) = 0;
%! betainc (1,1,b);
%!error <invalid value for TAIL> betainc (1,2,3, "foobar")