Mercurial > octave
view libinterp/corefcn/psi.cc @ 22197:e43d83253e28
refill multi-line macro definitions
Use the Emacs C++ mode style for line continuation markers in
multi-line macro definitions.
* make_int.cc, __dsearchn__.cc, __magick_read__.cc, besselj.cc,
bitfcns.cc, bsxfun.cc, cellfun.cc, data.cc, defun-dld.h, defun-int.h,
defun.h, det.cc, error.h, find.cc, gcd.cc, graphics.cc, interpreter.h,
jit-ir.h, jit-typeinfo.h, lookup.cc, ls-mat5.cc, max.cc, mexproto.h,
mxarray.in.h, oct-stream.cc, ordschur.cc, pr-output.cc, profiler.h,
psi.cc, regexp.cc, sparse-xdiv.cc, sparse-xpow.cc, tril.cc, txt-eng.h,
utils.cc, variables.cc, variables.h, xdiv.cc, xpow.cc, __glpk__.cc,
ov-base.cc, ov-base.h, ov-cell.cc, ov-ch-mat.cc, ov-classdef.cc,
ov-complex.cc, ov-cx-mat.cc, ov-cx-sparse.cc, ov-float.cc, ov-float.h,
ov-flt-complex.cc, ov-flt-cx-mat.cc, ov-flt-re-mat.cc,
ov-int-traits.h, ov-lazy-idx.h, ov-perm.cc, ov-re-mat.cc,
ov-re-sparse.cc, ov-scalar.cc, ov-scalar.h, ov-str-mat.cc,
ov-type-conv.h, ov.cc, ov.h, op-class.cc, op-int-conv.cc, op-int.h,
op-str-str.cc, ops.h, lex.ll, Array.cc, CMatrix.cc, CSparse.cc,
MArray.cc, MArray.h, MDiagArray2.cc, MDiagArray2.h, MSparse.h,
Sparse.cc, dMatrix.cc, dSparse.cc, fCMatrix.cc, fMatrix.cc,
idx-vector.cc, f77-fcn.h, quit.h, bsxfun-decl.h, bsxfun-defs.cc,
lo-specfun.cc, oct-convn.cc, oct-convn.h, oct-norm.cc, oct-norm.h,
oct-rand.cc, Sparse-op-decls.h, Sparse-op-defs.h, mx-inlines.cc,
mx-op-decl.h, mx-op-defs.h, mach-info.cc, oct-group.cc, oct-passwd.cc,
oct-syscalls.cc, oct-time.cc, data-conv.cc, kpse.cc, lo-ieee.h,
lo-macros.h, oct-cmplx.h, oct-glob.cc, oct-inttypes.cc,
oct-inttypes.h, oct-locbuf.h, oct-sparse.h, url-transfer.cc,
oct-conf-post.in.h, shared-fcns.h: Refill macro definitions.
author | John W. Eaton <jwe@octave.org> |
---|---|
date | Mon, 01 Aug 2016 12:40:18 -0400 |
parents | 112b20240c87 |
children | bac0d6f07a3e |
line wrap: on
line source
/* Copyright (C) 2015 Carnë Draug 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 <http://www.gnu.org/licenses/>. */ #if defined (HAVE_CONFIG_H) # include "config.h" #endif #include "ov.h" #include "defun.h" #include "error.h" #include "dNDArray.h" #include "fNDArray.h" #include "lo-specfun.h" DEFUN (psi, args, , doc: /* -*- texinfo -*- @deftypefn {} {} psi (@var{z}) @deftypefnx {} {} psi (@var{k}, @var{z}) Compute the psi (polygamma) function. The polygamma functions are the @var{k}th derivative of the logarithm of the gamma function. If unspecified, @var{k} defaults to zero. A value of zero computes the digamma function, a value of 1, the trigamma function, and so on. The digamma function is defined: @tex $$ \Psi (z) = {d (log (\Gamma (z))) \over dx} $$ @end tex @ifnottex @example @group psi (z) = d (log (gamma (z))) / dx @end group @end example @end ifnottex When computing the digamma function (when @var{k} equals zero), @var{z} can have any value real or complex value. However, for polygamma functions (@var{k} higher than 0), @var{z} must be real and non-negative. @seealso{gamma, gammainc, gammaln} @end deftypefn */) { int nargin = args.length (); if (nargin < 1 || nargin > 2) print_usage (); const octave_value oct_z = (nargin == 1) ? args(0) : args(1); const octave_idx_type k = (nargin == 1) ? 0 : args(0).idx_type_value ("psi: K must be an integer"); if (k < 0) error ("psi: K must be non-negative"); octave_value retval; if (k == 0) { #define FLOAT_BRANCH(T, A, M, E) \ if (oct_z.is_ ## T ##_type ()) \ { \ const A ## NDArray z = oct_z.M ## array_value (); \ A ## NDArray psi_z (z.dims ()); \ \ const E* zv = z.data (); \ E* psi_zv = psi_z.fortran_vec (); \ const octave_idx_type n = z.numel (); \ for (octave_idx_type i = 0; i < n; i++) \ *psi_zv++ = octave::math::psi (*zv++); \ \ retval = psi_z; \ } if (oct_z.is_complex_type ()) { FLOAT_BRANCH(double, Complex, complex_, Complex) else FLOAT_BRANCH(single, FloatComplex, float_complex_, FloatComplex) else error ("psi: Z must be a floating point"); } else { FLOAT_BRANCH(double, , , double) else FLOAT_BRANCH(single, Float, float_, float) else error ("psi: Z must be a floating point"); } #undef FLOAT_BRANCH } else { if (! oct_z.is_real_type ()) error ("psi: Z must be real value for polygamma (K > 0)"); #define FLOAT_BRANCH(T, A, M, E) \ if (oct_z.is_ ## T ##_type ()) \ { \ const A ## NDArray z = oct_z.M ## array_value (); \ A ## NDArray psi_z (z.dims ()); \ \ const E* zv = z.data (); \ E* psi_zv = psi_z.fortran_vec (); \ const octave_idx_type n = z.numel (); \ for (octave_idx_type i = 0; i < n; i++) \ { \ if (*zv < 0) \ error ("psi: Z must be non-negative for polygamma (K > 0)"); \ \ *psi_zv++ = octave::math::psi (k, *zv++); \ } \ retval = psi_z; \ } FLOAT_BRANCH(double, , , double) else FLOAT_BRANCH(single, Float, float_, float) else error ("psi: Z must be a floating point for polygamma (K > 0)"); #undef FLOAT_BRANCH } return retval; } /* %!shared em %! em = 0.577215664901532860606512090082402431042; # Euler-Mascheroni Constant %!assert (psi (ones (7, 3, 5)), repmat (-em, [7 3 5])) %!assert (psi ([0 1]), [-Inf -em]) %!assert (psi ([-20:1]), [repmat(-Inf, [1 21]) -em]) %!assert (psi (single ([0 1])), single ([-Inf -em])) ## Abramowitz and Stegun, page 258, eq 6.3.5 %!test %! z = [-100:-1 1:200] ./ 10; # drop the 0 %! assert (psi (z + 1), psi (z) + 1 ./ z, eps*1000); ## Abramowitz and Stegun, page 258, eq 6.3.2 %!assert (psi (1), -em) ## Abramowitz and Stegun, page 258, eq 6.3.3 %!assert (psi (1/2), -em - 2 * log (2)) ## The following tests are from Pascal Sebah and Xavier Gourdon (2002) ## "Introduction to the Gamma Function" ## Interesting identities of the digamma function, in section of 5.1.3 %!assert (psi (1/3), - em - (3/2) * log(3) - ((sqrt (3) / 6) * pi), eps*10) %!assert (psi (1/4), - em -3 * log (2) - pi/2, eps*10) %!assert (psi (1/6), - em -2 * log (2) - (3/2) * log (3) - ((sqrt (3) / 2) * pi), eps*10) ## First 6 zeros of the digamma function, in section of 5.1.5 (and also on ## Abramowitz and Stegun, page 258, eq 6.3.19) %!assert (psi ( 1.46163214496836234126265954232572132846819620400644), 0, eps) %!assert (psi (-0.504083008264455409258269304533302498955385182368579), 0, eps*2) %!assert (psi (-1.573498473162390458778286043690434612655040859116846), 0, eps*2) %!assert (psi (-2.610720868444144650001537715718724207951074010873480), 0, eps*10) %!assert (psi (-3.635293366436901097839181566946017713948423861193530), 0, eps*10) %!assert (psi (-4.653237761743142441714598151148207363719069416133868), 0, eps*100) ## Tests for complex values %!shared z %! z = [-100:-1 1:200] ./ 10; # drop the 0 ## Abramowitz and Stegun, page 259 eq 6.3.10 %!assert (real (psi (i*z)), real (psi (1 - i*z))) ## Abramowitz and Stegun, page 259 eq 6.3.11 %!assert (imag (psi (i*z)), 1/2 .* 1./z + 1/2 * pi * coth (pi * z), eps *10) ## Abramowitz and Stegun, page 259 eq 6.3.12 %!assert (imag (psi (1/2 + i*z)), 1/2 * pi * tanh (pi * z), eps*10) ## Abramowitz and Stegun, page 259 eq 6.3.13 %!assert (imag (psi (1 + i*z)), - 1./(2*z) + 1/2 * pi * coth (pi * z), eps*10) ## Abramowitz and Stegun, page 260 eq 6.4.5 %!test %! for z = 0:20 %! assert (psi (1, z + 0.5), %! 0.5 * (pi^2) - 4 * sum ((2*(1:z) -1) .^(-2)), %! eps*10); %! endfor ## Abramowitz and Stegun, page 260 eq 6.4.6 %!test %! z = 0.1:0.1:20; %! for n = 0:8 %! ## our precision goes down really quick when computing n is too high. %! assert (psi (n, z+1), %! psi (n, z) + ((-1)^n) * factorial (n) * (z.^(-n-1)), 0.1); %! endfor ## Test input validation %!error psi () %!error psi (1, 2, 3) %!error <Z must be> psi ("non numeric") %!error <conversion of 5.3 to int.* value failed> psi (5.3, 1) %!error <K must be non-negative> psi (-5, 1) %!error <Z must be non-negative for polygamma> psi (5, -1) %!error <Z must be a floating point> psi (5, uint8 (-1)) %!error <Z must be real value for polygamma> psi (5, 5i) */