GNU Octave  6.2.0
A high-level interpreted language, primarily intended for numerical computations, mostly compatible with Matlab
psi.cc
Go to the documentation of this file.
1 ////////////////////////////////////////////////////////////////////////
2 //
3 // Copyright (C) 2016-2021 The Octave Project Developers
4 //
5 // See the file COPYRIGHT.md in the top-level directory of this
6 // distribution or <https://octave.org/copyright/>.
7 //
8 // This file is part of Octave.
9 //
10 // Octave is free software: you can redistribute it and/or modify it
11 // under the terms of the GNU General Public License as published by
12 // the Free Software Foundation, either version 3 of the License, or
13 // (at your option) any later version.
14 //
15 // Octave is distributed in the hope that it will be useful, but
16 // WITHOUT ANY WARRANTY; without even the implied warranty of
17 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
18 // GNU General Public License for more details.
19 //
20 // You should have received a copy of the GNU General Public License
21 // along with Octave; see the file COPYING. If not, see
22 // <https://www.gnu.org/licenses/>.
23 //
24 ////////////////////////////////////////////////////////////////////////
25 
26 #if defined (HAVE_CONFIG_H)
27 # include "config.h"
28 #endif
29 
30 #include "ov.h"
31 #include "defun.h"
32 #include "error.h"
33 #include "dNDArray.h"
34 #include "fNDArray.h"
35 
36 #include "lo-specfun.h"
37 
38 DEFUN (psi, args, ,
39  doc: /* -*- texinfo -*-
40 @deftypefn {} {} psi (@var{z})
41 @deftypefnx {} {} psi (@var{k}, @var{z})
42 Compute the psi (polygamma) function.
43 
44 The polygamma functions are the @var{k}th derivative of the logarithm
45 of the gamma function. If unspecified, @var{k} defaults to zero. A value
46 of zero computes the digamma function, a value of 1, the trigamma function,
47 and so on.
48 
49 The digamma function is defined:
50 
51 @tex
52 $$
53 \Psi (z) = {d (log (\Gamma (z))) \over dx}
54 $$
55 @end tex
56 @ifnottex
57 
58 @example
59 @group
60 psi (z) = d (log (gamma (z))) / dx
61 @end group
62 @end example
63 
64 @end ifnottex
65 
66 When computing the digamma function (when @var{k} equals zero), @var{z}
67 can have any value real or complex value. However, for polygamma functions
68 (@var{k} higher than 0), @var{z} must be real and non-negative.
69 
70 @seealso{gamma, gammainc, gammaln}
71 @end deftypefn */)
72 {
73  int nargin = args.length ();
74 
75  if (nargin < 1 || nargin > 2)
76  print_usage ();
77 
78  const octave_value oct_z = (nargin == 1) ? args(0) : args(1);
79  const octave_idx_type k = (nargin == 1) ? 0 : args(0).xidx_type_value ("psi: K must be an integer");
80  if (k < 0)
81  error ("psi: K must be non-negative");
82 
84 
85  if (k == 0)
86  {
87 #define FLOAT_BRANCH(T, A, M, E) \
88  if (oct_z.is_ ## T ##_type ()) \
89  { \
90  const A ## NDArray z = oct_z.M ## array_value (); \
91  A ## NDArray psi_z (z.dims ()); \
92  \
93  const E *zv = z.data (); \
94  E *psi_zv = psi_z.fortran_vec (); \
95  const octave_idx_type n = z.numel (); \
96  for (octave_idx_type i = 0; i < n; i++) \
97  *psi_zv++ = octave::math::psi (*zv++); \
98  \
99  retval = psi_z; \
100  }
101 
102  if (oct_z.iscomplex ())
103  {
104  FLOAT_BRANCH(double, Complex, complex_, Complex)
105  else FLOAT_BRANCH(single, FloatComplex, float_complex_, FloatComplex)
106  else
107  error ("psi: Z must be a floating point");
108  }
109  else
110  {
111  FLOAT_BRANCH(double, , , double)
112  else FLOAT_BRANCH(single, Float, float_, float)
113  else
114  error ("psi: Z must be a floating point");
115  }
116 
117 #undef FLOAT_BRANCH
118  }
119  else
120  {
121  if (! oct_z.isreal ())
122  error ("psi: Z must be real value for polygamma (K > 0)");
123 
124 #define FLOAT_BRANCH(T, A, M, E) \
125  if (oct_z.is_ ## T ##_type ()) \
126  { \
127  const A ## NDArray z = oct_z.M ## array_value (); \
128  A ## NDArray psi_z (z.dims ()); \
129  \
130  const E *zv = z.data (); \
131  E *psi_zv = psi_z.fortran_vec (); \
132  const octave_idx_type n = z.numel (); \
133  for (octave_idx_type i = 0; i < n; i++) \
134  { \
135  if (*zv < 0) \
136  error ("psi: Z must be non-negative for polygamma (K > 0)"); \
137  \
138  *psi_zv++ = octave::math::psi (k, *zv++); \
139  } \
140  retval = psi_z; \
141  }
142 
143  FLOAT_BRANCH(double, , , double)
144  else FLOAT_BRANCH(single, Float, float_, float)
145  else
146  error ("psi: Z must be a floating point for polygamma (K > 0)");
147 
148 #undef FLOAT_BRANCH
149  }
150 
151  return retval;
152 }
153 
154 /*
155 %!shared em
156 %! em = 0.577215664901532860606512090082402431042; # Euler-Mascheroni Constant
157 
158 %!assert (psi (ones (7, 3, 5)), repmat (-em, [7 3 5]))
159 %!assert (psi ([0 1]), [-Inf -em])
160 %!assert (psi ([-20:1]), [repmat(-Inf, [1 21]) -em])
161 %!assert (psi (single ([0 1])), single ([-Inf -em]))
162 
163 ## Abramowitz and Stegun, page 258, eq 6.3.5
164 %!test
165 %! z = [-100:-1 1:200] ./ 10; # drop the 0
166 %! assert (psi (z + 1), psi (z) + 1 ./ z, eps*1000);
167 
168 ## Abramowitz and Stegun, page 258, eq 6.3.2
169 %!assert (psi (1), -em)
170 
171 ## Abramowitz and Stegun, page 258, eq 6.3.3
172 %!assert (psi (1/2), -em - 2 * log (2))
173 
174 ## The following tests are from Pascal Sebah and Xavier Gourdon (2002)
175 ## "Introduction to the Gamma Function"
176 
177 ## Interesting identities of the digamma function, in section of 5.1.3
178 %!assert (psi (1/3), - em - (3/2) * log(3) - ((sqrt (3) / 6) * pi), eps*10)
179 %!assert (psi (1/4), - em -3 * log (2) - pi/2, eps*10)
180 %!assert (psi (1/6), - em -2 * log (2) - (3/2) * log (3) - ((sqrt (3) / 2) * pi), eps*10)
181 
182 ## First 6 zeros of the digamma function, in section of 5.1.5 (and also on
183 ## Abramowitz and Stegun, page 258, eq 6.3.19)
184 %!assert (psi ( 1.46163214496836234126265954232572132846819620400644), 0, eps)
185 %!assert (psi (-0.504083008264455409258269304533302498955385182368579), 0, eps*2)
186 %!assert (psi (-1.573498473162390458778286043690434612655040859116846), 0, eps*2)
187 %!assert (psi (-2.610720868444144650001537715718724207951074010873480), 0, eps*10)
188 %!assert (psi (-3.635293366436901097839181566946017713948423861193530), 0, eps*10)
189 %!assert (psi (-4.653237761743142441714598151148207363719069416133868), 0, eps*100)
190 
191 ## Tests for complex values
192 %!shared z
193 %! z = [-100:-1 1:200] ./ 10; # drop the 0
194 
195 ## Abramowitz and Stegun, page 259 eq 6.3.10
196 %!assert (real (psi (i*z)), real (psi (1 - i*z)))
197 
198 ## Abramowitz and Stegun, page 259 eq 6.3.11
199 %!assert (imag (psi (i*z)), 1/2 .* 1./z + 1/2 * pi * coth (pi * z), eps *10)
200 
201 ## Abramowitz and Stegun, page 259 eq 6.3.12
202 %!assert (imag (psi (1/2 + i*z)), 1/2 * pi * tanh (pi * z), eps*10)
203 
204 ## Abramowitz and Stegun, page 259 eq 6.3.13
205 %!assert (imag (psi (1 + i*z)), - 1./(2*z) + 1/2 * pi * coth (pi * z), eps*10)
206 
207 ## Abramowitz and Stegun, page 260 eq 6.4.5
208 %!test
209 %! for z = 0:20
210 %! assert (psi (1, z + 0.5),
211 %! 0.5 * (pi^2) - 4 * sum ((2*(1:z) -1) .^(-2)),
212 %! eps*10);
213 %! endfor
214 
215 ## Abramowitz and Stegun, page 260 eq 6.4.6
216 %!test
217 %! z = 0.1:0.1:20;
218 %! for n = 0:8
219 %! ## our precision goes down really quick when computing n is too high.
220 %! assert (psi (n, z+1),
221 %! psi (n, z) + ((-1)^n) * factorial (n) * (z.^(-n-1)), 0.1);
222 %! endfor
223 
224 ## Test input validation
225 %!error psi ()
226 %!error psi (1, 2, 3)
227 %!error <Z must be> psi ("non numeric")
228 %!error <K must be an integer> psi ({5.3}, 1)
229 %!error <K must be non-negative> psi (-5, 1)
230 %!error <Z must be non-negative for polygamma> psi (5, -1)
231 %!error <Z must be a floating point> psi (5, uint8 (-1))
232 %!error <Z must be real value for polygamma> psi (5, 5i)
233 
234 */
bool isreal(void) const
Definition: ov.h:691
bool iscomplex(void) const
Definition: ov.h:694
OCTINTERP_API void print_usage(void)
Definition: defun.cc:53
#define DEFUN(name, args_name, nargout_name, doc)
Macro to define a builtin function.
Definition: defun.h:56
void error(const char *fmt,...)
Definition: error.cc:968
double psi(double z)
Definition: lo-specfun.cc:2099
std::complex< double > Complex
Definition: oct-cmplx.h:33
std::complex< float > FloatComplex
Definition: oct-cmplx.h:34
octave_value::octave_value(const Array< char > &chm, char type) return retval
Definition: ov.cc:811
#define FLOAT_BRANCH(T, A, M, E)