Object-oriented Scientific Computing Library: Version 0.910
gsl_inte_kronrod.h
00001 /*
00002   -------------------------------------------------------------------
00003   
00004   Copyright (C) 2006-2012, Jerry Gagelman
00005   and Andrew W. Steiner
00006   
00007   This file is part of O2scl.
00008   
00009   O2scl is free software; you can redistribute it and/or modify
00010   it under the terms of the GNU General Public License as published by
00011   the Free Software Foundation; either version 3 of the License, or
00012   (at your option) any later version.
00013   
00014   O2scl is distributed in the hope that it will be useful,
00015   but WITHOUT ANY WARRANTY; without even the implied warranty of
00016   MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00017   GNU General Public License for more details.
00018   
00019   You should have received a copy of the GNU General Public License
00020   along with O2scl. If not, see <http://www.gnu.org/licenses/>.
00021   
00022   -------------------------------------------------------------------
00023 */
00024 #ifndef O2SCL_GSL_INTE_KRONROD_H
00025 #define O2SCL_GSL_INTE_KRONROD_H
00026 
00027 #include <cmath>
00028 #include <gsl/gsl_integration.h>
00029 #include <o2scl/inte.h>
00030 #include <o2scl/gsl_inte.h>
00031 #include <o2scl/string_conv.h>
00032 
00033 /** \brief A namespace for the GSL adaptive Gauss-Kronrod integration
00034     coefficients
00035 
00036     The storage scheme follows that of QUADPACK: symmetry about the  
00037     origin permits storing only those abscissae in the interval 
00038     \f$ [0, 1). \f$ Similiarly for the weights. The odd-indexed
00039     abscissae \c xgk[1], \c xgk[3],... belong to the low-order Gauss 
00040     rule. For the full Gauss-Kronrod rule, the array \c wgk[] holds
00041     the weights corresponding to the abscissae \c xgk[]. For the
00042     low-order rule, \c wg[i] is the weight correpsonding to
00043     \c xgk[2*i+1].
00044  
00045     <b>Documentation from GSL</b>: \n    
00046     Gauss quadrature weights and kronrod quadrature abscissae and
00047     weights as evaluated with 80 decimal digit arithmetic by
00048     L. W. Fullerton, Bell Labs, Nov. 1981.
00049 */
00050 namespace o2scl_inte_gk_coeffs {
00051   
00052   /** \brief Abscissae of the 15-point Kronrod rule */
00053   static const double qk15_xgk[8]={
00054     0.991455371120812639206854697526329,
00055     0.949107912342758524526189684047851,
00056     0.864864423359769072789712788640926,
00057     0.741531185599394439863864773280788,
00058     0.586087235467691130294144838258730,
00059     0.405845151377397166906606412076961,
00060     0.207784955007898467600689403773245,
00061     0.000000000000000000000000000000000
00062   };
00063   
00064   /** \brief Weights of the 7-point Gauss rule */
00065   static const double qk15_wg[4] =     
00066     {
00067       0.129484966168869693270611432679082,
00068       0.279705391489276667901467771423780,
00069       0.381830050505118944950369775488975,
00070       0.417959183673469387755102040816327
00071     };
00072 
00073   /** \brief Weights of the 15-point Kronrod rule */
00074   static const double qk15_wgk[8] =    
00075     {
00076       0.022935322010529224963732008058970,
00077       0.063092092629978553290700663189204,
00078       0.104790010322250183839876322541518,
00079       0.140653259715525918745189590510238,
00080       0.169004726639267902826583426598550,
00081       0.190350578064785409913256402421014,
00082       0.204432940075298892414161999234649,
00083       0.209482141084727828012999174891714
00084     };
00085 
00086   /** \brief Abscissae of the 21-point Kronrod rule */
00087   static const double qk21_xgk[11] =   
00088     {
00089       0.995657163025808080735527280689003,
00090       0.973906528517171720077964012084452,
00091       0.930157491355708226001207180059508,
00092       0.865063366688984510732096688423493,
00093       0.780817726586416897063717578345042,
00094       0.679409568299024406234327365114874,
00095       0.562757134668604683339000099272694,
00096       0.433395394129247190799265943165784,
00097       0.294392862701460198131126603103866,
00098       0.148874338981631210884826001129720,
00099       0.000000000000000000000000000000000
00100     };
00101 
00102   /** \brief Weights of the 10-point Gauss rule */
00103   static const double qk21_wg[5] =     
00104     {
00105       0.066671344308688137593568809893332,
00106       0.149451349150580593145776339657697,
00107       0.219086362515982043995534934228163,
00108       0.269266719309996355091226921569469,
00109       0.295524224714752870173892994651338
00110     };
00111 
00112   /** \brief Weights of the 21-point Kronrod rule */
00113   static const double qk21_wgk[11] =   
00114     {
00115       0.011694638867371874278064396062192,
00116       0.032558162307964727478818972459390,
00117       0.054755896574351996031381300244580,
00118       0.075039674810919952767043140916190,
00119       0.093125454583697605535065465083366,
00120       0.109387158802297641899210590325805,
00121       0.123491976262065851077958109831074,
00122       0.134709217311473325928054001771707,
00123       0.142775938577060080797094273138717,
00124       0.147739104901338491374841515972068,
00125       0.149445554002916905664936468389821
00126     };
00127 
00128   /** \brief Abscissae of the 31-point Kronrod rule */
00129   static const double qk31_xgk[16] =   
00130     {
00131       0.998002298693397060285172840152271,
00132       0.987992518020485428489565718586613,
00133       0.967739075679139134257347978784337,
00134       0.937273392400705904307758947710209,
00135       0.897264532344081900882509656454496,
00136       0.848206583410427216200648320774217,
00137       0.790418501442465932967649294817947,
00138       0.724417731360170047416186054613938,
00139       0.650996741297416970533735895313275,
00140       0.570972172608538847537226737253911,
00141       0.485081863640239680693655740232351,
00142       0.394151347077563369897207370981045,
00143       0.299180007153168812166780024266389,
00144       0.201194093997434522300628303394596,
00145       0.101142066918717499027074231447392,
00146       0.000000000000000000000000000000000
00147     };
00148 
00149   /** \brief Weights of the 15-point Gauss rule */
00150   static const double qk31_wg[8] =     
00151     {
00152       0.030753241996117268354628393577204,
00153       0.070366047488108124709267416450667,
00154       0.107159220467171935011869546685869,
00155       0.139570677926154314447804794511028,
00156       0.166269205816993933553200860481209,
00157       0.186161000015562211026800561866423,
00158       0.198431485327111576456118326443839,
00159       0.202578241925561272880620199967519
00160     };
00161 
00162   /** \brief Weights of the 31-point Kronrod rule */
00163   static const double qk31_wgk[16] =   
00164     {
00165       0.005377479872923348987792051430128,
00166       0.015007947329316122538374763075807,
00167       0.025460847326715320186874001019653,
00168       0.035346360791375846222037948478360,
00169       0.044589751324764876608227299373280,
00170       0.053481524690928087265343147239430,
00171       0.062009567800670640285139230960803,
00172       0.069854121318728258709520077099147,
00173       0.076849680757720378894432777482659,
00174       0.083080502823133021038289247286104,
00175       0.088564443056211770647275443693774,
00176       0.093126598170825321225486872747346,
00177       0.096642726983623678505179907627589,
00178       0.099173598721791959332393173484603,
00179       0.100769845523875595044946662617570,
00180       0.101330007014791549017374792767493
00181     };
00182 
00183   /** \brief Abscissae of the 41-point Kronrod rule */
00184   static const double qk41_xgk[21] =   
00185     {
00186       0.998859031588277663838315576545863,
00187       0.993128599185094924786122388471320,
00188       0.981507877450250259193342994720217,
00189       0.963971927277913791267666131197277,
00190       0.940822633831754753519982722212443,
00191       0.912234428251325905867752441203298,
00192       0.878276811252281976077442995113078,
00193       0.839116971822218823394529061701521,
00194       0.795041428837551198350638833272788,
00195       0.746331906460150792614305070355642,
00196       0.693237656334751384805490711845932,
00197       0.636053680726515025452836696226286,
00198       0.575140446819710315342946036586425,
00199       0.510867001950827098004364050955251,
00200       0.443593175238725103199992213492640,
00201       0.373706088715419560672548177024927,
00202       0.301627868114913004320555356858592,
00203       0.227785851141645078080496195368575,
00204       0.152605465240922675505220241022678,
00205       0.076526521133497333754640409398838,
00206       0.000000000000000000000000000000000
00207     };
00208 
00209   /** \brief Weights of the 20-point Gauss rule */
00210   static const double qk41_wg[11] =    
00211     {
00212       0.017614007139152118311861962351853,
00213       0.040601429800386941331039952274932,
00214       0.062672048334109063569506535187042,
00215       0.083276741576704748724758143222046,
00216       0.101930119817240435036750135480350,
00217       0.118194531961518417312377377711382,
00218       0.131688638449176626898494499748163,
00219       0.142096109318382051329298325067165,
00220       0.149172986472603746787828737001969,
00221       0.152753387130725850698084331955098
00222     };
00223 
00224   /** \brief Weights of the 41-point Kronrod rule */
00225   static const double qk41_wgk[21] =   
00226     {
00227       0.003073583718520531501218293246031,
00228       0.008600269855642942198661787950102,
00229       0.014626169256971252983787960308868,
00230       0.020388373461266523598010231432755,
00231       0.025882133604951158834505067096153,
00232       0.031287306777032798958543119323801,
00233       0.036600169758200798030557240707211,
00234       0.041668873327973686263788305936895,
00235       0.046434821867497674720231880926108,
00236       0.050944573923728691932707670050345,
00237       0.055195105348285994744832372419777,
00238       0.059111400880639572374967220648594,
00239       0.062653237554781168025870122174255,
00240       0.065834597133618422111563556969398,
00241       0.068648672928521619345623411885368,
00242       0.071054423553444068305790361723210,
00243       0.073030690332786667495189417658913,
00244       0.074582875400499188986581418362488,
00245       0.075704497684556674659542775376617,
00246       0.076377867672080736705502835038061,
00247       0.076600711917999656445049901530102
00248     };
00249 
00250   /** \brief Abscissae of the 51-point Kronrod rule */
00251   static const double qk51_xgk[26] =   
00252     {
00253       0.999262104992609834193457486540341,
00254       0.995556969790498097908784946893902,
00255       0.988035794534077247637331014577406,
00256       0.976663921459517511498315386479594,
00257       0.961614986425842512418130033660167,
00258       0.942974571228974339414011169658471,
00259       0.920747115281701561746346084546331,
00260       0.894991997878275368851042006782805,
00261       0.865847065293275595448996969588340,
00262       0.833442628760834001421021108693570,
00263       0.797873797998500059410410904994307,
00264       0.759259263037357630577282865204361,
00265       0.717766406813084388186654079773298,
00266       0.673566368473468364485120633247622,
00267       0.626810099010317412788122681624518,
00268       0.577662930241222967723689841612654,
00269       0.526325284334719182599623778158010,
00270       0.473002731445714960522182115009192,
00271       0.417885382193037748851814394594572,
00272       0.361172305809387837735821730127641,
00273       0.303089538931107830167478909980339,
00274       0.243866883720988432045190362797452,
00275       0.183718939421048892015969888759528,
00276       0.122864692610710396387359818808037,
00277       0.061544483005685078886546392366797,
00278       0.000000000000000000000000000000000
00279     };
00280 
00281   /** \brief Weights of the 25-point Gauss rule */
00282   static const double qk51_wg[13] =    
00283     {
00284       0.011393798501026287947902964113235,
00285       0.026354986615032137261901815295299,
00286       0.040939156701306312655623487711646,
00287       0.054904695975835191925936891540473,
00288       0.068038333812356917207187185656708,
00289       0.080140700335001018013234959669111,
00290       0.091028261982963649811497220702892,
00291       0.100535949067050644202206890392686,
00292       0.108519624474263653116093957050117,
00293       0.114858259145711648339325545869556,
00294       0.119455763535784772228178126512901,
00295       0.122242442990310041688959518945852,
00296       0.123176053726715451203902873079050
00297     };
00298 
00299   /** \brief Weights of the 51-point Kronrod rule */
00300   static const double qk51_wgk[26] =   
00301     {
00302       0.001987383892330315926507851882843,
00303       0.005561932135356713758040236901066,
00304       0.009473973386174151607207710523655,
00305       0.013236229195571674813656405846976,
00306       0.016847817709128298231516667536336,
00307       0.020435371145882835456568292235939,
00308       0.024009945606953216220092489164881,
00309       0.027475317587851737802948455517811,
00310       0.030792300167387488891109020215229,
00311       0.034002130274329337836748795229551,
00312       0.037116271483415543560330625367620,
00313       0.040083825504032382074839284467076,
00314       0.042872845020170049476895792439495,
00315       0.045502913049921788909870584752660,
00316       0.047982537138836713906392255756915,
00317       0.050277679080715671963325259433440,
00318       0.052362885806407475864366712137873,
00319       0.054251129888545490144543370459876,
00320       0.055950811220412317308240686382747,
00321       0.057437116361567832853582693939506,
00322       0.058689680022394207961974175856788,
00323       0.059720340324174059979099291932562,
00324       0.060539455376045862945360267517565,
00325       0.061128509717053048305859030416293,
00326       0.061471189871425316661544131965264,
00327       0.061580818067832935078759824240066
00328     };
00329 
00330   /** \brief Abscissae of the 61-point Kronrod rule */
00331   static const double qk61_xgk[31] =   
00332     {
00333       0.999484410050490637571325895705811,
00334       0.996893484074649540271630050918695,
00335       0.991630996870404594858628366109486,
00336       0.983668123279747209970032581605663,
00337       0.973116322501126268374693868423707,
00338       0.960021864968307512216871025581798,
00339       0.944374444748559979415831324037439,
00340       0.926200047429274325879324277080474,
00341       0.905573307699907798546522558925958,
00342       0.882560535792052681543116462530226,
00343       0.857205233546061098958658510658944,
00344       0.829565762382768397442898119732502,
00345       0.799727835821839083013668942322683,
00346       0.767777432104826194917977340974503,
00347       0.733790062453226804726171131369528,
00348       0.697850494793315796932292388026640,
00349       0.660061064126626961370053668149271,
00350       0.620526182989242861140477556431189,
00351       0.579345235826361691756024932172540,
00352       0.536624148142019899264169793311073,
00353       0.492480467861778574993693061207709,
00354       0.447033769538089176780609900322854,
00355       0.400401254830394392535476211542661,
00356       0.352704725530878113471037207089374,
00357       0.304073202273625077372677107199257,
00358       0.254636926167889846439805129817805,
00359       0.204525116682309891438957671002025,
00360       0.153869913608583546963794672743256,
00361       0.102806937966737030147096751318001,
00362       0.051471842555317695833025213166723,
00363       0.000000000000000000000000000000000
00364     };
00365 
00366   /** \brief Weights of the 30-point Gauss rule */
00367   static const double qk61_wg[15] =    
00368     {
00369       0.007968192496166605615465883474674,
00370       0.018466468311090959142302131912047,
00371       0.028784707883323369349719179611292,
00372       0.038799192569627049596801936446348,
00373       0.048402672830594052902938140422808,
00374       0.057493156217619066481721689402056,
00375       0.065974229882180495128128515115962,
00376       0.073755974737705206268243850022191,
00377       0.080755895229420215354694938460530,
00378       0.086899787201082979802387530715126,
00379       0.092122522237786128717632707087619,
00380       0.096368737174644259639468626351810,
00381       0.099593420586795267062780282103569,
00382       0.101762389748405504596428952168554,
00383       0.102852652893558840341285636705415
00384     };
00385 
00386   /** \brief Weights of the 61-point Kronrod rule */
00387   static const double qk61_wgk[31] =   
00388     {
00389       0.001389013698677007624551591226760,
00390       0.003890461127099884051267201844516,
00391       0.006630703915931292173319826369750,
00392       0.009273279659517763428441146892024,
00393       0.011823015253496341742232898853251,
00394       0.014369729507045804812451432443580,
00395       0.016920889189053272627572289420322,
00396       0.019414141193942381173408951050128,
00397       0.021828035821609192297167485738339,
00398       0.024191162078080601365686370725232,
00399       0.026509954882333101610601709335075,
00400       0.028754048765041292843978785354334,
00401       0.030907257562387762472884252943092,
00402       0.032981447057483726031814191016854,
00403       0.034979338028060024137499670731468,
00404       0.036882364651821229223911065617136,
00405       0.038678945624727592950348651532281,
00406       0.040374538951535959111995279752468,
00407       0.041969810215164246147147541285970,
00408       0.043452539701356069316831728117073,
00409       0.044814800133162663192355551616723,
00410       0.046059238271006988116271735559374,
00411       0.047185546569299153945261478181099,
00412       0.048185861757087129140779492298305,
00413       0.049055434555029778887528165367238,
00414       0.049795683427074206357811569379942,
00415       0.050405921402782346840893085653585,
00416       0.050881795898749606492297473049805,
00417       0.051221547849258772170656282604944,
00418       0.051426128537459025933862879215781,
00419       0.051494729429451567558340433647099
00420     };
00421   
00422 }
00423 
00424 #ifndef DOXYGENP
00425 namespace o2scl {
00426 #endif
00427 
00428   /** \brief Integration workspace for the GSL integrators
00429 
00430       This is a simple rewrite of \c gsl_integration_workspace
00431       into a class.
00432 
00433       QUADPACK workspace documentation:
00434       \verbatim
00435       c    parameters (meaning at output)
00436       c       limit  - integer
00437       c                maximum number of error estimates the list
00438       c                can contain
00439       c          last   - integer
00440       c                number of error estimates currently in the list
00441       c          maxerr - integer
00442       c                maxerr points to the nrmax-th largest error
00443       c                estimate currently in the list
00444       c          ermax  - double precision
00445       c                nrmax-th largest error estimate
00446       c                ermax = elist(maxerr)
00447       c          elist  - double precision
00448       c                vector of dimension last containing
00449       c                the error estimates
00450       c          iord   - integer
00451       c                vector of dimension last, the first k elements
00452       c                of which contain pointers to the error
00453       c                estimates, such that
00454       c                elist(iord(1)),...,  elist(iord(k))
00455       c                form a decreasing sequence, with
00456       c                k = last if last.le.(limit/2+2), and
00457       c                k = limit+1-last otherwise
00458       c          nrmax  - integer
00459       c      maxerr = iord(nrmax)
00460       \endverbatim
00461       \verbatim
00462       c     alist   - real
00463       c               vector of dimension at least limit, the first
00464       c                last  elements of which are the left
00465       c               end points of the subintervals in the partition
00466       c               of the given integration range (a,b)
00467       c            blist   - real
00468       c               vector of dimension at least limit, the first
00469       c                last  elements of which are the right
00470       c               end points of the subintervals in the partition
00471       c               of the given integration range (a,b)
00472       c            rlist   - real
00473       c               vector of dimension at least limit, the first
00474       c                last  elements of which are the
00475       c               integral approximations on the subintervals
00476       c            elist   - real
00477       c               vector of dimension at least limit, the first
00478       c                last  elements of which are the moduli of the
00479       c               absolute error estimates on the subintervals
00480       c            iord    - integer
00481       c               vector of dimension at least limit, the first k
00482       c               elements of which are pointers to the
00483       c               error estimates over the subintervals,
00484       c               such that elist(iord(1)), ...,
00485       c               elist(iord(k)) form a decreasing sequence,
00486       c               with k = last if last.le.(limit/2+2), and
00487       c               k = limit+1-last otherwise
00488       c            last    - integer
00489       c               number of subintervals actually produced in the
00490       c               subdivision process
00491       \endverbatim
00492       
00493   */
00494   class gsl_inte_workspace {
00495     
00496   public:
00497     
00498     gsl_inte_workspace();
00499     
00500     ~gsl_inte_workspace();
00501 
00502     /// Maximum number of subintervals allocated
00503     size_t limit;
00504     /// Current number of subintervals being used
00505     size_t size;                
00506     /// Counter for extrapolation routine
00507     size_t nrmax;
00508     /// Index of current subinterval
00509     size_t i;
00510     /// Depth of subdivisions reached
00511     size_t maximum_level;
00512     /// Left endpoints of subintervals
00513     double *alist;
00514     /// Right endpoints of subintervals
00515     double *blist;
00516     /// Integral approximations on subintervals
00517     double *rlist;
00518     /// Integral error estimates
00519     double *elist;
00520     /// Linear ordering vector for sort routine
00521     size_t *order;      
00522     /// Numbers of subdivisions made
00523     size_t *level;      
00524     
00525     /// Allocate a workspace
00526     int allocate(size_t sz);
00527 
00528     /// Free allocated workspace memory
00529     int free();
00530     
00531     /** \brief Initialize the workspace for an integration with limits
00532         \c a and \c b.
00533     */
00534     int initialise(double a, double b);
00535 
00536     /** \brief Update the workspace with the result and error
00537         from the first integration
00538     */
00539     int set_initial_result(double result, double error);
00540 
00541     /** \brief Retrieve the ith result from the 
00542         workspace stack
00543         
00544         The workspace variable \c i is used to specify which 
00545         interval is requested.
00546     */
00547     int retrieve(double *a, double *b, double *r, double *e) const;
00548     
00549     /** \brief Sort the workspace stack
00550         
00551         This routine maintains the descending ordering in the list of
00552         the local error estimated resulting from the interval
00553         subdivision process. at each call two error estimates are
00554         inserted using the sequential search method, top-down for the
00555         largest error estimate and bottom-up for the smallest error
00556         estimate.
00557 
00558         Originally written in QUADPACK by R. Piessens and E. de
00559         Doncker, translated into C for GSL by Brian Gough, and then
00560         rewritten for \o2.
00561     */
00562     int qpsrt();
00563 
00564     /** \brief Determine which new subinterval to add to the workspace
00565         stack and perform update
00566     */
00567     int update(double a1, double b1, double area1, double error1,
00568                double a2, double b2, double area2, double error2);
00569     
00570     /// Add up all of the contributions to construct the final result
00571     double sum_results();
00572 
00573     /** \brief Test whether the proposed subdivision falls 
00574         before floating-point precision
00575     */
00576     int subinterval_too_small(double a1, double a2, double b2);
00577 
00578     /// Push a new interval to the workspace stack
00579     void append_interval(double a1, double b1, double area1, double error1);
00580 
00581   };
00582 
00583   /** \brief Basic Gauss-Kronrod integration class (GSL)
00584       
00585       This class provides the basic Gauss-Kronrod integration function
00586       and integration workspace for some of the GSL-based integration
00587       classes. 
00588 
00589       The main function of interest is \ref set_rule(), which
00590       sets the integration rule for the GSL integration classes
00591       which inherit from this base class. The argument to
00592       set rule should be selected from the following list:
00593       - 1: 15-point Gauss-Kronrod integration rule (default)
00594       - 2: 21-point rule
00595       - 3: 31-point rule
00596       - 4: 41-point rule
00597       - 5: 51-point rule
00598       - 6: 61-point rule
00599       Any value other than 1-6 forces the error handler to be
00600       called. All classes default to the 15-point rule,
00601       except for \ref gsl_inte_qags which defaults to
00602       the 21-point rule for singularities. 
00603 
00604       The integration coefficients for use with this class and
00605       associated children are stored in the \ref o2scl_inte_gk_coeffs
00606       namespace.
00607 
00608       \future Convert 'double *' objects to uvector
00609   */
00610   template<class func_t> class gsl_inte_kronrod : public gsl_inte, 
00611     public inte<func_t> {
00612       
00613   protected:
00614       
00615       /// The integration workspace
00616       gsl_inte_workspace *w;
00617       
00618       /// Size of Gauss-Kronrod arrays
00619       int n_gk;
00620 
00621       /// Gauss-Kronrod abscissae pointer
00622       const double *x_gk;
00623       
00624       /// Gauss weight pointer
00625       const double *w_g;
00626       
00627       /// Gauss-Kronrod weight pointer
00628       const double *w_gk;
00629       
00630       /// Scratch space
00631       double *f_v1;
00632 
00633       /// Scratch space
00634       double *f_v2;
00635 
00636   public:
00637 
00638       gsl_inte_kronrod() {
00639         w=new gsl_inte_workspace;
00640         w->allocate(1000);
00641         n_gk=0;
00642         set_rule(1);
00643       }
00644       
00645       ~gsl_inte_kronrod() {
00646         w->free();
00647         delete w;
00648         if (n_gk>0) {
00649           delete[] f_v1;
00650           delete[] f_v2;
00651         }
00652       }
00653 
00654       /** \brief Get the Gauss-Kronrod integration rule
00655 
00656           This returns the index of the GSL integration rule
00657           a number between 1 and 6 (inclusive)
00658       */
00659       int get_rule() {
00660         switch (n_gk) {
00661         case 8:
00662           return 1;
00663         case 11:
00664           return 2;
00665         case 16:
00666           return 3;
00667         case 21:
00668           return 4;
00669         case 26:
00670           return 5;
00671         case 31:
00672           return 6;
00673         }
00674         return 0;
00675       }
00676     
00677       /// Set the Gauss-Kronrod integration rule to be used
00678       void set_rule(int rule) {
00679         
00680         using namespace o2scl_inte_gk_coeffs;
00681         
00682         if (n_gk > 0) {
00683           delete[] f_v1;
00684           delete[] f_v2;
00685         }
00686 
00687         switch (rule) {
00688         case 1:
00689           n_gk=8;
00690           x_gk=qk15_xgk;
00691           w_g=qk15_wg;
00692           w_gk=qk15_wgk;
00693           break;
00694         case 2:
00695           n_gk=11;
00696           x_gk=qk21_xgk;
00697           w_g=qk21_wg;
00698           w_gk=qk21_wgk;
00699           break;
00700         case 3:
00701           n_gk=16;
00702           x_gk=qk31_xgk;
00703           w_g=qk31_wg;
00704           w_gk=qk31_wgk;
00705           break;
00706         case 4:
00707           n_gk=21;
00708           x_gk=qk41_xgk;
00709           w_g=qk41_wg;
00710           w_gk=qk41_wgk;
00711           break;
00712         case 5:
00713           n_gk=26;
00714           x_gk=qk51_xgk;
00715           w_g=qk51_wg;
00716           w_gk=qk51_wgk;
00717           break;
00718         case 6:
00719           n_gk=31;
00720           x_gk=qk61_xgk;
00721           w_g=qk61_wg;
00722           w_gk=qk61_wgk;
00723           break;
00724         default:
00725           O2SCL_ERR("Invalid rule.",gsl_einval);
00726           break;
00727         }
00728 
00729         f_v1=new double[n_gk];
00730         f_v2=new double[n_gk];
00731 
00732         return;
00733       }
00734 
00735       /** \brief Set the limit for the number of subdivisions of 
00736           the integration region (default 1000)
00737         
00738           If the value of \c size is zero, the error handler will
00739           be called.
00740       */
00741       int set_limit(size_t lim) {
00742         if (lim==0) {
00743           O2SCL_ERR2("Requested zero limit in ",
00744                      "gsl_inte_kronrod::set_limit().",gsl_einval);
00745         }
00746         w->free();
00747         w->allocate(lim);
00748         return 0;
00749       }
00750 
00751       /** \brief The base Gauss-Kronrod integration function template
00752           
00753           Given abcissas and weights, this performs the integration of
00754           \c func between \c a and \c b, providing a result with
00755           uncertainties.
00756         
00757           The Gauss-Kronrod rule uses \f$ 2m+1 \f$ abscissae \f$ x_1,
00758           \ldots, x_{2m+1} \in (a, b) \f$ to estimate the integral of a
00759           function as a linear combination of values,
00760           \f[ 
00761           \int_a^b f~dx \; \approx \; 
00762           Q_m^{GK}f \; \equiv \;
00763           \sum_{k=1}^{2m+1} \beta_k f(x_k),
00764           \f]
00765           where the weights \f$ \beta_1,\ldots, \beta_{2m+1} \f$ are
00766           intrinsic to the abscissae. The data are designed so that the
00767           even-indexed abscissae yield a coarser estimate,
00768           \f[
00769           Q_m^{G}f \; \equiv \;
00770           \sum_{j=1}^{m} \alpha_j f(x_{2j}),
00771           \f]
00772           and their difference \f$ |Q_m^Gf - Q_m^{GK}f| \f$ is the "raw"
00773           error estimate. The various quantities that the function
00774           computes are
00775           \f{eqnarray*}{
00776           \mathtt{result} &=& Q_m^{GK}f, \\
00777           \mathtt{resabs} &=& Q_m^{GK}|f|, \\
00778           \mathtt{resasc} &=& Q_m^{GK}(|f| - \mu),
00779           \quad \mu \equiv \frac{Q_m^{GK}f}{b - a}.
00780           \f}
00781           The "absolute" error \c abserr is computed from the raw error
00782           value using the function \ref gsl_inte::rescale_error.
00783         
00784           This function is designed for use with the values given in the
00785           o2scl_inte_gk_coeffs namespace.
00786         
00787           This function never calls the error handler.
00788 
00789           \future This function, in principle, is an integration
00790           object in and of itself, and could be implemented separately
00791           as an object of type \ref inte.
00792        */
00793       template<class func2_t> void gauss_kronrod_base
00794         (func2_t &func, double a, double b, double *result, 
00795          double *abserr, double *resabs, double *resasc) {
00796         
00797         const double center=0.5*(a+b);
00798         const double half_length=0.5*(b-a);
00799         const double abs_half_length=fabs (half_length);
00800           
00801         double f_center;
00802         f_center=func(center);
00803           
00804         double result_gauss=0.0;
00805         double result_kronrod=f_center*this->w_gk[this->n_gk-1];
00806  
00807         double result_abs=fabs (result_kronrod);
00808         double result_asc=0.0;
00809         double mean=0.0, err=0.0;
00810       
00811         if (this->n_gk % 2 == 0) {
00812           result_gauss=f_center*this->w_g[this->n_gk / 2-1];
00813         }
00814       
00815         // Sum up results from left half of interval
00816         for (int j=0; j < (this->n_gk-1) / 2; j++) {
00817 
00818           const int jtw=j*2+1;        
00819           const double abscissa=half_length*this->x_gk[jtw];
00820               
00821           double fval1, fval2;
00822           fval1=func(center-abscissa);
00823           fval2=func(center+abscissa);
00824               
00825           const double fsum=fval1+fval2;
00826           this->f_v1[jtw]=fval1;
00827           this->f_v2[jtw]=fval2;
00828           result_gauss+=this->w_g[j]*fsum;
00829           result_kronrod+=this->w_gk[jtw]*fsum;
00830           result_abs+=this->w_gk[jtw]*(fabs(fval1)+fabs(fval2));
00831         }
00832       
00833         // Sum up results from right half of interval
00834         for (int j=0; j < this->n_gk / 2; j++) {
00835 
00836           int jtwm1=j*2;
00837           const double abscissa=half_length*this->x_gk[jtwm1];
00838 
00839           double fval1, fval2;
00840           fval1=func(center-abscissa);
00841           fval2=func(center+abscissa);
00842 
00843           this->f_v1[jtwm1]=fval1;
00844           this->f_v2[jtwm1]=fval2;
00845           result_kronrod+=this->w_gk[jtwm1]*(fval1+fval2);
00846           result_abs+=this->w_gk[jtwm1]*(fabs(fval1)+fabs(fval2));
00847         }
00848       
00849         mean=result_kronrod*0.5;
00850       
00851         result_asc=this->w_gk[this->n_gk-1]*fabs(f_center-mean);
00852         
00853         for (int j=0;j<this->n_gk-1;j++) {
00854           result_asc+=this->w_gk[j]*(fabs(this->f_v1[j]-mean)+
00855                                      fabs(this->f_v2[j]-mean));
00856         }
00857       
00858         /* Scale by the width of the integration region */
00859       
00860         err=(result_kronrod-result_gauss)*half_length;
00861       
00862         result_kronrod*=half_length;
00863         result_abs*=abs_half_length;
00864         result_asc*=abs_half_length;
00865       
00866         *result=result_kronrod;
00867         *resabs=result_abs;
00868         *resasc=result_asc;
00869         *abserr=rescale_error(err,result_abs,result_asc);
00870 
00871         return;
00872       }
00873 
00874       /** \brief Integration wrapper for user-specified function type
00875       */
00876       virtual void gauss_kronrod
00877         (func_t &func, double a, double b,
00878          double *result, double *abserr, double *resabs, double *resasc) {
00879         return gauss_kronrod_base<func_t>
00880           (func,a,b,result,abserr,resabs,resasc);
00881                                           
00882       }
00883     
00884     };
00885   
00886 #ifndef DOXYGENP
00887 }
00888 #endif
00889 
00890 #endif
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines

Documentation generated with Doxygen. Provided under the GNU Free Documentation License (see License Information).

Get Object-oriented Scientific Computing
Lib at SourceForge.net. Fast, secure and Free Open Source software
downloads.