GNU Octave  9.1.0
A high-level interpreted language, primarily intended for numerical computations, mostly compatible with Matlab
lo-specfun.cc
Go to the documentation of this file.
1 ////////////////////////////////////////////////////////////////////////
2 //
3 // Copyright (C) 1996-2024 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 <cmath>
31 
32 #include <algorithm>
33 #include <limits>
34 #include <string>
35 
36 #include "CColVector.h"
37 #include "CMatrix.h"
38 #include "CNDArray.h"
39 #include "Faddeeva.hh"
40 #include "dMatrix.h"
41 #include "dNDArray.h"
42 #include "dRowVector.h"
43 #include "f77-fcn.h"
44 #include "fCColVector.h"
45 #include "fCMatrix.h"
46 #include "fCNDArray.h"
47 #include "fMatrix.h"
48 #include "fNDArray.h"
49 #include "fRowVector.h"
50 #include "lo-amos-proto.h"
51 #include "lo-error.h"
52 #include "lo-ieee.h"
53 #include "lo-mappers.h"
54 #include "lo-slatec-proto.h"
55 #include "lo-specfun.h"
56 #include "mx-inlines.cc"
57 
59 
61 
62 static inline Complex
63 bessel_return_value (const Complex& val, octave_idx_type ierr)
64 {
65  static const Complex inf_val
68 
69  static const Complex nan_val
72 
73  Complex retval;
74 
75  switch (ierr)
76  {
77  case 0:
78  case 3:
79  case 4:
80  retval = val;
81  break;
82 
83  case 2:
84  retval = inf_val;
85  break;
86 
87  default:
88  retval = nan_val;
89  break;
90  }
91 
92  return retval;
93 }
94 
95 static inline FloatComplex
96 bessel_return_value (const FloatComplex& val, octave_idx_type ierr)
97 {
98  static const FloatComplex inf_val
101 
102  static const FloatComplex nan_val
105 
106  FloatComplex retval;
107 
108  switch (ierr)
109  {
110  case 0:
111  case 3:
112  case 4:
113  retval = val;
114  break;
115 
116  case 2:
117  retval = inf_val;
118  break;
119 
120  default:
121  retval = nan_val;
122  break;
123  }
124 
125  return retval;
126 }
127 
128 Complex
129 airy (const Complex& z, bool deriv, bool scaled, octave_idx_type& ierr)
130 {
131  double ar = 0.0;
132  double ai = 0.0;
133 
134  double zr = z.real ();
135  double zi = z.imag ();
136 
137  F77_INT id = (deriv ? 1 : 0);
138  F77_INT nz, t_ierr;
139  F77_INT sc = (scaled ? 2 : 1);
140 
141  F77_FUNC (zairy, ZAIRY) (zr, zi, id, sc, ar, ai, nz, t_ierr);
142 
143  ierr = t_ierr;
144 
145  if (zi == 0.0 && (! scaled || zr >= 0.0))
146  ai = 0.0;
147 
148  return bessel_return_value (Complex (ar, ai), ierr);
149 }
150 
152 airy (const ComplexMatrix& z, bool deriv, bool scaled,
154 {
155  octave_idx_type nr = z.rows ();
156  octave_idx_type nc = z.cols ();
157 
158  ComplexMatrix retval (nr, nc);
159 
160  ierr.resize (dim_vector (nr, nc));
161 
162  for (octave_idx_type j = 0; j < nc; j++)
163  for (octave_idx_type i = 0; i < nr; i++)
164  retval(i, j) = airy (z(i, j), deriv, scaled, ierr(i, j));
165 
166  return retval;
167 }
168 
170 airy (const ComplexNDArray& z, bool deriv, bool scaled,
172 {
173  dim_vector dv = z.dims ();
174  octave_idx_type nel = dv.numel ();
175  ComplexNDArray retval (dv);
176 
177  ierr.resize (dv);
178 
179  for (octave_idx_type i = 0; i < nel; i++)
180  retval(i) = airy (z(i), deriv, scaled, ierr(i));
181 
182  return retval;
183 }
184 
186 airy (const FloatComplex& z, bool deriv, bool scaled,
188 {
189  FloatComplex a;
190 
191  F77_INT id = (deriv ? 1 : 0);
192  F77_INT nz, t_ierr;
193  F77_INT sc = (scaled ? 2 : 1);
194 
195  F77_FUNC (cairy, CAIRY) (F77_CONST_CMPLX_ARG (&z), id, sc,
196  F77_CMPLX_ARG (&a), nz, t_ierr);
197 
198  ierr = t_ierr;
199 
200  float ar = a.real ();
201  float ai = a.imag ();
202 
203  if (z.imag () == 0.0 && (! scaled || z.real () >= 0.0))
204  ai = 0.0;
205 
206  return bessel_return_value (FloatComplex (ar, ai), ierr);
207 }
208 
210 airy (const FloatComplexMatrix& z, bool deriv, bool scaled,
212 {
213  octave_idx_type nr = z.rows ();
214  octave_idx_type nc = z.cols ();
215 
216  FloatComplexMatrix retval (nr, nc);
217 
218  ierr.resize (dim_vector (nr, nc));
219 
220  for (octave_idx_type j = 0; j < nc; j++)
221  for (octave_idx_type i = 0; i < nr; i++)
222  retval(i, j) = airy (z(i, j), deriv, scaled, ierr(i, j));
223 
224  return retval;
225 }
226 
228 airy (const FloatComplexNDArray& z, bool deriv, bool scaled,
230 {
231  dim_vector dv = z.dims ();
232  octave_idx_type nel = dv.numel ();
233  FloatComplexNDArray retval (dv);
234 
235  ierr.resize (dv);
236 
237  for (octave_idx_type i = 0; i < nel; i++)
238  retval(i) = airy (z(i), deriv, scaled, ierr(i));
239 
240  return retval;
241 }
242 
243 static inline bool
244 is_integer_value (double x)
245 {
246  return x == static_cast<double> (static_cast<long> (x));
247 }
248 
249 static inline Complex
250 zbesj (const Complex& z, double alpha, int kode, octave_idx_type& ierr);
251 
252 static inline Complex
253 zbesy (const Complex& z, double alpha, int kode, octave_idx_type& ierr);
254 
255 static inline Complex
256 zbesi (const Complex& z, double alpha, int kode, octave_idx_type& ierr);
257 
258 static inline Complex
259 zbesk (const Complex& z, double alpha, int kode, octave_idx_type& ierr);
260 
261 static inline Complex
262 zbesh1 (const Complex& z, double alpha, int kode, octave_idx_type& ierr);
263 
264 static inline Complex
265 zbesh2 (const Complex& z, double alpha, int kode, octave_idx_type& ierr);
266 
267 static inline Complex
268 zbesj (const Complex& z, double alpha, int kode, octave_idx_type& ierr)
269 {
270  Complex retval;
271 
272  if (alpha >= 0.0)
273  {
274  double yr = 0.0;
275  double yi = 0.0;
276 
277  F77_INT nz, t_ierr;
278 
279  double zr = z.real ();
280  double zi = z.imag ();
281 
282  F77_FUNC (zbesj, ZBESJ) (zr, zi, alpha, kode, 1, &yr, &yi, nz, t_ierr);
283 
284  ierr = t_ierr;
285 
286  if (zi == 0.0 && zr >= 0.0)
287  yi = 0.0;
288 
289  retval = bessel_return_value (Complex (yr, yi), ierr);
290  }
291  else if (is_integer_value (alpha))
292  {
293  // zbesy can overflow as z->0, and cause troubles for generic case below
294  alpha = -alpha;
295  Complex tmp = zbesj (z, alpha, kode, ierr);
296  if ((static_cast<long> (alpha)) & 1)
297  tmp = - tmp;
298  retval = bessel_return_value (tmp, ierr);
299  }
300  else
301  {
302  alpha = -alpha;
303 
304  Complex tmp = cos (M_PI * alpha) * zbesj (z, alpha, kode, ierr);
305 
306  if (ierr == 0 || ierr == 3)
307  {
308  tmp -= sin (M_PI * alpha) * zbesy (z, alpha, kode, ierr);
309 
310  retval = bessel_return_value (tmp, ierr);
311  }
312  else
315  }
316 
317  return retval;
318 }
319 
320 static inline Complex
321 zbesy (const Complex& z, double alpha, int kode, octave_idx_type& ierr)
322 {
323  Complex retval;
324 
325  if (alpha >= 0.0)
326  {
327  double yr = 0.0;
328  double yi = 0.0;
329 
330  F77_INT nz, t_ierr;
331 
332  double wr, wi;
333 
334  double zr = z.real ();
335  double zi = z.imag ();
336 
337  ierr = 0;
338 
339  if (zr == 0.0 && zi == 0.0)
340  {
342  yi = 0.0;
343  }
344  else
345  {
346  F77_FUNC (zbesy, ZBESY) (zr, zi, alpha, kode, 1, &yr, &yi, nz,
347  &wr, &wi, t_ierr);
348 
349  ierr = t_ierr;
350 
351  if (zi == 0.0 && zr >= 0.0)
352  yi = 0.0;
353  }
354 
355  return bessel_return_value (Complex (yr, yi), ierr);
356  }
357  else if (is_integer_value (alpha - 0.5))
358  {
359  // zbesy can overflow as z->0, and cause troubles for generic case below
360  alpha = -alpha;
361  Complex tmp = zbesj (z, alpha, kode, ierr);
362  if ((static_cast<long> (alpha - 0.5)) & 1)
363  tmp = - tmp;
364  retval = bessel_return_value (tmp, ierr);
365  }
366  else
367  {
368  alpha = -alpha;
369 
370  Complex tmp = cos (M_PI * alpha) * zbesy (z, alpha, kode, ierr);
371 
372  if (ierr == 0 || ierr == 3)
373  {
374  tmp += sin (M_PI * alpha) * zbesj (z, alpha, kode, ierr);
375 
376  retval = bessel_return_value (tmp, ierr);
377  }
378  else
381  }
382 
383  return retval;
384 }
385 
386 static inline Complex
387 zbesi (const Complex& z, double alpha, int kode, octave_idx_type& ierr)
388 {
389  Complex retval;
390 
391  if (alpha >= 0.0)
392  {
393  double yr = 0.0;
394  double yi = 0.0;
395 
396  F77_INT nz, t_ierr;
397 
398  double zr = z.real ();
399  double zi = z.imag ();
400 
401  F77_FUNC (zbesi, ZBESI) (zr, zi, alpha, kode, 1, &yr, &yi, nz, t_ierr);
402 
403  ierr = t_ierr;
404 
405  if (zi == 0.0 && zr >= 0.0)
406  yi = 0.0;
407 
408  retval = bessel_return_value (Complex (yr, yi), ierr);
409  }
410  else if (is_integer_value (alpha))
411  {
412  // zbesi can overflow as z->0, and cause troubles for generic case below
413  alpha = -alpha;
414  Complex tmp = zbesi (z, alpha, kode, ierr);
415  retval = bessel_return_value (tmp, ierr);
416  }
417  else
418  {
419  alpha = -alpha;
420 
421  Complex tmp = zbesi (z, alpha, kode, ierr);
422 
423  if (ierr == 0 || ierr == 3)
424  {
425  Complex tmp2 = (2.0 / M_PI) * sin (M_PI * alpha)
426  * zbesk (z, alpha, kode, ierr);
427 
428  if (kode == 2)
429  {
430  // Compensate for different scaling factor of besk.
431  tmp2 *= exp (-z - std::abs (z.real ()));
432  }
433 
434  tmp += tmp2;
435 
436  retval = bessel_return_value (tmp, ierr);
437  }
438  else
441  }
442 
443  return retval;
444 }
445 
446 static inline Complex
447 zbesk (const Complex& z, double alpha, int kode, octave_idx_type& ierr)
448 {
449  Complex retval;
450 
451  if (alpha >= 0.0)
452  {
453  double yr = 0.0;
454  double yi = 0.0;
455 
456  F77_INT nz, t_ierr;
457 
458  double zr = z.real ();
459  double zi = z.imag ();
460 
461  ierr = 0;
462 
463  if (zr == 0.0 && zi == 0.0)
464  {
466  yi = 0.0;
467  }
468  else
469  {
470  F77_FUNC (zbesk, ZBESK) (zr, zi, alpha, kode, 1, &yr, &yi, nz,
471  t_ierr);
472 
473  ierr = t_ierr;
474 
475  if (zi == 0.0 && zr >= 0.0)
476  yi = 0.0;
477  }
478 
479  retval = bessel_return_value (Complex (yr, yi), ierr);
480  }
481  else
482  {
483  Complex tmp = zbesk (z, -alpha, kode, ierr);
484 
485  retval = bessel_return_value (tmp, ierr);
486  }
487 
488  return retval;
489 }
490 
491 static inline Complex
492 zbesh1 (const Complex& z, double alpha, int kode, octave_idx_type& ierr)
493 {
494  Complex retval;
495 
496  if (alpha >= 0.0)
497  {
498  double yr = 0.0;
499  double yi = 0.0;
500 
501  F77_INT nz, t_ierr;
502 
503  double zr = z.real ();
504  double zi = z.imag ();
505 
506  F77_FUNC (zbesh, ZBESH) (zr, zi, alpha, kode, 1, 1, &yr, &yi, nz,
507  t_ierr);
508 
509  ierr = t_ierr;
510 
511  retval = bessel_return_value (Complex (yr, yi), ierr);
512  }
513  else
514  {
515  alpha = -alpha;
516 
517  static const Complex eye = Complex (0.0, 1.0);
518 
519  Complex tmp = exp (M_PI * alpha * eye) * zbesh1 (z, alpha, kode, ierr);
520 
521  retval = bessel_return_value (tmp, ierr);
522  }
523 
524  return retval;
525 }
526 
527 static inline Complex
528 zbesh2 (const Complex& z, double alpha, int kode, octave_idx_type& ierr)
529 {
530  Complex retval;
531 
532  if (alpha >= 0.0)
533  {
534  double yr = 0.0;
535  double yi = 0.0;
536 
537  F77_INT nz, t_ierr;
538 
539  double zr = z.real ();
540  double zi = z.imag ();
541 
542  F77_FUNC (zbesh, ZBESH) (zr, zi, alpha, kode, 2, 1, &yr, &yi, nz,
543  t_ierr);
544 
545  ierr = t_ierr;
546 
547  retval = bessel_return_value (Complex (yr, yi), ierr);
548  }
549  else
550  {
551  alpha = -alpha;
552 
553  static const Complex eye = Complex (0.0, 1.0);
554 
555  Complex tmp = exp (-M_PI * alpha * eye) * zbesh2 (z, alpha, kode, ierr);
556 
557  retval = bessel_return_value (tmp, ierr);
558  }
559 
560  return retval;
561 }
562 
563 typedef Complex (*dptr) (const Complex&, double, int, octave_idx_type&);
564 
565 static inline Complex
566 do_bessel (dptr f, const char *, double alpha, const Complex& x,
567  bool scaled, octave_idx_type& ierr)
568 {
569  Complex retval;
570 
571  retval = f (x, alpha, (scaled ? 2 : 1), ierr);
572 
573  return retval;
574 }
575 
576 static inline ComplexMatrix
577 do_bessel (dptr f, const char *, double alpha, const ComplexMatrix& x,
578  bool scaled, Array<octave_idx_type>& ierr)
579 {
580  octave_idx_type nr = x.rows ();
581  octave_idx_type nc = x.cols ();
582 
583  ComplexMatrix retval (nr, nc);
584 
585  ierr.resize (dim_vector (nr, nc));
586 
587  for (octave_idx_type j = 0; j < nc; j++)
588  for (octave_idx_type i = 0; i < nr; i++)
589  retval(i, j) = f (x(i, j), alpha, (scaled ? 2 : 1), ierr(i, j));
590 
591  return retval;
592 }
593 
594 static inline ComplexMatrix
595 do_bessel (dptr f, const char *, const Matrix& alpha, const Complex& x,
596  bool scaled, Array<octave_idx_type>& ierr)
597 {
598  octave_idx_type nr = alpha.rows ();
599  octave_idx_type nc = alpha.cols ();
600 
601  ComplexMatrix retval (nr, nc);
602 
603  ierr.resize (dim_vector (nr, nc));
604 
605  for (octave_idx_type j = 0; j < nc; j++)
606  for (octave_idx_type i = 0; i < nr; i++)
607  retval(i, j) = f (x, alpha(i, j), (scaled ? 2 : 1), ierr(i, j));
608 
609  return retval;
610 }
611 
612 static inline ComplexMatrix
613 do_bessel (dptr f, const char *fn, const Matrix& alpha,
614  const ComplexMatrix& x, bool scaled, Array<octave_idx_type>& ierr)
615 {
616  ComplexMatrix retval;
617 
618  octave_idx_type x_nr = x.rows ();
619  octave_idx_type x_nc = x.cols ();
620 
621  octave_idx_type alpha_nr = alpha.rows ();
622  octave_idx_type alpha_nc = alpha.cols ();
623 
624  if (x_nr != alpha_nr || x_nc != alpha_nc)
625  (*current_liboctave_error_handler)
626  ("%s: the sizes of alpha and x must conform", fn);
627 
628  octave_idx_type nr = x_nr;
629  octave_idx_type nc = x_nc;
630 
631  retval.resize (nr, nc);
632 
633  ierr.resize (dim_vector (nr, nc));
634 
635  for (octave_idx_type j = 0; j < nc; j++)
636  for (octave_idx_type i = 0; i < nr; i++)
637  retval(i, j) = f (x(i, j), alpha(i, j), (scaled ? 2 : 1), ierr(i, j));
638 
639  return retval;
640 }
641 
642 static inline ComplexNDArray
643 do_bessel (dptr f, const char *, double alpha, const ComplexNDArray& x,
644  bool scaled, Array<octave_idx_type>& ierr)
645 {
646  dim_vector dv = x.dims ();
647  octave_idx_type nel = dv.numel ();
648  ComplexNDArray retval (dv);
649 
650  ierr.resize (dv);
651 
652  for (octave_idx_type i = 0; i < nel; i++)
653  retval(i) = f (x(i), alpha, (scaled ? 2 : 1), ierr(i));
654 
655  return retval;
656 }
657 
658 static inline ComplexNDArray
659 do_bessel (dptr f, const char *, const NDArray& alpha, const Complex& x,
660  bool scaled, Array<octave_idx_type>& ierr)
661 {
662  dim_vector dv = alpha.dims ();
663  octave_idx_type nel = dv.numel ();
664  ComplexNDArray retval (dv);
665 
666  ierr.resize (dv);
667 
668  for (octave_idx_type i = 0; i < nel; i++)
669  retval(i) = f (x, alpha(i), (scaled ? 2 : 1), ierr(i));
670 
671  return retval;
672 }
673 
674 static inline ComplexNDArray
675 do_bessel (dptr f, const char *fn, const NDArray& alpha,
676  const ComplexNDArray& x, bool scaled, Array<octave_idx_type>& ierr)
677 {
678  dim_vector dv = x.dims ();
679  ComplexNDArray retval;
680 
681  if (dv != alpha.dims ())
682  (*current_liboctave_error_handler)
683  ("%s: the sizes of alpha and x must conform", fn);
684 
685  octave_idx_type nel = dv.numel ();
686 
687  retval.resize (dv);
688  ierr.resize (dv);
689 
690  for (octave_idx_type i = 0; i < nel; i++)
691  retval(i) = f (x(i), alpha(i), (scaled ? 2 : 1), ierr(i));
692 
693  return retval;
694 }
695 
696 static inline ComplexMatrix
697 do_bessel (dptr f, const char *, const RowVector& alpha,
698  const ComplexColumnVector& x, bool scaled,
700 {
701  octave_idx_type nr = x.numel ();
702  octave_idx_type nc = alpha.numel ();
703 
704  ComplexMatrix retval (nr, nc);
705 
706  ierr.resize (dim_vector (nr, nc));
707 
708  for (octave_idx_type j = 0; j < nc; j++)
709  for (octave_idx_type i = 0; i < nr; i++)
710  retval(i, j) = f (x(i), alpha(j), (scaled ? 2 : 1), ierr(i, j));
711 
712  return retval;
713 }
714 
715 #define SS_BESSEL(name, fcn) \
716  Complex \
717  name (double alpha, const Complex& x, bool scaled, octave_idx_type& ierr) \
718  { \
719  return do_bessel (fcn, #name, alpha, x, scaled, ierr); \
720  }
721 
722 #define SM_BESSEL(name, fcn) \
723  ComplexMatrix \
724  name (double alpha, const ComplexMatrix& x, bool scaled, \
725  Array<octave_idx_type>& ierr) \
726  { \
727  return do_bessel (fcn, #name, alpha, x, scaled, ierr); \
728  }
729 
730 #define MS_BESSEL(name, fcn) \
731  ComplexMatrix \
732  name (const Matrix& alpha, const Complex& x, bool scaled, \
733  Array<octave_idx_type>& ierr) \
734  { \
735  return do_bessel (fcn, #name, alpha, x, scaled, ierr); \
736  }
737 
738 #define MM_BESSEL(name, fcn) \
739  ComplexMatrix \
740  name (const Matrix& alpha, const ComplexMatrix& x, bool scaled, \
741  Array<octave_idx_type>& ierr) \
742  { \
743  return do_bessel (fcn, #name, alpha, x, scaled, ierr); \
744  }
745 
746 #define SN_BESSEL(name, fcn) \
747  ComplexNDArray \
748  name (double alpha, const ComplexNDArray& x, bool scaled, \
749  Array<octave_idx_type>& ierr) \
750  { \
751  return do_bessel (fcn, #name, alpha, x, scaled, ierr); \
752  }
753 
754 #define NS_BESSEL(name, fcn) \
755  ComplexNDArray \
756  name (const NDArray& alpha, const Complex& x, bool scaled, \
757  Array<octave_idx_type>& ierr) \
758  { \
759  return do_bessel (fcn, #name, alpha, x, scaled, ierr); \
760  }
761 
762 #define NN_BESSEL(name, fcn) \
763  ComplexNDArray \
764  name (const NDArray& alpha, const ComplexNDArray& x, bool scaled, \
765  Array<octave_idx_type>& ierr) \
766  { \
767  return do_bessel (fcn, #name, alpha, x, scaled, ierr); \
768  }
769 
770 #define RC_BESSEL(name, fcn) \
771  ComplexMatrix \
772  name (const RowVector& alpha, const ComplexColumnVector& x, bool scaled, \
773  Array<octave_idx_type>& ierr) \
774  { \
775  return do_bessel (fcn, #name, alpha, x, scaled, ierr); \
776  }
777 
778 #define ALL_BESSEL(name, fcn) \
779  SS_BESSEL (name, fcn) \
780  SM_BESSEL (name, fcn) \
781  MS_BESSEL (name, fcn) \
782  MM_BESSEL (name, fcn) \
783  SN_BESSEL (name, fcn) \
784  NS_BESSEL (name, fcn) \
785  NN_BESSEL (name, fcn) \
786  RC_BESSEL (name, fcn)
787 
794 
795 #undef ALL_BESSEL
796 #undef SS_BESSEL
797 #undef SM_BESSEL
798 #undef MS_BESSEL
799 #undef MM_BESSEL
800 #undef SN_BESSEL
801 #undef NS_BESSEL
802 #undef NN_BESSEL
803 #undef RC_BESSEL
804 
805 static inline FloatComplex
806 cbesj (const FloatComplex& z, float alpha, int kode, octave_idx_type& ierr);
807 
808 static inline FloatComplex
809 cbesy (const FloatComplex& z, float alpha, int kode, octave_idx_type& ierr);
810 
811 static inline FloatComplex
812 cbesi (const FloatComplex& z, float alpha, int kode, octave_idx_type& ierr);
813 
814 static inline FloatComplex
815 cbesk (const FloatComplex& z, float alpha, int kode, octave_idx_type& ierr);
816 
817 static inline FloatComplex
818 cbesh1 (const FloatComplex& z, float alpha, int kode, octave_idx_type& ierr);
819 
820 static inline FloatComplex
821 cbesh2 (const FloatComplex& z, float alpha, int kode, octave_idx_type& ierr);
822 
823 static inline bool
824 is_integer_value (float x)
825 {
826  return x == static_cast<float> (static_cast<long> (x));
827 }
828 
829 static inline FloatComplex
830 cbesj (const FloatComplex& z, float alpha, int kode, octave_idx_type& ierr)
831 {
832  FloatComplex retval;
833 
834  if (alpha >= 0.0)
835  {
836  FloatComplex y = 0.0;
837 
838  F77_INT nz, t_ierr;
839 
840  F77_FUNC (cbesj, CBESJ) (F77_CONST_CMPLX_ARG (&z), alpha, kode, 1,
841  F77_CMPLX_ARG (&y), nz, t_ierr);
842 
843  ierr = t_ierr;
844 
845  if (z.imag () == 0.0 && z.real () >= 0.0)
846  y = FloatComplex (y.real (), 0.0);
847 
848  retval = bessel_return_value (y, ierr);
849  }
850  else if (is_integer_value (alpha))
851  {
852  // zbesy can overflow as z->0, and cause troubles for generic case below
853  alpha = -alpha;
854  FloatComplex tmp = cbesj (z, alpha, kode, ierr);
855  if ((static_cast<long> (alpha)) & 1)
856  tmp = - tmp;
857  retval = bessel_return_value (tmp, ierr);
858  }
859  else
860  {
861  alpha = -alpha;
862 
863  FloatComplex tmp = cosf (static_cast<float> (M_PI) * alpha)
864  * cbesj (z, alpha, kode, ierr);
865 
866  if (ierr == 0 || ierr == 3)
867  {
868  tmp -= sinf (static_cast<float> (M_PI) * alpha)
869  * cbesy (z, alpha, kode, ierr);
870 
871  retval = bessel_return_value (tmp, ierr);
872  }
873  else
876  }
877 
878  return retval;
879 }
880 
881 static inline FloatComplex
882 cbesy (const FloatComplex& z, float alpha, int kode, octave_idx_type& ierr)
883 {
884  FloatComplex retval;
885 
886  if (alpha >= 0.0)
887  {
888  FloatComplex y = 0.0;
889 
890  F77_INT nz, t_ierr;
891 
892  FloatComplex w;
893 
894  ierr = 0;
895 
896  if (z.real () == 0.0 && z.imag () == 0.0)
897  {
899  }
900  else
901  {
902  F77_FUNC (cbesy, CBESY) (F77_CONST_CMPLX_ARG (&z), alpha, kode, 1,
903  F77_CMPLX_ARG (&y), nz,
904  F77_CMPLX_ARG (&w), t_ierr);
905 
906  ierr = t_ierr;
907 
908  if (z.imag () == 0.0 && z.real () >= 0.0)
909  y = FloatComplex (y.real (), 0.0);
910  }
911 
912  return bessel_return_value (y, ierr);
913  }
914  else if (is_integer_value (alpha - 0.5))
915  {
916  // zbesy can overflow as z->0, and cause troubles for generic case below
917  alpha = -alpha;
918  FloatComplex tmp = cbesj (z, alpha, kode, ierr);
919  if ((static_cast<long> (alpha - 0.5)) & 1)
920  tmp = - tmp;
921  retval = bessel_return_value (tmp, ierr);
922  }
923  else
924  {
925  alpha = -alpha;
926 
927  FloatComplex tmp = cosf (static_cast<float> (M_PI) * alpha)
928  * cbesy (z, alpha, kode, ierr);
929 
930  if (ierr == 0 || ierr == 3)
931  {
932  tmp += sinf (static_cast<float> (M_PI) * alpha)
933  * cbesj (z, alpha, kode, ierr);
934 
935  retval = bessel_return_value (tmp, ierr);
936  }
937  else
940  }
941 
942  return retval;
943 }
944 
945 static inline FloatComplex
946 cbesi (const FloatComplex& z, float alpha, int kode, octave_idx_type& ierr)
947 {
948  FloatComplex retval;
949 
950  if (alpha >= 0.0)
951  {
952  FloatComplex y = 0.0;
953 
954  F77_INT nz, t_ierr;
955 
956  F77_FUNC (cbesi, CBESI) (F77_CONST_CMPLX_ARG (&z), alpha, kode, 1,
957  F77_CMPLX_ARG (&y), nz, t_ierr);
958 
959  ierr = t_ierr;
960 
961  if (z.imag () == 0.0 && z.real () >= 0.0)
962  y = FloatComplex (y.real (), 0.0);
963 
964  retval = bessel_return_value (y, ierr);
965  }
966  else
967  {
968  alpha = -alpha;
969 
970  FloatComplex tmp = cbesi (z, alpha, kode, ierr);
971 
972  if (ierr == 0 || ierr == 3)
973  {
974  FloatComplex tmp2 = static_cast<float> (2.0 / M_PI)
975  * sinf (static_cast<float> (M_PI) * alpha)
976  * cbesk (z, alpha, kode, ierr);
977 
978  if (kode == 2)
979  {
980  // Compensate for different scaling factor of besk.
981  tmp2 *= exp (-z - std::abs (z.real ()));
982  }
983 
984  tmp += tmp2;
985 
986  retval = bessel_return_value (tmp, ierr);
987  }
988  else
991  }
992 
993  return retval;
994 }
995 
996 static inline FloatComplex
997 cbesk (const FloatComplex& z, float alpha, int kode, octave_idx_type& ierr)
998 {
999  FloatComplex retval;
1000 
1001  if (alpha >= 0.0)
1002  {
1003  FloatComplex y = 0.0;
1004 
1005  F77_INT nz, t_ierr;
1006 
1007  ierr = 0;
1008 
1009  if (z.real () == 0.0 && z.imag () == 0.0)
1010  {
1012  }
1013  else
1014  {
1015  F77_FUNC (cbesk, CBESK) (F77_CONST_CMPLX_ARG (&z), alpha, kode, 1,
1016  F77_CMPLX_ARG (&y), nz, t_ierr);
1017 
1018  ierr = t_ierr;
1019 
1020  if (z.imag () == 0.0 && z.real () >= 0.0)
1021  y = FloatComplex (y.real (), 0.0);
1022  }
1023 
1024  retval = bessel_return_value (y, ierr);
1025  }
1026  else
1027  {
1028  FloatComplex tmp = cbesk (z, -alpha, kode, ierr);
1029 
1030  retval = bessel_return_value (tmp, ierr);
1031  }
1032 
1033  return retval;
1034 }
1035 
1036 static inline FloatComplex
1037 cbesh1 (const FloatComplex& z, float alpha, int kode, octave_idx_type& ierr)
1038 {
1039  FloatComplex retval;
1040 
1041  if (alpha >= 0.0)
1042  {
1043  FloatComplex y = 0.0;
1044 
1045  F77_INT nz, t_ierr;
1046 
1047  F77_FUNC (cbesh, CBESH) (F77_CONST_CMPLX_ARG (&z), alpha, kode, 1, 1,
1048  F77_CMPLX_ARG (&y), nz, t_ierr);
1049 
1050  ierr = t_ierr;
1051 
1052  retval = bessel_return_value (y, ierr);
1053  }
1054  else
1055  {
1056  alpha = -alpha;
1057 
1058  static const FloatComplex eye = FloatComplex (0.0, 1.0);
1059 
1060  FloatComplex tmp = exp (static_cast<float> (M_PI) * alpha * eye)
1061  * cbesh1 (z, alpha, kode, ierr);
1062 
1063  retval = bessel_return_value (tmp, ierr);
1064  }
1065 
1066  return retval;
1067 }
1068 
1069 static inline FloatComplex
1070 cbesh2 (const FloatComplex& z, float alpha, int kode, octave_idx_type& ierr)
1071 {
1072  FloatComplex retval;
1073 
1074  if (alpha >= 0.0)
1075  {
1076  FloatComplex y = 0.0;
1077 
1078  F77_INT nz, t_ierr;
1079 
1080  F77_FUNC (cbesh, CBESH) (F77_CONST_CMPLX_ARG (&z), alpha, kode, 2, 1,
1081  F77_CMPLX_ARG (&y), nz, t_ierr);
1082 
1083  ierr = t_ierr;
1084 
1085  retval = bessel_return_value (y, ierr);
1086  }
1087  else
1088  {
1089  alpha = -alpha;
1090 
1091  static const FloatComplex eye = FloatComplex (0.0, 1.0);
1092 
1093  FloatComplex tmp = exp (-static_cast<float> (M_PI) * alpha * eye)
1094  * cbesh2 (z, alpha, kode, ierr);
1095 
1096  retval = bessel_return_value (tmp, ierr);
1097  }
1098 
1099  return retval;
1100 }
1101 
1102 typedef FloatComplex (*fptr) (const FloatComplex&, float, int,
1103  octave_idx_type&);
1104 
1105 static inline FloatComplex
1106 do_bessel (fptr f, const char *, float alpha, const FloatComplex& x,
1107  bool scaled, octave_idx_type& ierr)
1108 {
1109  FloatComplex retval;
1110 
1111  retval = f (x, alpha, (scaled ? 2 : 1), ierr);
1112 
1113  return retval;
1114 }
1115 
1116 static inline FloatComplexMatrix
1117 do_bessel (fptr f, const char *, float alpha, const FloatComplexMatrix& x,
1118  bool scaled, Array<octave_idx_type>& ierr)
1119 {
1120  octave_idx_type nr = x.rows ();
1121  octave_idx_type nc = x.cols ();
1122 
1123  FloatComplexMatrix retval (nr, nc);
1124 
1125  ierr.resize (dim_vector (nr, nc));
1126 
1127  for (octave_idx_type j = 0; j < nc; j++)
1128  for (octave_idx_type i = 0; i < nr; i++)
1129  retval(i, j) = f (x(i, j), alpha, (scaled ? 2 : 1), ierr(i, j));
1130 
1131  return retval;
1132 }
1133 
1134 static inline FloatComplexMatrix
1135 do_bessel (fptr f, const char *, const FloatMatrix& alpha,
1136  const FloatComplex& x,
1137  bool scaled, Array<octave_idx_type>& ierr)
1138 {
1139  octave_idx_type nr = alpha.rows ();
1140  octave_idx_type nc = alpha.cols ();
1141 
1142  FloatComplexMatrix retval (nr, nc);
1143 
1144  ierr.resize (dim_vector (nr, nc));
1145 
1146  for (octave_idx_type j = 0; j < nc; j++)
1147  for (octave_idx_type i = 0; i < nr; i++)
1148  retval(i, j) = f (x, alpha(i, j), (scaled ? 2 : 1), ierr(i, j));
1149 
1150  return retval;
1151 }
1152 
1153 static inline FloatComplexMatrix
1154 do_bessel (fptr f, const char *fn, const FloatMatrix& alpha,
1155  const FloatComplexMatrix& x, bool scaled,
1157 {
1158  FloatComplexMatrix retval;
1159 
1160  octave_idx_type x_nr = x.rows ();
1161  octave_idx_type x_nc = x.cols ();
1162 
1163  octave_idx_type alpha_nr = alpha.rows ();
1164  octave_idx_type alpha_nc = alpha.cols ();
1165 
1166  if (x_nr != alpha_nr || x_nc != alpha_nc)
1167  (*current_liboctave_error_handler)
1168  ("%s: the sizes of alpha and x must conform", fn);
1169 
1170  octave_idx_type nr = x_nr;
1171  octave_idx_type nc = x_nc;
1172 
1173  retval.resize (nr, nc);
1174 
1175  ierr.resize (dim_vector (nr, nc));
1176 
1177  for (octave_idx_type j = 0; j < nc; j++)
1178  for (octave_idx_type i = 0; i < nr; i++)
1179  retval(i, j) = f (x(i, j), alpha(i, j), (scaled ? 2 : 1), ierr(i, j));
1180 
1181  return retval;
1182 }
1183 
1184 static inline FloatComplexNDArray
1185 do_bessel (fptr f, const char *, float alpha, const FloatComplexNDArray& x,
1186  bool scaled, Array<octave_idx_type>& ierr)
1187 {
1188  dim_vector dv = x.dims ();
1189  octave_idx_type nel = dv.numel ();
1190  FloatComplexNDArray retval (dv);
1191 
1192  ierr.resize (dv);
1193 
1194  for (octave_idx_type i = 0; i < nel; i++)
1195  retval(i) = f (x(i), alpha, (scaled ? 2 : 1), ierr(i));
1196 
1197  return retval;
1198 }
1199 
1200 static inline FloatComplexNDArray
1201 do_bessel (fptr f, const char *, const FloatNDArray& alpha,
1202  const FloatComplex& x, bool scaled, Array<octave_idx_type>& ierr)
1203 {
1204  dim_vector dv = alpha.dims ();
1205  octave_idx_type nel = dv.numel ();
1206  FloatComplexNDArray retval (dv);
1207 
1208  ierr.resize (dv);
1209 
1210  for (octave_idx_type i = 0; i < nel; i++)
1211  retval(i) = f (x, alpha(i), (scaled ? 2 : 1), ierr(i));
1212 
1213  return retval;
1214 }
1215 
1216 static inline FloatComplexNDArray
1217 do_bessel (fptr f, const char *fn, const FloatNDArray& alpha,
1218  const FloatComplexNDArray& x, bool scaled,
1220 {
1221  dim_vector dv = x.dims ();
1222  FloatComplexNDArray retval;
1223 
1224  if (dv != alpha.dims ())
1225  (*current_liboctave_error_handler)
1226  ("%s: the sizes of alpha and x must conform", fn);
1227 
1228  octave_idx_type nel = dv.numel ();
1229 
1230  retval.resize (dv);
1231  ierr.resize (dv);
1232 
1233  for (octave_idx_type i = 0; i < nel; i++)
1234  retval(i) = f (x(i), alpha(i), (scaled ? 2 : 1), ierr(i));
1235 
1236  return retval;
1237 }
1238 
1239 static inline FloatComplexMatrix
1240 do_bessel (fptr f, const char *, const FloatRowVector& alpha,
1241  const FloatComplexColumnVector& x, bool scaled,
1243 {
1244  octave_idx_type nr = x.numel ();
1245  octave_idx_type nc = alpha.numel ();
1246 
1247  FloatComplexMatrix retval (nr, nc);
1248 
1249  ierr.resize (dim_vector (nr, nc));
1250 
1251  for (octave_idx_type j = 0; j < nc; j++)
1252  for (octave_idx_type i = 0; i < nr; i++)
1253  retval(i, j) = f (x(i), alpha(j), (scaled ? 2 : 1), ierr(i, j));
1254 
1255  return retval;
1256 }
1257 
1258 #define SS_BESSEL(name, fcn) \
1259  FloatComplex \
1260  name (float alpha, const FloatComplex& x, bool scaled, \
1261  octave_idx_type& ierr) \
1262  { \
1263  return do_bessel (fcn, #name, alpha, x, scaled, ierr); \
1264  }
1265 
1266 #define SM_BESSEL(name, fcn) \
1267  FloatComplexMatrix \
1268  name (float alpha, const FloatComplexMatrix& x, bool scaled, \
1269  Array<octave_idx_type>& ierr) \
1270  { \
1271  return do_bessel (fcn, #name, alpha, x, scaled, ierr); \
1272  }
1273 
1274 #define MS_BESSEL(name, fcn) \
1275  FloatComplexMatrix \
1276  name (const FloatMatrix& alpha, const FloatComplex& x, bool scaled, \
1277  Array<octave_idx_type>& ierr) \
1278  { \
1279  return do_bessel (fcn, #name, alpha, x, scaled, ierr); \
1280  }
1281 
1282 #define MM_BESSEL(name, fcn) \
1283  FloatComplexMatrix \
1284  name (const FloatMatrix& alpha, const FloatComplexMatrix& x, \
1285  bool scaled, Array<octave_idx_type>& ierr) \
1286  { \
1287  return do_bessel (fcn, #name, alpha, x, scaled, ierr); \
1288  }
1289 
1290 #define SN_BESSEL(name, fcn) \
1291  FloatComplexNDArray \
1292  name (float alpha, const FloatComplexNDArray& x, bool scaled, \
1293  Array<octave_idx_type>& ierr) \
1294  { \
1295  return do_bessel (fcn, #name, alpha, x, scaled, ierr); \
1296  }
1297 
1298 #define NS_BESSEL(name, fcn) \
1299  FloatComplexNDArray \
1300  name (const FloatNDArray& alpha, const FloatComplex& x, \
1301  bool scaled, Array<octave_idx_type>& ierr) \
1302  { \
1303  return do_bessel (fcn, #name, alpha, x, scaled, ierr); \
1304  }
1305 
1306 #define NN_BESSEL(name, fcn) \
1307  FloatComplexNDArray \
1308  name (const FloatNDArray& alpha, const FloatComplexNDArray& x, \
1309  bool scaled, Array<octave_idx_type>& ierr) \
1310  { \
1311  return do_bessel (fcn, #name, alpha, x, scaled, ierr); \
1312  }
1313 
1314 #define RC_BESSEL(name, fcn) \
1315  FloatComplexMatrix \
1316  name (const FloatRowVector& alpha, \
1317  const FloatComplexColumnVector& x, bool scaled, \
1318  Array<octave_idx_type>& ierr) \
1319  { \
1320  return do_bessel (fcn, #name, alpha, x, scaled, ierr); \
1321  }
1322 
1323 #define ALL_BESSEL(name, fcn) \
1324  SS_BESSEL (name, fcn) \
1325  SM_BESSEL (name, fcn) \
1326  MS_BESSEL (name, fcn) \
1327  MM_BESSEL (name, fcn) \
1328  SN_BESSEL (name, fcn) \
1329  NS_BESSEL (name, fcn) \
1330  NN_BESSEL (name, fcn) \
1331  RC_BESSEL (name, fcn)
1332 
1339 
1340 #undef ALL_BESSEL
1341 #undef SS_BESSEL
1342 #undef SM_BESSEL
1343 #undef MS_BESSEL
1344 #undef MM_BESSEL
1345 #undef SN_BESSEL
1346 #undef NS_BESSEL
1347 #undef NN_BESSEL
1348 #undef RC_BESSEL
1349 
1350 Complex
1351 biry (const Complex& z, bool deriv, bool scaled, octave_idx_type& ierr)
1352 {
1353  double ar = 0.0;
1354  double ai = 0.0;
1355 
1356  double zr = z.real ();
1357  double zi = z.imag ();
1358 
1359  F77_INT id = (deriv ? 1 : 0);
1360  F77_INT t_ierr;
1361  F77_INT sc = (scaled ? 2 : 1);
1362 
1363  F77_FUNC (zbiry, ZBIRY) (zr, zi, id, sc, ar, ai, t_ierr);
1364 
1365  ierr = t_ierr;
1366 
1367  if (zi == 0.0 && (! scaled || zr >= 0.0))
1368  ai = 0.0;
1369 
1370  return bessel_return_value (Complex (ar, ai), ierr);
1371 }
1372 
1374 biry (const ComplexMatrix& z, bool deriv, bool scaled,
1376 {
1377  octave_idx_type nr = z.rows ();
1378  octave_idx_type nc = z.cols ();
1379 
1380  ComplexMatrix retval (nr, nc);
1381 
1382  ierr.resize (dim_vector (nr, nc));
1383 
1384  for (octave_idx_type j = 0; j < nc; j++)
1385  for (octave_idx_type i = 0; i < nr; i++)
1386  retval(i, j) = biry (z(i, j), deriv, scaled, ierr(i, j));
1387 
1388  return retval;
1389 }
1390 
1392 biry (const ComplexNDArray& z, bool deriv, bool scaled,
1394 {
1395  dim_vector dv = z.dims ();
1396  octave_idx_type nel = dv.numel ();
1397  ComplexNDArray retval (dv);
1398 
1399  ierr.resize (dv);
1400 
1401  for (octave_idx_type i = 0; i < nel; i++)
1402  retval(i) = biry (z(i), deriv, scaled, ierr(i));
1403 
1404  return retval;
1405 }
1406 
1408 biry (const FloatComplex& z, bool deriv, bool scaled,
1410 {
1411  FloatComplex a;
1412 
1413  F77_INT id = (deriv ? 1 : 0);
1414  F77_INT t_ierr;
1415  F77_INT sc = (scaled ? 2 : 1);
1416 
1417  F77_FUNC (cbiry, CBIRY) (F77_CONST_CMPLX_ARG (&z), id, sc,
1418  F77_CMPLX_ARG (&a), t_ierr);
1419 
1420  ierr = t_ierr;
1421 
1422  float ar = a.real ();
1423  float ai = a.imag ();
1424 
1425  if (z.imag () == 0.0 && (! scaled || z.real () >= 0.0))
1426  ai = 0.0;
1427 
1428  return bessel_return_value (FloatComplex (ar, ai), ierr);
1429 }
1430 
1432 biry (const FloatComplexMatrix& z, bool deriv, bool scaled,
1434 {
1435  octave_idx_type nr = z.rows ();
1436  octave_idx_type nc = z.cols ();
1437 
1438  FloatComplexMatrix retval (nr, nc);
1439 
1440  ierr.resize (dim_vector (nr, nc));
1441 
1442  for (octave_idx_type j = 0; j < nc; j++)
1443  for (octave_idx_type i = 0; i < nr; i++)
1444  retval(i, j) = biry (z(i, j), deriv, scaled, ierr(i, j));
1445 
1446  return retval;
1447 }
1448 
1450 biry (const FloatComplexNDArray& z, bool deriv, bool scaled,
1452 {
1453  dim_vector dv = z.dims ();
1454  octave_idx_type nel = dv.numel ();
1455  FloatComplexNDArray retval (dv);
1456 
1457  ierr.resize (dv);
1458 
1459  for (octave_idx_type i = 0; i < nel; i++)
1460  retval(i) = biry (z(i), deriv, scaled, ierr(i));
1461 
1462  return retval;
1463 }
1464 
1465 // Real and complex Dawson function (= scaled erfi) from Faddeeva package
1466 double
1467 dawson (double x) { return Faddeeva::Dawson (x); }
1468 float
1469 dawson (float x) { return Faddeeva::Dawson (x); }
1470 
1471 Complex
1472 dawson (const Complex& x)
1473 {
1474  return Faddeeva::Dawson (x);
1475 }
1476 
1479 {
1480  Complex xd (x.real (), x.imag ());
1481  Complex ret;
1482  ret = Faddeeva::Dawson (xd, std::numeric_limits<float>::epsilon ());
1483  return FloatComplex (ret.real (), ret.imag ());
1484 }
1485 
1486 void
1487 ellipj (double u, double m, double& sn, double& cn, double& dn, double& err)
1488 {
1489  static const int Nmax = 16;
1490  double m1, t=0, si_u, co_u, se_u, ta_u, b, c[Nmax], a[Nmax], phi;
1491  int n, Nn, ii;
1492 
1493  if (m < 0 || m > 1)
1494  {
1495  (*current_liboctave_warning_with_id_handler)
1496  ("Octave:ellipj-invalid-m",
1497  "ellipj: invalid M value, required value 0 <= M <= 1");
1498 
1499  sn = cn = dn = lo_ieee_nan_value ();
1500 
1501  return;
1502  }
1503 
1504  double sqrt_eps = std::sqrt (std::numeric_limits<double>::epsilon ());
1505  if (m < sqrt_eps)
1506  {
1507  // For small m, (Abramowitz and Stegun, Section 16.13)
1508  si_u = sin (u);
1509  co_u = cos (u);
1510  t = 0.25*m*(u - si_u*co_u);
1511  sn = si_u - t * co_u;
1512  cn = co_u + t * si_u;
1513  dn = 1 - 0.5*m*si_u*si_u;
1514  }
1515  else if ((1 - m) < sqrt_eps)
1516  {
1517  // For m1 = (1-m) small (Abramowitz and Stegun, Section 16.15)
1518  m1 = 1 - m;
1519  si_u = sinh (u);
1520  co_u = cosh (u);
1521  ta_u = tanh (u);
1522  se_u = 1/co_u;
1523  sn = ta_u + 0.25*m1*(si_u*co_u - u)*se_u*se_u;
1524  cn = se_u - 0.25*m1*(si_u*co_u - u)*ta_u*se_u;
1525  dn = se_u + 0.25*m1*(si_u*co_u + u)*ta_u*se_u;
1526  }
1527  else
1528  {
1529  // Arithmetic-Geometric Mean (AGM) algorithm
1530  // (Abramowitz and Stegun, Section 16.4)
1531  a[0] = 1;
1532  b = std::sqrt (1 - m);
1533  c[0] = std::sqrt (m);
1534  for (n = 1; n < Nmax; ++n)
1535  {
1536  a[n] = (a[n - 1] + b)/2;
1537  c[n] = (a[n - 1] - b)/2;
1538  b = std::sqrt (a[n - 1]*b);
1539  if (c[n]/a[n] < std::numeric_limits<double>::epsilon ()) break;
1540  }
1541  if (n >= Nmax - 1)
1542  {
1543  err = 1;
1544  return;
1545  }
1546  Nn = n;
1547  for (ii = 1; n > 0; ii *= 2, --n) {} // ii = pow(2,Nn)
1548  phi = ii*a[Nn]*u;
1549  for (n = Nn; n > 0; --n)
1550  {
1551  phi = (std::asin ((c[n]/a[n])* sin (phi)) + phi)/2;
1552  }
1553  sn = sin (phi);
1554  cn = cos (phi);
1555  dn = std::sqrt (1 - m*sn*sn);
1556  }
1557 }
1558 
1559 void
1560 ellipj (const Complex& u, double m, Complex& sn, Complex& cn, Complex& dn,
1561  double& err)
1562 {
1563  double m1 = 1 - m, ss1, cc1, dd1;
1564 
1565  ellipj (u.imag (), m1, ss1, cc1, dd1, err);
1566  if (u.real () == 0)
1567  {
1568  // u is pure imag: Jacoby imag. transf.
1569  sn = Complex (0, ss1/cc1);
1570  cn = 1/cc1; // cn.imag = 0;
1571  dn = dd1/cc1; // dn.imag = 0;
1572  }
1573  else
1574  {
1575  // u is generic complex
1576  double ss, cc, dd, ddd;
1577 
1578  ellipj (u.real (), m, ss, cc, dd, err);
1579  ddd = cc1*cc1 + m*ss*ss*ss1*ss1;
1580  sn = Complex (ss*dd1/ddd, cc*dd*ss1*cc1/ddd);
1581  cn = Complex (cc*cc1/ddd, -ss*dd*ss1*dd1/ddd);
1582  dn = Complex (dd*cc1*dd1/ddd, -m*ss*cc*ss1/ddd);
1583  }
1584 }
1585 
1586 // Complex error function from the Faddeeva package
1587 Complex
1588 erf (const Complex& x)
1589 {
1590  return Faddeeva::erf (x);
1591 }
1592 
1595 {
1596  Complex xd (x.real (), x.imag ());
1597  Complex ret = Faddeeva::erf (xd, std::numeric_limits<float>::epsilon ());
1598  return FloatComplex (ret.real (), ret.imag ());
1599 }
1600 
1601 // Complex complementary error function from the Faddeeva package
1602 Complex
1603 erfc (const Complex& x)
1604 {
1605  return Faddeeva::erfc (x);
1606 }
1607 
1610 {
1611  Complex xd (x.real (), x.imag ());
1612  Complex ret = Faddeeva::erfc (xd, std::numeric_limits<float>::epsilon ());
1613  return FloatComplex (ret.real (), ret.imag ());
1614 }
1615 
1616 // The algorithm for erfcinv is an adaptation of the erfinv algorithm
1617 // above from P. J. Acklam. It has been modified to run over the
1618 // different input domain of erfcinv. See the notes for erfinv for an
1619 // explanation.
1620 
1621 static double
1622 do_erfcinv (double x, bool refine)
1623 {
1624  // Coefficients of rational approximation.
1625  static const double a[] =
1626  {
1627  -2.806989788730439e+01, 1.562324844726888e+02,
1628  -1.951109208597547e+02, 9.783370457507161e+01,
1629  -2.168328665628878e+01, 1.772453852905383e+00
1630  };
1631  static const double b[] =
1632  {
1633  -5.447609879822406e+01, 1.615858368580409e+02,
1634  -1.556989798598866e+02, 6.680131188771972e+01,
1635  -1.328068155288572e+01
1636  };
1637  static const double c[] =
1638  {
1639  -5.504751339936943e-03, -2.279687217114118e-01,
1640  -1.697592457770869e+00, -1.802933168781950e+00,
1641  3.093354679843505e+00, 2.077595676404383e+00
1642  };
1643  static const double d[] =
1644  {
1645  7.784695709041462e-03, 3.224671290700398e-01,
1646  2.445134137142996e+00, 3.754408661907416e+00
1647  };
1648 
1649  static const double spi2 = 8.862269254527579e-01; // sqrt(pi)/2.
1650  static const double pbreak_lo = 0.04850; // 1-pbreak
1651  static const double pbreak_hi = 1.95150; // 1+pbreak
1652  double y;
1653 
1654  // Select case.
1655  if (x >= pbreak_lo && x <= pbreak_hi)
1656  {
1657  // Middle region.
1658  const double q = 0.5*(1-x), r = q*q;
1659  const double yn = (((((a[0]*r + a[1])*r + a[2])*r + a[3])*r + a[4])*r + a[5])*q;
1660  const double yd = ((((b[0]*r + b[1])*r + b[2])*r + b[3])*r + b[4])*r + 1.0;
1661  y = yn / yd;
1662  }
1663  else if (x > 0.0 && x < 2.0)
1664  {
1665  // Tail region.
1666  const double q = (x < 1
1667  ? std::sqrt (-2*std::log (0.5*x))
1668  : std::sqrt (-2*std::log (0.5*(2-x))));
1669 
1670  const double yn = ((((c[0]*q + c[1])*q + c[2])*q + c[3])*q + c[4])*q + c[5];
1671 
1672  const double yd = (((d[0]*q + d[1])*q + d[2])*q + d[3])*q + 1.0;
1673 
1674  y = yn / yd;
1675 
1676  if (x < pbreak_lo)
1677  y = -y;
1678  }
1679  else if (x == 0.0)
1680  return numeric_limits<double>::Inf ();
1681  else if (x == 2.0)
1682  return -numeric_limits<double>::Inf ();
1683  else
1684  return numeric_limits<double>::NaN ();
1685 
1686  if (refine)
1687  {
1688  // One iteration of Halley's method gives full precision.
1689  double u = (erf (y) - (1-x)) * spi2 * exp (y*y);
1690  y -= u / (1 + y*u);
1691  }
1692 
1693  return y;
1694 }
1695 
1696 double
1697 erfcinv (double x)
1698 {
1699  return do_erfcinv (x, true);
1700 }
1701 
1702 float
1703 erfcinv (float x)
1704 {
1705  return do_erfcinv (x, false);
1706 }
1707 
1708 // Real and complex scaled complementary error function from Faddeeva pkg.
1709 double
1710 erfcx (double x) { return Faddeeva::erfcx (x); }
1711 float
1712 erfcx (float x) { return Faddeeva::erfcx (x); }
1713 
1714 Complex
1715 erfcx (const Complex& x)
1716 {
1717  return Faddeeva::erfcx (x);
1718 }
1719 
1722 {
1723  Complex xd (x.real (), x.imag ());
1724  Complex ret;
1725  ret = Faddeeva::erfcx (xd, std::numeric_limits<float>::epsilon ());
1726  return FloatComplex (ret.real (), ret.imag ());
1727 }
1728 
1729 // Real and complex imaginary error function from Faddeeva package
1730 double
1731 erfi (double x) { return Faddeeva::erfi (x); }
1732 float
1733 erfi (float x) { return Faddeeva::erfi (x); }
1734 
1735 Complex
1736 erfi (const Complex& x)
1737 {
1738  return Faddeeva::erfi (x);
1739 }
1740 
1743 {
1744  Complex xd (x.real (), x.imag ());
1745  Complex ret = Faddeeva::erfi (xd, std::numeric_limits<float>::epsilon ());
1746  return FloatComplex (ret.real (), ret.imag ());
1747 }
1748 
1749 // This algorithm is due to P. J. Acklam.
1750 //
1751 // See http://home.online.no/~pjacklam/notes/invnorm/
1752 //
1753 // The rational approximation has relative accuracy 1.15e-9 in the whole
1754 // region. For doubles, it is refined by a single step of Halley's 3rd
1755 // order method. For single precision, the accuracy is already OK, so
1756 // we skip it to get faster evaluation.
1757 
1758 static double
1759 do_erfinv (double x, bool refine)
1760 {
1761  // Coefficients of rational approximation.
1762  static const double a[] =
1763  {
1764  -2.806989788730439e+01, 1.562324844726888e+02,
1765  -1.951109208597547e+02, 9.783370457507161e+01,
1766  -2.168328665628878e+01, 1.772453852905383e+00
1767  };
1768  static const double b[] =
1769  {
1770  -5.447609879822406e+01, 1.615858368580409e+02,
1771  -1.556989798598866e+02, 6.680131188771972e+01,
1772  -1.328068155288572e+01
1773  };
1774  static const double c[] =
1775  {
1776  -5.504751339936943e-03, -2.279687217114118e-01,
1777  -1.697592457770869e+00, -1.802933168781950e+00,
1778  3.093354679843505e+00, 2.077595676404383e+00
1779  };
1780  static const double d[] =
1781  {
1782  7.784695709041462e-03, 3.224671290700398e-01,
1783  2.445134137142996e+00, 3.754408661907416e+00
1784  };
1785 
1786  static const double spi2 = 8.862269254527579e-01; // sqrt(pi)/2.
1787  static const double pbreak = 0.95150;
1788  double ax = fabs (x), y;
1789 
1790  // Select case.
1791  if (ax <= pbreak)
1792  {
1793  // Middle region.
1794  const double q = 0.5 * x, r = q*q;
1795  const double yn = (((((a[0]*r + a[1])*r + a[2])*r + a[3])*r + a[4])*r + a[5])*q;
1796  const double yd = ((((b[0]*r + b[1])*r + b[2])*r + b[3])*r + b[4])*r + 1.0;
1797  y = yn / yd;
1798  }
1799  else if (ax < 1.0)
1800  {
1801  // Tail region.
1802  const double q = std::sqrt (-2*std::log (0.5*(1-ax)));
1803  const double yn = ((((c[0]*q + c[1])*q + c[2])*q + c[3])*q + c[4])*q + c[5];
1804  const double yd = (((d[0]*q + d[1])*q + d[2])*q + d[3])*q + 1.0;
1805  y = yn / yd * math::signum (-x);
1806  }
1807  else if (ax == 1.0)
1809  else
1810  return numeric_limits<double>::NaN ();
1811 
1812  if (refine)
1813  {
1814  // One iteration of Halley's method gives full precision.
1815  double u = (erf (y) - x) * spi2 * exp (y*y);
1816  y -= u / (1 + y*u);
1817  }
1818 
1819  return y;
1820 }
1821 
1822 double
1823 erfinv (double x)
1824 {
1825  return do_erfinv (x, true);
1826 }
1827 
1828 float
1829 erfinv (float x)
1830 {
1831  return do_erfinv (x, false);
1832 }
1833 
1834 Complex
1835 expm1 (const Complex& x)
1836 {
1837  Complex retval;
1838 
1839  if (std::abs (x) < 1)
1840  {
1841  double im = x.imag ();
1842  double u = expm1 (x.real ());
1843  double v = sin (im/2);
1844  v = -2*v*v;
1845  retval = Complex (u*v + u + v, (u+1) * sin (im));
1846  }
1847  else
1848  retval = std::exp (x) - Complex (1);
1849 
1850  return retval;
1851 }
1852 
1855 {
1856  FloatComplex retval;
1857 
1858  if (std::abs (x) < 1)
1859  {
1860  float im = x.imag ();
1861  float u = expm1 (x.real ());
1862  float v = sin (im/2);
1863  v = -2*v*v;
1864  retval = FloatComplex (u*v + u + v, (u+1) * sin (im));
1865  }
1866  else
1867  retval = std::exp (x) - FloatComplex (1);
1868 
1869  return retval;
1870 }
1871 
1872 double
1873 gamma (double x)
1874 {
1875  double result;
1876 
1877  // Special cases for (near) compatibility with Matlab instead of tgamma.
1878  // Matlab does not have -0.
1879 
1880  if (x == 0)
1881  result = (math::negative_sign (x)
1884  else if ((x < 0 && math::x_nint (x) == x)
1885  || math::isinf (x))
1886  result = numeric_limits<double>::Inf ();
1887  else if (math::isnan (x))
1888  result = numeric_limits<double>::NaN ();
1889  else
1890  result = std::tgamma (x);
1891 
1892  return result;
1893 }
1894 
1895 float
1896 gamma (float x)
1897 {
1898  float result;
1899 
1900  // Special cases for (near) compatibility with Matlab instead of tgamma.
1901  // Matlab does not have -0.
1902 
1903  if (x == 0)
1904  result = (math::negative_sign (x)
1907  else if ((x < 0 && math::x_nint (x) == x)
1908  || math::isinf (x))
1909  result = numeric_limits<float>::Inf ();
1910  else if (math::isnan (x))
1911  result = numeric_limits<float>::NaN ();
1912  else
1913  result = std::tgammaf (x);
1914 
1915  return result;
1916 }
1917 
1918 Complex
1919 log1p (const Complex& x)
1920 {
1921  Complex retval;
1922 
1923  double r = x.real (), i = x.imag ();
1924 
1925  if (fabs (r) < 0.5 && fabs (i) < 0.5)
1926  {
1927  double u = 2*r + r*r + i*i;
1928  retval = Complex (log1p (u / (1+std::sqrt (u+1))),
1929  atan2 (i, 1 + r));
1930  }
1931  else
1932  retval = std::log (Complex (1) + x);
1933 
1934  return retval;
1935 }
1936 
1939 {
1940  FloatComplex retval;
1941 
1942  float r = x.real (), i = x.imag ();
1943 
1944  if (fabs (r) < 0.5 && fabs (i) < 0.5)
1945  {
1946  float u = 2*r + r*r + i*i;
1947  retval = FloatComplex (log1p (u / (1+std::sqrt (u+1))),
1948  atan2 (i, 1 + r));
1949  }
1950  else
1951  retval = std::log (FloatComplex (1) + x);
1952 
1953  return retval;
1954 }
1955 
1956 static const double pi = 3.14159265358979323846;
1957 
1958 template <typename T>
1959 static inline T
1960 xlog (const T& x)
1961 {
1962  return log (x);
1963 }
1964 
1965 template <>
1966 inline double
1967 xlog (const double& x)
1968 {
1969  return std::log (x);
1970 }
1971 
1972 template <>
1973 inline float
1974 xlog (const float& x)
1975 {
1976  return std::log (x);
1977 }
1978 
1979 template <typename T>
1980 static T
1981 lanczos_approximation_psi (const T zc)
1982 {
1983  // Coefficients for C.Lanczos expansion of psi function from XLiFE++
1984  // gammaFunctions psi_coef[k] = - (2k+1) * lg_coef[k] (see melina++
1985  // gamma functions -1/12, 3/360,-5/1260, 7/1680,-9/1188,
1986  // 11*691/360360,-13/156, 15*3617/122400, ? , ?
1987  static const T dg_coeff[10] =
1988  {
1989  -0.83333333333333333e-1, 0.83333333333333333e-2,
1990  -0.39682539682539683e-2, 0.41666666666666667e-2,
1991  -0.75757575757575758e-2, 0.21092796092796093e-1,
1992  -0.83333333333333333e-1, 0.4432598039215686,
1993  -0.3053954330270122e+1, 0.125318899521531e+2
1994  };
1995 
1996  T overz2 = T (1.0) / (zc * zc);
1997  T overz2k = overz2;
1998 
1999  T p = 0;
2000  for (octave_idx_type k = 0; k < 10; k++, overz2k *= overz2)
2001  p += dg_coeff[k] * overz2k;
2002  p += xlog (zc) - T (0.5) / zc;
2003  return p;
2004 }
2005 
2006 template <typename T>
2007 T
2008 xpsi (T z)
2009 {
2010  static const double euler_mascheroni
2011  = 0.577215664901532860606512090082402431042;
2012 
2013  const bool is_int = (std::floor (z) == z);
2014 
2015  T p = 0;
2016  if (z <= 0)
2017  {
2018  // limits - zeros of the gamma function
2019  if (is_int)
2020  p = -numeric_limits<T>::Inf (); // Matlab returns -Inf for psi (0)
2021  else
2022  // Abramowitz and Stegun, page 259, eq 6.3.7
2023  p = psi (1 - z) - (pi / tan (pi * z));
2024  }
2025  else if (is_int)
2026  {
2027  // Abramowitz and Stegun, page 258, eq 6.3.2
2028  p = - euler_mascheroni;
2029  for (octave_idx_type k = z - 1; k > 0; k--)
2030  p += 1.0 / k;
2031  }
2032  else if (std::floor (z + 0.5) == z + 0.5)
2033  {
2034  // Abramowitz and Stegun, page 258, eq 6.3.3 and 6.3.4
2035  for (octave_idx_type k = z; k > 0; k--)
2036  p += 1.0 / (2 * k - 1);
2037 
2038  p = - euler_mascheroni - 2 * std::log (2) + 2 * (p);
2039  }
2040  else
2041  {
2042  // adapted from XLiFE++ gammaFunctions
2043 
2044  T zc = z;
2045  // Use formula for derivative of LogGamma(z)
2046  if (z < 10)
2047  {
2048  const signed char n = 10 - z;
2049  for (signed char k = n - 1; k >= 0; k--)
2050  p -= 1.0 / (k + z);
2051  zc += n;
2052  }
2053  p += lanczos_approximation_psi (zc);
2054  }
2055 
2056  return p;
2057 }
2058 
2059 // explicit instantiations
2060 double
2061 psi (double z) { return xpsi (z); }
2062 float
2063 psi (float z) { return xpsi (z); }
2064 
2065 template <typename T>
2066 std::complex<T>
2067 xpsi (const std::complex<T>& z)
2068 {
2069  // adapted from XLiFE++ gammaFunctions
2070 
2071  typedef typename std::complex<T>::value_type P;
2072 
2073  P z_r = z.real ();
2074  P z_ra = z_r;
2075 
2076  std::complex<T> dgam (0.0, 0.0);
2077  if (z.imag () == 0)
2078  dgam = std::complex<T> (psi (z_r), 0.0);
2079  else if (z_r < 0)
2080  dgam = psi (P (1.0) - z)- (P (pi) / tan (P (pi) * z));
2081  else
2082  {
2083  // Use formula for derivative of LogGamma(z)
2084  std::complex<T> z_m = z;
2085  if (z_ra < 8)
2086  {
2087  unsigned char n = 8 - z_ra;
2088  z_m = z + std::complex<T> (n, 0.0);
2089 
2090  // Recurrence formula. For | Re(z) | < 8, use recursively
2091  //
2092  // DiGamma(z) = DiGamma(z+1) - 1/z
2093  std::complex<T> z_p = z + P (n - 1);
2094  for (unsigned char k = n; k > 0; k--, z_p -= 1.0)
2095  dgam -= P (1.0) / z_p;
2096  }
2097 
2098  // for | Re(z) | > 8, use derivative of C.Lanczos expansion for
2099  // LogGamma
2100  //
2101  // psi(z) = log(z) - 1/(2z) - 1/12z^2 + 3/360z^4 - 5/1260z^6
2102  // + 7/1680z^8 - 9/1188z^10 + ...
2103  //
2104  // (Abramowitz&Stegun, page 259, formula 6.3.18
2105  dgam += lanczos_approximation_psi (z_m);
2106  }
2107  return dgam;
2108 }
2109 
2110 // explicit instantiations
2111 Complex
2112 psi (const Complex& z) { return xpsi (z); }
2114 psi (const FloatComplex& z) { return xpsi (z); }
2115 
2116 template <typename T>
2117 static inline void
2118 fortran_psifn (T z, octave_idx_type n, T& ans, octave_idx_type& ierr);
2119 
2120 template <>
2121 inline void
2122 fortran_psifn<double> (double z, octave_idx_type n_arg,
2123  double& ans, octave_idx_type& ierr)
2124 {
2125  F77_INT n = to_f77_int (n_arg);
2126  F77_INT flag = 0;
2127  F77_INT t_ierr;
2128  F77_XFCN (dpsifn, DPSIFN, (z, n, 1, 1, ans, flag, t_ierr));
2129  ierr = t_ierr;
2130 }
2131 
2132 template <>
2133 inline void
2134 fortran_psifn<float> (float z, octave_idx_type n_arg,
2135  float& ans, octave_idx_type& ierr)
2136 {
2137  F77_INT n = to_f77_int (n_arg);
2138  F77_INT flag = 0;
2139  F77_INT t_ierr;
2140  F77_XFCN (psifn, PSIFN, (z, n, 1, 1, ans, flag, t_ierr));
2141  ierr = t_ierr;
2142 }
2143 
2144 template <typename T>
2145 T
2147 {
2148  T ans;
2149  octave_idx_type ierr = 0;
2150  fortran_psifn<T> (z, n, ans, ierr);
2151  if (ierr == 0)
2152  {
2153  // Remember that psifn and dpsifn return scales values
2154  // When n is 1: do nothing since ((-1)**(n+1)/gamma(n+1)) == 1
2155  // When n is 0: change sign since ((-1)**(n+1)/gamma(n+1)) == -1
2156  if (n > 1)
2157  // FIXME: xgamma here is a killer for our precision since it grows
2158  // way too fast.
2159  ans = ans / (std::pow (-1.0, n + 1) / gamma (double (n+1)));
2160  else if (n == 0)
2161  ans = -ans;
2162  }
2163  else if (ierr == 2)
2164  ans = - numeric_limits<T>::Inf ();
2165  else // we probably never get here
2166  ans = numeric_limits<T>::NaN ();
2167 
2168  return ans;
2169 }
2170 
2171 double
2172 psi (octave_idx_type n, double z) { return xpsi (n, z); }
2173 float
2174 psi (octave_idx_type n, float z) { return xpsi (n, z); }
2175 
2176 Complex
2177 rc_lgamma (double x)
2178 {
2179  double result;
2180 
2181 #if defined (HAVE_LGAMMA_R)
2182  int sgngam;
2183  result = lgamma_r (x, &sgngam);
2184 #else
2185  result = std::lgamma (x);
2186  int sgngam = signgam;
2187 #endif
2188 
2189  if (sgngam < 0)
2190  return result + Complex (0., M_PI);
2191  else
2192  return result;
2193 }
2194 
2196 rc_lgamma (float x)
2197 {
2198  float result;
2199 
2200 #if defined (HAVE_LGAMMAF_R)
2201  int sgngam;
2202  result = lgammaf_r (x, &sgngam);
2203 #else
2204  result = std::lgammaf (x);
2205  int sgngam = signgam;
2206 #endif
2207 
2208  if (sgngam < 0)
2209  return result + FloatComplex (0., M_PI);
2210  else
2211  return result;
2212 }
2213 
2214 Complex
2215 rc_log1p (double x)
2216 {
2217  return (x < -1.0
2218  ? Complex (std::log (-(1.0 + x)), M_PI)
2219  : Complex (log1p (x)));
2220 }
2221 
2223 rc_log1p (float x)
2224 {
2225  return (x < -1.0f
2226  ? FloatComplex (std::log (-(1.0f + x)), M_PI)
2227  : FloatComplex (log1p (x)));
2228 }
2229 
2230 OCTAVE_END_NAMESPACE(math)
2231 OCTAVE_END_NAMESPACE(octave)
#define Inf
Definition: Faddeeva.cc:260
#define NaN
Definition: Faddeeva.cc:261
octave_value_list do_bessel(enum bessel_type type, const char *fcn, const octave_value_list &args, int nargout)
Definition: besselj.cc:86
subroutine cairy(Z, ID, KODE, AI, NZ, IERR)
Definition: cairy.f:2
subroutine cbesh(Z, FNU, KODE, M, N, CY, NZ, IERR)
Definition: cbesh.f:2
subroutine cbesi(Z, FNU, KODE, N, CY, NZ, IERR)
Definition: cbesi.f:2
subroutine cbesj(Z, FNU, KODE, N, CY, NZ, IERR)
Definition: cbesj.f:2
subroutine cbesk(Z, FNU, KODE, N, CY, NZ, IERR)
Definition: cbesk.f:2
subroutine cbesy(Z, FNU, KODE, N, CY, NZ, CWRK, IERR)
Definition: cbesy.f:2
subroutine cbiry(Z, ID, KODE, BI, IERR)
Definition: cbiry.f:2
octave_idx_type rows() const
Definition: Array.h:459
void resize(const dim_vector &dv, const T &rfv)
Size of the specified dimension.
Definition: Array-base.cc:1023
octave_idx_type cols() const
Definition: Array.h:469
const dim_vector & dims() const
Return a const-reference so that dims ()(i) works efficiently.
Definition: Array.h:503
octave_idx_type numel() const
Number of elements in the array.
Definition: Array.h:414
void resize(octave_idx_type nr, octave_idx_type nc, const Complex &rfv=Complex(0))
Definition: CMatrix.h:193
void resize(octave_idx_type nr, octave_idx_type nc, const FloatComplex &rfv=FloatComplex(0))
Definition: fCMatrix.h:201
Definition: dMatrix.h:42
Vector representing the dimensions (size) of an Array.
Definition: dim-vector.h:94
octave_idx_type numel(int n=0) const
Number of elements that a matrix with this dimensions would have.
Definition: dim-vector.h:335
OCTAVE_BEGIN_NAMESPACE(octave) static octave_value daspk_fcn
subroutine dpsifn(X, N, KODE, M, ANS, NZ, IERR)
Definition: dpsifn.f:3
#define F77_CONST_CMPLX_ARG(x)
Definition: f77-fcn.h:313
#define F77_CMPLX_ARG(x)
Definition: f77-fcn.h:310
#define F77_XFCN(f, F, args)
Definition: f77-fcn.h:45
octave_f77_int_type F77_INT
Definition: f77-fcn.h:306
double lo_ieee_nan_value()
Definition: lo-ieee.cc:84
Complex asin(const Complex &x)
Definition: lo-mappers.cc:107
bool negative_sign(double x)
Definition: lo-mappers.cc:181
bool isinf(double x)
Definition: lo-mappers.h:203
double signum(double x)
Definition: lo-mappers.h:229
bool isnan(bool)
Definition: lo-mappers.h:178
std::complex< T > floor(const std::complex< T > &x)
Definition: lo-mappers.h:130
T x_nint(T x)
Definition: lo-mappers.h:269
F77_RET_T const F77_DBLE const F77_DBLE F77_DBLE * d
F77_RET_T const F77_DBLE const F77_DBLE F77_DBLE const F77_INT F77_INT & ierr
F77_RET_T const F77_DBLE * x
F77_RET_T const F77_DBLE const F77_DBLE * f
Complex(* dptr)(const Complex &, double, int, octave_idx_type &)
Definition: lo-specfun.cc:563
void fortran_psifn< float >(float z, octave_idx_type n_arg, float &ans, octave_idx_type &ierr)
Definition: lo-specfun.cc:2134
Complex besselj(double alpha, const Complex &x, bool scaled, octave_idx_type &ierr)
Definition: lo-specfun.cc:788
double erfinv(double x)
Definition: lo-specfun.cc:1823
double erfi(double x)
Definition: lo-specfun.cc:1731
Complex besselh2(double alpha, const Complex &x, bool scaled, octave_idx_type &ierr)
Definition: lo-specfun.cc:793
double gamma(double x)
Definition: lo-specfun.cc:1873
double dawson(double x)
Definition: lo-specfun.cc:1467
Complex erf(const Complex &x)
Definition: lo-specfun.cc:1588
void fortran_psifn< double >(double z, octave_idx_type n_arg, double &ans, octave_idx_type &ierr)
Definition: lo-specfun.cc:2122
double erfcinv(double x)
Definition: lo-specfun.cc:1697
Complex rc_log1p(double x)
Definition: lo-specfun.cc:2215
T xpsi(T z)
Definition: lo-specfun.cc:2008
Complex bessely(double alpha, const Complex &x, bool scaled, octave_idx_type &ierr)
Definition: lo-specfun.cc:789
Complex besselk(double alpha, const Complex &x, bool scaled, octave_idx_type &ierr)
Definition: lo-specfun.cc:791
Complex besselh1(double alpha, const Complex &x, bool scaled, octave_idx_type &ierr)
Definition: lo-specfun.cc:792
Complex rc_lgamma(double x)
Definition: lo-specfun.cc:2177
#define ALL_BESSEL(name, fcn)
Definition: lo-specfun.cc:1323
double xlog(const double &x)
Definition: lo-specfun.cc:1967
FloatComplex(* fptr)(const FloatComplex &, float, int, octave_idx_type &)
Definition: lo-specfun.cc:1102
Complex airy(const Complex &z, bool deriv, bool scaled, octave_idx_type &ierr)
Definition: lo-specfun.cc:129
double erfcx(double x)
Definition: lo-specfun.cc:1710
Complex expm1(const Complex &x)
Definition: lo-specfun.cc:1835
void ellipj(double u, double m, double &sn, double &cn, double &dn, double &err)
Definition: lo-specfun.cc:1487
Complex besseli(double alpha, const Complex &x, bool scaled, octave_idx_type &ierr)
Definition: lo-specfun.cc:790
Complex biry(const Complex &z, bool deriv, bool scaled, octave_idx_type &ierr)
Definition: lo-specfun.cc:1351
Complex log1p(const Complex &x)
Definition: lo-specfun.cc:1919
Complex erfc(const Complex &x)
Definition: lo-specfun.cc:1603
double psi(double z)
Definition: lo-specfun.cc:2061
double lgamma(double x)
Definition: lo-specfun.h:336
T octave_idx_type m
Definition: mx-inlines.cc:781
octave_idx_type n
Definition: mx-inlines.cc:761
T * r
Definition: mx-inlines.cc:781
std::complex< double > w(std::complex< double > z, double relerr=0)
std::complex< double > erfc(std::complex< double > z, double relerr=0)
std::complex< double > erfcx(std::complex< double > z, double relerr=0)
std::complex< double > erfi(std::complex< double > z, double relerr=0)
std::complex< double > erf(std::complex< double > z, double relerr=0)
std::complex< double > Dawson(std::complex< double > z, double relerr=0)
std::complex< double > Complex
Definition: oct-cmplx.h:33
std::complex< float > FloatComplex
Definition: oct-cmplx.h:34
octave_int< T > pow(const octave_int< T > &a, const octave_int< T > &b)
subroutine psifn(X, N, KODE, M, ANS, NZ, IERR)
Definition: psifn.f:3
F77_RET_T F77_FUNC(xerbla, XERBLA)(F77_CONST_CHAR_ARG_DEF(s_arg
subroutine zairy(ZR, ZI, ID, KODE, AIR, AII, NZ, IERR)
Definition: zairy.f:2
subroutine zbesh(ZR, ZI, FNU, KODE, M, N, CYR, CYI, NZ, IERR)
Definition: zbesh.f:2
subroutine zbesi(ZR, ZI, FNU, KODE, N, CYR, CYI, NZ, IERR)
Definition: zbesi.f:2
subroutine zbesj(ZR, ZI, FNU, KODE, N, CYR, CYI, NZ, IERR)
Definition: zbesj.f:2
subroutine zbesk(ZR, ZI, FNU, KODE, N, CYR, CYI, NZ, IERR)
Definition: zbesk.f:2
subroutine zbesy(ZR, ZI, FNU, KODE, N, CYR, CYI, NZ, CWRKR, CWRKI, IERR)
Definition: zbesy.f:3
subroutine zbiry(ZR, ZI, ID, KODE, BIR, BII, IERR)
Definition: zbiry.f:2