3155
|
1 /* |
|
2 |
|
3 Copyright (C) 1997 John W. Eaton |
|
4 |
|
5 This file is part of Octave. |
|
6 |
|
7 Octave is free software; you can redistribute it and/or modify it |
|
8 under the terms of the GNU General Public License as published by the |
|
9 Free Software Foundation; either version 2, or (at your option) any |
|
10 later version. |
|
11 |
|
12 Octave is distributed in the hope that it will be useful, but WITHOUT |
|
13 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or |
|
14 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License |
|
15 for more details. |
|
16 |
|
17 You should have received a copy of the GNU General Public License |
|
18 along with Octave; see the file COPYING. If not, write to the Free |
5307
|
19 Software Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA |
|
20 02110-1301, USA. |
3155
|
21 |
|
22 */ |
|
23 |
|
24 #ifdef HAVE_CONFIG_H |
|
25 #include <config.h> |
|
26 #endif |
|
27 |
|
28 #include "lo-specfun.h" |
4153
|
29 #include "quit.h" |
3155
|
30 |
|
31 #include "defun-dld.h" |
|
32 #include "error.h" |
|
33 #include "gripes.h" |
|
34 #include "oct-obj.h" |
|
35 #include "utils.h" |
|
36 |
3220
|
37 enum bessel_type |
|
38 { |
|
39 BESSEL_J, |
|
40 BESSEL_Y, |
|
41 BESSEL_I, |
|
42 BESSEL_K, |
|
43 BESSEL_H1, |
|
44 BESSEL_H2 |
|
45 }; |
|
46 |
|
47 #define DO_BESSEL(type, alpha, x, scaled, ierr, result) \ |
3155
|
48 do \ |
|
49 { \ |
|
50 switch (type) \ |
|
51 { \ |
3220
|
52 case BESSEL_J: \ |
|
53 result = besselj (alpha, x, scaled, ierr); \ |
|
54 break; \ |
|
55 \ |
|
56 case BESSEL_Y: \ |
|
57 result = bessely (alpha, x, scaled, ierr); \ |
3155
|
58 break; \ |
|
59 \ |
3220
|
60 case BESSEL_I: \ |
|
61 result = besseli (alpha, x, scaled, ierr); \ |
3155
|
62 break; \ |
|
63 \ |
3220
|
64 case BESSEL_K: \ |
|
65 result = besselk (alpha, x, scaled, ierr); \ |
3155
|
66 break; \ |
|
67 \ |
3220
|
68 case BESSEL_H1: \ |
|
69 result = besselh1 (alpha, x, scaled, ierr); \ |
|
70 break; \ |
|
71 \ |
|
72 case BESSEL_H2: \ |
|
73 result = besselh2 (alpha, x, scaled, ierr); \ |
3155
|
74 break; \ |
|
75 \ |
|
76 default: \ |
|
77 break; \ |
|
78 } \ |
|
79 } \ |
|
80 while (0) |
|
81 |
3220
|
82 static inline Matrix |
5275
|
83 int_array2_to_matrix (const Array2<octave_idx_type>& a) |
3220
|
84 { |
5275
|
85 octave_idx_type nr = a.rows (); |
|
86 octave_idx_type nc = a.cols (); |
3220
|
87 |
|
88 Matrix retval (nr, nc); |
|
89 |
5275
|
90 for (octave_idx_type j = 0; j < nc; j++) |
|
91 for (octave_idx_type i = 0; i < nr; i++) |
4153
|
92 { |
|
93 OCTAVE_QUIT; |
|
94 |
5760
|
95 retval(i,j) = static_cast<double> (a(i,j)); |
4153
|
96 } |
3220
|
97 |
|
98 return retval; |
|
99 } |
|
100 |
4844
|
101 static inline NDArray |
5760
|
102 int_arrayN_to_array (const ArrayN<octave_idx_type>& a) |
4844
|
103 { |
|
104 dim_vector dv = a.dims (); |
|
105 int nel = dv.numel (); |
|
106 |
|
107 NDArray retval (dv); |
|
108 |
|
109 for (int i = 0; i < nel; i++) |
|
110 { |
|
111 OCTAVE_QUIT; |
|
112 |
5760
|
113 retval(i) = static_cast<double> (a(i)); |
4844
|
114 } |
|
115 |
|
116 return retval; |
|
117 } |
|
118 |
3155
|
119 static void |
3220
|
120 gripe_bessel_arg (const char *fn, const char *arg) |
3155
|
121 { |
3220
|
122 error ("%s: expecting scalar or matrix as %s argument", fn, arg); |
3155
|
123 } |
|
124 |
|
125 octave_value_list |
3220
|
126 do_bessel (enum bessel_type type, const char *fn, |
|
127 const octave_value_list& args, int nargout) |
3155
|
128 { |
3220
|
129 octave_value_list retval; |
3155
|
130 |
|
131 int nargin = args.length (); |
|
132 |
3220
|
133 if (nargin == 2 || nargin == 3) |
3155
|
134 { |
3220
|
135 bool scaled = (nargin == 3); |
|
136 |
3155
|
137 octave_value alpha_arg = args(0); |
3220
|
138 octave_value x_arg = args(1); |
3155
|
139 |
|
140 if (alpha_arg.is_scalar_type ()) |
|
141 { |
3220
|
142 double alpha = args(0).double_value (); |
3155
|
143 |
|
144 if (! error_state) |
|
145 { |
3220
|
146 if (x_arg.is_scalar_type ()) |
|
147 { |
|
148 Complex x = x_arg.complex_value (); |
|
149 |
|
150 if (! error_state) |
|
151 { |
5275
|
152 octave_idx_type ierr; |
3220
|
153 octave_value result; |
|
154 |
|
155 DO_BESSEL (type, alpha, x, scaled, ierr, result); |
|
156 |
|
157 if (nargout > 1) |
5760
|
158 retval(1) = static_cast<double> (ierr); |
3155
|
159 |
3220
|
160 retval(0) = result; |
|
161 } |
|
162 else |
|
163 gripe_bessel_arg (fn, "second"); |
|
164 } |
3155
|
165 else |
3220
|
166 { |
4844
|
167 ComplexNDArray x = x_arg.complex_array_value (); |
3220
|
168 |
|
169 if (! error_state) |
|
170 { |
5275
|
171 ArrayN<octave_idx_type> ierr; |
3220
|
172 octave_value result; |
|
173 |
|
174 DO_BESSEL (type, alpha, x, scaled, ierr, result); |
|
175 |
|
176 if (nargout > 1) |
4844
|
177 retval(1) = int_arrayN_to_array (ierr); |
3220
|
178 |
|
179 retval(0) = result; |
|
180 } |
|
181 else |
|
182 gripe_bessel_arg (fn, "second"); |
|
183 } |
3155
|
184 } |
|
185 else |
3220
|
186 gripe_bessel_arg (fn, "first"); |
3155
|
187 } |
|
188 else |
|
189 { |
4844
|
190 dim_vector dv0 = args(0).dims (); |
|
191 dim_vector dv1 = args(1).dims (); |
|
192 |
|
193 bool args0_is_row_vector = (dv0 (1) == dv0.numel ()); |
|
194 bool args1_is_col_vector = (dv1 (0) == dv1.numel ()); |
3155
|
195 |
4844
|
196 if (args0_is_row_vector && args1_is_col_vector) |
3155
|
197 { |
4844
|
198 RowVector ralpha = args(0).row_vector_value (); |
|
199 |
|
200 if (! error_state) |
3220
|
201 { |
4844
|
202 ComplexColumnVector cx = |
|
203 x_arg.complex_column_vector_value (); |
3220
|
204 |
|
205 if (! error_state) |
|
206 { |
5275
|
207 Array2<octave_idx_type> ierr; |
3220
|
208 octave_value result; |
3155
|
209 |
4844
|
210 DO_BESSEL (type, ralpha, cx, scaled, ierr, result); |
|
211 |
3220
|
212 if (nargout > 1) |
|
213 retval(1) = int_array2_to_matrix (ierr); |
|
214 |
|
215 retval(0) = result; |
|
216 } |
|
217 else |
|
218 gripe_bessel_arg (fn, "second"); |
|
219 } |
|
220 else |
4844
|
221 gripe_bessel_arg (fn, "first"); |
|
222 } |
|
223 else |
|
224 { |
|
225 NDArray alpha = args(0).array_value (); |
3220
|
226 |
4844
|
227 if (! error_state) |
|
228 { |
|
229 if (x_arg.is_scalar_type ()) |
|
230 { |
|
231 Complex x = x_arg.complex_value (); |
3155
|
232 |
4844
|
233 if (! error_state) |
3155
|
234 { |
5275
|
235 ArrayN<octave_idx_type> ierr; |
3220
|
236 octave_value result; |
|
237 |
|
238 DO_BESSEL (type, alpha, x, scaled, ierr, result); |
|
239 |
|
240 if (nargout > 1) |
4844
|
241 retval(1) = int_arrayN_to_array (ierr); |
3220
|
242 |
|
243 retval(0) = result; |
3155
|
244 } |
4844
|
245 else |
|
246 gripe_bessel_arg (fn, "second"); |
3155
|
247 } |
3220
|
248 else |
4844
|
249 { |
|
250 ComplexNDArray x = x_arg.complex_array_value (); |
|
251 |
|
252 if (! error_state) |
|
253 { |
5275
|
254 ArrayN<octave_idx_type> ierr; |
4844
|
255 octave_value result; |
|
256 |
|
257 DO_BESSEL (type, alpha, x, scaled, ierr, result); |
|
258 |
|
259 if (nargout > 1) |
|
260 retval(1) = int_arrayN_to_array (ierr); |
|
261 |
|
262 retval(0) = result; |
|
263 } |
|
264 else |
|
265 gripe_bessel_arg (fn, "second"); |
|
266 } |
3155
|
267 } |
4844
|
268 else |
|
269 gripe_bessel_arg (fn, "first"); |
3155
|
270 } |
|
271 } |
|
272 } |
|
273 else |
5823
|
274 print_usage (); |
3155
|
275 |
|
276 return retval; |
|
277 } |
|
278 |
3220
|
279 DEFUN_DLD (besselj, args, nargout, |
3459
|
280 "-*- texinfo -*-\n\ |
|
281 @deftypefn {Loadable Function} {[@var{j}, @var{ierr}] =} besselj (@var{alpha}, @var{x}, @var{opt})\n\ |
|
282 @deftypefnx {Loadable Function} {[@var{y}, @var{ierr}] =} bessely (@var{alpha}, @var{x}, @var{opt})\n\ |
|
283 @deftypefnx {Loadable Function} {[@var{i}, @var{ierr}] =} besseli (@var{alpha}, @var{x}, @var{opt})\n\ |
|
284 @deftypefnx {Loadable Function} {[@var{k}, @var{ierr}] =} besselk (@var{alpha}, @var{x}, @var{opt})\n\ |
|
285 @deftypefnx {Loadable Function} {[@var{h}, @var{ierr}] =} besselh (@var{alpha}, @var{k}, @var{x}, @var{opt})\n\ |
|
286 Compute Bessel or Hankel functions of various kinds:\n\ |
3155
|
287 \n\ |
3459
|
288 @table @code\n\ |
|
289 @item besselj\n\ |
|
290 Bessel functions of the first kind.\n\ |
|
291 @item bessely\n\ |
|
292 Bessel functions of the second kind.\n\ |
|
293 @item besseli\n\ |
|
294 Modified Bessel functions of the first kind.\n\ |
|
295 @item besselk\n\ |
|
296 Modified Bessel functions of the second kind.\n\ |
|
297 @item besselh\n\ |
|
298 Compute Hankel functions of the first (@var{k} = 1) or second (@var{k}\n\ |
|
299 = 2) kind.\n\ |
|
300 @end table\n\ |
3220
|
301 \n\ |
3568
|
302 If the argument @var{opt} is supplied, the result is scaled by the\n\ |
3459
|
303 @code{exp (-I*@var{x})} for @var{k} = 1 or @code{exp (I*@var{x})} for\n\ |
3499
|
304 @var{k} = 2.\n\ |
3220
|
305 \n\ |
3459
|
306 If @var{alpha} is a scalar, the result is the same size as @var{x}.\n\ |
|
307 If @var{x} is a scalar, the result is the same size as @var{alpha}.\n\ |
|
308 If @var{alpha} is a row vector and @var{x} is a column vector, the\n\ |
|
309 result is a matrix with @code{length (@var{x})} rows and\n\ |
|
310 @code{length (@var{alpha})} columns. Otherwise, @var{alpha} and\n\ |
|
311 @var{x} must conform and the result will be the same size.\n\ |
3155
|
312 \n\ |
3459
|
313 The value of @var{alpha} must be real. The value of @var{x} may be\n\ |
|
314 complex.\n\ |
|
315 \n\ |
|
316 If requested, @var{ierr} contains the following status information\n\ |
|
317 and is the same size as the result.\n\ |
3548
|
318 \n\ |
3459
|
319 @enumerate 0\n\ |
|
320 @item\n\ |
|
321 Normal return.\n\ |
|
322 @item\n\ |
|
323 Input error, return @code{NaN}.\n\ |
|
324 @item\n\ |
|
325 Overflow, return @code{Inf}.\n\ |
|
326 @item\n\ |
|
327 Loss of significance by argument reduction results in less than\n\ |
|
328 half of machine accuracy.\n\ |
|
329 @item\n\ |
|
330 Complete loss of significance by argument reduction, return @code{NaN}.\n\ |
|
331 @item\n\ |
|
332 Error---no computation, algorithm termination condition not met,\n\ |
|
333 return @code{NaN}.\n\ |
|
334 @end enumerate\n\ |
|
335 @end deftypefn") |
3155
|
336 { |
3220
|
337 return do_bessel (BESSEL_J, "besselj", args, nargout); |
3155
|
338 } |
|
339 |
3220
|
340 DEFUN_DLD (bessely, args, nargout, |
3459
|
341 "-*- texinfo -*-\n\ |
|
342 @deftypefn {Loadable Function} {[@var{y}, @var{ierr}] =} bessely (@var{alpha}, @var{x}, @var{opt})\n\ |
|
343 See besselj.\n\ |
|
344 @end deftypefn") |
3155
|
345 { |
3220
|
346 return do_bessel (BESSEL_Y, "bessely", args, nargout); |
3155
|
347 } |
|
348 |
3220
|
349 DEFUN_DLD (besseli, args, nargout, |
3459
|
350 "-*- texinfo -*-\n\ |
|
351 @deftypefn {Loadable Function} {[@var{i}, @var{ierr}] =} besseli (@var{alpha}, @var{x}, @var{opt})\n\ |
|
352 See besselj.\n\ |
|
353 @end deftypefn") |
3155
|
354 { |
3220
|
355 return do_bessel (BESSEL_I, "besseli", args, nargout); |
3155
|
356 } |
|
357 |
3220
|
358 DEFUN_DLD (besselk, args, nargout, |
3459
|
359 "-*- texinfo -*-\n\ |
|
360 @deftypefn {Loadable Function} {[@var{k}, @var{ierr}] =} besselk (@var{alpha}, @var{x}, @var{opt})\n\ |
|
361 See besselj.\n\ |
|
362 @end deftypefn") |
3220
|
363 { |
|
364 return do_bessel (BESSEL_K, "besselk", args, nargout); |
|
365 } |
|
366 |
|
367 DEFUN_DLD (besselh, args, nargout, |
3459
|
368 "-*- texinfo -*-\n\ |
|
369 @deftypefn {Loadable Function} {[@var{h}, @var{ierr}] =} besselh (@var{alpha}, @var{k}, @var{x}, @var{opt})\n\ |
|
370 See besselj.\n\ |
|
371 @end deftypefn") |
3220
|
372 { |
|
373 octave_value_list retval; |
|
374 |
|
375 int nargin = args.length (); |
|
376 |
|
377 if (nargin == 2) |
|
378 { |
|
379 retval = do_bessel (BESSEL_H1, "besselh", args, nargout); |
|
380 } |
|
381 else if (nargin == 3 || nargin == 4) |
|
382 { |
5275
|
383 octave_idx_type kind = args(1).int_value (); |
3220
|
384 |
3810
|
385 if (! error_state) |
3220
|
386 { |
|
387 octave_value_list tmp_args; |
|
388 |
|
389 if (nargin == 4) |
|
390 tmp_args(2) = args(3); |
|
391 |
|
392 tmp_args(1) = args(2); |
|
393 tmp_args(0) = args(0); |
|
394 |
|
395 if (kind == 1) |
|
396 retval = do_bessel (BESSEL_H1, "besselh", tmp_args, nargout); |
|
397 else if (kind == 2) |
|
398 retval = do_bessel (BESSEL_H2, "besselh", tmp_args, nargout); |
|
399 else |
|
400 error ("besselh: expecting K = 1 or 2"); |
|
401 } |
|
402 else |
|
403 error ("besselh: invalid value of K"); |
|
404 } |
|
405 else |
5823
|
406 print_usage (); |
3220
|
407 |
|
408 return retval; |
|
409 } |
|
410 |
|
411 DEFUN_DLD (airy, args, nargout, |
3459
|
412 "-*- texinfo -*-\n\ |
|
413 @deftypefn {Loadable Function} {[@var{a}, @var{ierr}] =} airy (@var{k}, @var{z}, @var{opt})\n\ |
3220
|
414 Compute Airy functions of the first and second kind, and their\n\ |
|
415 derivatives.\n\ |
|
416 \n\ |
3459
|
417 @example\n\ |
3220
|
418 K Function Scale factor (if a third argument is supplied)\n\ |
|
419 --- -------- ----------------------------------------------\n\ |
|
420 0 Ai (Z) exp ((2/3) * Z * sqrt (Z))\n\ |
|
421 1 dAi(Z)/dZ exp ((2/3) * Z * sqrt (Z))\n\ |
|
422 2 Bi (Z) exp (-abs (real ((2/3) * Z *sqrt (Z))))\n\ |
|
423 3 dBi(Z)/dZ exp (-abs (real ((2/3) * Z *sqrt (Z))))\n\ |
3459
|
424 @end example\n\ |
3220
|
425 \n\ |
3549
|
426 The function call @code{airy (@var{z})} is equivalent to\n\ |
3459
|
427 @code{airy (0, @var{z})}.\n\ |
3155
|
428 \n\ |
3549
|
429 The result is the same size as @var{z}.\n\ |
3220
|
430 \n\ |
3459
|
431 If requested, @var{ierr} contains the following status information and\n\ |
|
432 is the same size as the result.\n\ |
3548
|
433 \n\ |
3459
|
434 @enumerate 0\n\ |
|
435 @item\n\ |
|
436 Normal return.\n\ |
|
437 @item\n\ |
|
438 Input error, return @code{NaN}.\n\ |
|
439 @item\n\ |
|
440 Overflow, return @code{Inf}.\n\ |
|
441 @item\n\ |
|
442 Loss of significance by argument reduction results in less than half\n\ |
|
443 of machine accuracy.\n\ |
|
444 @item\n\ |
|
445 Complete loss of significance by argument reduction, return @code{NaN}.\n\ |
|
446 @item\n\ |
|
447 Error---no computation, algorithm termination condition not met,\n\ |
5448
|
448 return @code{NaN}.\n\ |
3459
|
449 @end enumerate\n\ |
|
450 @end deftypefn") |
3155
|
451 { |
3220
|
452 octave_value_list retval; |
|
453 |
|
454 int nargin = args.length (); |
|
455 |
|
456 if (nargin > 0 && nargin < 4) |
|
457 { |
|
458 bool scale = (nargin == 3); |
|
459 |
|
460 int kind = 0; |
|
461 |
4844
|
462 ComplexNDArray z; |
3220
|
463 |
|
464 if (nargin > 1) |
|
465 { |
5760
|
466 kind = args(0).int_value (); |
3220
|
467 |
|
468 if (! error_state) |
|
469 { |
|
470 if (kind < 0 || kind > 3) |
|
471 error ("airy: expecting K = 0, 1, 2, or 3"); |
|
472 } |
|
473 else |
|
474 error ("airy: expecting integer value for K"); |
|
475 } |
|
476 |
|
477 if (! error_state) |
|
478 { |
4844
|
479 z = args(nargin == 1 ? 0 : 1).complex_array_value (); |
3220
|
480 |
|
481 if (! error_state) |
|
482 { |
5275
|
483 ArrayN<octave_idx_type> ierr; |
3220
|
484 octave_value result; |
|
485 |
|
486 if (kind > 1) |
|
487 result = biry (z, kind == 3, scale, ierr); |
|
488 else |
|
489 result = airy (z, kind == 1, scale, ierr); |
|
490 |
|
491 if (nargout > 1) |
4844
|
492 retval(1) = int_arrayN_to_array (ierr); |
3220
|
493 |
|
494 retval(0) = result; |
|
495 } |
|
496 else |
|
497 error ("airy: expecting complex matrix for Z"); |
|
498 } |
|
499 } |
|
500 else |
5823
|
501 print_usage (); |
3220
|
502 |
|
503 return retval; |
3155
|
504 } |
|
505 |
|
506 /* |
|
507 ;;; Local Variables: *** |
|
508 ;;; mode: C++ *** |
|
509 ;;; End: *** |
|
510 */ |