Mercurial > octave
view libinterp/corefcn/psi.cc @ 21966:112b20240c87
move docstrings in C++ files out of C strings and into comments
* __contourc__.cc, __dispatch__.cc, __dsearchn__.cc, __ichol__.cc,
__ilu__.cc, __lin_interpn__.cc, __luinc__.cc, __magick_read__.cc,
__pchip_deriv__.cc, __qp__.cc, balance.cc, besselj.cc, betainc.cc,
bitfcns.cc, bsxfun.cc, cellfun.cc, colloc.cc, conv2.cc, daspk.cc,
dasrt.cc, dassl.cc, data.cc, debug.cc, defaults.cc, det.cc, dirfns.cc,
dlmread.cc, dot.cc, eig.cc, ellipj.cc, error.cc, fft.cc, fft2.cc,
fftn.cc, file-io.cc, filter.cc, find.cc, gammainc.cc, gcd.cc,
getgrent.cc, getpwent.cc, getrusage.cc, givens.cc, graphics.cc,
hash.cc, help.cc, hess.cc, hex2num.cc, input.cc, inv.cc, kron.cc,
load-path.cc, load-save.cc, lookup.cc, ls-oct-text.cc, lsode.cc,
lu.cc, mappers.cc, matrix_type.cc, max.cc, mgorth.cc, nproc.cc,
oct-hist.cc, octave-link.cc, ordschur.cc, pager.cc, pinv.cc,
pr-output.cc, profiler.cc, psi.cc, pt-jit.cc, quad.cc, quadcc.cc,
qz.cc, rand.cc, rcond.cc, regexp.cc, schur.cc, sighandlers.cc,
sparse.cc, spparms.cc, sqrtm.cc, str2double.cc, strfind.cc, strfns.cc,
sub2ind.cc, svd.cc, sylvester.cc, symtab.cc, syscalls.cc, sysdep.cc,
time.cc, toplev.cc, tril.cc, tsearch.cc, typecast.cc, urlwrite.cc,
utils.cc, variables.cc, __delaunayn__.cc, __eigs__.cc,
__fltk_uigetfile__.cc, __glpk__.cc, __init_fltk__.cc,
__init_gnuplot__.cc, __osmesa_print__.cc, __voronoi__.cc, amd.cc,
audiodevinfo.cc, audioread.cc, ccolamd.cc, chol.cc, colamd.cc,
convhulln.cc, dmperm.cc, fftw.cc, qr.cc, symbfact.cc, symrcm.cc,
ov-base.cc, ov-bool-mat.cc, ov-cell.cc, ov-class.cc, ov-classdef.cc,
ov-fcn-handle.cc, ov-fcn-inline.cc, ov-flt-re-mat.cc, ov-int16.cc,
ov-int32.cc, ov-int64.cc, ov-int8.cc, ov-java.cc, ov-null-mat.cc,
ov-oncleanup.cc, ov-range.cc, ov-re-mat.cc, ov-struct.cc,
ov-typeinfo.cc, ov-uint16.cc, ov-uint32.cc, ov-uint64.cc, ov-uint8.cc,
ov-usr-fcn.cc, ov.cc, octave.cc, pt-arg-list.cc, pt-binop.cc,
pt-eval.cc, pt-mat.cc, lex.ll, oct-parse.in.yy:
Docstrings are now comments instead of C strings.
* build-aux/mk-opts.pl: Emit docstrings as comments instead of C
strings.
* DASPK-opts.in, LSODE-opts.in: Don't quote " in docstring fragments.
* builtins.h: Include builtin-defun-decls.h unconditionally.
* defun.h (DEFUN, DEFUNX, DEFCONSTFUN): Simply emit declaration.
(DEFALIAS): Always expand to nothing.
* defun-dld.h: No special macro expansions for MAKE_BUILTINS.
(DEFUN_DLD): Use FORWARD_DECLARE_FUN.
(DEFUNX_DLD): Use FORWARD_DECLARE_FUNX.
* defun-int.h: No special macro expansions for MAKE_BUILTINS.
(FORWARD_DECLARE_FUN, FORWARD_DECLARE_FUNX): New macros.
(DEFINE_FUN_INSTALLER_FUN): If compiling an Octave source file, pass
"external-doc" to DEFINE_FUNX_INSTALLER_FUN.
(DEFUN_INTERNAL, DEFCONSTFUN_INTERNAL, DEFUNX_INTERNAL,
DEFALIAS_INTERNAL): Delete.
* common.mk (move_if_change_rule): New macro.
(simple_move_if_change_rule): Define using move_if_change_rule.
* find-defun-files.sh (DEFUN_PATTERN): Update. Don't transform file
name extension to ".df".
* libinterp/mk-pkg-add, gendoc.pl: Operate directly on source files.
* mkbuiltins: New argument, SRCDIR. Operate directly on source files.
* mkdefs: Delete.
* libinterp/module.mk (BUILT_SOURCES): Update list to contain only
files included in other source files.
(GENERATED_MAKE_BUILTINS_INCS, DEF_FILES): Delete.
(LIBINTERP_BUILT_DISTFILES): Include $(OPT_HANDLERS) here.
(LIBINTERP_BUILT_NODISTFILES): Not here. Remove $(ALL_DEF_FILES from
the list.
(libinterp_EXTRA_DIST): Remove mkdefs from the list.
(FOUND_DEFUN_FILES): Rename from SRC_DEF_FILES.
(DLDFCN_DEFUN_FILES): Rename from DLDFCN_DEF_FILES.
(SRC_DEFUN_FILES): Rename from SRC_DEF_FILES.
(ALL_DEFUN_FILES): Rename from ALL_DEF_FILES.
(%.df: %.cc): Delete pattern rule.
(libinterp/build-env-features.cc, libinterp/builtins.cc,
libinterp/dldfcn/PKG_ADD): Use mv instead of move-if-change.
(libinterp/builtins.cc, libinterp/builtin-defun-decls.h):
Update mkbuiltins command.
($(srcdir)/libinterp/DOCSTRINGS): Update gendoc.pl command.
* liboctave/module.mk (BUILT_SOURCES): Don't include
liboctave-build-info.cc in the list.
author | John W. Eaton <jwe@octave.org> |
---|---|
date | Tue, 21 Jun 2016 16:07:51 -0400 |
parents | cbd8cf0a8a5c |
children | e43d83253e28 |
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) */