lo-specfun.cc

Go to the documentation of this file.
00001 /*
00002 
00003 Copyright (C) 1996-2012 John W. Eaton
00004 Copyright (C) 2010 Jaroslav Hajek
00005 Copyright (C) 2010 VZLU Prague
00006 
00007 This file is part of Octave.
00008 
00009 Octave is free software; you can redistribute it and/or modify it
00010 under the terms of the GNU General Public License as published by the
00011 Free Software Foundation; either version 3 of the License, or (at your
00012 option) any later version.
00013 
00014 Octave is distributed in the hope that it will be useful, but WITHOUT
00015 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00016 FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
00017 for more details.
00018 
00019 You should have received a copy of the GNU General Public License
00020 along with Octave; see the file COPYING.  If not, see
00021 <http://www.gnu.org/licenses/>.
00022 
00023 */
00024 
00025 #ifdef HAVE_CONFIG_H
00026 #include <config.h>
00027 #endif
00028 
00029 #include "Range.h"
00030 #include "CColVector.h"
00031 #include "CMatrix.h"
00032 #include "dRowVector.h"
00033 #include "dMatrix.h"
00034 #include "dNDArray.h"
00035 #include "CNDArray.h"
00036 #include "fCColVector.h"
00037 #include "fCMatrix.h"
00038 #include "fRowVector.h"
00039 #include "fMatrix.h"
00040 #include "fNDArray.h"
00041 #include "fCNDArray.h"
00042 #include "f77-fcn.h"
00043 #include "lo-error.h"
00044 #include "lo-ieee.h"
00045 #include "lo-specfun.h"
00046 #include "mx-inlines.cc"
00047 #include "lo-mappers.h"
00048 
00049 #ifndef M_PI
00050 #define M_PI 3.14159265358979323846
00051 #endif
00052 
00053 extern "C"
00054 {
00055   F77_RET_T
00056   F77_FUNC (zbesj, ZBESJ) (const double&, const double&, const double&,
00057                            const octave_idx_type&, const octave_idx_type&,
00058                            double*, double*, octave_idx_type&,
00059                            octave_idx_type&);
00060 
00061   F77_RET_T
00062   F77_FUNC (zbesy, ZBESY) (const double&, const double&, const double&,
00063                            const octave_idx_type&, const octave_idx_type&,
00064                            double*, double*, octave_idx_type&, double*,
00065                            double*, octave_idx_type&);
00066 
00067   F77_RET_T
00068   F77_FUNC (zbesi, ZBESI) (const double&, const double&, const double&,
00069                            const octave_idx_type&, const octave_idx_type&,
00070                            double*, double*, octave_idx_type&,
00071                            octave_idx_type&);
00072 
00073   F77_RET_T
00074   F77_FUNC (zbesk, ZBESK) (const double&, const double&, const double&,
00075                            const octave_idx_type&, const octave_idx_type&,
00076                            double*, double*, octave_idx_type&,
00077                            octave_idx_type&);
00078 
00079   F77_RET_T
00080   F77_FUNC (zbesh, ZBESH) (const double&, const double&, const double&,
00081                            const octave_idx_type&, const octave_idx_type&,
00082                            const octave_idx_type&, double*, double*,
00083                            octave_idx_type&, octave_idx_type&);
00084 
00085   F77_RET_T
00086   F77_FUNC (cbesj, cBESJ) (const FloatComplex&, const float&,
00087                            const octave_idx_type&, const octave_idx_type&,
00088                            FloatComplex*, octave_idx_type&, octave_idx_type&);
00089 
00090   F77_RET_T
00091   F77_FUNC (cbesy, CBESY) (const FloatComplex&, const float&,
00092                            const octave_idx_type&, const octave_idx_type&,
00093                            FloatComplex*, octave_idx_type&,
00094                            FloatComplex*, octave_idx_type&);
00095 
00096   F77_RET_T
00097   F77_FUNC (cbesi, CBESI) (const FloatComplex&, const float&,
00098                            const octave_idx_type&, const octave_idx_type&,
00099                            FloatComplex*, octave_idx_type&, octave_idx_type&);
00100 
00101   F77_RET_T
00102   F77_FUNC (cbesk, CBESK) (const FloatComplex&, const float&,
00103                            const octave_idx_type&, const octave_idx_type&,
00104                            FloatComplex*, octave_idx_type&, octave_idx_type&);
00105 
00106   F77_RET_T
00107   F77_FUNC (cbesh, CBESH) (const FloatComplex&, const float&,
00108                            const octave_idx_type&, const octave_idx_type&,
00109                            const octave_idx_type&, FloatComplex*,
00110                            octave_idx_type&, octave_idx_type&);
00111 
00112   F77_RET_T
00113   F77_FUNC (zairy, ZAIRY) (const double&, const double&,
00114                            const octave_idx_type&, const octave_idx_type&,
00115                            double&, double&, octave_idx_type&,
00116                            octave_idx_type&);
00117 
00118   F77_RET_T
00119   F77_FUNC (cairy, CAIRY) (const float&, const float&, const octave_idx_type&,
00120                            const octave_idx_type&, float&, float&,
00121                            octave_idx_type&, octave_idx_type&);
00122 
00123   F77_RET_T
00124   F77_FUNC (zbiry, ZBIRY) (const double&, const double&,
00125                            const octave_idx_type&, const octave_idx_type&,
00126                            double&, double&, octave_idx_type&);
00127 
00128   F77_RET_T
00129   F77_FUNC (cbiry, CBIRY) (const float&, const float&, const octave_idx_type&,
00130                            const octave_idx_type&, float&, float&,
00131                            octave_idx_type&);
00132 
00133   F77_RET_T
00134   F77_FUNC (xdacosh, XDACOSH) (const double&, double&);
00135 
00136   F77_RET_T
00137   F77_FUNC (xacosh, XACOSH) (const float&, float&);
00138 
00139   F77_RET_T
00140   F77_FUNC (xdasinh, XDASINH) (const double&, double&);
00141 
00142   F77_RET_T
00143   F77_FUNC (xasinh, XASINH) (const float&, float&);
00144 
00145   F77_RET_T
00146   F77_FUNC (xdatanh, XDATANH) (const double&, double&);
00147 
00148   F77_RET_T
00149   F77_FUNC (xatanh, XATANH) (const float&, float&);
00150 
00151   F77_RET_T
00152   F77_FUNC (xderf, XDERF) (const double&, double&);
00153 
00154   F77_RET_T
00155   F77_FUNC (xerf, XERF) (const float&, float&);
00156 
00157   F77_RET_T
00158   F77_FUNC (xderfc, XDERFC) (const double&, double&);
00159 
00160   F77_RET_T
00161   F77_FUNC (xerfc, XERFC) (const float&, float&);
00162 
00163   F77_RET_T
00164   F77_FUNC (xdbetai, XDBETAI) (const double&, const double&,
00165                                const double&, double&);
00166 
00167   F77_RET_T
00168   F77_FUNC (xbetai, XBETAI) (const float&, const float&,
00169                              const float&, float&);
00170 
00171   F77_RET_T
00172   F77_FUNC (xdgamma, XDGAMMA) (const double&, double&);
00173 
00174   F77_RET_T
00175   F77_FUNC (xgamma, XGAMMA) (const float&, float&);
00176 
00177   F77_RET_T
00178   F77_FUNC (xgammainc, XGAMMAINC) (const double&, const double&, double&);
00179 
00180   F77_RET_T
00181   F77_FUNC (xsgammainc, XSGAMMAINC) (const float&, const float&, float&);
00182 
00183   F77_RET_T
00184   F77_FUNC (dlgams, DLGAMS) (const double&, double&, double&);
00185 
00186   F77_RET_T
00187   F77_FUNC (algams, ALGAMS) (const float&, float&, float&);
00188 }
00189 
00190 #if !defined (HAVE_ACOSH)
00191 double
00192 acosh (double x)
00193 {
00194   double retval;
00195   F77_XFCN (xdacosh, XDACOSH, (x, retval));
00196   return retval;
00197 }
00198 #endif
00199 
00200 #if !defined (HAVE_ACOSHF)
00201 float
00202 acoshf (float x)
00203 {
00204   float retval;
00205   F77_XFCN (xacosh, XACOSH, (x, retval));
00206   return retval;
00207 }
00208 #endif
00209 
00210 #if !defined (HAVE_ASINH)
00211 double
00212 asinh (double x)
00213 {
00214   double retval;
00215   F77_XFCN (xdasinh, XDASINH, (x, retval));
00216   return retval;
00217 }
00218 #endif
00219 
00220 #if !defined (HAVE_ASINHF)
00221 float
00222 asinhf (float x)
00223 {
00224   float retval;
00225   F77_XFCN (xasinh, XASINH, (x, retval));
00226   return retval;
00227 }
00228 #endif
00229 
00230 #if !defined (HAVE_ATANH)
00231 double
00232 atanh (double x)
00233 {
00234   double retval;
00235   F77_XFCN (xdatanh, XDATANH, (x, retval));
00236   return retval;
00237 }
00238 #endif
00239 
00240 #if !defined (HAVE_ATANHF)
00241 float
00242 atanhf (float x)
00243 {
00244   float retval;
00245   F77_XFCN (xatanh, XATANH, (x, retval));
00246   return retval;
00247 }
00248 #endif
00249 
00250 #if !defined (HAVE_ERF)
00251 double
00252 erf (double x)
00253 {
00254   double retval;
00255   F77_XFCN (xderf, XDERF, (x, retval));
00256   return retval;
00257 }
00258 #endif
00259 
00260 #if !defined (HAVE_ERFF)
00261 float
00262 erff (float x)
00263 {
00264   float retval;
00265   F77_XFCN (xerf, XERF, (x, retval));
00266   return retval;
00267 }
00268 #endif
00269 
00270 #if !defined (HAVE_ERFC)
00271 double
00272 erfc (double x)
00273 {
00274   double retval;
00275   F77_XFCN (xderfc, XDERFC, (x, retval));
00276   return retval;
00277 }
00278 #endif
00279 
00280 #if !defined (HAVE_ERFCF)
00281 float
00282 erfcf (float x)
00283 {
00284   float retval;
00285   F77_XFCN (xerfc, XERFC, (x, retval));
00286   return retval;
00287 }
00288 #endif
00289 
00290 double
00291 xgamma (double x)
00292 {
00293   double result;
00294 
00295   if (xisnan (x))
00296     result = x;
00297   else if ((x <= 0 && D_NINT (x) == x) || xisinf (x))
00298     result = octave_Inf;
00299   else
00300 #if defined (HAVE_TGAMMA)
00301     result = tgamma (x);
00302 #else
00303     F77_XFCN (xdgamma, XDGAMMA, (x, result));
00304 #endif
00305 
00306   return result;
00307 }
00308 
00309 double
00310 xlgamma (double x)
00311 {
00312 #if defined (HAVE_LGAMMA)
00313   return lgamma (x);
00314 #else
00315   double result;
00316   double sgngam;
00317 
00318   if (xisnan (x))
00319     result = x;
00320   else if ((x <= 0 && D_NINT (x) == x) || xisinf (x))
00321     result = octave_Inf;
00322   else
00323     F77_XFCN (dlgams, DLGAMS, (x, result, sgngam));
00324 
00325   return result;
00326 #endif
00327 }
00328 
00329 Complex
00330 rc_lgamma (double x)
00331 {
00332   double result;
00333 
00334 #if defined (HAVE_LGAMMA_R)
00335   int sgngam;
00336   result = lgamma_r (x, &sgngam);
00337 #else
00338   double sgngam;
00339 
00340   if (xisnan (x))
00341     result = x;
00342   else if ((x <= 0 && D_NINT (x) == x) || xisinf (x))
00343     result = octave_Inf;
00344   else
00345     F77_XFCN (dlgams, DLGAMS, (x, result, sgngam));
00346 
00347 #endif
00348 
00349   if (sgngam < 0)
00350     return result + Complex (0., M_PI);
00351   else
00352     return result;
00353 }
00354 
00355 float
00356 xgamma (float x)
00357 {
00358   float result;
00359 
00360   if (xisnan (x))
00361     result = x;
00362   else if ((x <= 0 && D_NINT (x) == x) || xisinf (x))
00363     result = octave_Float_Inf;
00364   else
00365 #if defined (HAVE_TGAMMAF)
00366     result = tgammaf (x);
00367 #else
00368     F77_XFCN (xgamma, XGAMMA, (x, result));
00369 #endif
00370 
00371   return result;
00372 }
00373 
00374 float
00375 xlgamma (float x)
00376 {
00377 #if defined (HAVE_LGAMMAF)
00378   return lgammaf (x);
00379 #else
00380   float result;
00381   float sgngam;
00382 
00383   if (xisnan (x))
00384     result = x;
00385   else if ((x <= 0 && D_NINT (x) == x) || xisinf (x))
00386     result = octave_Float_Inf;
00387   else
00388     F77_XFCN (algams, ALGAMS, (x, result, sgngam));
00389 
00390   return result;
00391 #endif
00392 }
00393 
00394 FloatComplex
00395 rc_lgamma (float x)
00396 {
00397   float result;
00398 
00399 #if defined (HAVE_LGAMMAF_R)
00400   int sgngam;
00401   result = lgammaf_r (x, &sgngam);
00402 #else
00403   float sgngam;
00404 
00405   if (xisnan (x))
00406     result = x;
00407   else if ((x <= 0 && D_NINT (x) == x) || xisinf (x))
00408     result = octave_Float_Inf;
00409   else
00410     F77_XFCN (algams, ALGAMS, (x, result, sgngam));
00411 
00412 #endif
00413 
00414   if (sgngam < 0)
00415     return result + FloatComplex (0., M_PI);
00416   else
00417     return result;
00418 }
00419 
00420 #if !defined (HAVE_EXPM1)
00421 double
00422 expm1 (double x)
00423 {
00424   double retval;
00425 
00426   double ax = fabs (x);
00427 
00428   if (ax < 0.1)
00429     {
00430       ax /= 16;
00431 
00432       // use Taylor series to calculate exp(x)-1.
00433       double t = ax;
00434       double s = 0;
00435       for (int i = 2; i < 7; i++)
00436         s += (t *= ax/i);
00437       s += ax;
00438 
00439       // use the identity (a+1)^2-1 = a*(a+2)
00440       double e = s;
00441       for (int i = 0; i < 4; i++)
00442         {
00443           s *= e + 2;
00444           e *= e + 2;
00445         }
00446 
00447       retval = (x > 0) ? s : -s / (1+s);
00448     }
00449   else
00450     retval = exp (x) - 1;
00451 
00452   return retval;
00453 }
00454 #endif
00455 
00456 Complex
00457 expm1(const Complex& x)
00458 {
00459   Complex retval;
00460 
00461   if (std:: abs (x) < 1)
00462     {
00463       double im = x.imag();
00464       double u = expm1 (x.real ());
00465       double v = sin (im/2);
00466       v = -2*v*v;
00467       retval = Complex (u*v + u + v, (u+1) * sin (im));
00468     }
00469   else
00470     retval = std::exp (x) - Complex (1);
00471 
00472   return retval;
00473 }
00474 
00475 #if !defined (HAVE_EXPM1F)
00476 float
00477 expm1f (float x)
00478 {
00479   float retval;
00480 
00481   float ax = fabs (x);
00482 
00483   if (ax < 0.1)
00484     {
00485       ax /= 16;
00486 
00487       // use Taylor series to calculate exp(x)-1.
00488       float t = ax;
00489       float s = 0;
00490       for (int i = 2; i < 7; i++)
00491         s += (t *= ax/i);
00492       s += ax;
00493 
00494       // use the identity (a+1)^2-1 = a*(a+2)
00495       float e = s;
00496       for (int i = 0; i < 4; i++)
00497         {
00498           s *= e + 2;
00499           e *= e + 2;
00500         }
00501 
00502       retval = (x > 0) ? s : -s / (1+s);
00503     }
00504   else
00505     retval = exp (x) - 1;
00506 
00507   return retval;
00508 }
00509 #endif
00510 
00511 FloatComplex
00512 expm1(const FloatComplex& x)
00513 {
00514   FloatComplex retval;
00515 
00516   if (std:: abs (x) < 1)
00517     {
00518       float im = x.imag();
00519       float u = expm1 (x.real ());
00520       float v = sin (im/2);
00521       v = -2*v*v;
00522       retval = FloatComplex (u*v + u + v, (u+1) * sin (im));
00523     }
00524   else
00525     retval = std::exp (x) - FloatComplex (1);
00526 
00527   return retval;
00528 }
00529 
00530 #if !defined (HAVE_LOG1P)
00531 double
00532 log1p (double x)
00533 {
00534   double retval;
00535 
00536   double ax = fabs (x);
00537 
00538   if (ax < 0.2)
00539     {
00540       // use approximation log (1+x) ~ 2*sum ((x/(2+x)).^ii ./ ii), ii = 1:2:2n+1
00541       double u = x / (2 + x), t = 1, s = 0;
00542       for (int i = 2; i < 12; i += 2)
00543         s += (t *= u*u) / (i+1);
00544 
00545       retval = 2 * (s + 1) * u;
00546     }
00547   else
00548     retval = log (1 + x);
00549 
00550   return retval;
00551 }
00552 #endif
00553 
00554 Complex
00555 log1p (const Complex& x)
00556 {
00557   Complex retval;
00558 
00559   double r = x.real (), i = x.imag();
00560 
00561   if (fabs (r) < 0.5 && fabs (i) < 0.5)
00562     {
00563       double u = 2*r + r*r + i*i;
00564       retval = Complex (log1p (u / (1+sqrt (u+1))),
00565                         atan2 (1 + r, i));
00566     }
00567   else
00568     retval = std::log (Complex(1) + x);
00569 
00570   return retval;
00571 }
00572 
00573 #if !defined (HAVE_CBRT)
00574 double cbrt (double x)
00575 {
00576   static const double one_third = 0.3333333333333333333;
00577   if (xfinite (x))
00578     {
00579       // Use pow.
00580       double y = std::pow (std::abs (x), one_third) * signum (x);
00581       // Correct for better accuracy.
00582       return (x / (y*y) + y + y) / 3;
00583     }
00584   else
00585     return x;
00586 }
00587 #endif
00588 
00589 #if !defined (HAVE_LOG1PF)
00590 float
00591 log1pf (float x)
00592 {
00593   float retval;
00594 
00595   float ax = fabs (x);
00596 
00597   if (ax < 0.2)
00598     {
00599       // use approximation log (1+x) ~ 2*sum ((x/(2+x)).^ii ./ ii), ii = 1:2:2n+1
00600       float u = x / (2 + x), t = 1, s = 0;
00601       for (int i = 2; i < 12; i += 2)
00602         s += (t *= u*u) / (i+1);
00603 
00604       retval = 2 * (s + 1) * u;
00605     }
00606   else
00607     retval = log (1 + x);
00608 
00609   return retval;
00610 }
00611 #endif
00612 
00613 FloatComplex
00614 log1p (const FloatComplex& x)
00615 {
00616   FloatComplex retval;
00617 
00618   float r = x.real (), i = x.imag();
00619 
00620   if (fabs (r) < 0.5 && fabs (i) < 0.5)
00621     {
00622       float u = 2*r + r*r + i*i;
00623       retval = FloatComplex (log1p (u / (1+sqrt (u+1))),
00624                         atan2 (1 + r, i));
00625     }
00626   else
00627     retval = std::log (FloatComplex(1) + x);
00628 
00629   return retval;
00630 }
00631 
00632 #if !defined (HAVE_CBRTF)
00633 float cbrtf (float x)
00634 {
00635   static const float one_third = 0.3333333333333333333f;
00636   if (xfinite (x))
00637     {
00638       // Use pow.
00639       float y = std::pow (std::abs (x), one_third) * signum (x);
00640       // Correct for better accuracy.
00641       return (x / (y*y) + y + y) / 3;
00642     }
00643   else
00644     return x;
00645 }
00646 #endif
00647 
00648 static inline Complex
00649 zbesj (const Complex& z, double alpha, int kode, octave_idx_type& ierr);
00650 
00651 static inline Complex
00652 zbesy (const Complex& z, double alpha, int kode, octave_idx_type& ierr);
00653 
00654 static inline Complex
00655 zbesi (const Complex& z, double alpha, int kode, octave_idx_type& ierr);
00656 
00657 static inline Complex
00658 zbesk (const Complex& z, double alpha, int kode, octave_idx_type& ierr);
00659 
00660 static inline Complex
00661 zbesh1 (const Complex& z, double alpha, int kode, octave_idx_type& ierr);
00662 
00663 static inline Complex
00664 zbesh2 (const Complex& z, double alpha, int kode, octave_idx_type& ierr);
00665 
00666 static inline Complex
00667 bessel_return_value (const Complex& val, octave_idx_type ierr)
00668 {
00669   static const Complex inf_val = Complex (octave_Inf, octave_Inf);
00670   static const Complex nan_val = Complex (octave_NaN, octave_NaN);
00671 
00672   Complex retval;
00673 
00674   switch (ierr)
00675     {
00676     case 0:
00677     case 3:
00678       retval = val;
00679       break;
00680 
00681     case 2:
00682       retval = inf_val;
00683       break;
00684 
00685     default:
00686       retval = nan_val;
00687       break;
00688     }
00689 
00690   return retval;
00691 }
00692 
00693 static inline bool
00694 is_integer_value (double x)
00695 {
00696   return x == static_cast<double> (static_cast<long> (x));
00697 }
00698 
00699 static inline Complex
00700 zbesj (const Complex& z, double alpha, int kode, octave_idx_type& ierr)
00701 {
00702   Complex retval;
00703 
00704   if (alpha >= 0.0)
00705     {
00706       double yr = 0.0;
00707       double yi = 0.0;
00708 
00709       octave_idx_type nz;
00710 
00711       double zr = z.real ();
00712       double zi = z.imag ();
00713 
00714       F77_FUNC (zbesj, ZBESJ) (zr, zi, alpha, 2, 1, &yr, &yi, nz, ierr);
00715 
00716       if (kode != 2)
00717         {
00718           double expz = exp (std::abs (zi));
00719           yr *= expz;
00720           yi *= expz;
00721         }
00722 
00723       if (zi == 0.0 && zr >= 0.0)
00724         yi = 0.0;
00725 
00726       retval = bessel_return_value (Complex (yr, yi), ierr);
00727     }
00728   else if (is_integer_value (alpha))
00729     {
00730       // zbesy can overflow as z->0, and cause troubles for generic case below
00731       alpha = -alpha;
00732       Complex tmp = zbesj (z, alpha, kode, ierr);
00733       if ((static_cast <long> (alpha)) & 1)
00734         tmp = - tmp;
00735       retval = bessel_return_value (tmp, ierr);
00736     }
00737   else
00738     {
00739       alpha = -alpha;
00740 
00741       Complex tmp = cos (M_PI * alpha) * zbesj (z, alpha, kode, ierr);
00742 
00743       if (ierr == 0 || ierr == 3)
00744         {
00745           tmp -= sin (M_PI * alpha) * zbesy (z, alpha, kode, ierr);
00746 
00747           retval = bessel_return_value (tmp, ierr);
00748         }
00749       else
00750         retval = Complex (octave_NaN, octave_NaN);
00751     }
00752 
00753   return retval;
00754 }
00755 
00756 static inline Complex
00757 zbesy (const Complex& z, double alpha, int kode, octave_idx_type& ierr)
00758 {
00759   Complex retval;
00760 
00761   if (alpha >= 0.0)
00762     {
00763       double yr = 0.0;
00764       double yi = 0.0;
00765 
00766       octave_idx_type nz;
00767 
00768       double wr, wi;
00769 
00770       double zr = z.real ();
00771       double zi = z.imag ();
00772 
00773       ierr = 0;
00774 
00775       if (zr == 0.0 && zi == 0.0)
00776         {
00777           yr = -octave_Inf;
00778           yi = 0.0;
00779         }
00780       else
00781         {
00782           F77_FUNC (zbesy, ZBESY) (zr, zi, alpha, 2, 1, &yr, &yi, nz,
00783                                    &wr, &wi, ierr);
00784 
00785           if (kode != 2)
00786             {
00787               double expz = exp (std::abs (zi));
00788               yr *= expz;
00789               yi *= expz;
00790             }
00791 
00792           if (zi == 0.0 && zr >= 0.0)
00793             yi = 0.0;
00794         }
00795 
00796       return bessel_return_value (Complex (yr, yi), ierr);
00797     }
00798   else if (is_integer_value (alpha - 0.5))
00799     {
00800       // zbesy can overflow as z->0, and cause troubles for generic case below
00801       alpha = -alpha;
00802       Complex tmp = zbesj (z, alpha, kode, ierr);
00803       if ((static_cast <long> (alpha - 0.5)) & 1)
00804         tmp = - tmp;
00805       retval = bessel_return_value (tmp, ierr);
00806     }
00807   else
00808     {
00809       alpha = -alpha;
00810 
00811       Complex tmp = cos (M_PI * alpha) * zbesy (z, alpha, kode, ierr);
00812 
00813       if (ierr == 0 || ierr == 3)
00814         {
00815           tmp += sin (M_PI * alpha) * zbesj (z, alpha, kode, ierr);
00816 
00817           retval = bessel_return_value (tmp, ierr);
00818         }
00819       else
00820         retval = Complex (octave_NaN, octave_NaN);
00821     }
00822 
00823   return retval;
00824 }
00825 
00826 static inline Complex
00827 zbesi (const Complex& z, double alpha, int kode, octave_idx_type& ierr)
00828 {
00829   Complex retval;
00830 
00831   if (alpha >= 0.0)
00832     {
00833       double yr = 0.0;
00834       double yi = 0.0;
00835 
00836       octave_idx_type nz;
00837 
00838       double zr = z.real ();
00839       double zi = z.imag ();
00840 
00841       F77_FUNC (zbesi, ZBESI) (zr, zi, alpha, 2, 1, &yr, &yi, nz, ierr);
00842 
00843       if (kode != 2)
00844         {
00845           double expz = exp (std::abs (zr));
00846           yr *= expz;
00847           yi *= expz;
00848         }
00849 
00850       if (zi == 0.0 && zr >= 0.0)
00851         yi = 0.0;
00852 
00853       retval = bessel_return_value (Complex (yr, yi), ierr);
00854     }
00855   else if (is_integer_value (alpha))
00856     {
00857       // zbesi can overflow as z->0, and cause troubles for generic case below
00858       alpha = -alpha;
00859       Complex tmp = zbesi (z, alpha, kode, ierr);
00860       retval = bessel_return_value (tmp, ierr);
00861     }
00862   else
00863     {
00864       alpha = -alpha;
00865 
00866       Complex tmp = zbesi (z, alpha, kode, ierr);
00867 
00868       if (ierr == 0 || ierr == 3)
00869         {
00870           Complex tmp2 = (2.0 / M_PI) * sin (M_PI * alpha)
00871             * zbesk (z, alpha, kode, ierr);
00872 
00873           if (kode == 2)
00874             {
00875               // Compensate for different scaling factor of besk.
00876               tmp2 *= exp(-z - std::abs(z.real()));
00877             }
00878 
00879           tmp += tmp2;
00880 
00881           retval = bessel_return_value (tmp, ierr);
00882         }
00883       else
00884         retval = Complex (octave_NaN, octave_NaN);
00885     }
00886 
00887   return retval;
00888 }
00889 
00890 static inline Complex
00891 zbesk (const Complex& z, double alpha, int kode, octave_idx_type& ierr)
00892 {
00893   Complex retval;
00894 
00895   if (alpha >= 0.0)
00896     {
00897       double yr = 0.0;
00898       double yi = 0.0;
00899 
00900       octave_idx_type nz;
00901 
00902       double zr = z.real ();
00903       double zi = z.imag ();
00904 
00905       ierr = 0;
00906 
00907       if (zr == 0.0 && zi == 0.0)
00908         {
00909           yr = octave_Inf;
00910           yi = 0.0;
00911         }
00912       else
00913         {
00914           F77_FUNC (zbesk, ZBESK) (zr, zi, alpha, 2, 1, &yr, &yi, nz, ierr);
00915 
00916           if (kode != 2)
00917             {
00918               Complex expz = exp (-z);
00919 
00920               double rexpz = real (expz);
00921               double iexpz = imag (expz);
00922 
00923               double tmp = yr*rexpz - yi*iexpz;
00924 
00925               yi = yr*iexpz + yi*rexpz;
00926               yr = tmp;
00927             }
00928 
00929           if (zi == 0.0 && zr >= 0.0)
00930             yi = 0.0;
00931         }
00932 
00933       retval = bessel_return_value (Complex (yr, yi), ierr);
00934     }
00935   else
00936     {
00937       Complex tmp = zbesk (z, -alpha, kode, ierr);
00938 
00939       retval = bessel_return_value (tmp, ierr);
00940     }
00941 
00942   return retval;
00943 }
00944 
00945 static inline Complex
00946 zbesh1 (const Complex& z, double alpha, int kode, octave_idx_type& ierr)
00947 {
00948   Complex retval;
00949 
00950   if (alpha >= 0.0)
00951     {
00952       double yr = 0.0;
00953       double yi = 0.0;
00954 
00955       octave_idx_type nz;
00956 
00957       double zr = z.real ();
00958       double zi = z.imag ();
00959 
00960       F77_FUNC (zbesh, ZBESH) (zr, zi, alpha, 2, 1, 1, &yr, &yi, nz, ierr);
00961 
00962       if (kode != 2)
00963         {
00964           Complex expz = exp (Complex (0.0, 1.0) * z);
00965 
00966           double rexpz = real (expz);
00967           double iexpz = imag (expz);
00968 
00969           double tmp = yr*rexpz - yi*iexpz;
00970 
00971           yi = yr*iexpz + yi*rexpz;
00972           yr = tmp;
00973         }
00974 
00975       retval = bessel_return_value (Complex (yr, yi), ierr);
00976     }
00977   else
00978     {
00979       alpha = -alpha;
00980 
00981       static const Complex eye = Complex (0.0, 1.0);
00982 
00983       Complex tmp = exp (M_PI * alpha * eye) * zbesh1 (z, alpha, kode, ierr);
00984 
00985       retval = bessel_return_value (tmp, ierr);
00986     }
00987 
00988   return retval;
00989 }
00990 
00991 static inline Complex
00992 zbesh2 (const Complex& z, double alpha, int kode, octave_idx_type& ierr)
00993 {
00994   Complex retval;
00995 
00996   if (alpha >= 0.0)
00997     {
00998       double yr = 0.0;
00999       double yi = 0.0;
01000 
01001       octave_idx_type nz;
01002 
01003       double zr = z.real ();
01004       double zi = z.imag ();
01005 
01006       F77_FUNC (zbesh, ZBESH) (zr, zi, alpha, 2, 2, 1, &yr, &yi, nz, ierr);
01007 
01008       if (kode != 2)
01009         {
01010           Complex expz = exp (-Complex (0.0, 1.0) * z);
01011 
01012           double rexpz = real (expz);
01013           double iexpz = imag (expz);
01014 
01015           double tmp = yr*rexpz - yi*iexpz;
01016 
01017           yi = yr*iexpz + yi*rexpz;
01018           yr = tmp;
01019         }
01020 
01021       retval = bessel_return_value (Complex (yr, yi), ierr);
01022     }
01023   else
01024     {
01025       alpha = -alpha;
01026 
01027       static const Complex eye = Complex (0.0, 1.0);
01028 
01029       Complex tmp = exp (-M_PI * alpha * eye) * zbesh2 (z, alpha, kode, ierr);
01030 
01031       retval = bessel_return_value (tmp, ierr);
01032     }
01033 
01034   return retval;
01035 }
01036 
01037 typedef Complex (*dptr) (const Complex&, double, int, octave_idx_type&);
01038 
01039 static inline Complex
01040 do_bessel (dptr f, const char *, double alpha, const Complex& x,
01041            bool scaled, octave_idx_type& ierr)
01042 {
01043   Complex retval;
01044 
01045   retval = f (x, alpha, (scaled ? 2 : 1), ierr);
01046 
01047   return retval;
01048 }
01049 
01050 static inline ComplexMatrix
01051 do_bessel (dptr f, const char *, double alpha, const ComplexMatrix& x,
01052            bool scaled, Array<octave_idx_type>& ierr)
01053 {
01054   octave_idx_type nr = x.rows ();
01055   octave_idx_type nc = x.cols ();
01056 
01057   ComplexMatrix retval (nr, nc);
01058 
01059   ierr.resize (dim_vector (nr, nc));
01060 
01061   for (octave_idx_type j = 0; j < nc; j++)
01062     for (octave_idx_type i = 0; i < nr; i++)
01063       retval(i,j) = f (x(i,j), alpha, (scaled ? 2 : 1), ierr(i,j));
01064 
01065   return retval;
01066 }
01067 
01068 static inline ComplexMatrix
01069 do_bessel (dptr f, const char *, const Matrix& alpha, const Complex& x,
01070            bool scaled, Array<octave_idx_type>& ierr)
01071 {
01072   octave_idx_type nr = alpha.rows ();
01073   octave_idx_type nc = alpha.cols ();
01074 
01075   ComplexMatrix retval (nr, nc);
01076 
01077   ierr.resize (dim_vector (nr, nc));
01078 
01079   for (octave_idx_type j = 0; j < nc; j++)
01080     for (octave_idx_type i = 0; i < nr; i++)
01081       retval(i,j) = f (x, alpha(i,j), (scaled ? 2 : 1), ierr(i,j));
01082 
01083   return retval;
01084 }
01085 
01086 static inline ComplexMatrix
01087 do_bessel (dptr f, const char *fn, const Matrix& alpha,
01088            const ComplexMatrix& x, bool scaled, Array<octave_idx_type>& ierr)
01089 {
01090   ComplexMatrix retval;
01091 
01092   octave_idx_type x_nr = x.rows ();
01093   octave_idx_type x_nc = x.cols ();
01094 
01095   octave_idx_type alpha_nr = alpha.rows ();
01096   octave_idx_type alpha_nc = alpha.cols ();
01097 
01098   if (x_nr == alpha_nr && x_nc == alpha_nc)
01099     {
01100       octave_idx_type nr = x_nr;
01101       octave_idx_type nc = x_nc;
01102 
01103       retval.resize (nr, nc);
01104 
01105       ierr.resize (dim_vector (nr, nc));
01106 
01107       for (octave_idx_type j = 0; j < nc; j++)
01108         for (octave_idx_type i = 0; i < nr; i++)
01109           retval(i,j) = f (x(i,j), alpha(i,j), (scaled ? 2 : 1), ierr(i,j));
01110     }
01111   else
01112     (*current_liboctave_error_handler)
01113       ("%s: the sizes of alpha and x must conform", fn);
01114 
01115   return retval;
01116 }
01117 
01118 static inline ComplexNDArray
01119 do_bessel (dptr f, const char *, double alpha, const ComplexNDArray& x,
01120            bool scaled, Array<octave_idx_type>& ierr)
01121 {
01122   dim_vector dv = x.dims ();
01123   octave_idx_type nel = dv.numel ();
01124   ComplexNDArray retval (dv);
01125 
01126   ierr.resize (dv);
01127 
01128   for (octave_idx_type i = 0; i < nel; i++)
01129       retval(i) = f (x(i), alpha, (scaled ? 2 : 1), ierr(i));
01130 
01131   return retval;
01132 }
01133 
01134 static inline ComplexNDArray
01135 do_bessel (dptr f, const char *, const NDArray& alpha, const Complex& x,
01136            bool scaled, Array<octave_idx_type>& ierr)
01137 {
01138   dim_vector dv = alpha.dims ();
01139   octave_idx_type nel = dv.numel ();
01140   ComplexNDArray retval (dv);
01141 
01142   ierr.resize (dv);
01143 
01144   for (octave_idx_type i = 0; i < nel; i++)
01145     retval(i) = f (x, alpha(i), (scaled ? 2 : 1), ierr(i));
01146 
01147   return retval;
01148 }
01149 
01150 static inline ComplexNDArray
01151 do_bessel (dptr f, const char *fn, const NDArray& alpha,
01152            const ComplexNDArray& x, bool scaled, Array<octave_idx_type>& ierr)
01153 {
01154   dim_vector dv = x.dims ();
01155   ComplexNDArray retval;
01156 
01157   if (dv == alpha.dims ())
01158     {
01159       octave_idx_type nel = dv.numel ();
01160 
01161       retval.resize (dv);
01162       ierr.resize (dv);
01163 
01164       for (octave_idx_type i = 0; i < nel; i++)
01165         retval(i) = f (x(i), alpha(i), (scaled ? 2 : 1), ierr(i));
01166     }
01167   else
01168     (*current_liboctave_error_handler)
01169       ("%s: the sizes of alpha and x must conform", fn);
01170 
01171   return retval;
01172 }
01173 
01174 static inline ComplexMatrix
01175 do_bessel (dptr f, const char *, const RowVector& alpha,
01176            const ComplexColumnVector& x, bool scaled, Array<octave_idx_type>& ierr)
01177 {
01178   octave_idx_type nr = x.length ();
01179   octave_idx_type nc = alpha.length ();
01180 
01181   ComplexMatrix retval (nr, nc);
01182 
01183   ierr.resize (dim_vector (nr, nc));
01184 
01185   for (octave_idx_type j = 0; j < nc; j++)
01186     for (octave_idx_type i = 0; i < nr; i++)
01187       retval(i,j) = f (x(i), alpha(j), (scaled ? 2 : 1), ierr(i,j));
01188 
01189   return retval;
01190 }
01191 
01192 #define SS_BESSEL(name, fcn) \
01193   Complex \
01194   name (double alpha, const Complex& x, bool scaled, octave_idx_type& ierr) \
01195   { \
01196     return do_bessel (fcn, #name, alpha, x, scaled, ierr); \
01197   }
01198 
01199 #define SM_BESSEL(name, fcn) \
01200   ComplexMatrix \
01201   name (double alpha, const ComplexMatrix& x, bool scaled, \
01202         Array<octave_idx_type>& ierr) \
01203   { \
01204     return do_bessel (fcn, #name, alpha, x, scaled, ierr); \
01205   }
01206 
01207 #define MS_BESSEL(name, fcn) \
01208   ComplexMatrix \
01209   name (const Matrix& alpha, const Complex& x, bool scaled, \
01210         Array<octave_idx_type>& ierr) \
01211   { \
01212     return do_bessel (fcn, #name, alpha, x, scaled, ierr); \
01213   }
01214 
01215 #define MM_BESSEL(name, fcn) \
01216   ComplexMatrix \
01217   name (const Matrix& alpha, const ComplexMatrix& x, bool scaled, \
01218         Array<octave_idx_type>& ierr) \
01219   { \
01220     return do_bessel (fcn, #name, alpha, x, scaled, ierr); \
01221   }
01222 
01223 #define SN_BESSEL(name, fcn) \
01224   ComplexNDArray \
01225   name (double alpha, const ComplexNDArray& x, bool scaled, \
01226         Array<octave_idx_type>& ierr) \
01227   { \
01228     return do_bessel (fcn, #name, alpha, x, scaled, ierr); \
01229   }
01230 
01231 #define NS_BESSEL(name, fcn) \
01232   ComplexNDArray \
01233   name (const NDArray& alpha, const Complex& x, bool scaled, \
01234         Array<octave_idx_type>& ierr) \
01235   { \
01236     return do_bessel (fcn, #name, alpha, x, scaled, ierr); \
01237   }
01238 
01239 #define NN_BESSEL(name, fcn) \
01240   ComplexNDArray \
01241   name (const NDArray& alpha, const ComplexNDArray& x, bool scaled, \
01242         Array<octave_idx_type>& ierr) \
01243   { \
01244     return do_bessel (fcn, #name, alpha, x, scaled, ierr); \
01245   }
01246 
01247 #define RC_BESSEL(name, fcn) \
01248   ComplexMatrix \
01249   name (const RowVector& alpha, const ComplexColumnVector& x, bool scaled, \
01250         Array<octave_idx_type>& ierr) \
01251   { \
01252     return do_bessel (fcn, #name, alpha, x, scaled, ierr); \
01253   }
01254 
01255 #define ALL_BESSEL(name, fcn) \
01256   SS_BESSEL (name, fcn) \
01257   SM_BESSEL (name, fcn) \
01258   MS_BESSEL (name, fcn) \
01259   MM_BESSEL (name, fcn) \
01260   SN_BESSEL (name, fcn) \
01261   NS_BESSEL (name, fcn) \
01262   NN_BESSEL (name, fcn) \
01263   RC_BESSEL (name, fcn)
01264 
01265 ALL_BESSEL (besselj, zbesj)
01266 ALL_BESSEL (bessely, zbesy)
01267 ALL_BESSEL (besseli, zbesi)
01268 ALL_BESSEL (besselk, zbesk)
01269 ALL_BESSEL (besselh1, zbesh1)
01270 ALL_BESSEL (besselh2, zbesh2)
01271 
01272 #undef ALL_BESSEL
01273 #undef SS_BESSEL
01274 #undef SM_BESSEL
01275 #undef MS_BESSEL
01276 #undef MM_BESSEL
01277 #undef SN_BESSEL
01278 #undef NS_BESSEL
01279 #undef NN_BESSEL
01280 #undef RC_BESSEL
01281 
01282 static inline FloatComplex
01283 cbesj (const FloatComplex& z, float alpha, int kode, octave_idx_type& ierr);
01284 
01285 static inline FloatComplex
01286 cbesy (const FloatComplex& z, float alpha, int kode, octave_idx_type& ierr);
01287 
01288 static inline FloatComplex
01289 cbesi (const FloatComplex& z, float alpha, int kode, octave_idx_type& ierr);
01290 
01291 static inline FloatComplex
01292 cbesk (const FloatComplex& z, float alpha, int kode, octave_idx_type& ierr);
01293 
01294 static inline FloatComplex
01295 cbesh1 (const FloatComplex& z, float alpha, int kode, octave_idx_type& ierr);
01296 
01297 static inline FloatComplex
01298 cbesh2 (const FloatComplex& z, float alpha, int kode, octave_idx_type& ierr);
01299 
01300 static inline FloatComplex
01301 bessel_return_value (const FloatComplex& val, octave_idx_type ierr)
01302 {
01303   static const FloatComplex inf_val = FloatComplex (octave_Float_Inf, octave_Float_Inf);
01304   static const FloatComplex nan_val = FloatComplex (octave_Float_NaN, octave_Float_NaN);
01305 
01306   FloatComplex retval;
01307 
01308   switch (ierr)
01309     {
01310     case 0:
01311     case 3:
01312       retval = val;
01313       break;
01314 
01315     case 2:
01316       retval = inf_val;
01317       break;
01318 
01319     default:
01320       retval = nan_val;
01321       break;
01322     }
01323 
01324   return retval;
01325 }
01326 
01327 static inline bool
01328 is_integer_value (float x)
01329 {
01330   return x == static_cast<float> (static_cast<long> (x));
01331 }
01332 
01333 static inline FloatComplex
01334 cbesj (const FloatComplex& z, float alpha, int kode, octave_idx_type& ierr)
01335 {
01336   FloatComplex retval;
01337 
01338   if (alpha >= 0.0)
01339     {
01340       FloatComplex y = 0.0;
01341 
01342       octave_idx_type nz;
01343 
01344       F77_FUNC (cbesj, CBESJ) (z, alpha, 2, 1, &y, nz, ierr);
01345 
01346       if (kode != 2)
01347         {
01348           float expz = exp (std::abs (imag (z)));
01349           y *= expz;
01350         }
01351 
01352       if (imag (z) == 0.0 && real (z) >= 0.0)
01353         y = FloatComplex (y.real (), 0.0);
01354 
01355       retval = bessel_return_value (y, ierr);
01356     }
01357   else if (is_integer_value (alpha))
01358     {
01359       // zbesy can overflow as z->0, and cause troubles for generic case below
01360       alpha = -alpha;
01361       FloatComplex tmp = cbesj (z, alpha, kode, ierr);
01362       if ((static_cast <long> (alpha)) & 1)
01363         tmp = - tmp;
01364       retval = bessel_return_value (tmp, ierr);
01365     }
01366   else
01367     {
01368       alpha = -alpha;
01369 
01370       FloatComplex tmp = cosf (static_cast<float> (M_PI) * alpha) * cbesj (z, alpha, kode, ierr);
01371 
01372       if (ierr == 0 || ierr == 3)
01373         {
01374           tmp -= sinf (static_cast<float> (M_PI) * alpha) * cbesy (z, alpha, kode, ierr);
01375 
01376           retval = bessel_return_value (tmp, ierr);
01377         }
01378       else
01379         retval = FloatComplex (octave_Float_NaN, octave_Float_NaN);
01380     }
01381 
01382   return retval;
01383 }
01384 
01385 static inline FloatComplex
01386 cbesy (const FloatComplex& z, float alpha, int kode, octave_idx_type& ierr)
01387 {
01388   FloatComplex retval;
01389 
01390   if (alpha >= 0.0)
01391     {
01392       FloatComplex y = 0.0;
01393 
01394       octave_idx_type nz;
01395 
01396       FloatComplex w;
01397 
01398       ierr = 0;
01399 
01400       if (real (z) == 0.0 && imag (z) == 0.0)
01401         {
01402           y = FloatComplex (-octave_Float_Inf, 0.0);
01403         }
01404       else
01405         {
01406           F77_FUNC (cbesy, CBESY) (z, alpha, 2, 1, &y, nz, &w, ierr);
01407 
01408           if (kode != 2)
01409             {
01410               float expz = exp (std::abs (imag (z)));
01411               y *= expz;
01412             }
01413 
01414           if (imag (z) == 0.0 && real (z) >= 0.0)
01415             y = FloatComplex (y.real (), 0.0);
01416         }
01417 
01418       return bessel_return_value (y, ierr);
01419     }
01420   else if (is_integer_value (alpha - 0.5))
01421     {
01422       // zbesy can overflow as z->0, and cause troubles for generic case below
01423       alpha = -alpha;
01424       FloatComplex tmp = cbesj (z, alpha, kode, ierr);
01425       if ((static_cast <long> (alpha - 0.5)) & 1)
01426         tmp = - tmp;
01427       retval = bessel_return_value (tmp, ierr);
01428     }
01429   else
01430     {
01431       alpha = -alpha;
01432 
01433       FloatComplex tmp = cosf (static_cast<float> (M_PI) * alpha) * cbesy (z, alpha, kode, ierr);
01434 
01435       if (ierr == 0 || ierr == 3)
01436         {
01437           tmp += sinf (static_cast<float> (M_PI) * alpha) * cbesj (z, alpha, kode, ierr);
01438 
01439           retval = bessel_return_value (tmp, ierr);
01440         }
01441       else
01442         retval = FloatComplex (octave_Float_NaN, octave_Float_NaN);
01443     }
01444 
01445   return retval;
01446 }
01447 
01448 static inline FloatComplex
01449 cbesi (const FloatComplex& z, float alpha, int kode, octave_idx_type& ierr)
01450 {
01451   FloatComplex retval;
01452 
01453   if (alpha >= 0.0)
01454     {
01455       FloatComplex y = 0.0;
01456 
01457       octave_idx_type nz;
01458 
01459       F77_FUNC (cbesi, CBESI) (z, alpha, 2, 1, &y, nz, ierr);
01460 
01461       if (kode != 2)
01462         {
01463           float expz = exp (std::abs (real (z)));
01464           y *= expz;
01465         }
01466 
01467       if (imag (z) == 0.0 && real (z) >= 0.0)
01468         y = FloatComplex (y.real (), 0.0);
01469 
01470       retval = bessel_return_value (y, ierr);
01471     }
01472   else
01473     {
01474       alpha = -alpha;
01475 
01476       FloatComplex tmp = cbesi (z, alpha, kode, ierr);
01477 
01478       if (ierr == 0 || ierr == 3)
01479         {
01480           FloatComplex tmp2 = static_cast<float> (2.0 / M_PI) * sinf (static_cast<float> (M_PI) * alpha)
01481             * cbesk (z, alpha, kode, ierr);
01482 
01483           if (kode == 2)
01484             {
01485               // Compensate for different scaling factor of besk.
01486               tmp2 *= exp(-z - std::abs(z.real()));
01487             }
01488 
01489           tmp += tmp2;
01490 
01491           retval = bessel_return_value (tmp, ierr);
01492         }
01493       else
01494         retval = FloatComplex (octave_Float_NaN, octave_Float_NaN);
01495     }
01496 
01497   return retval;
01498 }
01499 
01500 static inline FloatComplex
01501 cbesk (const FloatComplex& z, float alpha, int kode, octave_idx_type& ierr)
01502 {
01503   FloatComplex retval;
01504 
01505   if (alpha >= 0.0)
01506     {
01507       FloatComplex y = 0.0;
01508 
01509       octave_idx_type nz;
01510 
01511       ierr = 0;
01512 
01513       if (real (z) == 0.0 && imag (z) == 0.0)
01514         {
01515           y = FloatComplex (octave_Float_Inf, 0.0);
01516         }
01517       else
01518         {
01519           F77_FUNC (cbesk, CBESK) (z, alpha, 2, 1, &y, nz, ierr);
01520 
01521           if (kode != 2)
01522             {
01523               FloatComplex expz = exp (-z);
01524 
01525               float rexpz = real (expz);
01526               float iexpz = imag (expz);
01527 
01528               float tmp_r = real (y) * rexpz - imag (y) * iexpz;
01529               float tmp_i = real (y) * iexpz + imag (y) * rexpz;
01530 
01531               y = FloatComplex (tmp_r, tmp_i);
01532             }
01533 
01534           if (imag (z) == 0.0 && real (z) >= 0.0)
01535             y = FloatComplex (y.real (), 0.0);
01536         }
01537 
01538       retval = bessel_return_value (y, ierr);
01539     }
01540   else
01541     {
01542       FloatComplex tmp = cbesk (z, -alpha, kode, ierr);
01543 
01544       retval = bessel_return_value (tmp, ierr);
01545     }
01546 
01547   return retval;
01548 }
01549 
01550 static inline FloatComplex
01551 cbesh1 (const FloatComplex& z, float alpha, int kode, octave_idx_type& ierr)
01552 {
01553   FloatComplex retval;
01554 
01555   if (alpha >= 0.0)
01556     {
01557       FloatComplex y = 0.0;
01558 
01559       octave_idx_type nz;
01560 
01561       F77_FUNC (cbesh, CBESH) (z, alpha, 2, 1, 1, &y, nz, ierr);
01562 
01563       if (kode != 2)
01564         {
01565           FloatComplex expz = exp (FloatComplex (0.0, 1.0) * z);
01566 
01567           float rexpz = real (expz);
01568           float iexpz = imag (expz);
01569 
01570           float tmp_r = real (y) * rexpz - imag (y) * iexpz;
01571           float tmp_i = real (y) * iexpz + imag (y) * rexpz;
01572 
01573           y = FloatComplex (tmp_r, tmp_i);
01574         }
01575 
01576       retval = bessel_return_value (y, ierr);
01577     }
01578   else
01579     {
01580       alpha = -alpha;
01581 
01582       static const FloatComplex eye = FloatComplex (0.0, 1.0);
01583 
01584       FloatComplex tmp = exp (static_cast<float> (M_PI) * alpha * eye) * cbesh1 (z, alpha, kode, ierr);
01585 
01586       retval = bessel_return_value (tmp, ierr);
01587     }
01588 
01589   return retval;
01590 }
01591 
01592 static inline FloatComplex
01593 cbesh2 (const FloatComplex& z, float alpha, int kode, octave_idx_type& ierr)
01594 {
01595   FloatComplex retval;
01596 
01597   if (alpha >= 0.0)
01598     {
01599       FloatComplex y = 0.0;
01600 
01601       octave_idx_type nz;
01602 
01603       F77_FUNC (cbesh, CBESH) (z, alpha, 2, 2, 1, &y, nz, ierr);
01604 
01605       if (kode != 2)
01606         {
01607           FloatComplex expz = exp (-FloatComplex (0.0, 1.0) * z);
01608 
01609           float rexpz = real (expz);
01610           float iexpz = imag (expz);
01611 
01612           float tmp_r = real (y) * rexpz - imag (y) * iexpz;
01613           float tmp_i = real (y) * iexpz + imag (y) * rexpz;
01614 
01615           y = FloatComplex (tmp_r, tmp_i);
01616         }
01617 
01618       retval = bessel_return_value (y, ierr);
01619     }
01620   else
01621     {
01622       alpha = -alpha;
01623 
01624       static const FloatComplex eye = FloatComplex (0.0, 1.0);
01625 
01626       FloatComplex tmp = exp (-static_cast<float> (M_PI) * alpha * eye) * cbesh2 (z, alpha, kode, ierr);
01627 
01628       retval = bessel_return_value (tmp, ierr);
01629     }
01630 
01631   return retval;
01632 }
01633 
01634 typedef FloatComplex (*fptr) (const FloatComplex&, float, int, octave_idx_type&);
01635 
01636 static inline FloatComplex
01637 do_bessel (fptr f, const char *, float alpha, const FloatComplex& x,
01638            bool scaled, octave_idx_type& ierr)
01639 {
01640   FloatComplex retval;
01641 
01642   retval = f (x, alpha, (scaled ? 2 : 1), ierr);
01643 
01644   return retval;
01645 }
01646 
01647 static inline FloatComplexMatrix
01648 do_bessel (fptr f, const char *, float alpha, const FloatComplexMatrix& x,
01649            bool scaled, Array<octave_idx_type>& ierr)
01650 {
01651   octave_idx_type nr = x.rows ();
01652   octave_idx_type nc = x.cols ();
01653 
01654   FloatComplexMatrix retval (nr, nc);
01655 
01656   ierr.resize (dim_vector (nr, nc));
01657 
01658   for (octave_idx_type j = 0; j < nc; j++)
01659     for (octave_idx_type i = 0; i < nr; i++)
01660       retval(i,j) = f (x(i,j), alpha, (scaled ? 2 : 1), ierr(i,j));
01661 
01662   return retval;
01663 }
01664 
01665 static inline FloatComplexMatrix
01666 do_bessel (fptr f, const char *, const FloatMatrix& alpha, const FloatComplex& x,
01667            bool scaled, Array<octave_idx_type>& ierr)
01668 {
01669   octave_idx_type nr = alpha.rows ();
01670   octave_idx_type nc = alpha.cols ();
01671 
01672   FloatComplexMatrix retval (nr, nc);
01673 
01674   ierr.resize (dim_vector (nr, nc));
01675 
01676   for (octave_idx_type j = 0; j < nc; j++)
01677     for (octave_idx_type i = 0; i < nr; i++)
01678       retval(i,j) = f (x, alpha(i,j), (scaled ? 2 : 1), ierr(i,j));
01679 
01680   return retval;
01681 }
01682 
01683 static inline FloatComplexMatrix
01684 do_bessel (fptr f, const char *fn, const FloatMatrix& alpha,
01685            const FloatComplexMatrix& x, bool scaled, Array<octave_idx_type>& ierr)
01686 {
01687   FloatComplexMatrix retval;
01688 
01689   octave_idx_type x_nr = x.rows ();
01690   octave_idx_type x_nc = x.cols ();
01691 
01692   octave_idx_type alpha_nr = alpha.rows ();
01693   octave_idx_type alpha_nc = alpha.cols ();
01694 
01695   if (x_nr == alpha_nr && x_nc == alpha_nc)
01696     {
01697       octave_idx_type nr = x_nr;
01698       octave_idx_type nc = x_nc;
01699 
01700       retval.resize (nr, nc);
01701 
01702       ierr.resize (dim_vector (nr, nc));
01703 
01704       for (octave_idx_type j = 0; j < nc; j++)
01705         for (octave_idx_type i = 0; i < nr; i++)
01706           retval(i,j) = f (x(i,j), alpha(i,j), (scaled ? 2 : 1), ierr(i,j));
01707     }
01708   else
01709     (*current_liboctave_error_handler)
01710       ("%s: the sizes of alpha and x must conform", fn);
01711 
01712   return retval;
01713 }
01714 
01715 static inline FloatComplexNDArray
01716 do_bessel (fptr f, const char *, float alpha, const FloatComplexNDArray& x,
01717            bool scaled, Array<octave_idx_type>& ierr)
01718 {
01719   dim_vector dv = x.dims ();
01720   octave_idx_type nel = dv.numel ();
01721   FloatComplexNDArray retval (dv);
01722 
01723   ierr.resize (dv);
01724 
01725   for (octave_idx_type i = 0; i < nel; i++)
01726       retval(i) = f (x(i), alpha, (scaled ? 2 : 1), ierr(i));
01727 
01728   return retval;
01729 }
01730 
01731 static inline FloatComplexNDArray
01732 do_bessel (fptr f, const char *, const FloatNDArray& alpha, const FloatComplex& x,
01733            bool scaled, Array<octave_idx_type>& ierr)
01734 {
01735   dim_vector dv = alpha.dims ();
01736   octave_idx_type nel = dv.numel ();
01737   FloatComplexNDArray retval (dv);
01738 
01739   ierr.resize (dv);
01740 
01741   for (octave_idx_type i = 0; i < nel; i++)
01742     retval(i) = f (x, alpha(i), (scaled ? 2 : 1), ierr(i));
01743 
01744   return retval;
01745 }
01746 
01747 static inline FloatComplexNDArray
01748 do_bessel (fptr f, const char *fn, const FloatNDArray& alpha,
01749            const FloatComplexNDArray& x, bool scaled, Array<octave_idx_type>& ierr)
01750 {
01751   dim_vector dv = x.dims ();
01752   FloatComplexNDArray retval;
01753 
01754   if (dv == alpha.dims ())
01755     {
01756       octave_idx_type nel = dv.numel ();
01757 
01758       retval.resize (dv);
01759       ierr.resize (dv);
01760 
01761       for (octave_idx_type i = 0; i < nel; i++)
01762         retval(i) = f (x(i), alpha(i), (scaled ? 2 : 1), ierr(i));
01763     }
01764   else
01765     (*current_liboctave_error_handler)
01766       ("%s: the sizes of alpha and x must conform", fn);
01767 
01768   return retval;
01769 }
01770 
01771 static inline FloatComplexMatrix
01772 do_bessel (fptr f, const char *, const FloatRowVector& alpha,
01773            const FloatComplexColumnVector& x, bool scaled, Array<octave_idx_type>& ierr)
01774 {
01775   octave_idx_type nr = x.length ();
01776   octave_idx_type nc = alpha.length ();
01777 
01778   FloatComplexMatrix retval (nr, nc);
01779 
01780   ierr.resize (dim_vector (nr, nc));
01781 
01782   for (octave_idx_type j = 0; j < nc; j++)
01783     for (octave_idx_type i = 0; i < nr; i++)
01784       retval(i,j) = f (x(i), alpha(j), (scaled ? 2 : 1), ierr(i,j));
01785 
01786   return retval;
01787 }
01788 
01789 #define SS_BESSEL(name, fcn) \
01790   FloatComplex \
01791   name (float alpha, const FloatComplex& x, bool scaled, octave_idx_type& ierr) \
01792   { \
01793     return do_bessel (fcn, #name, alpha, x, scaled, ierr); \
01794   }
01795 
01796 #define SM_BESSEL(name, fcn) \
01797   FloatComplexMatrix \
01798   name (float alpha, const FloatComplexMatrix& x, bool scaled, \
01799         Array<octave_idx_type>& ierr) \
01800   { \
01801     return do_bessel (fcn, #name, alpha, x, scaled, ierr); \
01802   }
01803 
01804 #define MS_BESSEL(name, fcn) \
01805   FloatComplexMatrix \
01806   name (const FloatMatrix& alpha, const FloatComplex& x, bool scaled, \
01807         Array<octave_idx_type>& ierr) \
01808   { \
01809     return do_bessel (fcn, #name, alpha, x, scaled, ierr); \
01810   }
01811 
01812 #define MM_BESSEL(name, fcn) \
01813   FloatComplexMatrix \
01814   name (const FloatMatrix& alpha, const FloatComplexMatrix& x, bool scaled, \
01815         Array<octave_idx_type>& ierr) \
01816   { \
01817     return do_bessel (fcn, #name, alpha, x, scaled, ierr); \
01818   }
01819 
01820 #define SN_BESSEL(name, fcn) \
01821   FloatComplexNDArray \
01822   name (float alpha, const FloatComplexNDArray& x, bool scaled, \
01823         Array<octave_idx_type>& ierr) \
01824   { \
01825     return do_bessel (fcn, #name, alpha, x, scaled, ierr); \
01826   }
01827 
01828 #define NS_BESSEL(name, fcn) \
01829   FloatComplexNDArray \
01830   name (const FloatNDArray& alpha, const FloatComplex& x, bool scaled, \
01831         Array<octave_idx_type>& ierr) \
01832   { \
01833     return do_bessel (fcn, #name, alpha, x, scaled, ierr); \
01834   }
01835 
01836 #define NN_BESSEL(name, fcn) \
01837   FloatComplexNDArray \
01838   name (const FloatNDArray& alpha, const FloatComplexNDArray& x, bool scaled, \
01839         Array<octave_idx_type>& ierr) \
01840   { \
01841     return do_bessel (fcn, #name, alpha, x, scaled, ierr); \
01842   }
01843 
01844 #define RC_BESSEL(name, fcn) \
01845   FloatComplexMatrix \
01846   name (const FloatRowVector& alpha, const FloatComplexColumnVector& x, bool scaled, \
01847         Array<octave_idx_type>& ierr) \
01848   { \
01849     return do_bessel (fcn, #name, alpha, x, scaled, ierr); \
01850   }
01851 
01852 #define ALL_BESSEL(name, fcn) \
01853   SS_BESSEL (name, fcn) \
01854   SM_BESSEL (name, fcn) \
01855   MS_BESSEL (name, fcn) \
01856   MM_BESSEL (name, fcn) \
01857   SN_BESSEL (name, fcn) \
01858   NS_BESSEL (name, fcn) \
01859   NN_BESSEL (name, fcn) \
01860   RC_BESSEL (name, fcn)
01861 
01862 ALL_BESSEL (besselj, cbesj)
01863 ALL_BESSEL (bessely, cbesy)
01864 ALL_BESSEL (besseli, cbesi)
01865 ALL_BESSEL (besselk, cbesk)
01866 ALL_BESSEL (besselh1, cbesh1)
01867 ALL_BESSEL (besselh2, cbesh2)
01868 
01869 #undef ALL_BESSEL
01870 #undef SS_BESSEL
01871 #undef SM_BESSEL
01872 #undef MS_BESSEL
01873 #undef MM_BESSEL
01874 #undef SN_BESSEL
01875 #undef NS_BESSEL
01876 #undef NN_BESSEL
01877 #undef RC_BESSEL
01878 
01879 Complex
01880 airy (const Complex& z, bool deriv, bool scaled, octave_idx_type& ierr)
01881 {
01882   double ar = 0.0;
01883   double ai = 0.0;
01884 
01885   octave_idx_type nz;
01886 
01887   double zr = z.real ();
01888   double zi = z.imag ();
01889 
01890   octave_idx_type id = deriv ? 1 : 0;
01891 
01892   F77_FUNC (zairy, ZAIRY) (zr, zi, id, 2, ar, ai, nz, ierr);
01893 
01894   if (! scaled)
01895     {
01896       Complex expz = exp (- 2.0 / 3.0 * z * sqrt(z));
01897 
01898       double rexpz = real (expz);
01899       double iexpz = imag (expz);
01900 
01901       double tmp = ar*rexpz - ai*iexpz;
01902 
01903       ai = ar*iexpz + ai*rexpz;
01904       ar = tmp;
01905     }
01906 
01907   if (zi == 0.0 && (! scaled || zr >= 0.0))
01908     ai = 0.0;
01909 
01910   return bessel_return_value (Complex (ar, ai), ierr);
01911 }
01912 
01913 Complex
01914 biry (const Complex& z, bool deriv, bool scaled, octave_idx_type& ierr)
01915 {
01916   double ar = 0.0;
01917   double ai = 0.0;
01918 
01919   double zr = z.real ();
01920   double zi = z.imag ();
01921 
01922   octave_idx_type id = deriv ? 1 : 0;
01923 
01924   F77_FUNC (zbiry, ZBIRY) (zr, zi, id, 2, ar, ai, ierr);
01925 
01926   if (! scaled)
01927     {
01928       Complex expz = exp (std::abs (real (2.0 / 3.0 * z * sqrt (z))));
01929 
01930       double rexpz = real (expz);
01931       double iexpz = imag (expz);
01932 
01933       double tmp = ar*rexpz - ai*iexpz;
01934 
01935       ai = ar*iexpz + ai*rexpz;
01936       ar = tmp;
01937     }
01938 
01939   if (zi == 0.0 && (! scaled || zr >= 0.0))
01940     ai = 0.0;
01941 
01942   return bessel_return_value (Complex (ar, ai), ierr);
01943 }
01944 
01945 ComplexMatrix
01946 airy (const ComplexMatrix& z, bool deriv, bool scaled, Array<octave_idx_type>& ierr)
01947 {
01948   octave_idx_type nr = z.rows ();
01949   octave_idx_type nc = z.cols ();
01950 
01951   ComplexMatrix retval (nr, nc);
01952 
01953   ierr.resize (dim_vector (nr, nc));
01954 
01955   for (octave_idx_type j = 0; j < nc; j++)
01956     for (octave_idx_type i = 0; i < nr; i++)
01957       retval(i,j) = airy (z(i,j), deriv, scaled, ierr(i,j));
01958 
01959   return retval;
01960 }
01961 
01962 ComplexMatrix
01963 biry (const ComplexMatrix& z, bool deriv, bool scaled, Array<octave_idx_type>& ierr)
01964 {
01965   octave_idx_type nr = z.rows ();
01966   octave_idx_type nc = z.cols ();
01967 
01968   ComplexMatrix retval (nr, nc);
01969 
01970   ierr.resize (dim_vector (nr, nc));
01971 
01972   for (octave_idx_type j = 0; j < nc; j++)
01973     for (octave_idx_type i = 0; i < nr; i++)
01974       retval(i,j) = biry (z(i,j), deriv, scaled, ierr(i,j));
01975 
01976   return retval;
01977 }
01978 
01979 ComplexNDArray
01980 airy (const ComplexNDArray& z, bool deriv, bool scaled, Array<octave_idx_type>& ierr)
01981 {
01982   dim_vector dv = z.dims ();
01983   octave_idx_type nel = dv.numel ();
01984   ComplexNDArray retval (dv);
01985 
01986   ierr.resize (dv);
01987 
01988   for (octave_idx_type i = 0; i < nel; i++)
01989     retval (i) = airy (z(i), deriv, scaled, ierr(i));
01990 
01991   return retval;
01992 }
01993 
01994 ComplexNDArray
01995 biry (const ComplexNDArray& z, bool deriv, bool scaled, Array<octave_idx_type>& ierr)
01996 {
01997   dim_vector dv = z.dims ();
01998   octave_idx_type nel = dv.numel ();
01999   ComplexNDArray retval (dv);
02000 
02001   ierr.resize (dv);
02002 
02003   for (octave_idx_type i = 0; i < nel; i++)
02004     retval (i) = biry (z(i), deriv, scaled, ierr(i));
02005 
02006   return retval;
02007 }
02008 
02009 FloatComplex
02010 airy (const FloatComplex& z, bool deriv, bool scaled, octave_idx_type& ierr)
02011 {
02012   float ar = 0.0;
02013   float ai = 0.0;
02014 
02015   octave_idx_type nz;
02016 
02017   float zr = z.real ();
02018   float zi = z.imag ();
02019 
02020   octave_idx_type id = deriv ? 1 : 0;
02021 
02022   F77_FUNC (cairy, CAIRY) (zr, zi, id, 2, ar, ai, nz, ierr);
02023 
02024   if (! scaled)
02025     {
02026       FloatComplex expz = exp (- static_cast<float> (2.0 / 3.0) * z * sqrt(z));
02027 
02028       float rexpz = real (expz);
02029       float iexpz = imag (expz);
02030 
02031       float tmp = ar*rexpz - ai*iexpz;
02032 
02033       ai = ar*iexpz + ai*rexpz;
02034       ar = tmp;
02035     }
02036 
02037   if (zi == 0.0 && (! scaled || zr >= 0.0))
02038     ai = 0.0;
02039 
02040   return bessel_return_value (FloatComplex (ar, ai), ierr);
02041 }
02042 
02043 FloatComplex
02044 biry (const FloatComplex& z, bool deriv, bool scaled, octave_idx_type& ierr)
02045 {
02046   float ar = 0.0;
02047   float ai = 0.0;
02048 
02049   float zr = z.real ();
02050   float zi = z.imag ();
02051 
02052   octave_idx_type id = deriv ? 1 : 0;
02053 
02054   F77_FUNC (cbiry, CBIRY) (zr, zi, id, 2, ar, ai, ierr);
02055 
02056   if (! scaled)
02057     {
02058       FloatComplex expz = exp (std::abs (real (static_cast<float> (2.0 / 3.0) * z * sqrt (z))));
02059 
02060       float rexpz = real (expz);
02061       float iexpz = imag (expz);
02062 
02063       float tmp = ar*rexpz - ai*iexpz;
02064 
02065       ai = ar*iexpz + ai*rexpz;
02066       ar = tmp;
02067     }
02068 
02069   if (zi == 0.0 && (! scaled || zr >= 0.0))
02070     ai = 0.0;
02071 
02072   return bessel_return_value (FloatComplex (ar, ai), ierr);
02073 }
02074 
02075 FloatComplexMatrix
02076 airy (const FloatComplexMatrix& z, bool deriv, bool scaled, Array<octave_idx_type>& ierr)
02077 {
02078   octave_idx_type nr = z.rows ();
02079   octave_idx_type nc = z.cols ();
02080 
02081   FloatComplexMatrix retval (nr, nc);
02082 
02083   ierr.resize (dim_vector (nr, nc));
02084 
02085   for (octave_idx_type j = 0; j < nc; j++)
02086     for (octave_idx_type i = 0; i < nr; i++)
02087       retval(i,j) = airy (z(i,j), deriv, scaled, ierr(i,j));
02088 
02089   return retval;
02090 }
02091 
02092 FloatComplexMatrix
02093 biry (const FloatComplexMatrix& z, bool deriv, bool scaled, Array<octave_idx_type>& ierr)
02094 {
02095   octave_idx_type nr = z.rows ();
02096   octave_idx_type nc = z.cols ();
02097 
02098   FloatComplexMatrix retval (nr, nc);
02099 
02100   ierr.resize (dim_vector (nr, nc));
02101 
02102   for (octave_idx_type j = 0; j < nc; j++)
02103     for (octave_idx_type i = 0; i < nr; i++)
02104       retval(i,j) = biry (z(i,j), deriv, scaled, ierr(i,j));
02105 
02106   return retval;
02107 }
02108 
02109 FloatComplexNDArray
02110 airy (const FloatComplexNDArray& z, bool deriv, bool scaled, Array<octave_idx_type>& ierr)
02111 {
02112   dim_vector dv = z.dims ();
02113   octave_idx_type nel = dv.numel ();
02114   FloatComplexNDArray retval (dv);
02115 
02116   ierr.resize (dv);
02117 
02118   for (octave_idx_type i = 0; i < nel; i++)
02119     retval (i) = airy (z(i), deriv, scaled, ierr(i));
02120 
02121   return retval;
02122 }
02123 
02124 FloatComplexNDArray
02125 biry (const FloatComplexNDArray& z, bool deriv, bool scaled, Array<octave_idx_type>& ierr)
02126 {
02127   dim_vector dv = z.dims ();
02128   octave_idx_type nel = dv.numel ();
02129   FloatComplexNDArray retval (dv);
02130 
02131   ierr.resize (dv);
02132 
02133   for (octave_idx_type i = 0; i < nel; i++)
02134     retval (i) = biry (z(i), deriv, scaled, ierr(i));
02135 
02136   return retval;
02137 }
02138 
02139 static void
02140 gripe_betainc_nonconformant (octave_idx_type r1, octave_idx_type c1, octave_idx_type r2, octave_idx_type c2, octave_idx_type r3,
02141                              octave_idx_type c3)
02142 {
02143   (*current_liboctave_error_handler)
02144    ("betainc: nonconformant arguments (x is %dx%d, a is %dx%d, b is %dx%d)",
02145      r1, c1, r2, c2, r3, c3);
02146 }
02147 
02148 static void
02149 gripe_betainc_nonconformant (const dim_vector& d1, const dim_vector& d2,
02150                              const dim_vector& d3)
02151 {
02152   std::string d1_str = d1.str ();
02153   std::string d2_str = d2.str ();
02154   std::string d3_str = d3.str ();
02155 
02156   (*current_liboctave_error_handler)
02157   ("betainc: nonconformant arguments (x is %s, a is %s, b is %s)",
02158    d1_str.c_str (), d2_str.c_str (), d3_str.c_str ());
02159 }
02160 
02161 double
02162 betainc (double x, double a, double b)
02163 {
02164   double retval;
02165   F77_XFCN (xdbetai, XDBETAI, (x, a, b, retval));
02166   return retval;
02167 }
02168 
02169 Matrix
02170 betainc (double x, double a, const Matrix& b)
02171 {
02172   octave_idx_type nr = b.rows ();
02173   octave_idx_type nc = b.cols ();
02174 
02175   Matrix retval (nr, nc);
02176 
02177   for (octave_idx_type j = 0; j < nc; j++)
02178     for (octave_idx_type i = 0; i < nr; i++)
02179       retval(i,j) = betainc (x, a, b(i,j));
02180 
02181   return retval;
02182 }
02183 
02184 Matrix
02185 betainc (double x, const Matrix& a, double b)
02186 {
02187   octave_idx_type nr = a.rows ();
02188   octave_idx_type nc = a.cols ();
02189 
02190   Matrix retval (nr, nc);
02191 
02192   for (octave_idx_type j = 0; j < nc; j++)
02193     for (octave_idx_type i = 0; i < nr; i++)
02194       retval(i,j) = betainc (x, a(i,j), b);
02195 
02196   return retval;
02197 }
02198 
02199 Matrix
02200 betainc (double x, const Matrix& a, const Matrix& b)
02201 {
02202   Matrix retval;
02203 
02204   octave_idx_type a_nr = a.rows ();
02205   octave_idx_type a_nc = a.cols ();
02206 
02207   octave_idx_type b_nr = b.rows ();
02208   octave_idx_type b_nc = b.cols ();
02209 
02210   if (a_nr == b_nr && a_nc == b_nc)
02211     {
02212       retval.resize (a_nr, a_nc);
02213 
02214       for (octave_idx_type j = 0; j < a_nc; j++)
02215         for (octave_idx_type i = 0; i < a_nr; i++)
02216           retval(i,j) = betainc (x, a(i,j), b(i,j));
02217     }
02218   else
02219     gripe_betainc_nonconformant (1, 1, a_nr, a_nc, b_nr, b_nc);
02220 
02221   return retval;
02222 }
02223 
02224 NDArray
02225 betainc (double x, double a, const NDArray& b)
02226 {
02227   dim_vector dv = b.dims ();
02228   octave_idx_type nel = dv.numel ();
02229 
02230   NDArray retval (dv);
02231 
02232   for (octave_idx_type i = 0; i < nel; i++)
02233     retval (i) = betainc (x, a, b(i));
02234 
02235   return retval;
02236 }
02237 
02238 NDArray
02239 betainc (double x, const NDArray& a, double b)
02240 {
02241   dim_vector dv = a.dims ();
02242   octave_idx_type nel = dv.numel ();
02243 
02244   NDArray retval (dv);
02245 
02246   for (octave_idx_type i = 0; i < nel; i++)
02247     retval (i) = betainc (x, a(i), b);
02248 
02249   return retval;
02250 }
02251 
02252 NDArray
02253 betainc (double x, const NDArray& a, const NDArray& b)
02254 {
02255   NDArray retval;
02256   dim_vector dv = a.dims ();
02257 
02258   if (dv == b.dims ())
02259     {
02260       octave_idx_type nel = dv.numel ();
02261 
02262       retval.resize (dv);
02263 
02264       for (octave_idx_type i = 0; i < nel; i++)
02265         retval (i) = betainc (x, a(i), b(i));
02266     }
02267   else
02268     gripe_betainc_nonconformant (dim_vector (0, 0), dv, b.dims ());
02269 
02270   return retval;
02271 }
02272 
02273 
02274 Matrix
02275 betainc (const Matrix& x, double a, double b)
02276 {
02277   octave_idx_type nr = x.rows ();
02278   octave_idx_type nc = x.cols ();
02279 
02280   Matrix retval (nr, nc);
02281 
02282   for (octave_idx_type j = 0; j < nc; j++)
02283     for (octave_idx_type i = 0; i < nr; i++)
02284       retval(i,j) = betainc (x(i,j), a, b);
02285 
02286   return retval;
02287 }
02288 
02289 Matrix
02290 betainc (const Matrix& x, double a, const Matrix& b)
02291 {
02292   Matrix retval;
02293 
02294   octave_idx_type nr = x.rows ();
02295   octave_idx_type nc = x.cols ();
02296 
02297   octave_idx_type b_nr = b.rows ();
02298   octave_idx_type b_nc = b.cols ();
02299 
02300   if (nr == b_nr && nc == b_nc)
02301     {
02302       retval.resize (nr, nc);
02303 
02304       for (octave_idx_type j = 0; j < nc; j++)
02305         for (octave_idx_type i = 0; i < nr; i++)
02306           retval(i,j) = betainc (x(i,j), a, b(i,j));
02307     }
02308   else
02309     gripe_betainc_nonconformant (nr, nc, 1, 1, b_nr, b_nc);
02310 
02311   return retval;
02312 }
02313 
02314 Matrix
02315 betainc (const Matrix& x, const Matrix& a, double b)
02316 {
02317   Matrix retval;
02318 
02319   octave_idx_type nr = x.rows ();
02320   octave_idx_type nc = x.cols ();
02321 
02322   octave_idx_type a_nr = a.rows ();
02323   octave_idx_type a_nc = a.cols ();
02324 
02325   if (nr == a_nr && nc == a_nc)
02326     {
02327       retval.resize (nr, nc);
02328 
02329       for (octave_idx_type j = 0; j < nc; j++)
02330         for (octave_idx_type i = 0; i < nr; i++)
02331           retval(i,j) = betainc (x(i,j), a(i,j), b);
02332     }
02333   else
02334     gripe_betainc_nonconformant (nr, nc, a_nr, a_nc, 1, 1);
02335 
02336   return retval;
02337 }
02338 
02339 Matrix
02340 betainc (const Matrix& x, const Matrix& a, const Matrix& b)
02341 {
02342   Matrix retval;
02343 
02344   octave_idx_type nr = x.rows ();
02345   octave_idx_type nc = x.cols ();
02346 
02347   octave_idx_type a_nr = a.rows ();
02348   octave_idx_type a_nc = a.cols ();
02349 
02350   octave_idx_type b_nr = b.rows ();
02351   octave_idx_type b_nc = b.cols ();
02352 
02353   if (nr == a_nr && nr == b_nr && nc == a_nc && nc == b_nc)
02354     {
02355       retval.resize (nr, nc);
02356 
02357       for (octave_idx_type j = 0; j < nc; j++)
02358         for (octave_idx_type i = 0; i < nr; i++)
02359           retval(i,j) = betainc (x(i,j), a(i,j), b(i,j));
02360     }
02361   else
02362     gripe_betainc_nonconformant (nr, nc, a_nr, a_nc, b_nr, b_nc);
02363 
02364   return retval;
02365 }
02366 
02367 NDArray
02368 betainc (const NDArray& x, double a, double b)
02369 {
02370   dim_vector dv = x.dims ();
02371   octave_idx_type nel = dv.numel ();
02372 
02373   NDArray retval (dv);
02374 
02375   for (octave_idx_type i = 0; i < nel; i++)
02376     retval (i) = betainc (x(i), a, b);
02377 
02378   return retval;
02379 }
02380 
02381 NDArray
02382 betainc (const NDArray& x, double a, const NDArray& b)
02383 {
02384   NDArray retval;
02385   dim_vector dv = x.dims ();
02386 
02387   if (dv == b.dims ())
02388     {
02389       octave_idx_type nel = dv.numel ();
02390 
02391       retval.resize (dv);
02392 
02393       for (octave_idx_type i = 0; i < nel; i++)
02394         retval (i) = betainc (x(i), a, b(i));
02395     }
02396   else
02397     gripe_betainc_nonconformant (dv, dim_vector (0, 0), b.dims ());
02398 
02399   return retval;
02400 }
02401 
02402 NDArray
02403 betainc (const NDArray& x, const NDArray& a, double b)
02404 {
02405   NDArray retval;
02406   dim_vector dv = x.dims ();
02407 
02408   if (dv == a.dims ())
02409     {
02410       octave_idx_type nel = dv.numel ();
02411 
02412       retval.resize (dv);
02413 
02414       for (octave_idx_type i = 0; i < nel; i++)
02415         retval (i) = betainc (x(i), a(i), b);
02416     }
02417   else
02418     gripe_betainc_nonconformant (dv, a.dims (), dim_vector (0, 0));
02419 
02420   return retval;
02421 }
02422 
02423 NDArray
02424 betainc (const NDArray& x, const NDArray& a, const NDArray& b)
02425 {
02426   NDArray retval;
02427   dim_vector dv = x.dims ();
02428 
02429   if (dv == a.dims () && dv == b.dims ())
02430     {
02431       octave_idx_type nel = dv.numel ();
02432 
02433       retval.resize (dv);
02434 
02435       for (octave_idx_type i = 0; i < nel; i++)
02436         retval (i) = betainc (x(i), a(i), b(i));
02437     }
02438   else
02439     gripe_betainc_nonconformant (dv, a.dims (), b.dims ());
02440 
02441   return retval;
02442 }
02443 
02444 float
02445 betainc (float x, float a, float b)
02446 {
02447   float retval;
02448   F77_XFCN (xbetai, XBETAI, (x, a, b, retval));
02449   return retval;
02450 }
02451 
02452 FloatMatrix
02453 betainc (float x, float a, const FloatMatrix& b)
02454 {
02455   octave_idx_type nr = b.rows ();
02456   octave_idx_type nc = b.cols ();
02457 
02458   FloatMatrix retval (nr, nc);
02459 
02460   for (octave_idx_type j = 0; j < nc; j++)
02461     for (octave_idx_type i = 0; i < nr; i++)
02462       retval(i,j) = betainc (x, a, b(i,j));
02463 
02464   return retval;
02465 }
02466 
02467 FloatMatrix
02468 betainc (float x, const FloatMatrix& a, float b)
02469 {
02470   octave_idx_type nr = a.rows ();
02471   octave_idx_type nc = a.cols ();
02472 
02473   FloatMatrix retval (nr, nc);
02474 
02475   for (octave_idx_type j = 0; j < nc; j++)
02476     for (octave_idx_type i = 0; i < nr; i++)
02477       retval(i,j) = betainc (x, a(i,j), b);
02478 
02479   return retval;
02480 }
02481 
02482 FloatMatrix
02483 betainc (float x, const FloatMatrix& a, const FloatMatrix& b)
02484 {
02485   FloatMatrix retval;
02486 
02487   octave_idx_type a_nr = a.rows ();
02488   octave_idx_type a_nc = a.cols ();
02489 
02490   octave_idx_type b_nr = b.rows ();
02491   octave_idx_type b_nc = b.cols ();
02492 
02493   if (a_nr == b_nr && a_nc == b_nc)
02494     {
02495       retval.resize (a_nr, a_nc);
02496 
02497       for (octave_idx_type j = 0; j < a_nc; j++)
02498         for (octave_idx_type i = 0; i < a_nr; i++)
02499           retval(i,j) = betainc (x, a(i,j), b(i,j));
02500     }
02501   else
02502     gripe_betainc_nonconformant (1, 1, a_nr, a_nc, b_nr, b_nc);
02503 
02504   return retval;
02505 }
02506 
02507 FloatNDArray
02508 betainc (float x, float a, const FloatNDArray& b)
02509 {
02510   dim_vector dv = b.dims ();
02511   octave_idx_type nel = dv.numel ();
02512 
02513   FloatNDArray retval (dv);
02514 
02515   for (octave_idx_type i = 0; i < nel; i++)
02516     retval (i) = betainc (x, a, b(i));
02517 
02518   return retval;
02519 }
02520 
02521 FloatNDArray
02522 betainc (float x, const FloatNDArray& a, float b)
02523 {
02524   dim_vector dv = a.dims ();
02525   octave_idx_type nel = dv.numel ();
02526 
02527   FloatNDArray retval (dv);
02528 
02529   for (octave_idx_type i = 0; i < nel; i++)
02530     retval (i) = betainc (x, a(i), b);
02531 
02532   return retval;
02533 }
02534 
02535 FloatNDArray
02536 betainc (float x, const FloatNDArray& a, const FloatNDArray& b)
02537 {
02538   FloatNDArray retval;
02539   dim_vector dv = a.dims ();
02540 
02541   if (dv == b.dims ())
02542     {
02543       octave_idx_type nel = dv.numel ();
02544 
02545       retval.resize (dv);
02546 
02547       for (octave_idx_type i = 0; i < nel; i++)
02548         retval (i) = betainc (x, a(i), b(i));
02549     }
02550   else
02551     gripe_betainc_nonconformant (dim_vector (0, 0), dv, b.dims ());
02552 
02553   return retval;
02554 }
02555 
02556 
02557 FloatMatrix
02558 betainc (const FloatMatrix& x, float a, float b)
02559 {
02560   octave_idx_type nr = x.rows ();
02561   octave_idx_type nc = x.cols ();
02562 
02563   FloatMatrix retval (nr, nc);
02564 
02565   for (octave_idx_type j = 0; j < nc; j++)
02566     for (octave_idx_type i = 0; i < nr; i++)
02567       retval(i,j) = betainc (x(i,j), a, b);
02568 
02569   return retval;
02570 }
02571 
02572 FloatMatrix
02573 betainc (const FloatMatrix& x, float a, const FloatMatrix& b)
02574 {
02575   FloatMatrix retval;
02576 
02577   octave_idx_type nr = x.rows ();
02578   octave_idx_type nc = x.cols ();
02579 
02580   octave_idx_type b_nr = b.rows ();
02581   octave_idx_type b_nc = b.cols ();
02582 
02583   if (nr == b_nr && nc == b_nc)
02584     {
02585       retval.resize (nr, nc);
02586 
02587       for (octave_idx_type j = 0; j < nc; j++)
02588         for (octave_idx_type i = 0; i < nr; i++)
02589           retval(i,j) = betainc (x(i,j), a, b(i,j));
02590     }
02591   else
02592     gripe_betainc_nonconformant (nr, nc, 1, 1, b_nr, b_nc);
02593 
02594   return retval;
02595 }
02596 
02597 FloatMatrix
02598 betainc (const FloatMatrix& x, const FloatMatrix& a, float b)
02599 {
02600   FloatMatrix retval;
02601 
02602   octave_idx_type nr = x.rows ();
02603   octave_idx_type nc = x.cols ();
02604 
02605   octave_idx_type a_nr = a.rows ();
02606   octave_idx_type a_nc = a.cols ();
02607 
02608   if (nr == a_nr && nc == a_nc)
02609     {
02610       retval.resize (nr, nc);
02611 
02612       for (octave_idx_type j = 0; j < nc; j++)
02613         for (octave_idx_type i = 0; i < nr; i++)
02614           retval(i,j) = betainc (x(i,j), a(i,j), b);
02615     }
02616   else
02617     gripe_betainc_nonconformant (nr, nc, a_nr, a_nc, 1, 1);
02618 
02619   return retval;
02620 }
02621 
02622 FloatMatrix
02623 betainc (const FloatMatrix& x, const FloatMatrix& a, const FloatMatrix& b)
02624 {
02625   FloatMatrix retval;
02626 
02627   octave_idx_type nr = x.rows ();
02628   octave_idx_type nc = x.cols ();
02629 
02630   octave_idx_type a_nr = a.rows ();
02631   octave_idx_type a_nc = a.cols ();
02632 
02633   octave_idx_type b_nr = b.rows ();
02634   octave_idx_type b_nc = b.cols ();
02635 
02636   if (nr == a_nr && nr == b_nr && nc == a_nc && nc == b_nc)
02637     {
02638       retval.resize (nr, nc);
02639 
02640       for (octave_idx_type j = 0; j < nc; j++)
02641         for (octave_idx_type i = 0; i < nr; i++)
02642           retval(i,j) = betainc (x(i,j), a(i,j), b(i,j));
02643     }
02644   else
02645     gripe_betainc_nonconformant (nr, nc, a_nr, a_nc, b_nr, b_nc);
02646 
02647   return retval;
02648 }
02649 
02650 FloatNDArray
02651 betainc (const FloatNDArray& x, float a, float b)
02652 {
02653   dim_vector dv = x.dims ();
02654   octave_idx_type nel = dv.numel ();
02655 
02656   FloatNDArray retval (dv);
02657 
02658   for (octave_idx_type i = 0; i < nel; i++)
02659     retval (i) = betainc (x(i), a, b);
02660 
02661   return retval;
02662 }
02663 
02664 FloatNDArray
02665 betainc (const FloatNDArray& x, float a, const FloatNDArray& b)
02666 {
02667   FloatNDArray retval;
02668   dim_vector dv = x.dims ();
02669 
02670   if (dv == b.dims ())
02671     {
02672       octave_idx_type nel = dv.numel ();
02673 
02674       retval.resize (dv);
02675 
02676       for (octave_idx_type i = 0; i < nel; i++)
02677         retval (i) = betainc (x(i), a, b(i));
02678     }
02679   else
02680     gripe_betainc_nonconformant (dv, dim_vector (0, 0), b.dims ());
02681 
02682   return retval;
02683 }
02684 
02685 FloatNDArray
02686 betainc (const FloatNDArray& x, const FloatNDArray& a, float b)
02687 {
02688   FloatNDArray retval;
02689   dim_vector dv = x.dims ();
02690 
02691   if (dv == a.dims ())
02692     {
02693       octave_idx_type nel = dv.numel ();
02694 
02695       retval.resize (dv);
02696 
02697       for (octave_idx_type i = 0; i < nel; i++)
02698         retval (i) = betainc (x(i), a(i), b);
02699     }
02700   else
02701     gripe_betainc_nonconformant (dv, a.dims (), dim_vector (0, 0));
02702 
02703   return retval;
02704 }
02705 
02706 FloatNDArray
02707 betainc (const FloatNDArray& x, const FloatNDArray& a, const FloatNDArray& b)
02708 {
02709   FloatNDArray retval;
02710   dim_vector dv = x.dims ();
02711 
02712   if (dv == a.dims () && dv == b.dims ())
02713     {
02714       octave_idx_type nel = dv.numel ();
02715 
02716       retval.resize (dv);
02717 
02718       for (octave_idx_type i = 0; i < nel; i++)
02719         retval (i) = betainc (x(i), a(i), b(i));
02720     }
02721   else
02722     gripe_betainc_nonconformant (dv, a.dims (), b.dims ());
02723 
02724   return retval;
02725 }
02726 
02727 // FIXME -- there is still room for improvement here...
02728 
02729 double
02730 gammainc (double x, double a, bool& err)
02731 {
02732   double retval;
02733 
02734   err = false;
02735 
02736   if (a < 0.0 || x < 0.0)
02737     {
02738       (*current_liboctave_error_handler)
02739         ("gammainc: A and X must be non-negative");
02740 
02741       err = true;
02742     }
02743   else
02744     F77_XFCN (xgammainc, XGAMMAINC, (a, x, retval));
02745 
02746   return retval;
02747 }
02748 
02749 Matrix
02750 gammainc (double x, const Matrix& a)
02751 {
02752   octave_idx_type nr = a.rows ();
02753   octave_idx_type nc = a.cols ();
02754 
02755   Matrix result (nr, nc);
02756   Matrix retval;
02757 
02758   bool err;
02759 
02760   for (octave_idx_type j = 0; j < nc; j++)
02761     for (octave_idx_type i = 0; i < nr; i++)
02762       {
02763         result(i,j) = gammainc (x, a(i,j), err);
02764 
02765         if (err)
02766           goto done;
02767       }
02768 
02769   retval = result;
02770 
02771  done:
02772 
02773   return retval;
02774 }
02775 
02776 Matrix
02777 gammainc (const Matrix& x, double a)
02778 {
02779   octave_idx_type nr = x.rows ();
02780   octave_idx_type nc = x.cols ();
02781 
02782   Matrix result (nr, nc);
02783   Matrix retval;
02784 
02785   bool err;
02786 
02787   for (octave_idx_type j = 0; j < nc; j++)
02788     for (octave_idx_type i = 0; i < nr; i++)
02789       {
02790         result(i,j) = gammainc (x(i,j), a, err);
02791 
02792         if (err)
02793           goto done;
02794       }
02795 
02796   retval = result;
02797 
02798  done:
02799 
02800   return retval;
02801 }
02802 
02803 Matrix
02804 gammainc (const Matrix& x, const Matrix& a)
02805 {
02806   Matrix result;
02807   Matrix retval;
02808 
02809   octave_idx_type nr = x.rows ();
02810   octave_idx_type nc = x.cols ();
02811 
02812   octave_idx_type a_nr = a.rows ();
02813   octave_idx_type a_nc = a.cols ();
02814 
02815   if (nr == a_nr && nc == a_nc)
02816     {
02817       result.resize (nr, nc);
02818 
02819       bool err;
02820 
02821       for (octave_idx_type j = 0; j < nc; j++)
02822         for (octave_idx_type i = 0; i < nr; i++)
02823           {
02824             result(i,j) = gammainc (x(i,j), a(i,j), err);
02825 
02826             if (err)
02827               goto done;
02828           }
02829 
02830       retval = result;
02831     }
02832   else
02833     (*current_liboctave_error_handler)
02834       ("gammainc: nonconformant arguments (arg 1 is %dx%d, arg 2 is %dx%d)",
02835        nr, nc, a_nr, a_nc);
02836 
02837  done:
02838 
02839   return retval;
02840 }
02841 
02842 NDArray
02843 gammainc (double x, const NDArray& a)
02844 {
02845   dim_vector dv = a.dims ();
02846   octave_idx_type nel = dv.numel ();
02847 
02848   NDArray retval;
02849   NDArray result (dv);
02850 
02851   bool err;
02852 
02853   for (octave_idx_type i = 0; i < nel; i++)
02854     {
02855       result (i) = gammainc (x, a(i), err);
02856 
02857       if (err)
02858         goto done;
02859     }
02860 
02861   retval = result;
02862 
02863  done:
02864 
02865   return retval;
02866 }
02867 
02868 NDArray
02869 gammainc (const NDArray& x, double a)
02870 {
02871   dim_vector dv = x.dims ();
02872   octave_idx_type nel = dv.numel ();
02873 
02874   NDArray retval;
02875   NDArray result (dv);
02876 
02877   bool err;
02878 
02879   for (octave_idx_type i = 0; i < nel; i++)
02880     {
02881       result (i) = gammainc (x(i), a, err);
02882 
02883       if (err)
02884         goto done;
02885     }
02886 
02887   retval = result;
02888 
02889  done:
02890 
02891   return retval;
02892 }
02893 
02894 NDArray
02895 gammainc (const NDArray& x, const NDArray& a)
02896 {
02897   dim_vector dv = x.dims ();
02898   octave_idx_type nel = dv.numel ();
02899 
02900   NDArray retval;
02901   NDArray result;
02902 
02903   if (dv == a.dims ())
02904     {
02905       result.resize (dv);
02906 
02907       bool err;
02908 
02909       for (octave_idx_type i = 0; i < nel; i++)
02910         {
02911           result (i) = gammainc (x(i), a(i), err);
02912 
02913           if (err)
02914             goto done;
02915         }
02916 
02917       retval = result;
02918     }
02919   else
02920     {
02921       std::string x_str = dv.str ();
02922       std::string a_str = a.dims ().str ();
02923 
02924       (*current_liboctave_error_handler)
02925         ("gammainc: nonconformant arguments (arg 1 is %s, arg 2 is %s)",
02926          x_str.c_str (), a_str. c_str ());
02927     }
02928 
02929  done:
02930 
02931   return retval;
02932 }
02933 
02934 float
02935 gammainc (float x, float a, bool& err)
02936 {
02937   float retval;
02938 
02939   err = false;
02940 
02941   if (a < 0.0 || x < 0.0)
02942     {
02943       (*current_liboctave_error_handler)
02944         ("gammainc: A and X must be non-negative");
02945 
02946       err = true;
02947     }
02948   else
02949     F77_XFCN (xsgammainc, XSGAMMAINC, (a, x, retval));
02950 
02951   return retval;
02952 }
02953 
02954 FloatMatrix
02955 gammainc (float x, const FloatMatrix& a)
02956 {
02957   octave_idx_type nr = a.rows ();
02958   octave_idx_type nc = a.cols ();
02959 
02960   FloatMatrix result (nr, nc);
02961   FloatMatrix retval;
02962 
02963   bool err;
02964 
02965   for (octave_idx_type j = 0; j < nc; j++)
02966     for (octave_idx_type i = 0; i < nr; i++)
02967       {
02968         result(i,j) = gammainc (x, a(i,j), err);
02969 
02970         if (err)
02971           goto done;
02972       }
02973 
02974   retval = result;
02975 
02976  done:
02977 
02978   return retval;
02979 }
02980 
02981 FloatMatrix
02982 gammainc (const FloatMatrix& x, float a)
02983 {
02984   octave_idx_type nr = x.rows ();
02985   octave_idx_type nc = x.cols ();
02986 
02987   FloatMatrix result (nr, nc);
02988   FloatMatrix retval;
02989 
02990   bool err;
02991 
02992   for (octave_idx_type j = 0; j < nc; j++)
02993     for (octave_idx_type i = 0; i < nr; i++)
02994       {
02995         result(i,j) = gammainc (x(i,j), a, err);
02996 
02997         if (err)
02998           goto done;
02999       }
03000 
03001   retval = result;
03002 
03003  done:
03004 
03005   return retval;
03006 }
03007 
03008 FloatMatrix
03009 gammainc (const FloatMatrix& x, const FloatMatrix& a)
03010 {
03011   FloatMatrix result;
03012   FloatMatrix retval;
03013 
03014   octave_idx_type nr = x.rows ();
03015   octave_idx_type nc = x.cols ();
03016 
03017   octave_idx_type a_nr = a.rows ();
03018   octave_idx_type a_nc = a.cols ();
03019 
03020   if (nr == a_nr && nc == a_nc)
03021     {
03022       result.resize (nr, nc);
03023 
03024       bool err;
03025 
03026       for (octave_idx_type j = 0; j < nc; j++)
03027         for (octave_idx_type i = 0; i < nr; i++)
03028           {
03029             result(i,j) = gammainc (x(i,j), a(i,j), err);
03030 
03031             if (err)
03032               goto done;
03033           }
03034 
03035       retval = result;
03036     }
03037   else
03038     (*current_liboctave_error_handler)
03039       ("gammainc: nonconformant arguments (arg 1 is %dx%d, arg 2 is %dx%d)",
03040        nr, nc, a_nr, a_nc);
03041 
03042  done:
03043 
03044   return retval;
03045 }
03046 
03047 FloatNDArray
03048 gammainc (float x, const FloatNDArray& a)
03049 {
03050   dim_vector dv = a.dims ();
03051   octave_idx_type nel = dv.numel ();
03052 
03053   FloatNDArray retval;
03054   FloatNDArray result (dv);
03055 
03056   bool err;
03057 
03058   for (octave_idx_type i = 0; i < nel; i++)
03059     {
03060       result (i) = gammainc (x, a(i), err);
03061 
03062       if (err)
03063         goto done;
03064     }
03065 
03066   retval = result;
03067 
03068  done:
03069 
03070   return retval;
03071 }
03072 
03073 FloatNDArray
03074 gammainc (const FloatNDArray& x, float a)
03075 {
03076   dim_vector dv = x.dims ();
03077   octave_idx_type nel = dv.numel ();
03078 
03079   FloatNDArray retval;
03080   FloatNDArray result (dv);
03081 
03082   bool err;
03083 
03084   for (octave_idx_type i = 0; i < nel; i++)
03085     {
03086       result (i) = gammainc (x(i), a, err);
03087 
03088       if (err)
03089         goto done;
03090     }
03091 
03092   retval = result;
03093 
03094  done:
03095 
03096   return retval;
03097 }
03098 
03099 FloatNDArray
03100 gammainc (const FloatNDArray& x, const FloatNDArray& a)
03101 {
03102   dim_vector dv = x.dims ();
03103   octave_idx_type nel = dv.numel ();
03104 
03105   FloatNDArray retval;
03106   FloatNDArray result;
03107 
03108   if (dv == a.dims ())
03109     {
03110       result.resize (dv);
03111 
03112       bool err;
03113 
03114       for (octave_idx_type i = 0; i < nel; i++)
03115         {
03116           result (i) = gammainc (x(i), a(i), err);
03117 
03118           if (err)
03119             goto done;
03120         }
03121 
03122       retval = result;
03123     }
03124   else
03125     {
03126       std::string x_str = dv.str ();
03127       std::string a_str = a.dims ().str ();
03128 
03129       (*current_liboctave_error_handler)
03130         ("gammainc: nonconformant arguments (arg 1 is %s, arg 2 is %s)",
03131          x_str.c_str (), a_str. c_str ());
03132     }
03133 
03134  done:
03135 
03136   return retval;
03137 }
03138 
03139 
03140 Complex rc_log1p (double x)
03141 {
03142   const double pi = 3.14159265358979323846;
03143   return x < -1.0 ? Complex (log (-(1.0 + x)), pi) : Complex (log1p (x));
03144 }
03145 
03146 FloatComplex rc_log1p (float x)
03147 {
03148   const float pi = 3.14159265358979323846f;
03149   return x < -1.0f ? FloatComplex (logf (-(1.0f + x)), pi) : FloatComplex (log1pf (x));
03150 }
03151 
03152 // This algorithm is due to P. J. Acklam.
03153 // See http://home.online.no/~pjacklam/notes/invnorm/
03154 // The rational approximation has relative accuracy 1.15e-9 in the whole region.
03155 // For doubles, it is refined by a single step of Higham's 3rd order method.
03156 // For single precision, the accuracy is already OK, so we skip it to get
03157 // faster evaluation.
03158 
03159 static double do_erfinv (double x, bool refine)
03160 {
03161   // Coefficients of rational approximation.
03162   static const double a[] =
03163     { -2.806989788730439e+01,  1.562324844726888e+02,
03164       -1.951109208597547e+02,  9.783370457507161e+01,
03165       -2.168328665628878e+01,  1.772453852905383e+00 };
03166   static const double b[] =
03167     { -5.447609879822406e+01,  1.615858368580409e+02,
03168       -1.556989798598866e+02,  6.680131188771972e+01,
03169       -1.328068155288572e+01 };
03170   static const double c[] =
03171     { -5.504751339936943e-03, -2.279687217114118e-01,
03172       -1.697592457770869e+00, -1.802933168781950e+00,
03173        3.093354679843505e+00,  2.077595676404383e+00 };
03174   static const double d[] =
03175     {  7.784695709041462e-03,  3.224671290700398e-01,
03176        2.445134137142996e+00,  3.754408661907416e+00 };
03177 
03178   static const double spi2 =  8.862269254527579e-01; // sqrt(pi)/2.
03179   static const double pbreak = 0.95150;
03180   double ax = fabs (x), y;
03181 
03182   // Select case.
03183   if (ax <= pbreak)
03184     {
03185       // Middle region.
03186       const double q = 0.5 * x, r = q*q;
03187       const double yn = (((((a[0]*r + a[1])*r + a[2])*r + a[3])*r + a[4])*r + a[5])*q;
03188       const double yd = ((((b[0]*r + b[1])*r + b[2])*r + b[3])*r + b[4])*r + 1.0;
03189       y = yn / yd;
03190     }
03191   else if (ax < 1.0)
03192     {
03193       // Tail region.
03194       const double q = sqrt (-2*log (0.5*(1-ax)));
03195       const double yn = ((((c[0]*q + c[1])*q + c[2])*q + c[3])*q + c[4])*q + c[5];
03196       const double yd = (((d[0]*q + d[1])*q + d[2])*q + d[3])*q + 1.0;
03197       y = yn / yd * signum (-x);
03198     }
03199   else if (ax == 1.0)
03200     return octave_Inf * signum (x);
03201   else
03202     return octave_NaN;
03203 
03204   if (refine)
03205     {
03206       // One iteration of Halley's method gives full precision.
03207       double u = (erf(y) - x) * spi2 * exp (y*y);
03208       y -= u / (1 + y*u);
03209     }
03210 
03211   return y;
03212 }
03213 
03214 double erfinv (double x)
03215 {
03216   return do_erfinv (x, true);
03217 }
03218 
03219 float erfinv (float x)
03220 {
03221   return do_erfinv (x, false);
03222 }
03223 
03224 // Implementation based on the Fortran code by W.J.Cody
03225 // see http://www.netlib.org/specfun/erf.
03226 // Templatized and simplified workflow.
03227 
03228 // FIXME: Maybe this should be globally visible.
03229 static inline float erfc (float x) { return erfcf (x); }
03230 
03231 template <class T>
03232 static T
03233 erfcx_impl (T x)
03234 {
03235   static const T c[] =
03236     {
03237       5.64188496988670089e-1,8.88314979438837594,
03238       6.61191906371416295e+1,2.98635138197400131e+2,
03239       8.81952221241769090e+2,1.71204761263407058e+3,
03240       2.05107837782607147e+3,1.23033935479799725e+3,
03241       2.15311535474403846e-8
03242     };
03243 
03244   static const T d[] =
03245     {
03246       1.57449261107098347e+1,1.17693950891312499e+2,
03247       5.37181101862009858e+2,1.62138957456669019e+3,
03248       3.29079923573345963e+3,4.36261909014324716e+3,
03249       3.43936767414372164e+3,1.23033935480374942e+3
03250     };
03251 
03252   static const T p[] =
03253     {
03254       3.05326634961232344e-1,3.60344899949804439e-1,
03255       1.25781726111229246e-1,1.60837851487422766e-2,
03256       6.58749161529837803e-4,1.63153871373020978e-2
03257     };
03258 
03259   static const T q[] =
03260     {
03261       2.56852019228982242,1.87295284992346047,
03262       5.27905102951428412e-1,6.05183413124413191e-2,
03263       2.33520497626869185e-3
03264     };
03265 
03266   static const T sqrpi = 5.6418958354775628695e-1;
03267   static const T xhuge = sqrt (1.0 / std::numeric_limits<T>::epsilon ());
03268   static const T xneg = -sqrt (log (std::numeric_limits<T>::max ()/2.0));
03269 
03270   double y = fabs (x), result;
03271   if (x < xneg)
03272     result = octave_Inf;
03273   else if (y <= 0.46875)
03274     result = std::exp (x*x) * erfc (x);
03275   else
03276     {
03277       if (y <= 4.0)
03278         {
03279           double xnum = c[8]*y, xden = y;
03280           for (int i = 0; i < 7; i++)
03281             {
03282               xnum = (xnum + c[i]) * y;
03283               xden = (xden + d[i]) * y;
03284             }
03285 
03286           result = (xnum + c[7]) / (xden + d[7]);
03287         }
03288       else if (y <= xhuge)
03289         {
03290           double y2 = 1/(y*y), xnum = p[5]*y2, xden = y2;
03291           for (int i = 0; i < 4; i++)
03292             {
03293               xnum = (xnum + p[i]) * y2;
03294               xden = (xden + q[i]) * y2;
03295             }
03296 
03297           result = y2 * (xnum + p[4]) / (xden + q[4]);
03298           result = (sqrpi - result) / y;
03299         }
03300       else
03301         result = sqrpi / y;
03302 
03303       // Fix up negative argument.
03304       if (x < 0)
03305         {
03306           double y2 = ceil (x / 16.0) * 16.0, del = (x-y2)*(x+y2);
03307           result = 2*(std::exp(y2*y2) * std::exp(del)) - result;
03308         }
03309     }
03310 
03311   return result;
03312 }
03313 
03314 double erfcx (double x)
03315 {
03316   return erfcx_impl (x);
03317 }
03318 
03319 float erfcx (float x)
03320 {
03321   return erfcx_impl (x);
03322 }
 All Classes Files Functions Variables Typedefs Enumerations Enumerator Friends Defines