3146
|
1 /* |
|
2 |
|
3 Copyright (C) 1996 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 |
|
19 Software Foundation, 59 Temple Place - Suite 330, Boston, MA 02111-1307, USA. |
|
20 |
|
21 */ |
|
22 |
|
23 #ifdef HAVE_CONFIG_H |
|
24 #include <config.h> |
|
25 #endif |
|
26 |
|
27 #include "Range.h" |
3220
|
28 #include "CColVector.h" |
|
29 #include "CMatrix.h" |
|
30 #include "dRowVector.h" |
3146
|
31 #include "dMatrix.h" |
|
32 #include "f77-fcn.h" |
|
33 #include "lo-error.h" |
3220
|
34 #include "lo-ieee.h" |
|
35 #include "lo-specfun.h" |
3146
|
36 #include "mx-inlines.cc" |
|
37 |
4064
|
38 #ifndef M_PI |
|
39 #define M_PI 3.14159265358979323846 |
|
40 #endif |
|
41 |
3146
|
42 extern "C" |
|
43 { |
3887
|
44 int F77_FUNC (zbesj, ZBESJ) (const double&, const double&, |
4062
|
45 const double&, const int&, const int&, |
|
46 double*, double*, int&, int&); |
3146
|
47 |
3887
|
48 int F77_FUNC (zbesy, ZBESY) (const double&, const double&, |
4062
|
49 const double&, const int&, const int&, |
|
50 double*, double*, int&, |
|
51 double*, double*, int&); |
3220
|
52 |
3887
|
53 int F77_FUNC (zbesi, ZBESI) (const double&, const double&, |
4062
|
54 const double&, const int&, const int&, |
|
55 double*, double*, int&, int&); |
3146
|
56 |
3887
|
57 int F77_FUNC (zbesk, ZBESK) (const double&, const double&, |
4062
|
58 const double&, const int&, const int&, |
|
59 double*, double*, int&, int&); |
3220
|
60 |
3887
|
61 int F77_FUNC (zbesh, ZBESH) (const double&, const double&, |
4062
|
62 const double&, const int&, const int&, |
|
63 const int&, double*, double*, int&, int&); |
3146
|
64 |
3887
|
65 int F77_FUNC (zairy, ZAIRY) (const double&, const double&, |
4062
|
66 const int&, const int&, |
|
67 double&, double&, int&, int&); |
3220
|
68 |
3887
|
69 int F77_FUNC (zbiry, ZBIRY) (const double&, const double&, |
4062
|
70 const int&, const int&, |
|
71 double&, double&, int&); |
3146
|
72 |
3887
|
73 int F77_FUNC (xdacosh, XDACOSH) (const double&, double&); |
3146
|
74 |
3887
|
75 int F77_FUNC (xdasinh, XDASINH) (const double&, double&); |
3146
|
76 |
3887
|
77 int F77_FUNC (xdatanh, XDATANH) (const double&, double&); |
3146
|
78 |
3887
|
79 int F77_FUNC (xderf, XDERF) (const double&, double&); |
3146
|
80 |
3887
|
81 int F77_FUNC (xderfc, XDERFC) (const double&, double&); |
3146
|
82 |
3887
|
83 int F77_FUNC (xdbetai, XDBETAI) (const double&, const double&, |
4062
|
84 const double&, double&); |
3146
|
85 |
3887
|
86 int F77_FUNC (xdgamma, XDGAMMA) (const double&, double&); |
3146
|
87 |
4004
|
88 int F77_FUNC (xgammainc, XGAMMAINC) (const double&, const double&, double&); |
3146
|
89 |
3887
|
90 int F77_FUNC (dlgams, DLGAMS) (const double&, double&, double&); |
3146
|
91 } |
|
92 |
|
93 #if !defined (HAVE_ACOSH) |
|
94 double |
|
95 acosh (double x) |
|
96 { |
|
97 double retval; |
4180
|
98 F77_FUNC (xdacosh, XDACOSH) (x, retval); |
3146
|
99 return retval; |
|
100 } |
|
101 #endif |
|
102 |
|
103 #if !defined (HAVE_ASINH) |
|
104 double |
|
105 asinh (double x) |
|
106 { |
|
107 double retval; |
4180
|
108 F77_FUNC (xdasinh, XDASINH) (x, retval); |
3146
|
109 return retval; |
|
110 } |
|
111 #endif |
|
112 |
|
113 #if !defined (HAVE_ATANH) |
|
114 double |
|
115 atanh (double x) |
|
116 { |
|
117 double retval; |
4180
|
118 F77_FUNC (xdatanh, XDATANH) (x, retval); |
3146
|
119 return retval; |
|
120 } |
|
121 #endif |
|
122 |
|
123 #if !defined (HAVE_ERF) |
|
124 double |
|
125 erf (double x) |
|
126 { |
|
127 double retval; |
4180
|
128 F77_FUNC (xderf, XDERF) (x, retval); |
3146
|
129 return retval; |
|
130 } |
|
131 #endif |
|
132 |
|
133 #if !defined (HAVE_ERFC) |
|
134 double |
|
135 erfc (double x) |
|
136 { |
|
137 double retval; |
4180
|
138 F77_FUNC (xderfc, XDERFC) (x, retval); |
3146
|
139 return retval; |
|
140 } |
|
141 #endif |
|
142 |
|
143 double |
3156
|
144 xgamma (double x) |
3146
|
145 { |
3156
|
146 double result; |
4180
|
147 F77_FUNC (xdgamma, XDGAMMA) (x, result); |
3156
|
148 return result; |
3146
|
149 } |
|
150 |
|
151 double |
3156
|
152 xlgamma (double x) |
3146
|
153 { |
3156
|
154 double result; |
3146
|
155 double sgngam; |
4497
|
156 |
|
157 if (x < 0) |
|
158 (*current_liboctave_error_handler) |
|
159 ("xlgamma: argument must be nonnegative"); |
|
160 |
4180
|
161 F77_FUNC (dlgams, DLGAMS) (x, result, sgngam); |
4497
|
162 |
3156
|
163 return result; |
3146
|
164 } |
|
165 |
3220
|
166 static inline Complex |
|
167 zbesj (const Complex& z, double alpha, int kode, int& ierr); |
|
168 |
|
169 static inline Complex |
|
170 zbesy (const Complex& z, double alpha, int kode, int& ierr); |
|
171 |
|
172 static inline Complex |
|
173 zbesi (const Complex& z, double alpha, int kode, int& ierr); |
|
174 |
|
175 static inline Complex |
|
176 zbesk (const Complex& z, double alpha, int kode, int& ierr); |
|
177 |
|
178 static inline Complex |
|
179 zbesh1 (const Complex& z, double alpha, int kode, int& ierr); |
|
180 |
|
181 static inline Complex |
|
182 zbesh2 (const Complex& z, double alpha, int kode, int& ierr); |
|
183 |
|
184 static inline Complex |
|
185 bessel_return_value (const Complex& val, int ierr) |
3146
|
186 { |
3220
|
187 static const Complex inf_val = Complex (octave_Inf, octave_Inf); |
|
188 static const Complex nan_val = Complex (octave_NaN, octave_NaN); |
|
189 |
|
190 Complex retval; |
|
191 |
|
192 switch (ierr) |
|
193 { |
|
194 case 0: |
|
195 case 3: |
|
196 retval = val; |
|
197 break; |
|
198 |
|
199 case 2: |
|
200 retval = inf_val; |
|
201 break; |
|
202 |
|
203 default: |
|
204 retval = nan_val; |
|
205 break; |
|
206 } |
|
207 |
3146
|
208 return retval; |
|
209 } |
|
210 |
3220
|
211 static inline Complex |
|
212 zbesj (const Complex& z, double alpha, int kode, int& ierr) |
3146
|
213 { |
3220
|
214 Complex retval; |
|
215 |
|
216 if (alpha >= 0.0) |
|
217 { |
|
218 double yr = 0.0; |
|
219 double yi = 0.0; |
|
220 |
|
221 int nz; |
|
222 |
|
223 double zr = z.real (); |
|
224 double zi = z.imag (); |
|
225 |
4506
|
226 F77_FUNC (zbesj, ZBESJ) (zr, zi, alpha, 2, 1, &yr, &yi, nz, ierr); |
|
227 |
|
228 if (kode != 2) |
|
229 { |
|
230 double expz = exp (std::abs (zi)); |
|
231 yr *= expz; |
|
232 yi *= expz; |
|
233 } |
3220
|
234 |
4490
|
235 if (zi == 0.0 && zr >= 0.0) |
3220
|
236 yi = 0.0; |
|
237 |
|
238 retval = bessel_return_value (Complex (yr, yi), ierr); |
|
239 } |
|
240 else |
|
241 { |
|
242 alpha = -alpha; |
|
243 |
|
244 Complex tmp = cos (M_PI * alpha) * zbesj (z, alpha, kode, ierr); |
|
245 |
|
246 if (ierr == 0 || ierr == 3) |
|
247 { |
|
248 tmp -= sin (M_PI * alpha) * zbesy (z, alpha, kode, ierr); |
|
249 |
|
250 retval = bessel_return_value (tmp, ierr); |
|
251 } |
|
252 else |
|
253 retval = Complex (octave_NaN, octave_NaN); |
|
254 } |
|
255 |
3146
|
256 return retval; |
|
257 } |
|
258 |
3220
|
259 static inline Complex |
|
260 zbesy (const Complex& z, double alpha, int kode, int& ierr) |
3146
|
261 { |
3220
|
262 Complex retval; |
3146
|
263 |
|
264 if (alpha >= 0.0) |
|
265 { |
3220
|
266 double yr = 0.0; |
|
267 double yi = 0.0; |
|
268 |
|
269 int nz; |
|
270 |
|
271 double wr, wi; |
3146
|
272 |
3220
|
273 double zr = z.real (); |
|
274 double zi = z.imag (); |
|
275 |
|
276 ierr = 0; |
|
277 |
|
278 if (zr == 0.0 && zi == 0.0) |
3146
|
279 { |
3220
|
280 yr = -octave_Inf; |
|
281 yi = 0.0; |
|
282 } |
|
283 else |
|
284 { |
4506
|
285 F77_FUNC (zbesy, ZBESY) (zr, zi, alpha, 2, 1, &yr, &yi, nz, |
|
286 &wr, &wi, ierr); |
|
287 |
|
288 if (kode != 2) |
|
289 { |
|
290 double expz = exp (std::abs (zi)); |
|
291 yr *= expz; |
|
292 yi *= expz; |
|
293 } |
3220
|
294 |
4490
|
295 if (zi == 0.0 && zr >= 0.0) |
3220
|
296 yi = 0.0; |
3146
|
297 } |
|
298 |
3220
|
299 return bessel_return_value (Complex (yr, yi), ierr); |
|
300 } |
|
301 else |
|
302 { |
|
303 alpha = -alpha; |
3146
|
304 |
3220
|
305 Complex tmp = cos (M_PI * alpha) * zbesy (z, alpha, kode, ierr); |
3146
|
306 |
3220
|
307 if (ierr == 0 || ierr == 3) |
3146
|
308 { |
3220
|
309 tmp += sin (M_PI * alpha) * zbesj (z, alpha, kode, ierr); |
|
310 |
|
311 retval = bessel_return_value (tmp, ierr); |
|
312 } |
|
313 else |
|
314 retval = Complex (octave_NaN, octave_NaN); |
|
315 } |
|
316 |
|
317 return retval; |
|
318 } |
|
319 |
|
320 static inline Complex |
|
321 zbesi (const Complex& z, double alpha, int kode, int& ierr) |
|
322 { |
|
323 Complex retval; |
3146
|
324 |
3220
|
325 if (alpha >= 0.0) |
|
326 { |
|
327 double yr = 0.0; |
|
328 double yi = 0.0; |
|
329 |
|
330 int nz; |
3146
|
331 |
3220
|
332 double zr = z.real (); |
|
333 double zi = z.imag (); |
|
334 |
4506
|
335 F77_FUNC (zbesi, ZBESI) (zr, zi, alpha, 2, 1, &yr, &yi, nz, ierr); |
|
336 |
|
337 if (kode != 2) |
|
338 { |
|
339 double expz = exp (std::abs (zr)); |
|
340 yr *= expz; |
|
341 yi *= expz; |
|
342 } |
3146
|
343 |
4490
|
344 if (zi == 0.0 && zr >= 0.0) |
3220
|
345 yi = 0.0; |
|
346 |
|
347 retval = bessel_return_value (Complex (yr, yi), ierr); |
3146
|
348 } |
|
349 else |
3220
|
350 { |
|
351 alpha = -alpha; |
|
352 |
|
353 Complex tmp = zbesi (z, alpha, kode, ierr); |
|
354 |
|
355 if (ierr == 0 || ierr == 3) |
|
356 { |
|
357 tmp += (2.0 / M_PI) * sin (M_PI * alpha) |
|
358 * zbesk (z, alpha, kode, ierr); |
|
359 |
|
360 retval = bessel_return_value (tmp, ierr); |
|
361 } |
|
362 else |
|
363 retval = Complex (octave_NaN, octave_NaN); |
|
364 } |
|
365 |
|
366 return retval; |
|
367 } |
|
368 |
|
369 static inline Complex |
|
370 zbesk (const Complex& z, double alpha, int kode, int& ierr) |
|
371 { |
|
372 Complex retval; |
|
373 |
|
374 if (alpha >= 0.0) |
|
375 { |
|
376 double yr = 0.0; |
|
377 double yi = 0.0; |
|
378 |
|
379 int nz; |
|
380 |
|
381 double zr = z.real (); |
|
382 double zi = z.imag (); |
|
383 |
|
384 ierr = 0; |
|
385 |
|
386 if (zr == 0.0 && zi == 0.0) |
|
387 { |
|
388 yr = octave_Inf; |
|
389 yi = 0.0; |
|
390 } |
|
391 else |
|
392 { |
4506
|
393 F77_FUNC (zbesk, ZBESK) (zr, zi, alpha, 2, 1, &yr, &yi, nz, ierr); |
|
394 |
|
395 if (kode != 2) |
|
396 { |
|
397 Complex expz = exp (-z); |
|
398 |
|
399 double rexpz = real (expz); |
|
400 double iexpz = imag (expz); |
|
401 |
|
402 double tmp = yr*rexpz - yi*iexpz; |
|
403 |
|
404 yi = yr*iexpz + yi*rexpz; |
|
405 yr = tmp; |
|
406 } |
3220
|
407 |
4490
|
408 if (zi == 0.0 && zr >= 0.0) |
3220
|
409 yi = 0.0; |
|
410 } |
|
411 |
|
412 retval = bessel_return_value (Complex (yr, yi), ierr); |
|
413 } |
|
414 else |
|
415 { |
|
416 Complex tmp = zbesk (z, -alpha, kode, ierr); |
|
417 |
|
418 retval = bessel_return_value (tmp, ierr); |
|
419 } |
3146
|
420 |
|
421 return retval; |
|
422 } |
|
423 |
3220
|
424 static inline Complex |
|
425 zbesh1 (const Complex& z, double alpha, int kode, int& ierr) |
3146
|
426 { |
3220
|
427 Complex retval; |
3146
|
428 |
3220
|
429 if (alpha >= 0.0) |
3146
|
430 { |
3220
|
431 double yr = 0.0; |
|
432 double yi = 0.0; |
|
433 |
|
434 int nz; |
|
435 |
|
436 double zr = z.real (); |
|
437 double zi = z.imag (); |
3146
|
438 |
4506
|
439 F77_FUNC (zbesh, ZBESH) (zr, zi, alpha, 2, 1, 1, &yr, &yi, nz, ierr); |
|
440 |
|
441 if (kode != 2) |
|
442 { |
|
443 Complex expz = exp (Complex (0.0, 1.0) * z); |
|
444 |
|
445 double rexpz = real (expz); |
|
446 double iexpz = imag (expz); |
|
447 |
|
448 double tmp = yr*rexpz - yi*iexpz; |
|
449 |
|
450 yi = yr*iexpz + yi*rexpz; |
|
451 yr = tmp; |
|
452 } |
3146
|
453 |
3220
|
454 retval = bessel_return_value (Complex (yr, yi), ierr); |
|
455 } |
|
456 else |
|
457 { |
|
458 alpha = -alpha; |
|
459 |
|
460 static const Complex eye = Complex (0.0, 1.0); |
3146
|
461 |
3220
|
462 Complex tmp = exp (M_PI * alpha * eye) * zbesh1 (z, alpha, kode, ierr); |
3146
|
463 |
3220
|
464 retval = bessel_return_value (tmp, ierr); |
|
465 } |
3146
|
466 |
3220
|
467 return retval; |
|
468 } |
3146
|
469 |
3220
|
470 static inline Complex |
|
471 zbesh2 (const Complex& z, double alpha, int kode, int& ierr) |
|
472 { |
|
473 Complex retval; |
3146
|
474 |
3220
|
475 if (alpha >= 0.0) |
|
476 { |
|
477 double yr = 0.0; |
|
478 double yi = 0.0; |
|
479 |
|
480 int nz; |
3146
|
481 |
3220
|
482 double zr = z.real (); |
|
483 double zi = z.imag (); |
3146
|
484 |
4506
|
485 F77_FUNC (zbesh, ZBESH) (zr, zi, alpha, 2, 2, 1, &yr, &yi, nz, ierr); |
|
486 |
|
487 if (kode != 2) |
|
488 { |
|
489 Complex expz = exp (-Complex (0.0, 1.0) * z); |
|
490 |
|
491 double rexpz = real (expz); |
|
492 double iexpz = imag (expz); |
|
493 |
|
494 double tmp = yr*rexpz - yi*iexpz; |
|
495 |
|
496 yi = yr*iexpz + yi*rexpz; |
|
497 yr = tmp; |
|
498 } |
3220
|
499 |
|
500 retval = bessel_return_value (Complex (yr, yi), ierr); |
3146
|
501 } |
|
502 else |
3220
|
503 { |
|
504 alpha = -alpha; |
|
505 |
|
506 static const Complex eye = Complex (0.0, 1.0); |
|
507 |
|
508 Complex tmp = exp (-M_PI * alpha * eye) * zbesh2 (z, alpha, kode, ierr); |
|
509 |
|
510 retval = bessel_return_value (tmp, ierr); |
|
511 } |
|
512 |
|
513 return retval; |
|
514 } |
|
515 |
|
516 typedef Complex (*fptr) (const Complex&, double, int, int&); |
|
517 |
|
518 static inline Complex |
|
519 do_bessel (fptr f, const char *, double alpha, const Complex& x, |
|
520 bool scaled, int& ierr) |
|
521 { |
|
522 Complex retval; |
|
523 |
|
524 retval = f (x, alpha, (scaled ? 2 : 1), ierr); |
|
525 |
|
526 return retval; |
|
527 } |
|
528 |
|
529 static inline ComplexMatrix |
|
530 do_bessel (fptr f, const char *, double alpha, const ComplexMatrix& x, |
|
531 bool scaled, Array2<int>& ierr) |
|
532 { |
|
533 int nr = x.rows (); |
|
534 int nc = x.cols (); |
|
535 |
|
536 ComplexMatrix retval (nr, nc); |
|
537 |
|
538 ierr.resize (nr, nc); |
|
539 |
|
540 for (int j = 0; j < nc; j++) |
|
541 for (int i = 0; i < nr; i++) |
|
542 retval(i,j) = f (x(i,j), alpha, (scaled ? 2 : 1), ierr(i,j)); |
|
543 |
|
544 return retval; |
|
545 } |
|
546 |
|
547 static inline ComplexMatrix |
|
548 do_bessel (fptr f, const char *, const Matrix& alpha, const Complex& x, |
|
549 bool scaled, Array2<int>& ierr) |
|
550 { |
|
551 int nr = alpha.rows (); |
|
552 int nc = alpha.cols (); |
|
553 |
|
554 ComplexMatrix retval (nr, nc); |
|
555 |
|
556 ierr.resize (nr, nc); |
|
557 |
|
558 for (int j = 0; j < nc; j++) |
|
559 for (int i = 0; i < nr; i++) |
|
560 retval(i,j) = f (x, alpha(i,j), (scaled ? 2 : 1), ierr(i,j)); |
3146
|
561 |
|
562 return retval; |
|
563 } |
|
564 |
3220
|
565 static inline ComplexMatrix |
|
566 do_bessel (fptr f, const char *fn, const Matrix& alpha, |
|
567 const ComplexMatrix& x, bool scaled, Array2<int>& ierr) |
3146
|
568 { |
3220
|
569 ComplexMatrix retval; |
|
570 |
|
571 int x_nr = x.rows (); |
|
572 int x_nc = x.cols (); |
|
573 |
|
574 int alpha_nr = alpha.rows (); |
|
575 int alpha_nc = alpha.cols (); |
|
576 |
|
577 if (x_nr == alpha_nr && x_nc == alpha_nc) |
|
578 { |
|
579 int nr = x_nr; |
|
580 int nc = x_nc; |
|
581 |
|
582 retval.resize (nr, nc); |
|
583 |
|
584 ierr.resize (nr, nc); |
|
585 |
|
586 for (int j = 0; j < nc; j++) |
|
587 for (int i = 0; i < nr; i++) |
|
588 retval(i,j) = f (x(i,j), alpha(i,j), (scaled ? 2 : 1), ierr(i,j)); |
|
589 } |
|
590 else |
|
591 (*current_liboctave_error_handler) |
|
592 ("%s: the sizes of alpha and x must conform", fn); |
|
593 |
|
594 return retval; |
3146
|
595 } |
|
596 |
3220
|
597 static inline ComplexMatrix |
|
598 do_bessel (fptr f, const char *, const RowVector& alpha, |
|
599 const ComplexColumnVector& x, bool scaled, Array2<int>& ierr) |
3146
|
600 { |
3220
|
601 int nr = x.length (); |
|
602 int nc = alpha.length (); |
|
603 |
|
604 ComplexMatrix retval (nr, nc); |
3146
|
605 |
3220
|
606 ierr.resize (nr, nc); |
|
607 |
|
608 for (int j = 0; j < nc; j++) |
|
609 for (int i = 0; i < nr; i++) |
|
610 retval(i,j) = f (x(i), alpha(j), (scaled ? 2 : 1), ierr(i,j)); |
|
611 |
|
612 return retval; |
3146
|
613 } |
|
614 |
3220
|
615 #define SS_BESSEL(name, fcn) \ |
|
616 Complex \ |
|
617 name (double alpha, const Complex& x, bool scaled, int& ierr) \ |
|
618 { \ |
|
619 return do_bessel (fcn, #name, alpha, x, scaled, ierr); \ |
|
620 } |
|
621 |
|
622 #define SM_BESSEL(name, fcn) \ |
|
623 ComplexMatrix \ |
|
624 name (double alpha, const ComplexMatrix& x, bool scaled, \ |
|
625 Array2<int>& ierr) \ |
|
626 { \ |
|
627 return do_bessel (fcn, #name, alpha, x, scaled, ierr); \ |
|
628 } |
|
629 |
|
630 #define MS_BESSEL(name, fcn) \ |
|
631 ComplexMatrix \ |
|
632 name (const Matrix& alpha, const Complex& x, bool scaled, \ |
|
633 Array2<int>& ierr) \ |
|
634 { \ |
|
635 return do_bessel (fcn, #name, alpha, x, scaled, ierr); \ |
|
636 } |
|
637 |
|
638 #define MM_BESSEL(name, fcn) \ |
|
639 ComplexMatrix \ |
|
640 name (const Matrix& alpha, const ComplexMatrix& x, bool scaled, \ |
|
641 Array2<int>& ierr) \ |
|
642 { \ |
|
643 return do_bessel (fcn, #name, alpha, x, scaled, ierr); \ |
|
644 } |
|
645 |
|
646 #define RC_BESSEL(name, fcn) \ |
|
647 ComplexMatrix \ |
|
648 name (const RowVector& alpha, const ComplexColumnVector& x, bool scaled, \ |
|
649 Array2<int>& ierr) \ |
|
650 { \ |
|
651 return do_bessel (fcn, #name, alpha, x, scaled, ierr); \ |
|
652 } |
|
653 |
|
654 #define ALL_BESSEL(name, fcn) \ |
|
655 SS_BESSEL (name, fcn) \ |
|
656 SM_BESSEL (name, fcn) \ |
|
657 MS_BESSEL (name, fcn) \ |
|
658 MM_BESSEL (name, fcn) \ |
|
659 RC_BESSEL (name, fcn) |
|
660 |
|
661 ALL_BESSEL (besselj, zbesj) |
|
662 ALL_BESSEL (bessely, zbesy) |
|
663 ALL_BESSEL (besseli, zbesi) |
|
664 ALL_BESSEL (besselk, zbesk) |
|
665 ALL_BESSEL (besselh1, zbesh1) |
|
666 ALL_BESSEL (besselh2, zbesh2) |
|
667 |
|
668 Complex |
|
669 airy (const Complex& z, bool deriv, bool scaled, int& ierr) |
3146
|
670 { |
3220
|
671 double ar = 0.0; |
|
672 double ai = 0.0; |
|
673 |
|
674 int nz; |
|
675 |
|
676 double zr = z.real (); |
|
677 double zi = z.imag (); |
3146
|
678 |
3220
|
679 int id = deriv ? 1 : 0; |
|
680 |
4506
|
681 F77_FUNC (zairy, ZAIRY) (zr, zi, id, 2, ar, ai, nz, ierr); |
|
682 |
|
683 if (! scaled) |
|
684 { |
|
685 Complex expz = exp (- 2.0 / 3.0 * z * sqrt(z)); |
3220
|
686 |
4506
|
687 double rexpz = real (expz); |
|
688 double iexpz = imag (expz); |
|
689 |
|
690 double tmp = ar*rexpz - ai*iexpz; |
|
691 |
|
692 ai = ar*iexpz + ai*rexpz; |
|
693 ar = tmp; |
|
694 } |
3220
|
695 |
4490
|
696 if (zi == 0.0 && (! scaled || zr >= 0.0)) |
3225
|
697 ai = 0.0; |
|
698 |
3220
|
699 return bessel_return_value (Complex (ar, ai), ierr); |
3146
|
700 } |
|
701 |
3220
|
702 Complex |
|
703 biry (const Complex& z, bool deriv, bool scaled, int& ierr) |
3146
|
704 { |
3220
|
705 double ar = 0.0; |
|
706 double ai = 0.0; |
|
707 |
|
708 double zr = z.real (); |
|
709 double zi = z.imag (); |
|
710 |
|
711 int id = deriv ? 1 : 0; |
|
712 |
4506
|
713 F77_FUNC (zbiry, ZBIRY) (zr, zi, id, 2, ar, ai, ierr); |
|
714 |
|
715 if (! scaled) |
|
716 { |
|
717 Complex expz = exp (std::abs (real (2.0 / 3.0 * z * sqrt (z)))); |
3220
|
718 |
4506
|
719 double rexpz = real (expz); |
|
720 double iexpz = imag (expz); |
|
721 |
|
722 double tmp = ar*rexpz - ai*iexpz; |
|
723 |
|
724 ai = ar*iexpz + ai*rexpz; |
|
725 ar = tmp; |
|
726 } |
3220
|
727 |
4490
|
728 if (zi == 0.0 && (! scaled || zr >= 0.0)) |
3225
|
729 ai = 0.0; |
|
730 |
3220
|
731 return bessel_return_value (Complex (ar, ai), ierr); |
3146
|
732 } |
|
733 |
3220
|
734 ComplexMatrix |
|
735 airy (const ComplexMatrix& z, bool deriv, bool scaled, Array2<int>& ierr) |
3146
|
736 { |
3220
|
737 int nr = z.rows (); |
|
738 int nc = z.cols (); |
|
739 |
|
740 ComplexMatrix retval (nr, nc); |
|
741 |
|
742 ierr.resize (nr, nc); |
|
743 |
|
744 for (int j = 0; j < nc; j++) |
|
745 for (int i = 0; i < nr; i++) |
|
746 retval(i,j) = airy (z(i,j), deriv, scaled, ierr(i,j)); |
|
747 |
|
748 return retval; |
3146
|
749 } |
|
750 |
3220
|
751 ComplexMatrix |
|
752 biry (const ComplexMatrix& z, bool deriv, bool scaled, Array2<int>& ierr) |
3146
|
753 { |
3220
|
754 int nr = z.rows (); |
|
755 int nc = z.cols (); |
|
756 |
|
757 ComplexMatrix retval (nr, nc); |
|
758 |
|
759 ierr.resize (nr, nc); |
|
760 |
|
761 for (int j = 0; j < nc; j++) |
|
762 for (int i = 0; i < nr; i++) |
|
763 retval(i,j) = biry (z(i,j), deriv, scaled, ierr(i,j)); |
|
764 |
|
765 return retval; |
3146
|
766 } |
|
767 |
|
768 static void |
|
769 gripe_betainc_nonconformant (int r1, int c1, int r2, int c2, int r3, |
|
770 int c3) |
|
771 { |
|
772 (*current_liboctave_error_handler) |
|
773 ("betainc: nonconformant arguments (x is %dx%d, a is %dx%d, b is %dx%d)", |
|
774 r1, c1, r2, c2, r3, c3); |
|
775 } |
|
776 |
|
777 double |
|
778 betainc (double x, double a, double b) |
|
779 { |
|
780 double retval; |
4180
|
781 F77_FUNC (xdbetai, XDBETAI) (x, a, b, retval); |
3146
|
782 return retval; |
|
783 } |
|
784 |
|
785 Matrix |
|
786 betainc (double x, double a, const Matrix& b) |
|
787 { |
|
788 int nr = b.rows (); |
|
789 int nc = b.cols (); |
|
790 |
|
791 Matrix retval (nr, nc); |
|
792 |
|
793 for (int j = 0; j < nc; j++) |
|
794 for (int i = 0; i < nr; i++) |
|
795 retval(i,j) = betainc (x, a, b(i,j)); |
|
796 |
|
797 return retval; |
|
798 } |
|
799 |
|
800 Matrix |
|
801 betainc (double x, const Matrix& a, double b) |
|
802 { |
|
803 int nr = a.rows (); |
|
804 int nc = a.cols (); |
|
805 |
|
806 Matrix retval (nr, nc); |
|
807 |
|
808 for (int j = 0; j < nc; j++) |
|
809 for (int i = 0; i < nr; i++) |
|
810 retval(i,j) = betainc (x, a(i,j), b); |
|
811 |
|
812 return retval; |
|
813 } |
|
814 |
|
815 Matrix |
|
816 betainc (double x, const Matrix& a, const Matrix& b) |
|
817 { |
|
818 Matrix retval; |
|
819 |
|
820 int a_nr = a.rows (); |
|
821 int a_nc = a.cols (); |
|
822 |
|
823 int b_nr = b.rows (); |
|
824 int b_nc = b.cols (); |
|
825 |
|
826 if (a_nr == b_nr && a_nc == b_nc) |
|
827 { |
|
828 retval.resize (a_nr, a_nc); |
|
829 |
|
830 for (int j = 0; j < a_nc; j++) |
|
831 for (int i = 0; i < a_nr; i++) |
|
832 retval(i,j) = betainc (x, a(i,j), b(i,j)); |
|
833 } |
|
834 else |
|
835 gripe_betainc_nonconformant (1, 1, a_nr, a_nc, b_nr, b_nc); |
|
836 |
|
837 return retval; |
|
838 } |
|
839 |
|
840 Matrix |
|
841 betainc (const Matrix& x, double a, double b) |
|
842 { |
|
843 int nr = x.rows (); |
|
844 int nc = x.cols (); |
|
845 |
|
846 Matrix retval (nr, nc); |
|
847 |
|
848 for (int j = 0; j < nc; j++) |
|
849 for (int i = 0; i < nr; i++) |
|
850 retval(i,j) = betainc (x(i,j), a, b); |
|
851 |
|
852 return retval; |
|
853 } |
|
854 |
|
855 Matrix |
|
856 betainc (const Matrix& x, double a, const Matrix& b) |
|
857 { |
|
858 Matrix retval; |
|
859 |
|
860 int nr = x.rows (); |
|
861 int nc = x.cols (); |
|
862 |
|
863 int b_nr = b.rows (); |
|
864 int b_nc = b.cols (); |
|
865 |
|
866 if (nr == b_nr && nc == b_nc) |
|
867 { |
|
868 retval.resize (nr, nc); |
|
869 |
|
870 for (int j = 0; j < nc; j++) |
|
871 for (int i = 0; i < nr; i++) |
|
872 retval(i,j) = betainc (x(i,j), a, b(i,j)); |
|
873 } |
|
874 else |
|
875 gripe_betainc_nonconformant (nr, nc, 1, 1, b_nr, b_nc); |
|
876 |
|
877 return retval; |
|
878 } |
|
879 |
|
880 Matrix |
|
881 betainc (const Matrix& x, const Matrix& a, double b) |
|
882 { |
|
883 Matrix retval; |
|
884 |
|
885 int nr = x.rows (); |
|
886 int nc = x.cols (); |
|
887 |
|
888 int a_nr = a.rows (); |
|
889 int a_nc = a.cols (); |
|
890 |
|
891 if (nr == a_nr && nc == a_nc) |
|
892 { |
|
893 retval.resize (nr, nc); |
|
894 |
|
895 for (int j = 0; j < nc; j++) |
|
896 for (int i = 0; i < nr; i++) |
|
897 retval(i,j) = betainc (x(i,j), a(i,j), b); |
|
898 } |
|
899 else |
|
900 gripe_betainc_nonconformant (nr, nc, a_nr, a_nc, 1, 1); |
|
901 |
|
902 return retval; |
|
903 } |
|
904 |
|
905 Matrix |
|
906 betainc (const Matrix& x, const Matrix& a, const Matrix& b) |
|
907 { |
|
908 Matrix retval; |
|
909 |
|
910 int nr = x.rows (); |
|
911 int nc = x.cols (); |
|
912 |
|
913 int a_nr = a.rows (); |
|
914 int a_nc = a.cols (); |
|
915 |
|
916 int b_nr = b.rows (); |
|
917 int b_nc = b.cols (); |
|
918 |
|
919 if (nr == a_nr && nr == b_nr && nc == a_nc && nc == b_nc) |
|
920 { |
|
921 retval.resize (nr, nc); |
|
922 |
|
923 for (int j = 0; j < nc; j++) |
|
924 for (int i = 0; i < nr; i++) |
|
925 retval(i,j) = betainc (x(i,j), a(i,j), b(i,j)); |
|
926 } |
|
927 else |
|
928 gripe_betainc_nonconformant (nr, nc, a_nr, a_nc, b_nr, b_nc); |
|
929 |
|
930 return retval; |
|
931 } |
|
932 |
3164
|
933 // XXX FIXME XXX -- there is still room for improvement here... |
|
934 |
3146
|
935 double |
4004
|
936 gammainc (double x, double a, bool& err) |
3146
|
937 { |
|
938 double retval; |
3164
|
939 |
4004
|
940 err = false; |
3164
|
941 |
4004
|
942 if (a < 0.0 || x < 0.0) |
|
943 { |
|
944 (*current_liboctave_error_handler) |
|
945 ("gammainc: A and X must be non-negative"); |
|
946 |
|
947 err = true; |
|
948 } |
|
949 else |
4180
|
950 F77_FUNC (xgammainc, XGAMMAINC) (a, x, retval); |
3164
|
951 |
3146
|
952 return retval; |
|
953 } |
|
954 |
|
955 Matrix |
|
956 gammainc (double x, const Matrix& a) |
|
957 { |
|
958 int nr = a.rows (); |
|
959 int nc = a.cols (); |
|
960 |
4004
|
961 Matrix result (nr, nc); |
|
962 Matrix retval; |
|
963 |
|
964 bool err; |
3146
|
965 |
|
966 for (int j = 0; j < nc; j++) |
|
967 for (int i = 0; i < nr; i++) |
4004
|
968 { |
|
969 result(i,j) = gammainc (x, a(i,j), err); |
|
970 |
|
971 if (err) |
|
972 goto done; |
|
973 } |
|
974 |
|
975 retval = result; |
|
976 |
|
977 done: |
3146
|
978 |
|
979 return retval; |
|
980 } |
|
981 |
|
982 Matrix |
|
983 gammainc (const Matrix& x, double a) |
|
984 { |
|
985 int nr = x.rows (); |
|
986 int nc = x.cols (); |
|
987 |
4004
|
988 Matrix result (nr, nc); |
|
989 Matrix retval; |
|
990 |
|
991 bool err; |
3146
|
992 |
|
993 for (int j = 0; j < nc; j++) |
|
994 for (int i = 0; i < nr; i++) |
4004
|
995 { |
|
996 result(i,j) = gammainc (x(i,j), a, err); |
|
997 |
|
998 if (err) |
|
999 goto done; |
|
1000 } |
|
1001 |
|
1002 retval = result; |
|
1003 |
|
1004 done: |
3146
|
1005 |
|
1006 return retval; |
|
1007 } |
|
1008 |
|
1009 Matrix |
|
1010 gammainc (const Matrix& x, const Matrix& a) |
|
1011 { |
4004
|
1012 Matrix result; |
3146
|
1013 Matrix retval; |
|
1014 |
|
1015 int nr = x.rows (); |
|
1016 int nc = x.cols (); |
|
1017 |
|
1018 int a_nr = a.rows (); |
|
1019 int a_nc = a.cols (); |
|
1020 |
|
1021 if (nr == a_nr && nc == a_nc) |
|
1022 { |
4004
|
1023 result.resize (nr, nc); |
|
1024 |
|
1025 bool err; |
3146
|
1026 |
|
1027 for (int j = 0; j < nc; j++) |
|
1028 for (int i = 0; i < nr; i++) |
4004
|
1029 { |
|
1030 result(i,j) = gammainc (x(i,j), a(i,j), err); |
|
1031 |
|
1032 if (err) |
|
1033 goto done; |
|
1034 } |
|
1035 |
|
1036 retval = result; |
3146
|
1037 } |
|
1038 else |
|
1039 (*current_liboctave_error_handler) |
|
1040 ("gammainc: nonconformant arguments (arg 1 is %dx%d, arg 2 is %dx%d)", |
|
1041 nr, nc, a_nr, a_nc); |
|
1042 |
4004
|
1043 done: |
|
1044 |
3146
|
1045 return retval; |
|
1046 } |
|
1047 |
|
1048 /* |
|
1049 ;;; Local Variables: *** |
|
1050 ;;; mode: C++ *** |
|
1051 ;;; End: *** |
|
1052 */ |