SpecialFunctions.cpp

Go to the documentation of this file.
00001 #pragma ident "$Id: SpecialFunctions.cpp 3140 2012-06-18 15:03:02Z susancummins $"
00002 
00009 //============================================================================
00010 //
00011 //  This file is part of GPSTk, the GPS Toolkit.
00012 //
00013 //  The GPSTk is free software; you can redistribute it and/or modify
00014 //  it under the terms of the GNU Lesser General Public License as published
00015 //  by the Free Software Foundation; either version 2.1 of the License, or
00016 //  any later version.
00017 //
00018 //  The GPSTk is distributed in the hope that it will be useful,
00019 //  but WITHOUT ANY WARRANTY; without even the implied warranty of
00020 //  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00021 //  GNU Lesser General Public License for more details.
00022 //
00023 //  You should have received a copy of the GNU Lesser General Public
00024 //  License along with GPSTk; if not, write to the Free Software Foundation,
00025 //  Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110, USA
00026 //
00027 //  Dagoberto Salazar - gAGE ( http://www.gage.es ). 2008
00028 //
00029 //============================================================================
00030 
00031 
00032 #include "SpecialFunctions.hpp"
00033 
00034 
00035 using namespace std;
00036 
00037 namespace gpstk
00038 {
00039 
00040 
00041       /* Computes the Gamma function using a simple Lanczos approximation.
00042        *
00043        * This implementation typically gives 15 correct decimal places, and
00044        * it is adapted from free Python code found in:
00045        *
00046        * http://en.wikipedia.org/wiki/Lanczos_approximation
00047        *
00048        * \warning Be aware that Gamma function is not defined for 0, -1, -2,...
00049        */
00050    double gamma(double val)
00051    {
00052 
00053       double inf( 9.0e+99 );
00054 
00055       if( val == 0.0 )
00056       {
00057          return inf;
00058       }
00059 
00060       if ( ( val < 0.0 ) &&
00061            ( floor(val) == val ) )
00062       {
00063          return inf;
00064       }
00065 
00066          // Set the number of coefficients being used
00067       int g(7);
00068       const double lanczos_coef[] = { 0.99999999999980993,
00069                                       676.5203681218851,
00070                                       -1259.1392167224028,
00071                                       771.32342877765313,
00072                                       -176.61502916214059,
00073                                       12.507343278686905,
00074                                       -0.13857109526572012,
00075                                       9.9843695780195716e-6,
00076                                       1.5056327351493116e-7 };
00077 
00078       if(val < 0.5)
00079       {
00080          return ( PI / (sin(PI*val)*gamma(1.0-val)) );
00081       }
00082       else
00083       {
00084          val -= 1.0;
00085 
00086          double x( lanczos_coef[0] );
00087 
00088          for(int i = 1; i<g+2; i++)
00089          {
00090             x += lanczos_coef[i]/(val+(double)i);
00091          }
00092 
00093          double t (val + static_cast<double>(g) + 0.5);
00094 
00095          return ( 2.5066282746310002 * pow( t, (val+0.5) ) * exp(-t) * x );
00096 
00097       }
00098 
00099    }  // End of function 'gamma()'
00100 
00101 
00102 
00103       /* Computes the natural logarithm of Gamma function
00104        * using the Lanczos approximation.
00105        *
00106        * \warning This version does not work for values <= 0.0
00107        */
00108    double lngamma(double val)
00109    {
00110 
00111       double inf( 9.0e+99 );
00112 
00113       if( val <= 0.0 )
00114       {
00115          return inf;
00116       }
00117 
00118          // Set the number of coefficients being used
00119       int g(7);
00120       const double lanczos_coef[] = { 0.99999999999980993,
00121                                       676.5203681218851,
00122                                       -1259.1392167224028,
00123                                       771.32342877765313,
00124                                       -176.61502916214059,
00125                                       12.507343278686905,
00126                                       -0.13857109526572012,
00127                                       9.9843695780195716e-6,
00128                                       1.5056327351493116e-7 };
00129 
00130       if(val < 0.5)
00131       {
00132          return ( 1.1447298858494002 - (log(sin(PI*val)) + lngamma(1.0-val)) );
00133       }
00134       else
00135       {
00136          val -= 1.0;
00137 
00138          double x( lanczos_coef[0] );
00139 
00140          for(int i = 1; i<g+2; i++)
00141          {
00142             x += lanczos_coef[i]/(val+(double)i);
00143          }
00144 
00145          double t (val + static_cast<double>(g) + 0.5);
00146 
00147          return ( 0.918938533204672741781 + (val+0.5)*log(t) + (-t) + log(x) );
00148       }
00149 
00150 
00151    }  // End of function 'lngamma()'
00152 
00153 
00154 
00155       // Auxiliar Kummer function.
00156    double kummerFunc(const double a, const double z);
00157 
00158 
00159 
00160       // We compute the lower incomplete gamma function g(a,z) using the
00161       // formula:
00162       //
00163       //        g(a,z) = z**a * exp(-z) * S(a,z) / a
00164       //
00165       // where:
00166       //
00167       //                  oo
00168       //                 ___            k
00169       //                \              z
00170       //   S(a,z) = 1 +  )     ------------------.
00171       //                /___   (a+1)(a+2)...(a+k)
00172       //                k = 1
00173       //
00174       // S(a,z) is computed with the Kummer function "kummerFunc()".
00175    double kummerFunc(const double aval, const double zval)
00176    {
00177 
00178       double eps(1.0e-15);    // Small threshold controling precision
00179 
00180       double z( abs(zval) );    // We only allow positive values of 'z'
00181       double a( abs(aval) );    // We only allow positive values of 'a'
00182 
00183       double den(a);          // Variable to store denominator
00184 
00185       double s(1.0);          // Result will be stored here
00186 
00187       double coef(1.0);       // Initialize coefficient
00188 
00189       while( abs(coef) > eps)
00190       {
00191          coef = coef*z;   // Compute numerator
00192          den += 1.0;
00193          coef = coef/den;  // Compute coefficient
00194          s += coef;        // Add new coefficient to result
00195       }
00196 
00197       return s;
00198 
00199    }  // End of function 'kummerFunc()'
00200 
00201 
00202 
00203       // Lower incomplete gamma function.
00204    double lower_gamma(const double a, const double z)
00205    {
00206 
00207       double zp( abs(z) );    // We only allow positive values of 'z'
00208       double ap( abs(a) );    // We only allow positive values of 'a'
00209 
00210       double s( kummerFunc(ap, zp) );
00211 
00212       return exp(log(zp)*ap) * exp(-zp) * s / ap;
00213 
00214    }  // End of function 'lower_gamma()'
00215 
00216 
00217 
00218       // Upper incomplete gamma function.
00219    double upper_gamma(const double a, const double z)
00220    {
00221 
00222       return ( gamma(a) - lower_gamma(a, z) );
00223 
00224    }  // End of function 'upper_gamma()'
00225 
00226 
00227 
00228       // Lower incomplete regularized gamma function P(a,z).
00229    double gammaP(const double a, const double z)
00230    {
00231       return ( lower_gamma(a,z) / gamma(a) );
00232    }
00233 
00234 
00235 
00236       // Upper incomplete regularized gamma function Q(a,z).
00237    double gammaQ(const double a, const double z)
00238    {
00239       return ( 1.0 - gammaP(a,z) );
00240    }
00241 
00242 
00243 
00244       /* Computes factorial of integer number n.
00245        *
00246        * This implementation typically gives 15 correct decimal places, and
00247        * returns the result as double.
00248        */
00249    double factorial(const int n)
00250    {
00251 
00252          // Return 0 if n<0
00253       if( n < 0 ) return 0.0;
00254 
00255       double facttable[] = { 1.0, 1.0, 2.0, 6.0, 24.0, 120.0, 720.0, 5040.0,
00256                              40320.0, 362880.0, 3628800.0, 39916800.0,
00257                              479001600.0, 6227020800.0, 87178291200.0,
00258                              1307674368000.0 };
00259 
00260          // Use a look-up table for small values of n
00261       if( (n >= 0) && (n <= 15) )
00262       {
00263          return facttable[n];
00264       }
00265 
00266          // Use gamma function properties for big values of n
00267       return gamma(n+1.0);
00268 
00269    }
00270 
00271             /* Computes factorial of double number d.
00272        */
00273    double factorial(const double d)
00274    {
00275         // copy the input
00276       double n(d);
00277 
00278       double returnValue(0.0);
00279 
00280       if (n < 0.0) 
00281       {
00282          return returnValue;
00283 
00284          // Exception e("Illegal value: ");
00285          // GPSTK_THROW(e);
00286       } 
00287       else if ((n == 0.0)||(n == 1.0)) 
00288       {
00289          returnValue = 1.0;
00290 
00291       } else    
00292       {
00293          returnValue = n * factorial(n - 1.0);
00294       }
00295 
00296       return returnValue;
00297    }
00298 
00299 
00300 
00301       // Auxiliar error function #1. erf(x) for x in [ 0, 0.84375 ]
00302    double erf1(const double x);
00303 
00304 
00305       // Auxiliar error function #2. erf(x) for x in [ 0.84375, 1.25 ]
00306    double erf2(const double x);
00307 
00308 
00309       // Auxiliar error function #3. erf(x) for x in [ 1.25, 2.857142 ]
00310    double erf3(const double x);
00311 
00312 
00313       // Auxiliar error function #4. erf(x) for x in [ 2.857142, 6.0 ]
00314    double erf4(const double x);
00315 
00316 
00317       // Auxiliar error function #5. erf(x) for x in [ 6.0, inf ]
00318    double erf5(const double x);
00319 
00320 
00321 
00322       /* Error function.
00323        *
00324        * This is a C++ implementation of the free Python code found in:
00325        *
00326        *   http://code.activestate.com/recipes/576391/
00327        *
00328        * Such code was based in a C code base with OpenBSD license from:
00329        *
00330        * ====================================================
00331        * Copyright (C) 1993 by Sun Microsystems, Inc. All rights reserved.
00332        *
00333        * Developed at SunPro, a Sun Microsystems, Inc. business.
00334        * Permission to use, copy, modify, and distribute this
00335        * software is freely granted, provided that this notice
00336        * is preserved.
00337        * ====================================================
00338        */
00339    double erf(const double x)
00340    {
00341 
00342       /*
00343        *                            |x
00344        *                    2       |\
00345        *     erf(x)  =  ---------   |   exp(-t*t)dt
00346        *                 sqrt(pi)  \|
00347        *                            |0
00348        *
00349        *     erfc(x) =  1-erf(x)
00350        *  Note that
00351        *              erf(-x) = -erf(x)
00352        *             erfc(-x) = 2 - erfc(x)
00353        *
00354        * Method:
00355        *      1. For |x| in [0, 0.84375]
00356        *          erf(x)  = x + x*R(x^2)
00357        *          erfc(x) = 1 - erf(x)           if x in [-.84375,0.25]
00358        *                  = 0.5 + ((0.5-x)-x*R)  if x in [0.25,0.84375]
00359        *
00360        *         where R = P/Q where P is an odd poly of degree 8 and
00361        *         Q is an odd poly of degree 10.
00362        *
00363        *                      | R - (erf(x)-x)/x | <= 2^(-57.90)
00364        *
00365        *
00366        *         Remark. The formula is derived by noting
00367        *            erf(x) = (2/sqrt(pi))*(x - x^3/3 + x^5/10 - x^7/42 + ....)
00368        *         and that
00369        *            2/sqrt(pi) = 1.128379167095512573896158903121545171688
00370        *         is close to one. The interval is chosen because the fix
00371        *         point of erf(x) is near 0.6174 (i.e., erf(x)=x when x is
00372        *         near 0.6174), and by some experiment, 0.84375 is chosen to
00373        *         guarantee the error is less than one ulp for erf.
00374        *
00375        *      2. For |x| in [0.84375,1.25], let s = |x| - 1, and
00376        *         c = 0.84506291151 rounded to single (24 bits)
00377        *              erf(x)  = sign(x) * (c  + P1(s)/Q1(s))
00378        *              erfc(x) = (1-c)  - P1(s)/Q1(s) if x > 0
00379        *                        1+(c+P1(s)/Q1(s))    if x < 0
00380        *              |P1/Q1 - (erf(|x|)-c)| <= 2^(-59.06)
00381        *
00382        *         Remark: here we use the taylor series expansion at x=1.
00383        *              erf(1+s) = erf(1) + s*Poly(s)
00384        *                       = 0.845.. + P1(s)/Q1(s)
00385        *         That is, we use rational approximation to approximate
00386        *              erf(1+s) - (c = (single)0.84506291151)
00387        *         Note that |P1/Q1|< 0.078 for x in [0.84375,1.25]
00388        *         where
00389        *              P1(s) = degree 6 poly in s
00390        *              Q1(s) = degree 6 poly in s
00391        *
00392        *      3. For x in [1.25,1/0.35(~2.857143)],
00393        *              erfc(x) = (1/x)*exp(-x*x-0.5625+R1/S1)
00394        *              erf(x)  = 1 - erfc(x)
00395        *         where
00396        *              R1(z) = degree 7 poly in z, (z=1/x^2)
00397        *              S1(z) = degree 8 poly in z
00398        *
00399        *      4. For x in [1/0.35,28]
00400        *              erfc(x) = (1/x)*exp(-x*x-0.5625+R2/S2) if x > 0
00401        *                      = 2.0 - (1/x)*exp(-x*x-0.5625+R2/S2) if -6<x<0
00402        *                      = 2.0 - tiny      (if x <= -6)
00403        *              erf(x)  = sign(x)*(1.0 - erfc(x)) if x < 6, else
00404        *              erf(x)  = sign(x)*(1.0 - tiny)
00405        *         where
00406        *              R2(z) = degree 6 poly in z, (z=1/x^2)
00407        *              S2(z) = degree 7 poly in z
00408        *
00409        *         Note1:
00410        *            To compute exp(-x*x-0.5625+R/S), let s be a single
00411        *            precision number and s := x; then
00412        *              -x*x = -s*s + (s-x)*(s+x)
00413        *              exp(-x*x-0.5626+R/S) =
00414        *                               exp(-s*s-0.5625)*exp((s-x)*(s+x)+R/S);
00415        *
00416        *         Note2:
00417        *            Here 4 and 5 make use of the asymptotic series
00418        *                         exp(-x*x)
00419        *              erfc(x) ~ ---------- * ( 1 + Poly(1/x^2) )
00420        *                        x*sqrt(pi)
00421        *            We use rational approximation to approximate
00422        *              g(s)=f(1/x^2) = log(erfc(x)*x) - x*x + 0.5625
00423        *            Here is the error bound for R1/S1 and R2/S2
00424        *              |R1/S1 - f(x)|  < 2**(-62.57)
00425        *              |R2/S2 - f(x)|  < 2**(-61.52)
00426        *
00427        *      5. For inf > x >= 28
00428        *              erf(x)  = sign(x) *(1 - tiny)  (raise inexact)
00429        *              erfc(x) = tiny*tiny (raise underflow) if x > 0
00430        *                      = 2 - tiny if x<0
00431        *
00432        *      7. Special case:
00433        *              erf(0)  = 0, erf(inf)  = 1, erf(-inf) = -1,
00434        *              erfc(0) = 1, erfc(inf) = 0, erfc(-inf) = 2,
00435        *              erfc/erf(NaN) is NaN
00436        */
00437 
00438       double inf( 9e+99 );
00439 
00440        if (x >= inf)
00441       {
00442          return 1.0;
00443       }
00444 
00445       if (x <= -inf)
00446       {
00447          return -1.0;
00448       }
00449 
00450       if ( fabs(x) < 0.84375 )
00451       {
00452          return erf1(x);
00453       }
00454       else if ( ( 0.84375 <= fabs(x) ) &&
00455                 ( fabs(x) < 1.25 ) )
00456       {
00457          return erf2(x);
00458       }
00459       else if ( ( 1.25 <= fabs(x) ) &&
00460                 ( fabs(x) < 2.857142) )
00461       {
00462          return erf3(x);
00463       }
00464       else if ( ( 2.857142 <= fabs(x) ) &&
00465                 ( fabs(x) < 6.0 ) )
00466       {
00467          return erf4(x);
00468       }
00469 
00470          // If we got here, it means that ( fabs(x) >= 6.0 )
00471       return erf5(x);
00472 
00473    }  // End of function 'erf()'
00474 
00475 
00476 
00477       // Auxiliar error function #1. erf(x) for x in [0,0.84375].
00478    double erf1(const double x)
00479    {
00480 
00481       int i;
00482 
00483          // Get the base-2 exponent of input 'x'
00484       frexp (x , &i);
00485 
00486       if ( fabs( static_cast<double>(i) ) > 28.0 )
00487       {
00488          if( fabs( static_cast<double>(i) ) > 57.0 )
00489          {
00490             return ( 1.128379167095512586316e+0 * x );
00491          }
00492 
00493          return ( x * ( 2.28379167095512586316e-01 ) );
00494       }
00495 
00496       double z( x*x );
00497       double r( 1.28379167095512558561e-01
00498                + z * ( -3.25042107247001499370e-01
00499                + z * ( -2.84817495755985104766e-02
00500                + z * ( -5.77027029648944159157e-03
00501                + z * -2.37630166566501626084e-05 ) ) ) );
00502 
00503 
00504       double s( 1.0 + z * ( 3.97917223959155352819e-01
00505                     + z * ( 6.50222499887672944485e-02
00506                     + z * ( 5.08130628187576562776e-03
00507                     + z * ( 1.32494738004321644526e-04
00508                     + z * -3.96022827877536812320e-06 ) ) ) ) );
00509 
00510       return ( x * (1.0 + r/s ) );
00511 
00512    }  // End of function 'erf1()'
00513 
00514 
00515 
00516       // Auxiliar error function #2. erf(x) for x in [0.84375,1.25]
00517    double erf2(const double x)
00518    {
00519 
00520       double s( fabs(x) - 1.0 );
00521 
00522       double P( -2.36211856075265944077e-03
00523                + s * ( 4.14856118683748331666e-01
00524                + s * ( -3.72207876035701323847e-01
00525                + s * ( 3.18346619901161753674e-01
00526                + s * ( -1.10894694282396677476e-01
00527                + s * ( 3.54783043256182359371e-02
00528                + s * -2.16637559486879084300e-03 ) ) ) ) ) );
00529 
00530       double Q( 1.0 + s * ( 1.06420880400844228286e-01
00531                     + s * ( 5.40397917702171048937e-01
00532                     + s * ( 7.18286544141962662868e-02
00533                     + s * ( 1.26171219808761642112e-01
00534                     + s * ( 1.36370839120290507362e-02
00535                     + s * 1.19844998467991074170e-02 ) ) ) ) ) );
00536 
00537       if( x >= 0.0 )
00538       {
00539          return ( 8.45062911510467529297e-01 + P/Q );
00540       }
00541       else
00542       {
00543          return ( -8.45062911510467529297e-01 - P/Q );
00544       }
00545 
00546    }  // End of function 'erf2()'
00547 
00548 
00549 
00550       // Auxiliar error function #3. erf(x) for x in [1.25,2.857142]
00551    double erf3(const double xval)
00552    {
00553 
00554       double x0(xval);
00555       double x( fabs(xval) );
00556       double s( 1.0/(x*x) );
00557 
00558       double R( -9.86494403484714822705e-03
00559                + s * ( -6.93858572707181764372e-01
00560                + s * ( -1.05586262253232909814e+01
00561                + s * ( -6.23753324503260060396e+01
00562                + s * ( -1.62396669462573470355e+02
00563                + s * ( -1.84605092906711035994e+02
00564                + s * ( -8.12874355063065934246e+01
00565                + s * -9.81432934416914548592e+00 ) ) ) ) ) ) );
00566 
00567       double S( 1.0 + s * ( 1.96512716674392571292e+01
00568                     + s * ( 1.37657754143519042600e+02
00569                     + s * ( 4.34565877475229228821e+02
00570                     + s * ( 6.45387271733267880336e+02
00571                     + s * ( 4.29008140027567833386e+02
00572                     + s * ( 1.08635005541779435134e+02
00573                     + s * ( 6.57024977031928170135e+00
00574                     + s * -6.04244152148580987438e-02 ) ) ) ) ) ) ) );
00575 
00576       double r( exp(-x0*x0-0.5625) * exp( (x0-x)*(x0+x)+R/S) );
00577 
00578       if( x0 >= 0.0 )
00579       {
00580          return ( 1.0 - r/x );
00581       }
00582       else
00583       {
00584          return ( r/x - 1.0 );
00585       }
00586 
00587    }  // End of function 'erf3()'
00588 
00589 
00590 
00591       // Auxiliar error function #4. erf(x) for x in [ 2.857142, 6.0 ]
00592    double erf4(const double xval)
00593    {
00594 
00595       double x0(xval);
00596       double x( fabs(xval) );
00597       double s( 1.0/(x*x) );
00598 
00599       double R( -9.86494292470009928597e-03
00600                + s * ( -7.99283237680523006574e-01
00601                + s * ( -1.77579549177547519889e+01
00602                + s * ( -1.60636384855821916062e+02
00603                + s * ( -6.37566443368389627722e+02
00604                + s * ( -1.02509513161107724954e+03
00605                + s * -4.83519191608651397019e+02 ) ) ) ) ) );
00606 
00607       double S( 1.0 + s * ( 3.03380607434824582924e+01
00608                     + s * ( 3.25792512996573918826e+02
00609                     + s * ( 1.53672958608443695994e+03
00610                     + s * ( 3.19985821950859553908e+03
00611                     + s * ( 2.55305040643316442583e+03
00612                     + s * ( 4.74528541206955367215e+02
00613                     + s * -2.24409524465858183362e+01 ) ) ) ) ) ) );
00614 
00615       double r( exp( -x0 * x0 - 0.5625 ) * exp( (x0-x)*(x0+x) + R/S ) );
00616 
00617       if( x0 >= 0.0 )
00618       {
00619          return ( 1.0 - r/x );
00620       }
00621       else
00622       {
00623          return ( r/x - 1.0 );
00624       }
00625 
00626    }  // End of function 'erf4()'
00627 
00628 
00629 
00630       // Auxiliar error function #5. erf(x) for x in [ 6.0, inf ]
00631    double erf5(const double x)
00632    {
00633 
00634       double tiny( 1e-99 );
00635 
00636       if ( x > 0.0 )
00637       {
00638          return ( 1.0 - tiny );
00639       }
00640       else
00641       {
00642          return ( tiny - 1.0 );
00643       }
00644 
00645    }  // End of function 'erf5()'
00646 
00647 
00648 
00649       // Complementary error function.
00650    double erfc(const double x)
00651    {
00652 
00653       return ( 1.0 - erf(x) );
00654 
00655    }  // End of function 'erfc()'
00656 
00657 
00658 
00659       /* Inverse of error function.
00660        *
00661        * \ warning Value "z" must be in the range (-1, 1)
00662        */
00663    double inverf(const double z)
00664    {
00665 
00666       double inf( 9.0e+99 );
00667 
00668          // Check limits
00669       if( z >= 1.0 ) return inf;
00670       if( z <= -1.0 ) return -inf;
00671 
00672       double z2PI( z*z*PI );
00673 
00674          // Compute the approximation found in:
00675          //    http://en.wikipedia.org/wiki/Error_function
00676          // The module of this approximation is always under fabs(inverf())
00677       double x( 0.88622692545275794 * z
00678               * ( 1.0 + z2PI * ( 0.083333333333333333
00679                       + z2PI * ( 0.014583333333333334
00680                       + z2PI * ( 0.0031498015873015874
00681                       + z2PI * ( 0.00075248704805996468
00682                       + z2PI * ( 0.0001907475361251403
00683                       + z2PI * ( 1.8780048076923078e-5
00684                       + z2PI * ( 1.3623642420578133e-5 ) ) ) ) ) ) ) ) );
00685 
00686       double delta(1.0);
00687       double threshold( 1.0e-10 );
00688       int iter( 0 );
00689 
00690       while ( ( fabs(delta) > threshold ) &&
00691               ( iter < 100 ) )
00692       {
00693 
00694             // Use the Newton-Raphson method. Denominator is d(erf(x))/dx
00695          delta = (z - gpstk::erf(x)) / (1.1283791670955126*std::exp(-x*x));
00696 
00697          x = x + delta;
00698 
00699          ++iter;
00700 
00701       }
00702 
00703       return x;
00704 
00705    }  // End of function 'inverf()'
00706 
00707 
00708 
00709       /* Beta function.
00710        *
00711        * \warning This version may not work for values > 130.0
00712        */
00713    double beta(const double x, const double y)
00714    {
00715 
00716       return ( gamma(x) * gamma(y) / gamma( x + y ) );
00717 
00718    }  // End of function 'beta()'
00719 
00720 
00721 
00722       /* Computes the natural logarithm of Beta function
00723        *
00724        * \warning This version does not work for values <= 0.0
00725        */
00726    double lnbeta(double x, double y)
00727    {
00728 
00729       x = fabs(x);
00730       y = fabs(y);
00731 
00732       return ( lngamma(x) + lngamma(y) - lngamma( x + y ) );
00733 
00734    }  // End of function 'lbeta()'
00735 
00736 
00737 
00738       /* Auxiliar function to compute incomplete beta function.
00739        * Power series for incomplete beta integral.
00740        * Use when b*x is small and x not too close to 1.
00741        */
00742    double incompletebetaps(const double x, const double a, const double b)
00743    {
00744 
00745          // Declare working variables
00746       double epsilon(1.0e-30);
00747       double big(1.0e99);
00748       double small(1.0e-99);
00749       double maxgam(171.0);   // Approximate maximum input for gamma function
00750       double ai( 1.0/a );
00751       double u( (1.0-b) * x );
00752       double v( u / (a+1.0) );
00753       double t1(v);
00754       double t(u);
00755       double n(2.0);
00756       double s(0.0);
00757       double z( epsilon*ai );
00758 
00759       while( std::fabs(v) > z )
00760       {
00761          u = (n-b) * x / n;
00762          t = t * u;
00763          v = t/(a+n);
00764          s = s+v;
00765          n = n+1.0;
00766       }
00767 
00768       s = s + t1 + ai;
00769       u = a * std::log(x);
00770 
00771          // Check if we can compute gamma values
00772       if( (a+b < maxgam) &&
00773           ( std::fabs(u) < std::log(big) ) )
00774       {
00775 
00776          t = gamma(a+b) / ( gamma(a) * gamma(b) );
00777          s = s * t * pow(x, a);
00778 
00779       }
00780       else
00781       {
00782 
00783          t = lngamma(a+b) - lngamma(a)- lngamma(b) + u + std::log(s);
00784 
00785          if( t < std::log(small) )
00786          {
00787             s = 0.0;
00788          }
00789          else
00790          {
00791             s = std::exp(t);
00792          }
00793 
00794       }  // End of block 'if( (a+b < maxgam) && ( std::fabs(u) ...'
00795 
00796       return s;
00797 
00798    }  // End of function 'incompletebetaps()'
00799 
00800 
00801 
00802       // Auxiliar function to compute incomplete beta function.
00803    double incompletebetafe(const double x, const double a, const double b)
00804    {
00805 
00806          // Declare auxiliar parameters
00807       double epsilon(1.0e-30);
00808       double big( 1.0e16 );
00809       double biginv( 1.0/big );
00810 
00811       double k1( a );
00812       double k2( a+b );
00813       double k3( a );
00814       double k4( a+1.0 );
00815       double k5( 1.0 );
00816       double k6( b-1.0 );
00817       double k7( k4 );
00818       double k8( a+2.0 );
00819 
00820       double pkm1( 1.0 );
00821       double pkm2( 0.0 );
00822       double qkm1( 1.0 );
00823       double qkm2( 1.0 );
00824 
00825       double ans( 1.0 );
00826       double r( 1.0 );
00827       double thresh( 3.0 * epsilon );
00828 
00829       int n( 0 );
00830 
00831          // Start iteration
00832       do
00833       {
00834 
00835          double xk( -x * k1 * k2 / (k3*k4) );
00836          double pk( pkm1 + pkm2 * xk );
00837          double qk( qkm1 + qkm2 * xk );
00838 
00839          pkm2 = pkm1;
00840          pkm1 = pk;
00841          qkm2 = qkm1;
00842          qkm1 = qk;
00843 
00844          xk = x * k5 * k6 / (k7*k8);
00845          pk = pkm1 + pkm2 * xk;
00846          qk = qkm1 + qkm2 * xk;
00847 
00848          pkm2 = pkm1;
00849          pkm1 = pk;
00850          qkm2 = qkm1;
00851          qkm1 = qk;
00852 
00853          if( qk != 0.0 )
00854          {
00855             r = pk/qk;
00856          }
00857 
00858          double t;
00859 
00860          if( r != 0.0 )
00861          {
00862             t = std::fabs( (ans-r)/r );
00863             ans = r;
00864          }
00865          else
00866          {
00867             t = 1.0;
00868          }
00869 
00870             // Exit if we are under threshold
00871          if( t < thresh )
00872          {
00873             break;
00874          }
00875 
00876             // Recompute auxiliar parameters
00877          k1 = k1 + 1.0;
00878          k2 = k2 + 1.0;
00879          k3 = k3 + 2.0;
00880          k4 = k4 + 2.0;
00881          k5 = k5 + 1.0;
00882          k6 = k6 - 1.0;
00883          k7 = k7 + 2.0;
00884          k8 = k8 + 2.0;
00885 
00886             // Check if qk, pk are too big
00887          if( (std::fabs(qk) + std::fabs(pk)) > big )
00888          {
00889             pkm2 = pkm2 * biginv;
00890             pkm1 = pkm1 * biginv;
00891             qkm2 = qkm2 * biginv;
00892             qkm1 = qkm1 * biginv;
00893          }
00894 
00895             // Check if qk, pk are too small
00896          if( (std::fabs(qk) < biginv) ||
00897              (std::fabs(pk) < biginv ) )
00898          {
00899             pkm2 = pkm2 * big;
00900             pkm1 = pkm1 * big;
00901             qkm2 = qkm2 * big;
00902             qkm1 = qkm1 * big;
00903          }
00904 
00905             // Increase n
00906          n = n+1;
00907 
00908       }
00909       while( n < 300 );
00910 
00911          // Return result
00912       return ans;
00913 
00914    }  // End of function 'incompletebetafe()'
00915 
00916 
00917 
00918       // Auxiliar function to compute incomplete beta function.
00919    double incompletebetafe2(const double x, const double a, const double b)
00920    {
00921 
00922          // Declare auxiliar parameters
00923       double epsilon(1.0e-30);
00924       double big( 1.0e16 );
00925       double biginv( 1.0/big );
00926 
00927       double k1( a );
00928       double k2( b-1.0 );
00929       double k3( a );
00930       double k4( a+1.0 );
00931       double k5( 1.0 );
00932       double k6( a+b );
00933       double k7( a+1.0 );
00934       double k8( a+2.0 );
00935 
00936       double pkm1( 1.0 );
00937       double pkm2( 0.0 );
00938       double qkm1( 1.0 );
00939       double qkm2( 1.0 );
00940 
00941       double z( x / (1.0-x) );
00942       double ans( 1.0 );
00943       double r( 1.0 );
00944       double thresh( 3.0 * epsilon );
00945 
00946       int n( 0 );
00947 
00948          // Start iteration
00949       do
00950       {
00951 
00952          double xk( -z * k1 * k2 / (k3*k4) );
00953          double pk( pkm1 + pkm2 * xk );
00954          double qk( qkm1 + qkm2 * xk );
00955 
00956          pkm2 = pkm1;
00957          pkm1 = pk;
00958          qkm2 = qkm1;
00959          qkm1 = qk;
00960 
00961          xk = z * k5 * k6 / (k7*k8);
00962          pk = pkm1 + pkm2 * xk;
00963          qk = qkm1 + qkm2 * xk;
00964 
00965          pkm2 = pkm1;
00966          pkm1 = pk;
00967          qkm2 = qkm1;
00968          qkm1 = qk;
00969 
00970          if( qk != 0.0 )
00971          {
00972             r = pk/qk;
00973          }
00974 
00975          double t;
00976 
00977          if( r != 0.0 )
00978          {
00979             t = std::fabs( (ans-r) / r );
00980             ans = r;
00981          }
00982          else
00983          {
00984             t = 1.0;
00985          }
00986 
00987             // Exit if we are under threshold
00988          if( t < thresh )
00989          {
00990             break;
00991          }
00992 
00993 
00994             // Recompute auxiliar parameters
00995          k1 = k1 + 1.0;
00996          k2 = k2 - 1.0;
00997          k3 = k3 + 2.0;
00998          k4 = k4 + 2.0;
00999          k5 = k5 + 1.0;
01000          k6 = k6 + 1.0;
01001          k7 = k7 + 2.0;
01002          k8 = k8 + 2.0;
01003 
01004             // Check if qk, pk are too big
01005          if( (std::fabs(qk) + std::fabs(pk)) > big )
01006          {
01007             pkm2 = pkm2 * biginv;
01008             pkm1 = pkm1 * biginv;
01009             qkm2 = qkm2 * biginv;
01010             qkm1 = qkm1 * biginv;
01011          }
01012 
01013             // Check if qk, pk are too small
01014          if( (std::fabs(qk) < biginv) ||
01015              (std::fabs(pk) < biginv ) )
01016          {
01017             pkm2 = pkm2 * big;
01018             pkm1 = pkm1 * big;
01019             qkm2 = qkm2 * big;
01020             qkm1 = qkm1 * big;
01021          }
01022 
01023             // Increase n
01024          n = n+1;
01025 
01026       }
01027       while( n < 300 );
01028 
01029          // Return result
01030       return ans;
01031 
01032    }  // End of function 'incompletebetafe2()'
01033 
01034 
01035 
01036       /* Computes the regularized incomplete Beta function Ix(a,b).
01037        *
01038        * This code is a C++ implementation and adaptation from code found
01039        * in Cephes Math Library Release 2.8, copyright by Stephen L. Moshier,
01040        * released under a BSD license.
01041        */
01042    double regIncompleteBeta(const double x, const double a, const double b)
01043       throw(InvalidParameter)
01044    {
01045 
01046          // Let's do some basic checks
01047       if( (a <= 0.0) || (b <= 0.0) )
01048       {
01049          InvalidParameter e("Function 'regIncompleteBeta()': 'a' and 'b' must \
01050 be greater than zero." );
01051          GPSTK_THROW(e);
01052       }
01053 
01054       if( (x < 0.0) || (x > 1.0) )
01055       {
01056          InvalidParameter e("Function 'regIncompleteBeta()': 'x' must be \
01057 within the interval [0,1]." );
01058          GPSTK_THROW(e);
01059       }
01060 
01061       if (x == 0.0) return 0.0;
01062       if (x == 1.0) return 1.0;
01063 
01064          // Declare working parameters
01065       int flag(0);
01066       double epsilon(1.0e-30);
01067       double xx(x);
01068       double aa(a);
01069       double bb(b);
01070 
01071 
01072          // Check if bb*xx is small and xx not too close to 1.
01073       if( (bb*xx <= 1.0) &&
01074           (xx <= 0.95) )
01075       {
01076          return incompletebetaps(xx, aa, bb);
01077       }
01078 
01079       double w( 1.0-xx );
01080       double t, xc;
01081 
01082       if( xx > (aa/(aa+bb)) )
01083       {
01084          flag = 1;
01085          t = aa;
01086          aa = bb;
01087          bb = t;
01088          xc = xx;
01089          xx = w;
01090       }
01091       else
01092       {
01093          xc = w;
01094       }
01095 
01096       double result(0.0);
01097 
01098       if( (flag==1)      &&
01099           (bb*xx <= 1.0) &&
01100           (xx <= 0.95) )
01101       {
01102 
01103             // bb*xx is small and xx not too close to 1.
01104          t = incompletebetaps(xx, aa, bb);
01105 
01106          if( t <= epsilon )
01107          {
01108             result = 1.0 - epsilon;
01109          }
01110          else
01111          {
01112             result = 1.0 - t;
01113          }
01114 
01115          return result;
01116 
01117       }
01118 
01119       double y( xx * (aa+bb-2.0) - (aa-1.0) );
01120 
01121       if( y < 0.0 )
01122       {
01123          w = incompletebetafe(xx, aa, bb);
01124       }
01125       else
01126       {
01127          w = incompletebetafe2(xx, aa, bb) / xc;
01128       }
01129 
01130       t = std::pow(xc, bb);
01131       t = t * std::pow(xx, aa);
01132       t = t/aa;
01133       t = t*w;
01134       t = t * ( gamma(aa+bb) / ( gamma(aa) * gamma(bb) ) );
01135 
01136       if( flag==1 )
01137       {
01138          if( t <= epsilon )
01139          {
01140             result = 1.0 - epsilon;
01141          }
01142          else
01143          {
01144             result = 1.0 - t;
01145          }
01146       }
01147       else
01148       {
01149          result = t;
01150       }
01151 
01152       return result;
01153 
01154    }  // End of function 'incompleteBeta()'
01155 
01156 
01157 
01158 }  // End of namespace gpstk

Generated on Sun May 19 03:31:14 2013 for GPS ToolKit Software Library by  doxygen 1.3.9.1