annotate src/DLD-FUNCTIONS/besselj.cc @ 11904:059fadc0cbc3 release-3-0-x

fix scaling factor for negative alpha in zbesi,cbesi add bessel function tests add all scaling factors to bessel documentation
author Brian Gough <bjg@gnu.org>
date Mon, 12 Jan 2009 10:42:03 +0100
parents 120f3135952f
children
Ignore whitespace changes - Everywhere: Within whitespace: At end of lines:
rev   line source
3155
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
1 /*
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
2
7017
a1dbe9d80eee [project @ 2007-10-12 21:27:11 by jwe]
jwe
parents: 7016
diff changeset
3 Copyright (C) 1997, 1998, 1999, 2000, 2001, 2002, 2004, 2005, 2006,
a1dbe9d80eee [project @ 2007-10-12 21:27:11 by jwe]
jwe
parents: 7016
diff changeset
4 2007 John W. Eaton
3155
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
5
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
6 This file is part of Octave.
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
7
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
8 Octave is free software; you can redistribute it and/or modify it
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
9 under the terms of the GNU General Public License as published by the
7016
93c65f2a5668 [project @ 2007-10-12 06:40:56 by jwe]
jwe
parents: 5823
diff changeset
10 Free Software Foundation; either version 3 of the License, or (at your
93c65f2a5668 [project @ 2007-10-12 06:40:56 by jwe]
jwe
parents: 5823
diff changeset
11 option) any later version.
3155
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
12
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
13 Octave is distributed in the hope that it will be useful, but WITHOUT
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
14 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
15 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
16 for more details.
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
17
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
18 You should have received a copy of the GNU General Public License
7016
93c65f2a5668 [project @ 2007-10-12 06:40:56 by jwe]
jwe
parents: 5823
diff changeset
19 along with Octave; see the file COPYING. If not, see
93c65f2a5668 [project @ 2007-10-12 06:40:56 by jwe]
jwe
parents: 5823
diff changeset
20 <http://www.gnu.org/licenses/>.
3155
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
21
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
22 */
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
23
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
24 #ifdef HAVE_CONFIG_H
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
25 #include <config.h>
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
26 #endif
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
27
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
28 #include "lo-specfun.h"
4153
6b96ce9f5743 [project @ 2002-11-06 20:38:49 by jwe]
jwe
parents: 3810
diff changeset
29 #include "quit.h"
3155
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
30
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
31 #include "defun-dld.h"
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
32 #include "error.h"
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
33 #include "gripes.h"
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
34 #include "oct-obj.h"
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
35 #include "utils.h"
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
36
3220
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
37 enum bessel_type
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
38 {
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
39 BESSEL_J,
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
40 BESSEL_Y,
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
41 BESSEL_I,
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
42 BESSEL_K,
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
43 BESSEL_H1,
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
44 BESSEL_H2
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
45 };
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
46
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
47 #define DO_BESSEL(type, alpha, x, scaled, ierr, result) \
3155
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
48 do \
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
49 { \
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
50 switch (type) \
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
51 { \
3220
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
52 case BESSEL_J: \
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
53 result = besselj (alpha, x, scaled, ierr); \
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
54 break; \
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
55 \
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
56 case BESSEL_Y: \
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
57 result = bessely (alpha, x, scaled, ierr); \
3155
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
58 break; \
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
59 \
3220
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
60 case BESSEL_I: \
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
61 result = besseli (alpha, x, scaled, ierr); \
3155
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
62 break; \
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
63 \
3220
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
64 case BESSEL_K: \
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
65 result = besselk (alpha, x, scaled, ierr); \
3155
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
66 break; \
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
67 \
3220
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
68 case BESSEL_H1: \
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
69 result = besselh1 (alpha, x, scaled, ierr); \
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
70 break; \
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
71 \
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
72 case BESSEL_H2: \
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
73 result = besselh2 (alpha, x, scaled, ierr); \
3155
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
74 break; \
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
75 \
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
76 default: \
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
77 break; \
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
78 } \
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
79 } \
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
80 while (0)
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
81
3220
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
82 static inline Matrix
5275
23b37da9fd5b [project @ 2005-04-08 16:07:35 by jwe]
jwe
parents: 4844
diff changeset
83 int_array2_to_matrix (const Array2<octave_idx_type>& a)
3220
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
84 {
5275
23b37da9fd5b [project @ 2005-04-08 16:07:35 by jwe]
jwe
parents: 4844
diff changeset
85 octave_idx_type nr = a.rows ();
23b37da9fd5b [project @ 2005-04-08 16:07:35 by jwe]
jwe
parents: 4844
diff changeset
86 octave_idx_type nc = a.cols ();
3220
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
87
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
88 Matrix retval (nr, nc);
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
89
5275
23b37da9fd5b [project @ 2005-04-08 16:07:35 by jwe]
jwe
parents: 4844
diff changeset
90 for (octave_idx_type j = 0; j < nc; j++)
23b37da9fd5b [project @ 2005-04-08 16:07:35 by jwe]
jwe
parents: 4844
diff changeset
91 for (octave_idx_type i = 0; i < nr; i++)
4153
6b96ce9f5743 [project @ 2002-11-06 20:38:49 by jwe]
jwe
parents: 3810
diff changeset
92 {
6b96ce9f5743 [project @ 2002-11-06 20:38:49 by jwe]
jwe
parents: 3810
diff changeset
93 OCTAVE_QUIT;
6b96ce9f5743 [project @ 2002-11-06 20:38:49 by jwe]
jwe
parents: 3810
diff changeset
94
5760
8d7162924bd3 [project @ 2006-04-14 04:01:37 by jwe]
jwe
parents: 5448
diff changeset
95 retval(i,j) = static_cast<double> (a(i,j));
4153
6b96ce9f5743 [project @ 2002-11-06 20:38:49 by jwe]
jwe
parents: 3810
diff changeset
96 }
3220
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
97
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
98 return retval;
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
99 }
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
100
4844
9f7ef92b50b0 [project @ 2004-04-02 17:26:53 by jwe]
jwe
parents: 4153
diff changeset
101 static inline NDArray
5760
8d7162924bd3 [project @ 2006-04-14 04:01:37 by jwe]
jwe
parents: 5448
diff changeset
102 int_arrayN_to_array (const ArrayN<octave_idx_type>& a)
4844
9f7ef92b50b0 [project @ 2004-04-02 17:26:53 by jwe]
jwe
parents: 4153
diff changeset
103 {
9f7ef92b50b0 [project @ 2004-04-02 17:26:53 by jwe]
jwe
parents: 4153
diff changeset
104 dim_vector dv = a.dims ();
9f7ef92b50b0 [project @ 2004-04-02 17:26:53 by jwe]
jwe
parents: 4153
diff changeset
105 int nel = dv.numel ();
9f7ef92b50b0 [project @ 2004-04-02 17:26:53 by jwe]
jwe
parents: 4153
diff changeset
106
9f7ef92b50b0 [project @ 2004-04-02 17:26:53 by jwe]
jwe
parents: 4153
diff changeset
107 NDArray retval (dv);
9f7ef92b50b0 [project @ 2004-04-02 17:26:53 by jwe]
jwe
parents: 4153
diff changeset
108
9f7ef92b50b0 [project @ 2004-04-02 17:26:53 by jwe]
jwe
parents: 4153
diff changeset
109 for (int i = 0; i < nel; i++)
9f7ef92b50b0 [project @ 2004-04-02 17:26:53 by jwe]
jwe
parents: 4153
diff changeset
110 {
9f7ef92b50b0 [project @ 2004-04-02 17:26:53 by jwe]
jwe
parents: 4153
diff changeset
111 OCTAVE_QUIT;
9f7ef92b50b0 [project @ 2004-04-02 17:26:53 by jwe]
jwe
parents: 4153
diff changeset
112
5760
8d7162924bd3 [project @ 2006-04-14 04:01:37 by jwe]
jwe
parents: 5448
diff changeset
113 retval(i) = static_cast<double> (a(i));
4844
9f7ef92b50b0 [project @ 2004-04-02 17:26:53 by jwe]
jwe
parents: 4153
diff changeset
114 }
9f7ef92b50b0 [project @ 2004-04-02 17:26:53 by jwe]
jwe
parents: 4153
diff changeset
115
9f7ef92b50b0 [project @ 2004-04-02 17:26:53 by jwe]
jwe
parents: 4153
diff changeset
116 return retval;
9f7ef92b50b0 [project @ 2004-04-02 17:26:53 by jwe]
jwe
parents: 4153
diff changeset
117 }
9f7ef92b50b0 [project @ 2004-04-02 17:26:53 by jwe]
jwe
parents: 4153
diff changeset
118
3155
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
119 static void
3220
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
120 gripe_bessel_arg (const char *fn, const char *arg)
3155
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
121 {
3220
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
122 error ("%s: expecting scalar or matrix as %s argument", fn, arg);
3155
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
123 }
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
124
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
125 octave_value_list
3220
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
126 do_bessel (enum bessel_type type, const char *fn,
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
127 const octave_value_list& args, int nargout)
3155
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
128 {
3220
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
129 octave_value_list retval;
3155
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
130
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
131 int nargin = args.length ();
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
132
3220
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
133 if (nargin == 2 || nargin == 3)
3155
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
134 {
3220
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
135 bool scaled = (nargin == 3);
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
136
3155
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
137 octave_value alpha_arg = args(0);
3220
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
138 octave_value x_arg = args(1);
3155
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
139
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
140 if (alpha_arg.is_scalar_type ())
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
141 {
3220
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
142 double alpha = args(0).double_value ();
3155
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
143
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
144 if (! error_state)
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
145 {
3220
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
146 if (x_arg.is_scalar_type ())
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
147 {
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
148 Complex x = x_arg.complex_value ();
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
149
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
150 if (! error_state)
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
151 {
5275
23b37da9fd5b [project @ 2005-04-08 16:07:35 by jwe]
jwe
parents: 4844
diff changeset
152 octave_idx_type ierr;
3220
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
153 octave_value result;
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
154
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
155 DO_BESSEL (type, alpha, x, scaled, ierr, result);
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
156
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
157 if (nargout > 1)
5760
8d7162924bd3 [project @ 2006-04-14 04:01:37 by jwe]
jwe
parents: 5448
diff changeset
158 retval(1) = static_cast<double> (ierr);
3155
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
159
3220
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
160 retval(0) = result;
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
161 }
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
162 else
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
163 gripe_bessel_arg (fn, "second");
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
164 }
3155
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
165 else
3220
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
166 {
4844
9f7ef92b50b0 [project @ 2004-04-02 17:26:53 by jwe]
jwe
parents: 4153
diff changeset
167 ComplexNDArray x = x_arg.complex_array_value ();
3220
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
168
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
169 if (! error_state)
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
170 {
5275
23b37da9fd5b [project @ 2005-04-08 16:07:35 by jwe]
jwe
parents: 4844
diff changeset
171 ArrayN<octave_idx_type> ierr;
3220
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
172 octave_value result;
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
173
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
174 DO_BESSEL (type, alpha, x, scaled, ierr, result);
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
175
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
176 if (nargout > 1)
4844
9f7ef92b50b0 [project @ 2004-04-02 17:26:53 by jwe]
jwe
parents: 4153
diff changeset
177 retval(1) = int_arrayN_to_array (ierr);
3220
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
178
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
179 retval(0) = result;
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
180 }
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
181 else
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
182 gripe_bessel_arg (fn, "second");
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
183 }
3155
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
184 }
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
185 else
3220
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
186 gripe_bessel_arg (fn, "first");
3155
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
187 }
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
188 else
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
189 {
4844
9f7ef92b50b0 [project @ 2004-04-02 17:26:53 by jwe]
jwe
parents: 4153
diff changeset
190 dim_vector dv0 = args(0).dims ();
9f7ef92b50b0 [project @ 2004-04-02 17:26:53 by jwe]
jwe
parents: 4153
diff changeset
191 dim_vector dv1 = args(1).dims ();
9f7ef92b50b0 [project @ 2004-04-02 17:26:53 by jwe]
jwe
parents: 4153
diff changeset
192
9f7ef92b50b0 [project @ 2004-04-02 17:26:53 by jwe]
jwe
parents: 4153
diff changeset
193 bool args0_is_row_vector = (dv0 (1) == dv0.numel ());
9f7ef92b50b0 [project @ 2004-04-02 17:26:53 by jwe]
jwe
parents: 4153
diff changeset
194 bool args1_is_col_vector = (dv1 (0) == dv1.numel ());
3155
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
195
4844
9f7ef92b50b0 [project @ 2004-04-02 17:26:53 by jwe]
jwe
parents: 4153
diff changeset
196 if (args0_is_row_vector && args1_is_col_vector)
3155
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
197 {
4844
9f7ef92b50b0 [project @ 2004-04-02 17:26:53 by jwe]
jwe
parents: 4153
diff changeset
198 RowVector ralpha = args(0).row_vector_value ();
9f7ef92b50b0 [project @ 2004-04-02 17:26:53 by jwe]
jwe
parents: 4153
diff changeset
199
9f7ef92b50b0 [project @ 2004-04-02 17:26:53 by jwe]
jwe
parents: 4153
diff changeset
200 if (! error_state)
3220
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
201 {
4844
9f7ef92b50b0 [project @ 2004-04-02 17:26:53 by jwe]
jwe
parents: 4153
diff changeset
202 ComplexColumnVector cx =
9f7ef92b50b0 [project @ 2004-04-02 17:26:53 by jwe]
jwe
parents: 4153
diff changeset
203 x_arg.complex_column_vector_value ();
3220
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
204
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
205 if (! error_state)
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
206 {
5275
23b37da9fd5b [project @ 2005-04-08 16:07:35 by jwe]
jwe
parents: 4844
diff changeset
207 Array2<octave_idx_type> ierr;
3220
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
208 octave_value result;
3155
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
209
4844
9f7ef92b50b0 [project @ 2004-04-02 17:26:53 by jwe]
jwe
parents: 4153
diff changeset
210 DO_BESSEL (type, ralpha, cx, scaled, ierr, result);
9f7ef92b50b0 [project @ 2004-04-02 17:26:53 by jwe]
jwe
parents: 4153
diff changeset
211
3220
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
212 if (nargout > 1)
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
213 retval(1) = int_array2_to_matrix (ierr);
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
214
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
215 retval(0) = result;
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
216 }
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
217 else
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
218 gripe_bessel_arg (fn, "second");
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
219 }
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
220 else
4844
9f7ef92b50b0 [project @ 2004-04-02 17:26:53 by jwe]
jwe
parents: 4153
diff changeset
221 gripe_bessel_arg (fn, "first");
9f7ef92b50b0 [project @ 2004-04-02 17:26:53 by jwe]
jwe
parents: 4153
diff changeset
222 }
9f7ef92b50b0 [project @ 2004-04-02 17:26:53 by jwe]
jwe
parents: 4153
diff changeset
223 else
9f7ef92b50b0 [project @ 2004-04-02 17:26:53 by jwe]
jwe
parents: 4153
diff changeset
224 {
9f7ef92b50b0 [project @ 2004-04-02 17:26:53 by jwe]
jwe
parents: 4153
diff changeset
225 NDArray alpha = args(0).array_value ();
3220
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
226
4844
9f7ef92b50b0 [project @ 2004-04-02 17:26:53 by jwe]
jwe
parents: 4153
diff changeset
227 if (! error_state)
9f7ef92b50b0 [project @ 2004-04-02 17:26:53 by jwe]
jwe
parents: 4153
diff changeset
228 {
9f7ef92b50b0 [project @ 2004-04-02 17:26:53 by jwe]
jwe
parents: 4153
diff changeset
229 if (x_arg.is_scalar_type ())
9f7ef92b50b0 [project @ 2004-04-02 17:26:53 by jwe]
jwe
parents: 4153
diff changeset
230 {
9f7ef92b50b0 [project @ 2004-04-02 17:26:53 by jwe]
jwe
parents: 4153
diff changeset
231 Complex x = x_arg.complex_value ();
3155
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
232
4844
9f7ef92b50b0 [project @ 2004-04-02 17:26:53 by jwe]
jwe
parents: 4153
diff changeset
233 if (! error_state)
3155
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
234 {
5275
23b37da9fd5b [project @ 2005-04-08 16:07:35 by jwe]
jwe
parents: 4844
diff changeset
235 ArrayN<octave_idx_type> ierr;
3220
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
236 octave_value result;
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
237
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
238 DO_BESSEL (type, alpha, x, scaled, ierr, result);
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
239
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
240 if (nargout > 1)
4844
9f7ef92b50b0 [project @ 2004-04-02 17:26:53 by jwe]
jwe
parents: 4153
diff changeset
241 retval(1) = int_arrayN_to_array (ierr);
3220
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
242
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
243 retval(0) = result;
3155
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
244 }
4844
9f7ef92b50b0 [project @ 2004-04-02 17:26:53 by jwe]
jwe
parents: 4153
diff changeset
245 else
9f7ef92b50b0 [project @ 2004-04-02 17:26:53 by jwe]
jwe
parents: 4153
diff changeset
246 gripe_bessel_arg (fn, "second");
3155
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
247 }
3220
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
248 else
4844
9f7ef92b50b0 [project @ 2004-04-02 17:26:53 by jwe]
jwe
parents: 4153
diff changeset
249 {
9f7ef92b50b0 [project @ 2004-04-02 17:26:53 by jwe]
jwe
parents: 4153
diff changeset
250 ComplexNDArray x = x_arg.complex_array_value ();
9f7ef92b50b0 [project @ 2004-04-02 17:26:53 by jwe]
jwe
parents: 4153
diff changeset
251
9f7ef92b50b0 [project @ 2004-04-02 17:26:53 by jwe]
jwe
parents: 4153
diff changeset
252 if (! error_state)
9f7ef92b50b0 [project @ 2004-04-02 17:26:53 by jwe]
jwe
parents: 4153
diff changeset
253 {
5275
23b37da9fd5b [project @ 2005-04-08 16:07:35 by jwe]
jwe
parents: 4844
diff changeset
254 ArrayN<octave_idx_type> ierr;
4844
9f7ef92b50b0 [project @ 2004-04-02 17:26:53 by jwe]
jwe
parents: 4153
diff changeset
255 octave_value result;
9f7ef92b50b0 [project @ 2004-04-02 17:26:53 by jwe]
jwe
parents: 4153
diff changeset
256
9f7ef92b50b0 [project @ 2004-04-02 17:26:53 by jwe]
jwe
parents: 4153
diff changeset
257 DO_BESSEL (type, alpha, x, scaled, ierr, result);
9f7ef92b50b0 [project @ 2004-04-02 17:26:53 by jwe]
jwe
parents: 4153
diff changeset
258
9f7ef92b50b0 [project @ 2004-04-02 17:26:53 by jwe]
jwe
parents: 4153
diff changeset
259 if (nargout > 1)
9f7ef92b50b0 [project @ 2004-04-02 17:26:53 by jwe]
jwe
parents: 4153
diff changeset
260 retval(1) = int_arrayN_to_array (ierr);
9f7ef92b50b0 [project @ 2004-04-02 17:26:53 by jwe]
jwe
parents: 4153
diff changeset
261
9f7ef92b50b0 [project @ 2004-04-02 17:26:53 by jwe]
jwe
parents: 4153
diff changeset
262 retval(0) = result;
9f7ef92b50b0 [project @ 2004-04-02 17:26:53 by jwe]
jwe
parents: 4153
diff changeset
263 }
9f7ef92b50b0 [project @ 2004-04-02 17:26:53 by jwe]
jwe
parents: 4153
diff changeset
264 else
9f7ef92b50b0 [project @ 2004-04-02 17:26:53 by jwe]
jwe
parents: 4153
diff changeset
265 gripe_bessel_arg (fn, "second");
9f7ef92b50b0 [project @ 2004-04-02 17:26:53 by jwe]
jwe
parents: 4153
diff changeset
266 }
3155
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
267 }
4844
9f7ef92b50b0 [project @ 2004-04-02 17:26:53 by jwe]
jwe
parents: 4153
diff changeset
268 else
9f7ef92b50b0 [project @ 2004-04-02 17:26:53 by jwe]
jwe
parents: 4153
diff changeset
269 gripe_bessel_arg (fn, "first");
3155
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
270 }
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
271 }
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
272 }
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
273 else
5823
080c08b192d8 [project @ 2006-05-19 05:32:17 by jwe]
jwe
parents: 5760
diff changeset
274 print_usage ();
3155
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
275
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
276 return retval;
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
277 }
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
278
3220
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
279 DEFUN_DLD (besselj, args, nargout,
3459
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
280 "-*- texinfo -*-\n\
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
281 @deftypefn {Loadable Function} {[@var{j}, @var{ierr}] =} besselj (@var{alpha}, @var{x}, @var{opt})\n\
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
282 @deftypefnx {Loadable Function} {[@var{y}, @var{ierr}] =} bessely (@var{alpha}, @var{x}, @var{opt})\n\
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
283 @deftypefnx {Loadable Function} {[@var{i}, @var{ierr}] =} besseli (@var{alpha}, @var{x}, @var{opt})\n\
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
284 @deftypefnx {Loadable Function} {[@var{k}, @var{ierr}] =} besselk (@var{alpha}, @var{x}, @var{opt})\n\
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
285 @deftypefnx {Loadable Function} {[@var{h}, @var{ierr}] =} besselh (@var{alpha}, @var{k}, @var{x}, @var{opt})\n\
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
286 Compute Bessel or Hankel functions of various kinds:\n\
3155
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
287 \n\
3459
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
288 @table @code\n\
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
289 @item besselj\n\
11904
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
290 Bessel functions of the first kind. If the argument @var{opt} is supplied, \n\
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
291 the result is multiplied by @code{exp(-abs(imag(x)))}.\n\
3459
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
292 @item bessely\n\
11904
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
293 Bessel functions of the second kind. If the argument @var{opt} is supplied,\n\
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
294 the result is multiplied by @code{exp(-abs(imag(x)))}.\n\
3459
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
295 @item besseli\n\
11904
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
296 Modified Bessel functions of the first kind. If the argument @var{opt} is supplied,\n\
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
297 the result is multiplied by @code{exp(-abs(real(x)))}.\n\
3459
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
298 @item besselk\n\
11904
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
299 Modified Bessel functions of the second kind. If the argument @var{opt} is supplied,\n\
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
300 the result is multiplied by @code{exp(x)}.\n\
3459
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
301 @item besselh\n\
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
302 Compute Hankel functions of the first (@var{k} = 1) or second (@var{k}\n\
11904
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
303 = 2) kind. If the argument @var{opt} is supplied, the result is multiplied by\n\
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
304 @code{exp (-I*@var{x})} for @var{k} = 1 or @code{exp (I*@var{x})} for\n\
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
305 @var{k} = 2.\n\
3459
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
306 @end table\n\
3220
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
307 \n\
3459
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
308 If @var{alpha} is a scalar, the result is the same size as @var{x}.\n\
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
309 If @var{x} is a scalar, the result is the same size as @var{alpha}.\n\
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
310 If @var{alpha} is a row vector and @var{x} is a column vector, the\n\
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
311 result is a matrix with @code{length (@var{x})} rows and\n\
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
312 @code{length (@var{alpha})} columns. Otherwise, @var{alpha} and\n\
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
313 @var{x} must conform and the result will be the same size.\n\
3155
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
314 \n\
3459
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
315 The value of @var{alpha} must be real. The value of @var{x} may be\n\
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
316 complex.\n\
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
317 \n\
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
318 If requested, @var{ierr} contains the following status information\n\
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
319 and is the same size as the result.\n\
3548
ab7fa5a8f23f [project @ 2000-02-03 01:17:15 by jwe]
jwe
parents: 3499
diff changeset
320 \n\
3459
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
321 @enumerate 0\n\
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
322 @item\n\
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
323 Normal return.\n\
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
324 @item\n\
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
325 Input error, return @code{NaN}.\n\
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
326 @item\n\
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
327 Overflow, return @code{Inf}.\n\
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
328 @item\n\
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
329 Loss of significance by argument reduction results in less than\n\
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
330 half of machine accuracy.\n\
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
331 @item\n\
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
332 Complete loss of significance by argument reduction, return @code{NaN}.\n\
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
333 @item\n\
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
334 Error---no computation, algorithm termination condition not met,\n\
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
335 return @code{NaN}.\n\
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
336 @end enumerate\n\
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
337 @end deftypefn")
3155
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
338 {
3220
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
339 return do_bessel (BESSEL_J, "besselj", args, nargout);
3155
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
340 }
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
341
3220
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
342 DEFUN_DLD (bessely, args, nargout,
3459
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
343 "-*- texinfo -*-\n\
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
344 @deftypefn {Loadable Function} {[@var{y}, @var{ierr}] =} bessely (@var{alpha}, @var{x}, @var{opt})\n\
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
345 See besselj.\n\
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
346 @end deftypefn")
3155
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
347 {
3220
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
348 return do_bessel (BESSEL_Y, "bessely", args, nargout);
3155
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
349 }
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
350
3220
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
351 DEFUN_DLD (besseli, args, nargout,
3459
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
352 "-*- texinfo -*-\n\
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
353 @deftypefn {Loadable Function} {[@var{i}, @var{ierr}] =} besseli (@var{alpha}, @var{x}, @var{opt})\n\
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
354 See besselj.\n\
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
355 @end deftypefn")
3155
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
356 {
3220
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
357 return do_bessel (BESSEL_I, "besseli", args, nargout);
3155
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
358 }
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
359
3220
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
360 DEFUN_DLD (besselk, args, nargout,
3459
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
361 "-*- texinfo -*-\n\
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
362 @deftypefn {Loadable Function} {[@var{k}, @var{ierr}] =} besselk (@var{alpha}, @var{x}, @var{opt})\n\
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
363 See besselj.\n\
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
364 @end deftypefn")
3220
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
365 {
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
366 return do_bessel (BESSEL_K, "besselk", args, nargout);
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
367 }
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
368
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
369 DEFUN_DLD (besselh, args, nargout,
3459
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
370 "-*- texinfo -*-\n\
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
371 @deftypefn {Loadable Function} {[@var{h}, @var{ierr}] =} besselh (@var{alpha}, @var{k}, @var{x}, @var{opt})\n\
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
372 See besselj.\n\
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
373 @end deftypefn")
3220
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
374 {
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
375 octave_value_list retval;
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
376
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
377 int nargin = args.length ();
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
378
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
379 if (nargin == 2)
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
380 {
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
381 retval = do_bessel (BESSEL_H1, "besselh", args, nargout);
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
382 }
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
383 else if (nargin == 3 || nargin == 4)
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
384 {
5275
23b37da9fd5b [project @ 2005-04-08 16:07:35 by jwe]
jwe
parents: 4844
diff changeset
385 octave_idx_type kind = args(1).int_value ();
3220
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
386
3810
f19f00339363 [project @ 2001-03-29 21:03:01 by jwe]
jwe
parents: 3568
diff changeset
387 if (! error_state)
3220
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
388 {
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
389 octave_value_list tmp_args;
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
390
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
391 if (nargin == 4)
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
392 tmp_args(2) = args(3);
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
393
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
394 tmp_args(1) = args(2);
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
395 tmp_args(0) = args(0);
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
396
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
397 if (kind == 1)
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
398 retval = do_bessel (BESSEL_H1, "besselh", tmp_args, nargout);
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
399 else if (kind == 2)
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
400 retval = do_bessel (BESSEL_H2, "besselh", tmp_args, nargout);
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
401 else
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
402 error ("besselh: expecting K = 1 or 2");
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
403 }
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
404 else
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
405 error ("besselh: invalid value of K");
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
406 }
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
407 else
5823
080c08b192d8 [project @ 2006-05-19 05:32:17 by jwe]
jwe
parents: 5760
diff changeset
408 print_usage ();
3220
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
409
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
410 return retval;
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
411 }
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
412
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
413 DEFUN_DLD (airy, args, nargout,
3459
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
414 "-*- texinfo -*-\n\
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
415 @deftypefn {Loadable Function} {[@var{a}, @var{ierr}] =} airy (@var{k}, @var{z}, @var{opt})\n\
3220
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
416 Compute Airy functions of the first and second kind, and their\n\
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
417 derivatives.\n\
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
418 \n\
3459
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
419 @example\n\
7031
120f3135952f [project @ 2007-10-15 15:30:03 by jwe]
jwe
parents: 7017
diff changeset
420 K Function Scale factor (if 'opt' is supplied)\n\
120f3135952f [project @ 2007-10-15 15:30:03 by jwe]
jwe
parents: 7017
diff changeset
421 --- -------- ---------------------------------------\n\
120f3135952f [project @ 2007-10-15 15:30:03 by jwe]
jwe
parents: 7017
diff changeset
422 0 Ai (Z) exp ((2/3) * Z * sqrt (Z))\n\
120f3135952f [project @ 2007-10-15 15:30:03 by jwe]
jwe
parents: 7017
diff changeset
423 1 dAi(Z)/dZ exp ((2/3) * Z * sqrt (Z))\n\
120f3135952f [project @ 2007-10-15 15:30:03 by jwe]
jwe
parents: 7017
diff changeset
424 2 Bi (Z) exp (-abs (real ((2/3) * Z *sqrt (Z))))\n\
120f3135952f [project @ 2007-10-15 15:30:03 by jwe]
jwe
parents: 7017
diff changeset
425 3 dBi(Z)/dZ exp (-abs (real ((2/3) * Z *sqrt (Z))))\n\
3459
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
426 @end example\n\
3220
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
427 \n\
3549
03025e79d8b9 [project @ 2000-02-03 01:28:24 by jwe]
jwe
parents: 3548
diff changeset
428 The function call @code{airy (@var{z})} is equivalent to\n\
3459
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
429 @code{airy (0, @var{z})}.\n\
3155
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
430 \n\
3549
03025e79d8b9 [project @ 2000-02-03 01:28:24 by jwe]
jwe
parents: 3548
diff changeset
431 The result is the same size as @var{z}.\n\
3220
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
432 \n\
3459
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
433 If requested, @var{ierr} contains the following status information and\n\
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
434 is the same size as the result.\n\
3548
ab7fa5a8f23f [project @ 2000-02-03 01:17:15 by jwe]
jwe
parents: 3499
diff changeset
435 \n\
3459
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
436 @enumerate 0\n\
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
437 @item\n\
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
438 Normal return.\n\
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
439 @item\n\
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
440 Input error, return @code{NaN}.\n\
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
441 @item\n\
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
442 Overflow, return @code{Inf}.\n\
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
443 @item\n\
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
444 Loss of significance by argument reduction results in less than half\n\
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
445 of machine accuracy.\n\
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
446 @item\n\
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
447 Complete loss of significance by argument reduction, return @code{NaN}.\n\
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
448 @item\n\
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
449 Error---no computation, algorithm termination condition not met,\n\
5448
ebe5d7d15522 [project @ 2005-09-14 18:55:04 by jwe]
jwe
parents: 5307
diff changeset
450 return @code{NaN}.\n\
3459
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
451 @end enumerate\n\
8e36c45e3a61 [project @ 2000-01-19 10:26:18 by jwe]
jwe
parents: 3325
diff changeset
452 @end deftypefn")
3155
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
453 {
3220
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
454 octave_value_list retval;
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
455
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
456 int nargin = args.length ();
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
457
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
458 if (nargin > 0 && nargin < 4)
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
459 {
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
460 bool scale = (nargin == 3);
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
461
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
462 int kind = 0;
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
463
4844
9f7ef92b50b0 [project @ 2004-04-02 17:26:53 by jwe]
jwe
parents: 4153
diff changeset
464 ComplexNDArray z;
3220
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
465
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
466 if (nargin > 1)
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
467 {
5760
8d7162924bd3 [project @ 2006-04-14 04:01:37 by jwe]
jwe
parents: 5448
diff changeset
468 kind = args(0).int_value ();
3220
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
469
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
470 if (! error_state)
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
471 {
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
472 if (kind < 0 || kind > 3)
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
473 error ("airy: expecting K = 0, 1, 2, or 3");
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
474 }
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
475 else
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
476 error ("airy: expecting integer value for K");
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
477 }
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
478
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
479 if (! error_state)
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
480 {
4844
9f7ef92b50b0 [project @ 2004-04-02 17:26:53 by jwe]
jwe
parents: 4153
diff changeset
481 z = args(nargin == 1 ? 0 : 1).complex_array_value ();
3220
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
482
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
483 if (! error_state)
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
484 {
5275
23b37da9fd5b [project @ 2005-04-08 16:07:35 by jwe]
jwe
parents: 4844
diff changeset
485 ArrayN<octave_idx_type> ierr;
3220
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
486 octave_value result;
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
487
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
488 if (kind > 1)
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
489 result = biry (z, kind == 3, scale, ierr);
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
490 else
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
491 result = airy (z, kind == 1, scale, ierr);
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
492
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
493 if (nargout > 1)
4844
9f7ef92b50b0 [project @ 2004-04-02 17:26:53 by jwe]
jwe
parents: 4153
diff changeset
494 retval(1) = int_arrayN_to_array (ierr);
3220
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
495
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
496 retval(0) = result;
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
497 }
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
498 else
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
499 error ("airy: expecting complex matrix for Z");
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
500 }
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
501 }
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
502 else
5823
080c08b192d8 [project @ 2006-05-19 05:32:17 by jwe]
jwe
parents: 5760
diff changeset
503 print_usage ();
3220
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
504
3deb1105fbc1 [project @ 1998-11-19 00:06:30 by jwe]
jwe
parents: 3187
diff changeset
505 return retval;
3155
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
506 }
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
507
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
508 /*
11904
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
509 %! # Test values computed with GP/PARI version 2.3.3
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
510 %!
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
511 %!shared alpha, x, jx, yx, ix, kx, nix
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
512 %!
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
513 %! # Bessel functions, even order, positive and negative x
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
514 %! alpha = 2; x = 1.25;
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
515 %! jx = 0.1710911312405234823613091417;
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
516 %! yx = -1.193199310178553861283790424;
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
517 %! ix = 0.2220184483766341752692212604;
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
518 %! kx = 0.9410016167388185767085460540;
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
519 %!
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
520 %!assert(besselj(alpha,x), jx, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
521 %!assert(bessely(alpha,x), yx, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
522 %!assert(besseli(alpha,x), ix, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
523 %!assert(besselk(alpha,x), kx, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
524 %!assert(besselh(alpha,1,x), jx + I*yx, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
525 %!assert(besselh(alpha,2,x), jx - I*yx, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
526 %!
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
527 %!assert(besselj(alpha,x,1), jx*exp(-abs(imag(x))), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
528 %!assert(bessely(alpha,x,1), yx*exp(-abs(imag(x))), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
529 %!assert(besseli(alpha,x,1), ix*exp(-abs(real(x))), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
530 %!assert(besselk(alpha,x,1), kx*exp(x), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
531 %!assert(besselh(alpha,1,x,1), (jx + I*yx)*exp(-I*x), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
532 %!assert(besselh(alpha,2,x,1), (jx - I*yx)*exp(I*x), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
533 %!
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
534 %!assert(besselj(-alpha,x), jx, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
535 %!assert(bessely(-alpha,x), yx, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
536 %!assert(besseli(-alpha,x), ix, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
537 %!assert(besselk(-alpha,x), kx, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
538 %!assert(besselh(-alpha,1,x), jx + I*yx, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
539 %!assert(besselh(-alpha,2,x), jx - I*yx, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
540 %!
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
541 %!assert(besselj(-alpha,x,1), jx*exp(-abs(imag(x))), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
542 %!assert(bessely(-alpha,x,1), yx*exp(-abs(imag(x))), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
543 %!assert(besseli(-alpha,x,1), ix*exp(-abs(real(x))), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
544 %!assert(besselk(-alpha,x,1), kx*exp(x), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
545 %!assert(besselh(-alpha,1,x,1), (jx + I*yx)*exp(-I*x), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
546 %!assert(besselh(-alpha,2,x,1), (jx - I*yx)*exp(I*x), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
547 %!
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
548 %! x *= -1;
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
549 %! yx = -1.193199310178553861283790424 + 0.3421822624810469647226182835*I;
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
550 %! kx = 0.9410016167388185767085460540 - 0.6974915263814386815610060884*I;
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
551 %!
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
552 %!assert(besselj(alpha,x), jx, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
553 %!assert(bessely(alpha,x), yx, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
554 %!assert(besseli(alpha,x), ix, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
555 %!assert(besselk(alpha,x), kx, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
556 %!assert(besselh(alpha,1,x), jx + I*yx, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
557 %!assert(besselh(alpha,2,x), jx - I*yx, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
558 %!
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
559 %!assert(besselj(alpha,x,1), jx*exp(-abs(imag(x))), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
560 %!assert(bessely(alpha,x,1), yx*exp(-abs(imag(x))), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
561 %!assert(besseli(alpha,x,1), ix*exp(-abs(real(x))), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
562 %!assert(besselk(alpha,x,1), kx*exp(x), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
563 %!assert(besselh(alpha,1,x,1), (jx + I*yx)*exp(-I*x), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
564 %!assert(besselh(alpha,2,x,1), (jx - I*yx)*exp(I*x), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
565 %!
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
566 %! # Bessel functions, odd order, positive and negative x
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
567 %! alpha = 3; x = 2.5;
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
568 %! jx = 0.2166003910391135247666890035;
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
569 %! yx = -0.7560554967536709968379029772;
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
570 %! ix = 0.4743704087780355895548240179;
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
571 %! kx = 0.2682271463934492027663765197;
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
572 %!
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
573 %!assert(besselj(alpha,x), jx, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
574 %!assert(bessely(alpha,x), yx, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
575 %!assert(besseli(alpha,x), ix, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
576 %!assert(besselk(alpha,x), kx, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
577 %!assert(besselh(alpha,1,x), jx + I*yx, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
578 %!assert(besselh(alpha,2,x), jx - I*yx, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
579 %!
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
580 %!assert(besselj(alpha,x,1), jx*exp(-abs(imag(x))), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
581 %!assert(bessely(alpha,x,1), yx*exp(-abs(imag(x))), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
582 %!assert(besseli(alpha,x,1), ix*exp(-abs(real(x))), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
583 %!assert(besselk(alpha,x,1), kx*exp(x), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
584 %!assert(besselh(alpha,1,x,1), (jx + I*yx)*exp(-I*x), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
585 %!assert(besselh(alpha,2,x,1), (jx - I*yx)*exp(I*x), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
586 %!
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
587 %!assert(besselj(-alpha,x), -jx, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
588 %!assert(bessely(-alpha,x), -yx, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
589 %!assert(besseli(-alpha,x), ix, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
590 %!assert(besselk(-alpha,x), kx, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
591 %!assert(besselh(-alpha,1,x), -(jx + I*yx), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
592 %!assert(besselh(-alpha,2,x), -(jx - I*yx), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
593 %!
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
594 %!assert(besselj(-alpha,x,1), -jx*exp(-abs(imag(x))), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
595 %!assert(bessely(-alpha,x,1), -yx*exp(-abs(imag(x))), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
596 %!assert(besseli(-alpha,x,1), ix*exp(-abs(real(x))), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
597 %!assert(besselk(-alpha,x,1), kx*exp(x), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
598 %!assert(besselh(-alpha,1,x,1), -(jx + I*yx)*exp(-I*x), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
599 %!assert(besselh(-alpha,2,x,1), -(jx - I*yx)*exp(I*x), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
600 %!
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
601 %! x *= -1;
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
602 %! jx = -jx;
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
603 %! yx = 0.7560554967536709968379029772 - 0.4332007820782270495333780070*I;
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
604 %! ix = -ix;
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
605 %! kx = -0.2682271463934492027663765197 - 1.490278591297463775542004240*I;
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
606 %!
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
607 %!assert(besselj(alpha,x), jx, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
608 %!assert(bessely(alpha,x), yx, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
609 %!assert(besseli(alpha,x), ix, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
610 %!assert(besselk(alpha,x), kx, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
611 %!assert(besselh(alpha,1,x), jx + I*yx, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
612 %!assert(besselh(alpha,2,x), jx - I*yx, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
613 %!
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
614 %!assert(besselj(alpha,x,1), jx*exp(-abs(imag(x))), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
615 %!assert(bessely(alpha,x,1), yx*exp(-abs(imag(x))), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
616 %!assert(besseli(alpha,x,1), ix*exp(-abs(real(x))), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
617 %!assert(besselk(alpha,x,1), kx*exp(x), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
618 %!assert(besselh(alpha,1,x,1), (jx + I*yx)*exp(-I*x), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
619 %!assert(besselh(alpha,2,x,1), (jx - I*yx)*exp(I*x), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
620 %!
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
621 %! # Bessel functions, fractional order, positive and negative x
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
622 %!
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
623 %! alpha = 3.5; x = 2.75;
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
624 %! jx = 0.1691636439842384154644784389;
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
625 %! yx = -0.8301381935499356070267953387;
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
626 %! ix = 0.3930540878794826310979363668;
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
627 %! kx = 0.2844099013460621170288192503;
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
628 %!
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
629 %!assert(besselj(alpha,x), jx, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
630 %!assert(bessely(alpha,x), yx, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
631 %!assert(besseli(alpha,x), ix, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
632 %!assert(besselk(alpha,x), kx, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
633 %!assert(besselh(alpha,1,x), jx + I*yx, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
634 %!assert(besselh(alpha,2,x), jx - I*yx, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
635 %!
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
636 %!assert(besselj(alpha,x,1), jx*exp(-abs(imag(x))), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
637 %!assert(bessely(alpha,x,1), yx*exp(-abs(imag(x))), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
638 %!assert(besseli(alpha,x,1), ix*exp(-abs(real(x))), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
639 %!assert(besselk(alpha,x,1), kx*exp(x), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
640 %!assert(besselh(alpha,1,x,1), (jx + I*yx)*exp(-I*x), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
641 %!assert(besselh(alpha,2,x,1), (jx - I*yx)*exp(I*x), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
642 %!
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
643 %! nix = 0.2119931212254662995364461998;
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
644 %!
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
645 %!assert(besselj(-alpha,x), yx, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
646 %!assert(bessely(-alpha,x), -jx, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
647 %!assert(besseli(-alpha,x), nix, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
648 %!assert(besselk(-alpha,x), kx, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
649 %!assert(besselh(-alpha,1,x), -I*(jx + I*yx), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
650 %!assert(besselh(-alpha,2,x), I*(jx - I*yx), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
651 %!
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
652 %!assert(besselj(-alpha,x,1), yx*exp(-abs(imag(x))), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
653 %!assert(bessely(-alpha,x,1), -jx*exp(-abs(imag(x))), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
654 %!assert(besseli(-alpha,x,1), nix*exp(-abs(real(x))), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
655 %!assert(besselk(-alpha,x,1), kx*exp(x), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
656 %!assert(besselh(-alpha,1,x,1), -I*(jx + I*yx)*exp(-I*x), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
657 %!assert(besselh(-alpha,2,x,1), I*(jx - I*yx)*exp(I*x), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
658 %!
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
659 %! x *= -1;
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
660 %! jx *= -I;
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
661 %! yx = -0.8301381935499356070267953387*I;
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
662 %! ix *= -I;
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
663 %! kx = -0.9504059335995575096509874508*I;
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
664 %!
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
665 %!assert(besselj(alpha,x), jx, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
666 %!assert(bessely(alpha,x), yx, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
667 %!assert(besseli(alpha,x), ix, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
668 %!assert(besselk(alpha,x), kx, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
669 %!assert(besselh(alpha,1,x), jx + I*yx, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
670 %!assert(besselh(alpha,2,x), jx - I*yx, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
671 %!
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
672 %!assert(besselj(alpha,x,1), jx*exp(-abs(imag(x))), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
673 %!assert(bessely(alpha,x,1), yx*exp(-abs(imag(x))), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
674 %!assert(besseli(alpha,x,1), ix*exp(-abs(real(x))), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
675 %!assert(besselk(alpha,x,1), kx*exp(x), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
676 %!assert(besselh(alpha,1,x,1), (jx + I*yx)*exp(-I*x), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
677 %!assert(besselh(alpha,2,x,1), (jx - I*yx)*exp(I*x), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
678 %!
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
679 %! # Bessel functions, even order, complex x
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
680 %!
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
681 %! alpha = 2; x = 1.25 + 3.625 * I;
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
682 %! jx = -1.299533366810794494030065917 + 4.370833116012278943267479589*I;
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
683 %! yx = -4.370357232383223896393056727 - 1.283083391453582032688834041*I;
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
684 %! ix = -0.6717801680341515541002273932 - 0.2314623443930774099910228553*I;
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
685 %! kx = -0.01108009888623253515463783379 + 0.2245218229358191588208084197*I;
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
686 %!
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
687 %!assert(besselj(alpha,x), jx, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
688 %!assert(bessely(alpha,x), yx, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
689 %!assert(besseli(alpha,x), ix, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
690 %!assert(besselk(alpha,x), kx, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
691 %!assert(besselh(alpha,1,x), jx + I*yx, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
692 %!assert(besselh(alpha,2,x), jx - I*yx, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
693 %!
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
694 %!assert(besselj(alpha,x,1), jx*exp(-abs(imag(x))), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
695 %!assert(bessely(alpha,x,1), yx*exp(-abs(imag(x))), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
696 %!assert(besseli(alpha,x,1), ix*exp(-abs(real(x))), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
697 %!assert(besselk(alpha,x,1), kx*exp(x), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
698 %!assert(besselh(alpha,1,x,1), (jx + I*yx)*exp(-I*x), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
699 %!assert(besselh(alpha,2,x,1), (jx - I*yx)*exp(I*x), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
700 %!
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
701 %!assert(besselj(-alpha,x), jx, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
702 %!assert(bessely(-alpha,x), yx, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
703 %!assert(besseli(-alpha,x), ix, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
704 %!assert(besselk(-alpha,x), kx, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
705 %!assert(besselh(-alpha,1,x), jx + I*yx, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
706 %!assert(besselh(-alpha,2,x), jx - I*yx, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
707 %!
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
708 %!assert(besselj(-alpha,x,1), jx*exp(-abs(imag(x))), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
709 %!assert(bessely(-alpha,x,1), yx*exp(-abs(imag(x))), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
710 %!assert(besseli(-alpha,x,1), ix*exp(-abs(real(x))), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
711 %!assert(besselk(-alpha,x,1), kx*exp(x), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
712 %!assert(besselh(-alpha,1,x,1), (jx + I*yx)*exp(-I*x), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
713 %!assert(besselh(-alpha,2,x,1), (jx - I*yx)*exp(I*x), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
714 %!
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
715 %! # Bessel functions, odd order, complex x
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
716 %!
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
717 %! alpha = 3; x = 2.5 + 1.875 * I;
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
718 %! jx = 0.1330721523048277493333458596 + 0.5386295217249660078754395597*I;
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
719 %! yx = -0.6485072392105829901122401551 + 0.2608129289785456797046996987*I;
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
720 %! ix = -0.6182064685486998097516365709 + 0.4677561094683470065767989920*I;
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
721 %! kx = -0.1568585587733540007867882337 - 0.05185853709490846050505141321*I;
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
722 %!
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
723 %!assert(besselj(alpha,x), jx, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
724 %!assert(bessely(alpha,x), yx, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
725 %!assert(besseli(alpha,x), ix, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
726 %!assert(besselk(alpha,x), kx, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
727 %!assert(besselh(alpha,1,x), jx + I*yx, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
728 %!assert(besselh(alpha,2,x), jx - I*yx, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
729 %!
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
730 %!assert(besselj(alpha,x,1), jx*exp(-abs(imag(x))), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
731 %!assert(bessely(alpha,x,1), yx*exp(-abs(imag(x))), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
732 %!assert(besseli(alpha,x,1), ix*exp(-abs(real(x))), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
733 %!assert(besselk(alpha,x,1), kx*exp(x), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
734 %!assert(besselh(alpha,1,x,1), (jx + I*yx)*exp(-I*x), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
735 %!assert(besselh(alpha,2,x,1), (jx - I*yx)*exp(I*x), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
736 %!
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
737 %!assert(besselj(-alpha,x), -jx, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
738 %!assert(bessely(-alpha,x), -yx, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
739 %!assert(besseli(-alpha,x), ix, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
740 %!assert(besselk(-alpha,x), kx, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
741 %!assert(besselh(-alpha,1,x), -(jx + I*yx), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
742 %!assert(besselh(-alpha,2,x), -(jx - I*yx), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
743 %!
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
744 %!assert(besselj(-alpha,x,1), -jx*exp(-abs(imag(x))), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
745 %!assert(bessely(-alpha,x,1), -yx*exp(-abs(imag(x))), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
746 %!assert(besseli(-alpha,x,1), ix*exp(-abs(real(x))), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
747 %!assert(besselk(-alpha,x,1), kx*exp(x), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
748 %!assert(besselh(-alpha,1,x,1), -(jx + I*yx)*exp(-I*x), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
749 %!assert(besselh(-alpha,2,x,1), -(jx - I*yx)*exp(I*x), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
750 %!
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
751 %! # Bessel functions, fractional order, complex x
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
752 %!
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
753 %! alpha = 3.5; x = 1.75 + 4.125 * I;
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
754 %! jx = -3.018566131370455929707009100 - 0.7585648436793900607704057611*I;
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
755 %! yx = 0.7772278839106298215614791107 - 3.018518722313849782683792010*I;
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
756 %! ix = 0.2100873577220057189038160913 - 0.6551765604618246531254970926*I;
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
757 %! kx = 0.1757147290513239935341488069 + 0.08772348296883849205562558311*I;
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
758 %!
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
759 %!assert(besselj(alpha,x), jx, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
760 %!assert(bessely(alpha,x), yx, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
761 %!assert(besseli(alpha,x), ix, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
762 %!assert(besselk(alpha,x), kx, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
763 %!assert(besselh(alpha,1,x), jx + I*yx, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
764 %!assert(besselh(alpha,2,x), jx - I*yx, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
765 %!
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
766 %!assert(besselj(alpha,x,1), jx*exp(-abs(imag(x))), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
767 %!assert(bessely(alpha,x,1), yx*exp(-abs(imag(x))), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
768 %!assert(besseli(alpha,x,1), ix*exp(-abs(real(x))), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
769 %!assert(besselk(alpha,x,1), kx*exp(x), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
770 %!assert(besselh(alpha,1,x,1), (jx + I*yx)*exp(-I*x), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
771 %!assert(besselh(alpha,2,x,1), (jx - I*yx)*exp(I*x), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
772 %!
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
773 %! nix = 0.09822388691172060573913739253 - 0.7110230642207380127317227407*I;
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
774 %!
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
775 %!assert(besselj(-alpha,x), yx, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
776 %!assert(bessely(-alpha,x), -jx, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
777 %!assert(besseli(-alpha,x), nix, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
778 %!assert(besselk(-alpha,x), kx, 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
779 %!assert(besselh(-alpha,1,x), -I*(jx + I*yx), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
780 %!assert(besselh(-alpha,2,x), I*(jx - I*yx), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
781 %!
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
782 %!assert(besselj(-alpha,x,1), yx*exp(-abs(imag(x))), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
783 %!assert(bessely(-alpha,x,1), -jx*exp(-abs(imag(x))), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
784 %!assert(besseli(-alpha,x,1), nix*exp(-abs(real(x))), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
785 %!assert(besselk(-alpha,x,1), kx*exp(x), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
786 %!assert(besselh(-alpha,1,x,1), -I*(jx + I*yx)*exp(-I*x), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
787 %!assert(besselh(-alpha,2,x,1), I*(jx - I*yx)*exp(I*x), 100*eps)
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
788 */
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
789
059fadc0cbc3 fix scaling factor for negative alpha in zbesi,cbesi
Brian Gough <bjg@gnu.org>
parents: 7031
diff changeset
790 /*
3155
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
791 ;;; Local Variables: ***
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
792 ;;; mode: C++ ***
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
793 ;;; End: ***
1016520a9d38 [project @ 1998-02-18 21:51:50 by jwe]
jwe
parents:
diff changeset
794 */