gsl_inte_qng.h

00001 /*
00002   -------------------------------------------------------------------
00003 
00004   Copyright (C) 2006, 2007, 2008, Andrew W. Steiner
00005 
00006   This file is part of O2scl.
00007 
00008   O2scl is free software; you can redistribute it and/or modify
00009   it under the terms of the GNU General Public License as published by
00010   the Free Software Foundation; either version 3 of the License, or
00011   (at your option) any later version.
00012 
00013   O2scl is distributed in the hope that it will be useful,
00014   but WITHOUT ANY WARRANTY; without even the implied warranty of
00015   MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00016   GNU General Public License for more details.
00017 
00018   You should have received a copy of the GNU General Public License
00019   along with O2scl. If not, see <http://www.gnu.org/licenses/>.
00020 
00021   -------------------------------------------------------------------
00022 */
00023 #ifndef O2SCL_GSL_INTE_QNG_H
00024 #define O2SCL_GSL_INTE_QNG_H
00025 
00026 #include <o2scl/inte.h>
00027 #include <o2scl/gsl_inte.h>
00028  
00029 /** \brief A namespace for the quadrature coefficients for 
00030     non-adaptive integration
00031     
00032     <b>Documentation from GSL</b>: \n
00033     Gauss-Kronrod-Patterson quadrature coefficients for use in
00034     quadpack routine qng. These coefficients were calculated with
00035     101 decimal digit arithmetic by L. W. Fullerton, Bell Labs, Nov
00036     1981. 
00037 
00038 */
00039 namespace o2scl_inte_qng_coeffs {
00040 
00041   /** x1, abscissae common to the 10-, 21-, 43- and 87-point rule */
00042   static const double x1[5] = {
00043     0.973906528517171720077964012084452,
00044     0.865063366688984510732096688423493,
00045     0.679409568299024406234327365114874,
00046     0.433395394129247190799265943165784,
00047     0.148874338981631210884826001129720
00048   } ;
00049 
00050   /** w10, weights of the 10-point formula */
00051   static const double w10[5] = {
00052     0.066671344308688137593568809893332,
00053     0.149451349150580593145776339657697,
00054     0.219086362515982043995534934228163,
00055     0.269266719309996355091226921569469,
00056     0.295524224714752870173892994651338
00057   } ;
00058 
00059   /** x2, abscissae common to the 21-, 43- and 87-point rule */
00060   static const double x2[5] = {
00061     0.995657163025808080735527280689003,
00062     0.930157491355708226001207180059508,
00063     0.780817726586416897063717578345042,
00064     0.562757134668604683339000099272694,
00065     0.294392862701460198131126603103866
00066   } ;
00067 
00068   /** w21a, weights of the 21-point formula for abscissae x1 */
00069   static const double w21a[5] = {
00070     0.032558162307964727478818972459390,
00071     0.075039674810919952767043140916190,
00072     0.109387158802297641899210590325805,
00073     0.134709217311473325928054001771707,
00074     0.147739104901338491374841515972068
00075   } ;
00076 
00077   /** w21b, weights of the 21-point formula for abscissae x2 */
00078   static const double w21b[6] = {
00079     0.011694638867371874278064396062192,
00080     0.054755896574351996031381300244580,
00081     0.093125454583697605535065465083366,
00082     0.123491976262065851077958109831074,
00083     0.142775938577060080797094273138717,
00084     0.149445554002916905664936468389821
00085   } ;
00086 
00087   /** x3, abscissae common to the 43- and 87-point rule */
00088   static const double x3[11] = {
00089     0.999333360901932081394099323919911,
00090     0.987433402908088869795961478381209,
00091     0.954807934814266299257919200290473,
00092     0.900148695748328293625099494069092,
00093     0.825198314983114150847066732588520,
00094     0.732148388989304982612354848755461,
00095     0.622847970537725238641159120344323,
00096     0.499479574071056499952214885499755,
00097     0.364901661346580768043989548502644,
00098     0.222254919776601296498260928066212,
00099     0.074650617461383322043914435796506
00100   } ;
00101 
00102   /** w43a, weights of the 43-point formula for abscissae x1, x3 */
00103   static const double w43a[10] = {
00104     0.016296734289666564924281974617663,
00105     0.037522876120869501461613795898115,
00106     0.054694902058255442147212685465005,
00107     0.067355414609478086075553166302174,
00108     0.073870199632393953432140695251367,
00109     0.005768556059769796184184327908655,
00110     0.027371890593248842081276069289151,
00111     0.046560826910428830743339154433824,
00112     0.061744995201442564496240336030883,
00113     0.071387267268693397768559114425516
00114   } ;
00115 
00116   /** w43b, weights of the 43-point formula for abscissae x3 */
00117   static const double w43b[12] = {
00118     0.001844477640212414100389106552965,
00119     0.010798689585891651740465406741293,
00120     0.021895363867795428102523123075149,
00121     0.032597463975345689443882222526137,
00122     0.042163137935191811847627924327955,
00123     0.050741939600184577780189020092084,
00124     0.058379395542619248375475369330206,
00125     0.064746404951445885544689259517511,
00126     0.069566197912356484528633315038405,
00127     0.072824441471833208150939535192842,
00128     0.074507751014175118273571813842889,
00129     0.074722147517403005594425168280423
00130   } ;
00131 
00132   /** x4, abscissae of the 87-point rule */
00133   static const double x4[22] = {
00134     0.999902977262729234490529830591582,
00135     0.997989895986678745427496322365960,
00136     0.992175497860687222808523352251425,
00137     0.981358163572712773571916941623894,
00138     0.965057623858384619128284110607926,
00139     0.943167613133670596816416634507426,
00140     0.915806414685507209591826430720050,
00141     0.883221657771316501372117548744163,
00142     0.845710748462415666605902011504855,
00143     0.803557658035230982788739474980964,
00144     0.757005730685495558328942793432020,
00145     0.706273209787321819824094274740840,
00146     0.651589466501177922534422205016736,
00147     0.593223374057961088875273770349144,
00148     0.531493605970831932285268948562671,
00149     0.466763623042022844871966781659270,
00150     0.399424847859218804732101665817923,
00151     0.329874877106188288265053371824597,
00152     0.258503559202161551802280975429025,
00153     0.185695396568346652015917141167606,
00154     0.111842213179907468172398359241362,
00155     0.037352123394619870814998165437704
00156   } ;
00157 
00158   /** w87a, weights of the 87-point formula for abscissae x1, x2, x3 */
00159   static const double w87a[21] = {
00160     0.008148377384149172900002878448190,
00161     0.018761438201562822243935059003794,
00162     0.027347451050052286161582829741283,
00163     0.033677707311637930046581056957588,
00164     0.036935099820427907614589586742499,
00165     0.002884872430211530501334156248695,
00166     0.013685946022712701888950035273128,
00167     0.023280413502888311123409291030404,
00168     0.030872497611713358675466394126442,
00169     0.035693633639418770719351355457044,
00170     0.000915283345202241360843392549948,
00171     0.005399280219300471367738743391053,
00172     0.010947679601118931134327826856808,
00173     0.016298731696787335262665703223280,
00174     0.021081568889203835112433060188190,
00175     0.025370969769253827243467999831710,
00176     0.029189697756475752501446154084920,
00177     0.032373202467202789685788194889595,
00178     0.034783098950365142750781997949596,
00179     0.036412220731351787562801163687577,
00180     0.037253875503047708539592001191226
00181   } ;
00182 
00183   /** w87b, weights of the 87-point formula for abscissae x4    */
00184   static const double w87b[23] = {
00185     0.000274145563762072350016527092881,
00186     0.001807124155057942948341311753254,
00187     0.004096869282759164864458070683480,
00188     0.006758290051847378699816577897424,
00189     0.009549957672201646536053581325377,
00190     0.012329447652244853694626639963780,
00191     0.015010447346388952376697286041943,
00192     0.017548967986243191099665352925900,
00193     0.019938037786440888202278192730714,
00194     0.022194935961012286796332102959499,
00195     0.024339147126000805470360647041454,
00196     0.026374505414839207241503786552615,
00197     0.028286910788771200659968002987960,
00198     0.030052581128092695322521110347341,
00199     0.031646751371439929404586051078883,
00200     0.033050413419978503290785944862689,
00201     0.034255099704226061787082821046821,
00202     0.035262412660156681033782717998428,
00203     0.036076989622888701185500318003895,
00204     0.036698604498456094498018047441094,
00205     0.037120549269832576114119958413599,
00206     0.037334228751935040321235449094698,
00207     0.037361073762679023410321241766599
00208   } ;
00209 
00210 }
00211 
00212 #ifndef DOXYGENP
00213 namespace o2scl {
00214 #endif
00215 
00216   /** 
00217       \brief Non-adaptive integration from a to b (GSL)
00218 
00219       integ() uses 10-point, 21-point, 43-point, and 87-point
00220       Gauss-Kronrod integration successively until the integral
00221       is returned within the accuracy specified by tolx and
00222       tolf.
00223 
00224       \future Compare directly with GSL as is done in gsl_inte_qag_ts.
00225   */
00226   template<class param_t, class func_t> class gsl_inte_qng :
00227   public inte<param_t,func_t>, public gsl_inte {
00228     
00229   public:
00230     
00231     /** 
00232         \brief The number of function evalutions for the last integration
00233 
00234         Set to either 0, 21, 43, or 87, depending on the number of
00235         function evaluations that were used. This variable is zero if
00236         an error occurs before any function evaluations were performed
00237         and is never equal 10, since in the 10-point method, the
00238         21-point result is used to estimate the error. If the function
00239         fails to achieve the desired precision, feval is set to 88. 
00240     */
00241     size_t feval;
00242     
00243     /** \brief Integrate function \c func from \c a to \c b.
00244      */
00245     virtual double integ(func_t &func, double a, double b, param_t &pa) {
00246       double res, err;
00247       integ_err(func,a,b,pa,res,err);
00248       this->interror=err;
00249       return res;
00250     }
00251 
00252     /** \brief Integrate function \c func from \c a to \c b
00253         giving result \c res and error \c err
00254     */
00255     virtual int integ_err(func_t &func, double a, double b, 
00256                           param_t &pa, double &res, double &err2) {
00257     
00258       double fv1[5], fv2[5], fv3[5], fv4[5];
00259 
00260       /* array of function values which have been computed */
00261       double savfun[21];  
00262     
00263       /* 10, 21, 43 and 87 point results */
00264       double res10, res21, res43, res87;    
00265 
00266       double result_kronrod, err;
00267     
00268       /* approximation to the integral of abs(f) */
00269       double resabs; 
00270     
00271       /* approximation to the integral of abs(f-i/(b-a)) */
00272       double resasc; 
00273     
00274       const double half_length =  0.5 * (b - a);
00275       const double abs_half_length = fabs (half_length);
00276       const double center = 0.5 * (b + a);
00277 
00278       double f_center;
00279       func(center,f_center,pa);
00280     
00281       int k;
00282     
00283       if (this->tolx <= 0 && (this->tolf < 50 * GSL_DBL_EPSILON || 
00284                           this->tolf < 0.5e-28)) {
00285         res = 0;
00286         err2 = 0;
00287         feval = 0;
00288 
00289         std::string estr="Tolerance cannot be achieved with given ";
00290         estr+="value of 'tolx' and 'tolf' in gsl_inte_qng::integ_err().";
00291         set_err_ret(estr.c_str(),gsl_ebadtol);
00292       };
00293       
00294       /* Compute the integral using the 10- and 21-point formula. */
00295     
00296       res10 = 0;
00297       res21 = o2scl_inte_qng_coeffs::w21b[5] * f_center;
00298       resabs = o2scl_inte_qng_coeffs::w21b[5] * fabs (f_center);
00299     
00300     
00301       for (k = 0; k < 5; k++) {
00302         const double abscissa = half_length * o2scl_inte_qng_coeffs::x1[k];
00303         double fval1, fval2;
00304         func(center+abscissa,fval1,pa);
00305         func(center-abscissa,fval2,pa);
00306         const double fval = fval1 + fval2;
00307         res10 += o2scl_inte_qng_coeffs::w10[k] * fval;
00308         res21 += o2scl_inte_qng_coeffs::w21a[k] * fval;
00309         resabs += o2scl_inte_qng_coeffs::w21a[k] * 
00310           (fabs (fval1) + fabs (fval2));
00311         savfun[k] = fval;
00312         fv1[k] = fval1;
00313         fv2[k] = fval2;
00314       }
00315     
00316       for (k = 0; k < 5; k++) {
00317         const double abscissa = half_length * o2scl_inte_qng_coeffs::x2[k];
00318         double fval1, fval2;
00319         func(center+abscissa,fval1,pa);
00320         func(center-abscissa,fval2,pa);
00321         const double fval = fval1 + fval2;
00322         res21 += o2scl_inte_qng_coeffs::w21b[k] * fval;
00323         resabs += o2scl_inte_qng_coeffs::w21b[k] * 
00324           (fabs (fval1) + fabs (fval2));
00325         savfun[k + 5] = fval;
00326         fv3[k] = fval1;
00327         fv4[k] = fval2;
00328       }
00329       
00330       resabs *= abs_half_length ;
00331     
00332       {
00333         const double mean = 0.5 * res21;
00334       
00335         resasc = o2scl_inte_qng_coeffs::w21b[5] * fabs (f_center - mean);
00336       
00337         for (k = 0; k < 5; k++) {
00338           resasc +=
00339             (o2scl_inte_qng_coeffs::w21a[k] * (fabs (fv1[k] - mean) + 
00340                                              fabs (fv2[k] - mean))
00341              + o2scl_inte_qng_coeffs::w21b[k] * (fabs (fv3[k] - mean) + 
00342                                                fabs (fv4[k] - mean)));
00343         }
00344         resasc *= abs_half_length ;
00345       }
00346       
00347       result_kronrod = res21 * half_length;
00348       
00349       err = rescale_error ((res21 - res10) * half_length, resabs, resasc) ;
00350     
00351       /* test for convergence. */
00352     
00353       if (err < this->tolx || err < this->tolf * fabs (result_kronrod)) {
00354         res = result_kronrod ;
00355         err2 = err ;
00356         feval = 21;
00357         return gsl_success;
00358       }
00359       
00360       /* compute the integral using the 43-point formula. */
00361     
00362       res43 = o2scl_inte_qng_coeffs::w43b[11] * f_center;
00363     
00364       for (k = 0; k < 10; k++) {
00365         res43 += savfun[k] * o2scl_inte_qng_coeffs::w43a[k];
00366       }
00367       
00368       for (k = 0; k < 11; k++) {
00369         const double abscissa = half_length * o2scl_inte_qng_coeffs::x3[k];
00370         double fval1, fval2;
00371         func(center+abscissa,fval1,pa);
00372         func(center-abscissa,fval2,pa);
00373         const double fval = fval1+fval2;
00374         res43 += fval * o2scl_inte_qng_coeffs::w43b[k];
00375         savfun[k + 10] = fval;
00376       }
00377       
00378       /*  test for convergence */
00379       result_kronrod = res43 * half_length;
00380       err = rescale_error ((res43 - res21) * half_length, resabs, resasc);
00381       
00382       if (err < this->tolx || err < this->tolf * fabs (result_kronrod)) {
00383         res = result_kronrod ;
00384         err2 = err ;
00385         feval = 43;
00386         return gsl_success;
00387       }
00388       
00389       /* compute the integral using the 87-point formula. */
00390     
00391       res87 = o2scl_inte_qng_coeffs::w87b[22] * f_center;
00392     
00393       for (k = 0; k < 21; k++) {
00394         res87 += savfun[k] * o2scl_inte_qng_coeffs::w87a[k];
00395       }
00396       
00397       for (k = 0; k < 22; k++) {
00398         const double abscissa = half_length * o2scl_inte_qng_coeffs::x4[k];
00399         double fval1, fval2;
00400         func(center+abscissa,fval1,pa);
00401         func(center-abscissa,fval2,pa);
00402         res87 += o2scl_inte_qng_coeffs::w87b[k] * (fval1+fval2);
00403       }
00404     
00405       /* test for convergence */
00406     
00407       result_kronrod = res87 * half_length ;
00408     
00409       err = rescale_error ((res87 - res43) * half_length, resabs, resasc);
00410     
00411       if (err < this->tolx || err < this->tolf * fabs (result_kronrod)) {
00412         res = result_kronrod ;
00413         err2 = err ;
00414         feval = 87;
00415         return gsl_success;
00416       }
00417       
00418       /* failed to converge */
00419       res = result_kronrod ;
00420       err2 = err ;
00421       feval = 88;
00422     
00423       std::string estr="Failed to reach tolerance ";
00424       estr+="in gsl_inte_qng::integ_err().";
00425       set_err_ret(estr.c_str(),gsl_etol);
00426     }
00427   
00428     /// Return string denoting type ("gsl_inte_qng")
00429     const char *type() { return "gsl_inte_qng"; }
00430 
00431   };
00432 
00433 #ifndef DOXYGENP
00434 }
00435 #endif
00436 
00437 #endif

Documentation generated with Doxygen and provided under the GNU Free Documentation License. See License Information for details.

Project hosting provided by SourceForge.net Logo, O2scl Sourceforge Project Page