besselj.cc

Go to the documentation of this file.
00001 /*
00002 
00003 Copyright (C) 1997-2012 John W. Eaton
00004 
00005 This file is part of Octave.
00006 
00007 Octave is free software; you can redistribute it and/or modify it
00008 under the terms of the GNU General Public License as published by the
00009 Free Software Foundation; either version 3 of the License, or (at your
00010 option) any later version.
00011 
00012 Octave is distributed in the hope that it will be useful, but WITHOUT
00013 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00014 FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
00015 for more details.
00016 
00017 You should have received a copy of the GNU General Public License
00018 along with Octave; see the file COPYING.  If not, see
00019 <http://www.gnu.org/licenses/>.
00020 
00021 */
00022 
00023 #ifdef HAVE_CONFIG_H
00024 #include <config.h>
00025 #endif
00026 
00027 #include "lo-specfun.h"
00028 #include "quit.h"
00029 
00030 #include "defun-dld.h"
00031 #include "error.h"
00032 #include "gripes.h"
00033 #include "oct-obj.h"
00034 #include "utils.h"
00035 
00036 enum bessel_type
00037 {
00038   BESSEL_J,
00039   BESSEL_Y,
00040   BESSEL_I,
00041   BESSEL_K,
00042   BESSEL_H1,
00043   BESSEL_H2
00044 };
00045 
00046 #define DO_BESSEL(type, alpha, x, scaled, ierr, result) \
00047   do \
00048     { \
00049       switch (type) \
00050         { \
00051           case BESSEL_J: \
00052             result = besselj (alpha, x, scaled, ierr); \
00053             break; \
00054  \
00055           case BESSEL_Y: \
00056             result = bessely (alpha, x, scaled, ierr); \
00057             break; \
00058  \
00059           case BESSEL_I: \
00060             result = besseli (alpha, x, scaled, ierr); \
00061             break; \
00062  \
00063           case BESSEL_K: \
00064             result = besselk (alpha, x, scaled, ierr); \
00065             break; \
00066  \
00067           case BESSEL_H1: \
00068             result = besselh1 (alpha, x, scaled, ierr); \
00069             break; \
00070  \
00071           case BESSEL_H2: \
00072             result = besselh2 (alpha, x, scaled, ierr); \
00073             break; \
00074  \
00075           default: \
00076             break; \
00077         } \
00078     } \
00079   while (0)
00080 
00081 static void
00082 gripe_bessel_arg (const char *fn, const char *arg)
00083 {
00084   error ("%s: expecting scalar or matrix as %s argument", fn, arg);
00085 }
00086 
00087 octave_value_list
00088 do_bessel (enum bessel_type type, const char *fn,
00089            const octave_value_list& args, int nargout)
00090 {
00091   octave_value_list retval;
00092 
00093   int nargin = args.length ();
00094 
00095   if (nargin == 2 || nargin == 3)
00096     {
00097       bool scaled = (nargin == 3);
00098 
00099       octave_value alpha_arg = args(0);
00100       octave_value x_arg = args(1);
00101 
00102       if (alpha_arg.is_single_type () || x_arg.is_single_type ())
00103         {
00104           if (alpha_arg.is_scalar_type ())
00105             {
00106               float alpha = args(0).float_value ();
00107 
00108               if (! error_state)
00109                 {
00110                   if (x_arg.is_scalar_type ())
00111                     {
00112                       FloatComplex x = x_arg.float_complex_value ();
00113 
00114                       if (! error_state)
00115                         {
00116                           octave_idx_type ierr;
00117                           octave_value result;
00118 
00119                           DO_BESSEL (type, alpha, x, scaled, ierr, result);
00120 
00121                           if (nargout > 1)
00122                             retval(1) = static_cast<float> (ierr);
00123 
00124                           retval(0) = result;
00125                         }
00126                       else
00127                         gripe_bessel_arg (fn, "second");
00128                     }
00129                   else
00130                     {
00131                       FloatComplexNDArray x = x_arg.float_complex_array_value ();
00132 
00133                       if (! error_state)
00134                         {
00135                           Array<octave_idx_type> ierr;
00136                           octave_value result;
00137 
00138                           DO_BESSEL (type, alpha, x, scaled, ierr, result);
00139 
00140                           if (nargout > 1)
00141                             retval(1) = NDArray (ierr);
00142 
00143                           retval(0) = result;
00144                         }
00145                       else
00146                         gripe_bessel_arg (fn, "second");
00147                     }
00148                 }
00149               else
00150                 gripe_bessel_arg (fn, "first");
00151             }
00152           else
00153             {
00154               dim_vector dv0 = args(0).dims ();
00155               dim_vector dv1 = args(1).dims ();
00156 
00157               bool args0_is_row_vector = (dv0 (1) == dv0.numel ());
00158               bool args1_is_col_vector = (dv1 (0) == dv1.numel ());
00159 
00160               if (args0_is_row_vector && args1_is_col_vector)
00161                 {
00162                   FloatRowVector ralpha = args(0).float_row_vector_value ();
00163 
00164                   if (! error_state)
00165                     {
00166                       FloatComplexColumnVector cx =
00167                         x_arg.float_complex_column_vector_value ();
00168 
00169                       if (! error_state)
00170                         {
00171                           Array<octave_idx_type> ierr;
00172                           octave_value result;
00173 
00174                           DO_BESSEL (type, ralpha, cx, scaled, ierr, result);
00175 
00176                           if (nargout > 1)
00177                             retval(1) = NDArray (ierr);
00178 
00179                           retval(0) = result;
00180                         }
00181                       else
00182                         gripe_bessel_arg (fn, "second");
00183                     }
00184                   else
00185                     gripe_bessel_arg (fn, "first");
00186                 }
00187               else
00188                 {
00189                   FloatNDArray alpha = args(0).float_array_value ();
00190 
00191                   if (! error_state)
00192                     {
00193                       if (x_arg.is_scalar_type ())
00194                         {
00195                           FloatComplex x = x_arg.float_complex_value ();
00196 
00197                           if (! error_state)
00198                             {
00199                               Array<octave_idx_type> ierr;
00200                               octave_value result;
00201 
00202                               DO_BESSEL (type, alpha, x, scaled, ierr, result);
00203 
00204                               if (nargout > 1)
00205                                 retval(1) = NDArray (ierr);
00206 
00207                               retval(0) = result;
00208                             }
00209                           else
00210                             gripe_bessel_arg (fn, "second");
00211                         }
00212                       else
00213                         {
00214                           FloatComplexNDArray x = x_arg.float_complex_array_value ();
00215 
00216                           if (! error_state)
00217                             {
00218                               Array<octave_idx_type> ierr;
00219                               octave_value result;
00220 
00221                               DO_BESSEL (type, alpha, x, scaled, ierr, result);
00222 
00223                               if (nargout > 1)
00224                                 retval(1) = NDArray (ierr);
00225 
00226                               retval(0) = result;
00227                             }
00228                           else
00229                             gripe_bessel_arg (fn, "second");
00230                         }
00231                     }
00232                   else
00233                     gripe_bessel_arg (fn, "first");
00234                 }
00235             }
00236         }
00237       else
00238         {
00239           if (alpha_arg.is_scalar_type ())
00240             {
00241               double alpha = args(0).double_value ();
00242 
00243               if (! error_state)
00244                 {
00245                   if (x_arg.is_scalar_type ())
00246                     {
00247                       Complex x = x_arg.complex_value ();
00248 
00249                       if (! error_state)
00250                         {
00251                           octave_idx_type ierr;
00252                           octave_value result;
00253 
00254                           DO_BESSEL (type, alpha, x, scaled, ierr, result);
00255 
00256                           if (nargout > 1)
00257                             retval(1) = static_cast<double> (ierr);
00258 
00259                           retval(0) = result;
00260                         }
00261                       else
00262                         gripe_bessel_arg (fn, "second");
00263                     }
00264                   else
00265                     {
00266                       ComplexNDArray x = x_arg.complex_array_value ();
00267 
00268                       if (! error_state)
00269                         {
00270                           Array<octave_idx_type> ierr;
00271                           octave_value result;
00272 
00273                           DO_BESSEL (type, alpha, x, scaled, ierr, result);
00274 
00275                           if (nargout > 1)
00276                             retval(1) = NDArray (ierr);
00277 
00278                           retval(0) = result;
00279                         }
00280                       else
00281                         gripe_bessel_arg (fn, "second");
00282                     }
00283                 }
00284               else
00285                 gripe_bessel_arg (fn, "first");
00286             }
00287           else
00288             {
00289               dim_vector dv0 = args(0).dims ();
00290               dim_vector dv1 = args(1).dims ();
00291 
00292               bool args0_is_row_vector = (dv0 (1) == dv0.numel ());
00293               bool args1_is_col_vector = (dv1 (0) == dv1.numel ());
00294 
00295               if (args0_is_row_vector && args1_is_col_vector)
00296                 {
00297                   RowVector ralpha = args(0).row_vector_value ();
00298 
00299                   if (! error_state)
00300                     {
00301                       ComplexColumnVector cx =
00302                         x_arg.complex_column_vector_value ();
00303 
00304                       if (! error_state)
00305                         {
00306                           Array<octave_idx_type> ierr;
00307                           octave_value result;
00308 
00309                           DO_BESSEL (type, ralpha, cx, scaled, ierr, result);
00310 
00311                           if (nargout > 1)
00312                             retval(1) = NDArray (ierr);
00313 
00314                           retval(0) = result;
00315                         }
00316                       else
00317                         gripe_bessel_arg (fn, "second");
00318                     }
00319                   else
00320                     gripe_bessel_arg (fn, "first");
00321                 }
00322               else
00323                 {
00324                   NDArray alpha = args(0).array_value ();
00325 
00326                   if (! error_state)
00327                     {
00328                       if (x_arg.is_scalar_type ())
00329                         {
00330                           Complex x = x_arg.complex_value ();
00331 
00332                           if (! error_state)
00333                             {
00334                               Array<octave_idx_type> ierr;
00335                               octave_value result;
00336 
00337                               DO_BESSEL (type, alpha, x, scaled, ierr, result);
00338 
00339                               if (nargout > 1)
00340                                 retval(1) = NDArray (ierr);
00341 
00342                               retval(0) = result;
00343                             }
00344                           else
00345                             gripe_bessel_arg (fn, "second");
00346                         }
00347                       else
00348                         {
00349                           ComplexNDArray x = x_arg.complex_array_value ();
00350 
00351                           if (! error_state)
00352                             {
00353                               Array<octave_idx_type> ierr;
00354                               octave_value result;
00355 
00356                               DO_BESSEL (type, alpha, x, scaled, ierr, result);
00357 
00358                               if (nargout > 1)
00359                                 retval(1) = NDArray (ierr);
00360 
00361                               retval(0) = result;
00362                             }
00363                           else
00364                             gripe_bessel_arg (fn, "second");
00365                         }
00366                     }
00367                   else
00368                     gripe_bessel_arg (fn, "first");
00369                 }
00370             }
00371         }
00372     }
00373   else
00374     print_usage ();
00375 
00376   return retval;
00377 }
00378 
00379 DEFUN_DLD (besselj, args, nargout,
00380   "-*- texinfo -*-\n\
00381 @deftypefn  {Loadable Function} {[@var{j}, @var{ierr}] =} besselj (@var{alpha}, @var{x}, @var{opt})\n\
00382 @deftypefnx {Loadable Function} {[@var{y}, @var{ierr}] =} bessely (@var{alpha}, @var{x}, @var{opt})\n\
00383 @deftypefnx {Loadable Function} {[@var{i}, @var{ierr}] =} besseli (@var{alpha}, @var{x}, @var{opt})\n\
00384 @deftypefnx {Loadable Function} {[@var{k}, @var{ierr}] =} besselk (@var{alpha}, @var{x}, @var{opt})\n\
00385 @deftypefnx {Loadable Function} {[@var{h}, @var{ierr}] =} besselh (@var{alpha}, @var{k}, @var{x}, @var{opt})\n\
00386 Compute Bessel or Hankel functions of various kinds:\n\
00387 \n\
00388 @table @code\n\
00389 @item besselj\n\
00390 Bessel functions of the first kind.  If the argument @var{opt} is supplied,\n\
00391 the result is multiplied by @code{exp(-abs(imag(@var{x})))}.\n\
00392 \n\
00393 @item bessely\n\
00394 Bessel functions of the second kind.  If the argument @var{opt} is supplied,\n\
00395 the result is multiplied by @code{exp(-abs(imag(@var{x})))}.\n\
00396 \n\
00397 @item besseli\n\
00398 \n\
00399 Modified Bessel functions of the first kind.  If the argument @var{opt} is\n\
00400 supplied, the result is multiplied by @code{exp(-abs(real(@var{x})))}.\n\
00401 \n\
00402 @item besselk\n\
00403 \n\
00404 Modified Bessel functions of the second kind.  If the argument @var{opt} is\n\
00405 supplied, the result is multiplied by @code{exp(@var{x})}.\n\
00406 \n\
00407 @item besselh\n\
00408 Compute Hankel functions of the first (@var{k} = 1) or second (@var{k}\n\
00409 = 2) kind.  If the argument @var{opt} is supplied, the result is multiplied\n\
00410 by @code{exp (-I*@var{x})} for @var{k} = 1 or @code{exp (I*@var{x})} for\n\
00411 @var{k} = 2.\n\
00412 @end table\n\
00413 \n\
00414 If @var{alpha} is a scalar, the result is the same size as @var{x}.\n\
00415 If @var{x} is a scalar, the result is the same size as @var{alpha}.\n\
00416 If @var{alpha} is a row vector and @var{x} is a column vector, the\n\
00417 result is a matrix with @code{length (@var{x})} rows and\n\
00418 @code{length (@var{alpha})} columns.  Otherwise, @var{alpha} and\n\
00419 @var{x} must conform and the result will be the same size.\n\
00420 \n\
00421 The value of @var{alpha} must be real.  The value of @var{x} may be\n\
00422 complex.\n\
00423 \n\
00424 If requested, @var{ierr} contains the following status information\n\
00425 and is the same size as the result.\n\
00426 \n\
00427 @enumerate 0\n\
00428 @item\n\
00429 Normal return.\n\
00430 \n\
00431 @item\n\
00432 Input error, return @code{NaN}.\n\
00433 \n\
00434 @item\n\
00435 Overflow, return @code{Inf}.\n\
00436 \n\
00437 @item\n\
00438 Loss of significance by argument reduction results in less than\n\
00439 half of machine accuracy.\n\
00440 \n\
00441 @item\n\
00442 Complete loss of significance by argument reduction, return @code{NaN}.\n\
00443 \n\
00444 @item\n\
00445 Error---no computation, algorithm termination condition not met,\n\
00446 return @code{NaN}.\n\
00447 @end enumerate\n\
00448 @end deftypefn")
00449 {
00450   return do_bessel (BESSEL_J, "besselj", args, nargout);
00451 }
00452 
00453 DEFUN_DLD (bessely, args, nargout,
00454   "-*- texinfo -*-\n\
00455 @deftypefn {Loadable Function} {[@var{y}, @var{ierr}] =} bessely (@var{alpha}, @var{x}, @var{opt})\n\
00456 See besselj.\n\
00457 @end deftypefn")
00458 {
00459   return do_bessel (BESSEL_Y, "bessely", args, nargout);
00460 }
00461 
00462 DEFUN_DLD (besseli, args, nargout,
00463   "-*- texinfo -*-\n\
00464 @deftypefn {Loadable Function} {[@var{i}, @var{ierr}] =} besseli (@var{alpha}, @var{x}, @var{opt})\n\
00465 See besselj.\n\
00466 @end deftypefn")
00467 {
00468   return do_bessel (BESSEL_I, "besseli", args, nargout);
00469 }
00470 
00471 DEFUN_DLD (besselk, args, nargout,
00472   "-*- texinfo -*-\n\
00473 @deftypefn {Loadable Function} {[@var{k}, @var{ierr}] =} besselk (@var{alpha}, @var{x}, @var{opt})\n\
00474 See besselj.\n\
00475 @end deftypefn")
00476 {
00477   return do_bessel (BESSEL_K, "besselk", args, nargout);
00478 }
00479 
00480 DEFUN_DLD (besselh, args, nargout,
00481   "-*- texinfo -*-\n\
00482 @deftypefn {Loadable Function} {[@var{h}, @var{ierr}] =} besselh (@var{alpha}, @var{k}, @var{x}, @var{opt})\n\
00483 See besselj.\n\
00484 @end deftypefn")
00485 {
00486   octave_value_list retval;
00487 
00488   int nargin = args.length ();
00489 
00490   if (nargin == 2)
00491     {
00492       retval = do_bessel (BESSEL_H1, "besselh", args, nargout);
00493     }
00494   else if (nargin == 3 || nargin == 4)
00495     {
00496       octave_idx_type kind = args(1).int_value ();
00497 
00498       if (! error_state)
00499         {
00500           octave_value_list tmp_args;
00501 
00502           if (nargin == 4)
00503             tmp_args(2) = args(3);
00504 
00505           tmp_args(1) = args(2);
00506           tmp_args(0) = args(0);
00507 
00508           if (kind == 1)
00509             retval = do_bessel (BESSEL_H1, "besselh", tmp_args, nargout);
00510           else if (kind == 2)
00511             retval = do_bessel (BESSEL_H2, "besselh", tmp_args, nargout);
00512           else
00513             error ("besselh: expecting K = 1 or 2");
00514         }
00515       else
00516         error ("besselh: invalid value of K");
00517     }
00518   else
00519     print_usage ();
00520 
00521   return retval;
00522 }
00523 
00524 DEFUN_DLD (airy, args, nargout,
00525   "-*- texinfo -*-\n\
00526 @deftypefn {Loadable Function} {[@var{a}, @var{ierr}] =} airy (@var{k}, @var{z}, @var{opt})\n\
00527 Compute Airy functions of the first and second kind, and their\n\
00528 derivatives.\n\
00529 \n\
00530 @example\n\
00531 @group\n\
00532  K   Function   Scale factor (if 'opt' is supplied)\n\
00533 ---  --------   ---------------------------------------\n\
00534  0   Ai (Z)     exp ((2/3) * Z * sqrt (Z))\n\
00535  1   dAi(Z)/dZ  exp ((2/3) * Z * sqrt (Z))\n\
00536  2   Bi (Z)     exp (-abs (real ((2/3) * Z *sqrt (Z))))\n\
00537  3   dBi(Z)/dZ  exp (-abs (real ((2/3) * Z *sqrt (Z))))\n\
00538 @end group\n\
00539 @end example\n\
00540 \n\
00541 The function call @code{airy (@var{z})} is equivalent to\n\
00542 @code{airy (0, @var{z})}.\n\
00543 \n\
00544 The result is the same size as @var{z}.\n\
00545 \n\
00546 If requested, @var{ierr} contains the following status information and\n\
00547 is the same size as the result.\n\
00548 \n\
00549 @enumerate 0\n\
00550 @item\n\
00551 Normal return.\n\
00552 \n\
00553 @item\n\
00554 Input error, return @code{NaN}.\n\
00555 \n\
00556 @item\n\
00557 Overflow, return @code{Inf}.\n\
00558 \n\
00559 @item\n\
00560 Loss of significance by argument reduction results in less than half\n\
00561  of machine accuracy.\n\
00562 \n\
00563 @item\n\
00564 Complete loss of significance by argument reduction, return @code{NaN}.\n\
00565 \n\
00566 @item\n\
00567 Error---no computation, algorithm termination condition not met,\n\
00568 return @code{NaN}.\n\
00569 @end enumerate\n\
00570 @end deftypefn")
00571 {
00572   octave_value_list retval;
00573 
00574   int nargin = args.length ();
00575 
00576   if (nargin > 0 && nargin < 4)
00577     {
00578       bool scale = (nargin == 3);
00579 
00580       int kind = 0;
00581 
00582       if (nargin > 1)
00583         {
00584           kind = args(0).int_value ();
00585 
00586           if (! error_state)
00587             {
00588               if (kind < 0 || kind > 3)
00589                 error ("airy: expecting K = 0, 1, 2, or 3");
00590             }
00591           else
00592             error ("airy: K must be an integer value");
00593         }
00594 
00595       if (! error_state)
00596         {
00597           int idx = nargin == 1 ? 0 : 1;
00598 
00599           if (args (idx).is_single_type ())
00600             {
00601               FloatComplexNDArray z = args(idx).float_complex_array_value ();
00602 
00603               if (! error_state)
00604                 {
00605                   Array<octave_idx_type> ierr;
00606                   octave_value result;
00607 
00608                   if (kind > 1)
00609                     result = biry (z, kind == 3, scale, ierr);
00610                   else
00611                     result = airy (z, kind == 1, scale, ierr);
00612 
00613                   if (nargout > 1)
00614                     retval(1) = NDArray (ierr);
00615 
00616                   retval(0) = result;
00617                 }
00618               else
00619                 error ("airy: Z must be a complex matrix");
00620             }
00621           else
00622             {
00623               ComplexNDArray z = args(idx).complex_array_value ();
00624 
00625               if (! error_state)
00626                 {
00627                   Array<octave_idx_type> ierr;
00628                   octave_value result;
00629 
00630                   if (kind > 1)
00631                     result = biry (z, kind == 3, scale, ierr);
00632                   else
00633                     result = airy (z, kind == 1, scale, ierr);
00634 
00635                   if (nargout > 1)
00636                     retval(1) = NDArray (ierr);
00637 
00638                   retval(0) = result;
00639                 }
00640               else
00641                 error ("airy: Z must be a complex matrix");
00642             }
00643         }
00644     }
00645   else
00646     print_usage ();
00647 
00648   return retval;
00649 }
00650 
00651 /*
00652 %! # Test values computed with GP/PARI version 2.3.3
00653 %!
00654 %!shared alpha, x, jx, yx, ix, kx, nix
00655 %!
00656 %! # Bessel functions, even order, positive and negative x
00657 %! alpha = 2; x = 1.25;
00658 %! jx = 0.1710911312405234823613091417;
00659 %! yx = -1.193199310178553861283790424;
00660 %! ix = 0.2220184483766341752692212604;
00661 %! kx = 0.9410016167388185767085460540;
00662 %!
00663 %!assert(besselj(alpha,x), jx, 100*eps)
00664 %!assert(bessely(alpha,x), yx, 100*eps)
00665 %!assert(besseli(alpha,x), ix, 100*eps)
00666 %!assert(besselk(alpha,x), kx, 100*eps)
00667 %!assert(besselh(alpha,1,x), jx + I*yx, 100*eps)
00668 %!assert(besselh(alpha,2,x), jx - I*yx, 100*eps)
00669 %!
00670 %!assert(besselj(alpha,x,1), jx*exp(-abs(imag(x))), 100*eps)
00671 %!assert(bessely(alpha,x,1), yx*exp(-abs(imag(x))), 100*eps)
00672 %!assert(besseli(alpha,x,1), ix*exp(-abs(real(x))), 100*eps)
00673 %!assert(besselk(alpha,x,1), kx*exp(x), 100*eps)
00674 %!assert(besselh(alpha,1,x,1), (jx + I*yx)*exp(-I*x), 100*eps)
00675 %!assert(besselh(alpha,2,x,1), (jx - I*yx)*exp(I*x), 100*eps)
00676 %!
00677 %!assert(besselj(-alpha,x), jx, 100*eps)
00678 %!assert(bessely(-alpha,x), yx, 100*eps)
00679 %!assert(besseli(-alpha,x), ix, 100*eps)
00680 %!assert(besselk(-alpha,x), kx, 100*eps)
00681 %!assert(besselh(-alpha,1,x), jx + I*yx, 100*eps)
00682 %!assert(besselh(-alpha,2,x), jx - I*yx, 100*eps)
00683 %!
00684 %!assert(besselj(-alpha,x,1), jx*exp(-abs(imag(x))), 100*eps)
00685 %!assert(bessely(-alpha,x,1), yx*exp(-abs(imag(x))), 100*eps)
00686 %!assert(besseli(-alpha,x,1), ix*exp(-abs(real(x))), 100*eps)
00687 %!assert(besselk(-alpha,x,1), kx*exp(x), 100*eps)
00688 %!assert(besselh(-alpha,1,x,1), (jx + I*yx)*exp(-I*x), 100*eps)
00689 %!assert(besselh(-alpha,2,x,1), (jx - I*yx)*exp(I*x), 100*eps)
00690 %!
00691 %! x *= -1;
00692 %! yx = -1.193199310178553861283790424 + 0.3421822624810469647226182835*I;
00693 %! kx = 0.9410016167388185767085460540 - 0.6974915263814386815610060884*I;
00694 %!
00695 %!assert(besselj(alpha,x), jx, 100*eps)
00696 %!assert(bessely(alpha,x), yx, 100*eps)
00697 %!assert(besseli(alpha,x), ix, 100*eps)
00698 %!assert(besselk(alpha,x), kx, 100*eps)
00699 %!assert(besselh(alpha,1,x), jx + I*yx, 100*eps)
00700 %!assert(besselh(alpha,2,x), jx - I*yx, 100*eps)
00701 %!
00702 %!assert(besselj(alpha,x,1), jx*exp(-abs(imag(x))), 100*eps)
00703 %!assert(bessely(alpha,x,1), yx*exp(-abs(imag(x))), 100*eps)
00704 %!assert(besseli(alpha,x,1), ix*exp(-abs(real(x))), 100*eps)
00705 %!assert(besselk(alpha,x,1), kx*exp(x), 100*eps)
00706 %!assert(besselh(alpha,1,x,1), (jx + I*yx)*exp(-I*x), 100*eps)
00707 %!assert(besselh(alpha,2,x,1), (jx - I*yx)*exp(I*x), 100*eps)
00708 %!
00709 %! # Bessel functions, odd order, positive and negative x
00710 %! alpha = 3; x = 2.5;
00711 %! jx = 0.2166003910391135247666890035;
00712 %! yx = -0.7560554967536709968379029772;
00713 %! ix = 0.4743704087780355895548240179;
00714 %! kx = 0.2682271463934492027663765197;
00715 %!
00716 %!assert(besselj(alpha,x), jx, 100*eps)
00717 %!assert(bessely(alpha,x), yx, 100*eps)
00718 %!assert(besseli(alpha,x), ix, 100*eps)
00719 %!assert(besselk(alpha,x), kx, 100*eps)
00720 %!assert(besselh(alpha,1,x), jx + I*yx, 100*eps)
00721 %!assert(besselh(alpha,2,x), jx - I*yx, 100*eps)
00722 %!
00723 %!assert(besselj(alpha,x,1), jx*exp(-abs(imag(x))), 100*eps)
00724 %!assert(bessely(alpha,x,1), yx*exp(-abs(imag(x))), 100*eps)
00725 %!assert(besseli(alpha,x,1), ix*exp(-abs(real(x))), 100*eps)
00726 %!assert(besselk(alpha,x,1), kx*exp(x), 100*eps)
00727 %!assert(besselh(alpha,1,x,1), (jx + I*yx)*exp(-I*x), 100*eps)
00728 %!assert(besselh(alpha,2,x,1), (jx - I*yx)*exp(I*x), 100*eps)
00729 %!
00730 %!assert(besselj(-alpha,x), -jx, 100*eps)
00731 %!assert(bessely(-alpha,x), -yx, 100*eps)
00732 %!assert(besseli(-alpha,x), ix, 100*eps)
00733 %!assert(besselk(-alpha,x), kx, 100*eps)
00734 %!assert(besselh(-alpha,1,x), -(jx + I*yx), 100*eps)
00735 %!assert(besselh(-alpha,2,x), -(jx - I*yx), 100*eps)
00736 %!
00737 %!assert(besselj(-alpha,x,1), -jx*exp(-abs(imag(x))), 100*eps)
00738 %!assert(bessely(-alpha,x,1), -yx*exp(-abs(imag(x))), 100*eps)
00739 %!assert(besseli(-alpha,x,1), ix*exp(-abs(real(x))), 100*eps)
00740 %!assert(besselk(-alpha,x,1), kx*exp(x), 100*eps)
00741 %!assert(besselh(-alpha,1,x,1), -(jx + I*yx)*exp(-I*x), 100*eps)
00742 %!assert(besselh(-alpha,2,x,1), -(jx - I*yx)*exp(I*x), 100*eps)
00743 %!
00744 %! x *= -1;
00745 %! jx = -jx;
00746 %! yx = 0.7560554967536709968379029772 - 0.4332007820782270495333780070*I;
00747 %! ix = -ix;
00748 %! kx = -0.2682271463934492027663765197 - 1.490278591297463775542004240*I;
00749 %!
00750 %!assert(besselj(alpha,x), jx, 100*eps)
00751 %!assert(bessely(alpha,x), yx, 100*eps)
00752 %!assert(besseli(alpha,x), ix, 100*eps)
00753 %!assert(besselk(alpha,x), kx, 100*eps)
00754 %!assert(besselh(alpha,1,x), jx + I*yx, 100*eps)
00755 %!assert(besselh(alpha,2,x), jx - I*yx, 100*eps)
00756 %!
00757 %!assert(besselj(alpha,x,1), jx*exp(-abs(imag(x))), 100*eps)
00758 %!assert(bessely(alpha,x,1), yx*exp(-abs(imag(x))), 100*eps)
00759 %!assert(besseli(alpha,x,1), ix*exp(-abs(real(x))), 100*eps)
00760 %!assert(besselk(alpha,x,1), kx*exp(x), 100*eps)
00761 %!assert(besselh(alpha,1,x,1), (jx + I*yx)*exp(-I*x), 100*eps)
00762 %!assert(besselh(alpha,2,x,1), (jx - I*yx)*exp(I*x), 100*eps)
00763 %!
00764 %! # Bessel functions, fractional order, positive and negative x
00765 %!
00766 %! alpha = 3.5; x = 2.75;
00767 %! jx = 0.1691636439842384154644784389;
00768 %! yx = -0.8301381935499356070267953387;
00769 %! ix = 0.3930540878794826310979363668;
00770 %! kx = 0.2844099013460621170288192503;
00771 %!
00772 %!assert(besselj(alpha,x), jx, 100*eps)
00773 %!assert(bessely(alpha,x), yx, 100*eps)
00774 %!assert(besseli(alpha,x), ix, 100*eps)
00775 %!assert(besselk(alpha,x), kx, 100*eps)
00776 %!assert(besselh(alpha,1,x), jx + I*yx, 100*eps)
00777 %!assert(besselh(alpha,2,x), jx - I*yx, 100*eps)
00778 %!
00779 %!assert(besselj(alpha,x,1), jx*exp(-abs(imag(x))), 100*eps)
00780 %!assert(bessely(alpha,x,1), yx*exp(-abs(imag(x))), 100*eps)
00781 %!assert(besseli(alpha,x,1), ix*exp(-abs(real(x))), 100*eps)
00782 %!assert(besselk(alpha,x,1), kx*exp(x), 100*eps)
00783 %!assert(besselh(alpha,1,x,1), (jx + I*yx)*exp(-I*x), 100*eps)
00784 %!assert(besselh(alpha,2,x,1), (jx - I*yx)*exp(I*x), 100*eps)
00785 %!
00786 %! nix = 0.2119931212254662995364461998;
00787 %!
00788 %!assert(besselj(-alpha,x), yx, 100*eps)
00789 %!assert(bessely(-alpha,x), -jx, 100*eps)
00790 %!assert(besseli(-alpha,x), nix, 100*eps)
00791 %!assert(besselk(-alpha,x), kx, 100*eps)
00792 %!assert(besselh(-alpha,1,x), -I*(jx + I*yx), 100*eps)
00793 %!assert(besselh(-alpha,2,x), I*(jx - I*yx), 100*eps)
00794 %!
00795 %!assert(besselj(-alpha,x,1), yx*exp(-abs(imag(x))), 100*eps)
00796 %!assert(bessely(-alpha,x,1), -jx*exp(-abs(imag(x))), 100*eps)
00797 %!assert(besseli(-alpha,x,1), nix*exp(-abs(real(x))), 100*eps)
00798 %!assert(besselk(-alpha,x,1), kx*exp(x), 100*eps)
00799 %!assert(besselh(-alpha,1,x,1), -I*(jx + I*yx)*exp(-I*x), 100*eps)
00800 %!assert(besselh(-alpha,2,x,1), I*(jx - I*yx)*exp(I*x), 100*eps)
00801 %!
00802 %! x *= -1;
00803 %! jx *= -I;
00804 %! yx = -0.8301381935499356070267953387*I;
00805 %! ix *= -I;
00806 %! kx = -0.9504059335995575096509874508*I;
00807 %!
00808 %!assert(besselj(alpha,x), jx, 100*eps)
00809 %!assert(bessely(alpha,x), yx, 100*eps)
00810 %!assert(besseli(alpha,x), ix, 100*eps)
00811 %!assert(besselk(alpha,x), kx, 100*eps)
00812 %!assert(besselh(alpha,1,x), jx + I*yx, 100*eps)
00813 %!assert(besselh(alpha,2,x), jx - I*yx, 100*eps)
00814 %!
00815 %!assert(besselj(alpha,x,1), jx*exp(-abs(imag(x))), 100*eps)
00816 %!assert(bessely(alpha,x,1), yx*exp(-abs(imag(x))), 100*eps)
00817 %!assert(besseli(alpha,x,1), ix*exp(-abs(real(x))), 100*eps)
00818 %!assert(besselk(alpha,x,1), kx*exp(x), 100*eps)
00819 %!assert(besselh(alpha,1,x,1), (jx + I*yx)*exp(-I*x), 100*eps)
00820 %!assert(besselh(alpha,2,x,1), (jx - I*yx)*exp(I*x), 100*eps)
00821 %!
00822 %! # Bessel functions, even order, complex x
00823 %!
00824 %! alpha = 2; x = 1.25 + 3.625 * I;
00825 %! jx = -1.299533366810794494030065917 + 4.370833116012278943267479589*I;
00826 %! yx = -4.370357232383223896393056727 - 1.283083391453582032688834041*I;
00827 %! ix = -0.6717801680341515541002273932 - 0.2314623443930774099910228553*I;
00828 %! kx = -0.01108009888623253515463783379 + 0.2245218229358191588208084197*I;
00829 %!
00830 %!assert(besselj(alpha,x), jx, 100*eps)
00831 %!assert(bessely(alpha,x), yx, 100*eps)
00832 %!assert(besseli(alpha,x), ix, 100*eps)
00833 %!assert(besselk(alpha,x), kx, 100*eps)
00834 %!assert(besselh(alpha,1,x), jx + I*yx, 100*eps)
00835 %!assert(besselh(alpha,2,x), jx - I*yx, 100*eps)
00836 %!
00837 %!assert(besselj(alpha,x,1), jx*exp(-abs(imag(x))), 100*eps)
00838 %!assert(bessely(alpha,x,1), yx*exp(-abs(imag(x))), 100*eps)
00839 %!assert(besseli(alpha,x,1), ix*exp(-abs(real(x))), 100*eps)
00840 %!assert(besselk(alpha,x,1), kx*exp(x), 100*eps)
00841 %!assert(besselh(alpha,1,x,1), (jx + I*yx)*exp(-I*x), 100*eps)
00842 %!assert(besselh(alpha,2,x,1), (jx - I*yx)*exp(I*x), 100*eps)
00843 %!
00844 %!assert(besselj(-alpha,x), jx, 100*eps)
00845 %!assert(bessely(-alpha,x), yx, 100*eps)
00846 %!assert(besseli(-alpha,x), ix, 100*eps)
00847 %!assert(besselk(-alpha,x), kx, 100*eps)
00848 %!assert(besselh(-alpha,1,x), jx + I*yx, 100*eps)
00849 %!assert(besselh(-alpha,2,x), jx - I*yx, 100*eps)
00850 %!
00851 %!assert(besselj(-alpha,x,1), jx*exp(-abs(imag(x))), 100*eps)
00852 %!assert(bessely(-alpha,x,1), yx*exp(-abs(imag(x))), 100*eps)
00853 %!assert(besseli(-alpha,x,1), ix*exp(-abs(real(x))), 100*eps)
00854 %!assert(besselk(-alpha,x,1), kx*exp(x), 100*eps)
00855 %!assert(besselh(-alpha,1,x,1), (jx + I*yx)*exp(-I*x), 100*eps)
00856 %!assert(besselh(-alpha,2,x,1), (jx - I*yx)*exp(I*x), 100*eps)
00857 %!
00858 %! # Bessel functions, odd order, complex x
00859 %!
00860 %! alpha = 3; x = 2.5 + 1.875 * I;
00861 %! jx = 0.1330721523048277493333458596 + 0.5386295217249660078754395597*I;
00862 %! yx = -0.6485072392105829901122401551 + 0.2608129289785456797046996987*I;
00863 %! ix = -0.6182064685486998097516365709 + 0.4677561094683470065767989920*I;
00864 %! kx = -0.1568585587733540007867882337 - 0.05185853709490846050505141321*I;
00865 %!
00866 %!assert(besselj(alpha,x), jx, 100*eps)
00867 %!assert(bessely(alpha,x), yx, 100*eps)
00868 %!assert(besseli(alpha,x), ix, 100*eps)
00869 %!assert(besselk(alpha,x), kx, 100*eps)
00870 %!assert(besselh(alpha,1,x), jx + I*yx, 100*eps)
00871 %!assert(besselh(alpha,2,x), jx - I*yx, 100*eps)
00872 %!
00873 %!assert(besselj(alpha,x,1), jx*exp(-abs(imag(x))), 100*eps)
00874 %!assert(bessely(alpha,x,1), yx*exp(-abs(imag(x))), 100*eps)
00875 %!assert(besseli(alpha,x,1), ix*exp(-abs(real(x))), 100*eps)
00876 %!assert(besselk(alpha,x,1), kx*exp(x), 100*eps)
00877 %!assert(besselh(alpha,1,x,1), (jx + I*yx)*exp(-I*x), 100*eps)
00878 %!assert(besselh(alpha,2,x,1), (jx - I*yx)*exp(I*x), 100*eps)
00879 %!
00880 %!assert(besselj(-alpha,x), -jx, 100*eps)
00881 %!assert(bessely(-alpha,x), -yx, 100*eps)
00882 %!assert(besseli(-alpha,x), ix, 100*eps)
00883 %!assert(besselk(-alpha,x), kx, 100*eps)
00884 %!assert(besselh(-alpha,1,x), -(jx + I*yx), 100*eps)
00885 %!assert(besselh(-alpha,2,x), -(jx - I*yx), 100*eps)
00886 %!
00887 %!assert(besselj(-alpha,x,1), -jx*exp(-abs(imag(x))), 100*eps)
00888 %!assert(bessely(-alpha,x,1), -yx*exp(-abs(imag(x))), 100*eps)
00889 %!assert(besseli(-alpha,x,1), ix*exp(-abs(real(x))), 100*eps)
00890 %!assert(besselk(-alpha,x,1), kx*exp(x), 100*eps)
00891 %!assert(besselh(-alpha,1,x,1), -(jx + I*yx)*exp(-I*x), 100*eps)
00892 %!assert(besselh(-alpha,2,x,1), -(jx - I*yx)*exp(I*x), 100*eps)
00893 %!
00894 %! # Bessel functions, fractional order, complex x
00895 %!
00896 %! alpha = 3.5; x = 1.75 + 4.125 * I;
00897 %! jx = -3.018566131370455929707009100 - 0.7585648436793900607704057611*I;
00898 %! yx = 0.7772278839106298215614791107 - 3.018518722313849782683792010*I;
00899 %! ix = 0.2100873577220057189038160913 - 0.6551765604618246531254970926*I;
00900 %! kx = 0.1757147290513239935341488069 + 0.08772348296883849205562558311*I;
00901 %!
00902 %!assert(besselj(alpha,x), jx, 100*eps)
00903 %!assert(bessely(alpha,x), yx, 100*eps)
00904 %!assert(besseli(alpha,x), ix, 100*eps)
00905 %!assert(besselk(alpha,x), kx, 100*eps)
00906 %!assert(besselh(alpha,1,x), jx + I*yx, 100*eps)
00907 %!assert(besselh(alpha,2,x), jx - I*yx, 100*eps)
00908 %!
00909 %!assert(besselj(alpha,x,1), jx*exp(-abs(imag(x))), 100*eps)
00910 %!assert(bessely(alpha,x,1), yx*exp(-abs(imag(x))), 100*eps)
00911 %!assert(besseli(alpha,x,1), ix*exp(-abs(real(x))), 100*eps)
00912 %!assert(besselk(alpha,x,1), kx*exp(x), 100*eps)
00913 %!assert(besselh(alpha,1,x,1), (jx + I*yx)*exp(-I*x), 100*eps)
00914 %!assert(besselh(alpha,2,x,1), (jx - I*yx)*exp(I*x), 100*eps)
00915 %!
00916 %!  nix = 0.09822388691172060573913739253 - 0.7110230642207380127317227407*I;
00917 %!
00918 %!assert(besselj(-alpha,x), yx, 100*eps)
00919 %!assert(bessely(-alpha,x), -jx, 100*eps)
00920 %!assert(besseli(-alpha,x), nix, 100*eps)
00921 %!assert(besselk(-alpha,x), kx, 100*eps)
00922 %!assert(besselh(-alpha,1,x), -I*(jx + I*yx), 100*eps)
00923 %!assert(besselh(-alpha,2,x), I*(jx - I*yx), 100*eps)
00924 %!
00925 %!assert(besselj(-alpha,x,1), yx*exp(-abs(imag(x))), 100*eps)
00926 %!assert(bessely(-alpha,x,1), -jx*exp(-abs(imag(x))), 100*eps)
00927 %!assert(besseli(-alpha,x,1), nix*exp(-abs(real(x))), 100*eps)
00928 %!assert(besselk(-alpha,x,1), kx*exp(x), 100*eps)
00929 %!assert(besselh(-alpha,1,x,1), -I*(jx + I*yx)*exp(-I*x), 100*eps)
00930 %!assert(besselh(-alpha,2,x,1), I*(jx - I*yx)*exp(I*x), 100*eps)
00931 */
 All Classes Files Functions Variables Typedefs Enumerations Enumerator Friends Defines