00001 /* 00002 ------------------------------------------------------------------- 00003 00004 Copyright (C) 2006, 2007, 2008, 2009, 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_QAG_B_H 00024 #define O2SCL_GSL_INTE_QAG_B_H 00025 00026 #include <cmath> 00027 #include <gsl/gsl_integration.h> 00028 #include <o2scl/inte.h> 00029 #include <o2scl/gsl_inte.h> 00030 00031 /** 00032 \brief A namespace for the GSL adaptive integration coefficients 00033 00034 In the abcissa arrays, the odd-indexed entries are the abcissae of 00035 the Gauss rule, while the even entries are the Kronrod abcissae 00036 which are added to the Gauss rule. 00037 00038 <b>Documentation from GSL</b>: \n 00039 Gauss quadrature weights and kronrod quadrature abscissae and 00040 weights as evaluated with 80 decimal digit arithmetic by 00041 L. W. Fullerton, Bell Labs, Nov. 1981. 00042 */ 00043 namespace o2scl_inte_qag_coeffs { 00044 00045 /** \brief Abscissae of the 15-point Kronrod rule */ 00046 static const double qk15_xgk[8] = 00047 { 00048 0.991455371120812639206854697526329, 00049 0.949107912342758524526189684047851, 00050 0.864864423359769072789712788640926, 00051 0.741531185599394439863864773280788, 00052 0.586087235467691130294144838258730, 00053 0.405845151377397166906606412076961, 00054 0.207784955007898467600689403773245, 00055 0.000000000000000000000000000000000 00056 }; 00057 00058 /** \brief Weights of the 7-point Gauss rule */ 00059 static const double qk15_wg[4] = 00060 { 00061 0.129484966168869693270611432679082, 00062 0.279705391489276667901467771423780, 00063 0.381830050505118944950369775488975, 00064 0.417959183673469387755102040816327 00065 }; 00066 00067 /** \brief Weights of the 15-point Kronrod rule */ 00068 static const double qk15_wgk[8] = 00069 { 00070 0.022935322010529224963732008058970, 00071 0.063092092629978553290700663189204, 00072 0.104790010322250183839876322541518, 00073 0.140653259715525918745189590510238, 00074 0.169004726639267902826583426598550, 00075 0.190350578064785409913256402421014, 00076 0.204432940075298892414161999234649, 00077 0.209482141084727828012999174891714 00078 }; 00079 00080 /** \brief Abscissae of the 21-point Kronrod rule */ 00081 static const double qk21_xgk[11] = 00082 { 00083 0.995657163025808080735527280689003, 00084 0.973906528517171720077964012084452, 00085 0.930157491355708226001207180059508, 00086 0.865063366688984510732096688423493, 00087 0.780817726586416897063717578345042, 00088 0.679409568299024406234327365114874, 00089 0.562757134668604683339000099272694, 00090 0.433395394129247190799265943165784, 00091 0.294392862701460198131126603103866, 00092 0.148874338981631210884826001129720, 00093 0.000000000000000000000000000000000 00094 }; 00095 00096 /** \brief Weights of the 10-point Gauss rule */ 00097 static const double qk21_wg[5] = 00098 { 00099 0.066671344308688137593568809893332, 00100 0.149451349150580593145776339657697, 00101 0.219086362515982043995534934228163, 00102 0.269266719309996355091226921569469, 00103 0.295524224714752870173892994651338 00104 }; 00105 00106 /** \brief Weights of the 21-point Kronrod rule */ 00107 static const double qk21_wgk[11] = 00108 { 00109 0.011694638867371874278064396062192, 00110 0.032558162307964727478818972459390, 00111 0.054755896574351996031381300244580, 00112 0.075039674810919952767043140916190, 00113 0.093125454583697605535065465083366, 00114 0.109387158802297641899210590325805, 00115 0.123491976262065851077958109831074, 00116 0.134709217311473325928054001771707, 00117 0.142775938577060080797094273138717, 00118 0.147739104901338491374841515972068, 00119 0.149445554002916905664936468389821 00120 }; 00121 00122 /** \brief Abscissae of the 31-point Kronrod rule */ 00123 static const double qk31_xgk[16] = 00124 { 00125 0.998002298693397060285172840152271, 00126 0.987992518020485428489565718586613, 00127 0.967739075679139134257347978784337, 00128 0.937273392400705904307758947710209, 00129 0.897264532344081900882509656454496, 00130 0.848206583410427216200648320774217, 00131 0.790418501442465932967649294817947, 00132 0.724417731360170047416186054613938, 00133 0.650996741297416970533735895313275, 00134 0.570972172608538847537226737253911, 00135 0.485081863640239680693655740232351, 00136 0.394151347077563369897207370981045, 00137 0.299180007153168812166780024266389, 00138 0.201194093997434522300628303394596, 00139 0.101142066918717499027074231447392, 00140 0.000000000000000000000000000000000 00141 }; 00142 00143 /** \brief Weights of the 15-point Gauss rule */ 00144 static const double qk31_wg[8] = 00145 { 00146 0.030753241996117268354628393577204, 00147 0.070366047488108124709267416450667, 00148 0.107159220467171935011869546685869, 00149 0.139570677926154314447804794511028, 00150 0.166269205816993933553200860481209, 00151 0.186161000015562211026800561866423, 00152 0.198431485327111576456118326443839, 00153 0.202578241925561272880620199967519 00154 }; 00155 00156 /** \brief Weights of the 31-point Kronrod rule */ 00157 static const double qk31_wgk[16] = 00158 { 00159 0.005377479872923348987792051430128, 00160 0.015007947329316122538374763075807, 00161 0.025460847326715320186874001019653, 00162 0.035346360791375846222037948478360, 00163 0.044589751324764876608227299373280, 00164 0.053481524690928087265343147239430, 00165 0.062009567800670640285139230960803, 00166 0.069854121318728258709520077099147, 00167 0.076849680757720378894432777482659, 00168 0.083080502823133021038289247286104, 00169 0.088564443056211770647275443693774, 00170 0.093126598170825321225486872747346, 00171 0.096642726983623678505179907627589, 00172 0.099173598721791959332393173484603, 00173 0.100769845523875595044946662617570, 00174 0.101330007014791549017374792767493 00175 }; 00176 00177 /** \brief Abscissae of the 41-point Kronrod rule */ 00178 static const double qk41_xgk[21] = 00179 { 00180 0.998859031588277663838315576545863, 00181 0.993128599185094924786122388471320, 00182 0.981507877450250259193342994720217, 00183 0.963971927277913791267666131197277, 00184 0.940822633831754753519982722212443, 00185 0.912234428251325905867752441203298, 00186 0.878276811252281976077442995113078, 00187 0.839116971822218823394529061701521, 00188 0.795041428837551198350638833272788, 00189 0.746331906460150792614305070355642, 00190 0.693237656334751384805490711845932, 00191 0.636053680726515025452836696226286, 00192 0.575140446819710315342946036586425, 00193 0.510867001950827098004364050955251, 00194 0.443593175238725103199992213492640, 00195 0.373706088715419560672548177024927, 00196 0.301627868114913004320555356858592, 00197 0.227785851141645078080496195368575, 00198 0.152605465240922675505220241022678, 00199 0.076526521133497333754640409398838, 00200 0.000000000000000000000000000000000 00201 }; 00202 00203 /** \brief Weights of the 20-point Gauss rule */ 00204 static const double qk41_wg[11] = 00205 { 00206 0.017614007139152118311861962351853, 00207 0.040601429800386941331039952274932, 00208 0.062672048334109063569506535187042, 00209 0.083276741576704748724758143222046, 00210 0.101930119817240435036750135480350, 00211 0.118194531961518417312377377711382, 00212 0.131688638449176626898494499748163, 00213 0.142096109318382051329298325067165, 00214 0.149172986472603746787828737001969, 00215 0.152753387130725850698084331955098 00216 }; 00217 00218 /** \brief Weights of the 41-point Kronrod rule */ 00219 static const double qk41_wgk[21] = 00220 { 00221 0.003073583718520531501218293246031, 00222 0.008600269855642942198661787950102, 00223 0.014626169256971252983787960308868, 00224 0.020388373461266523598010231432755, 00225 0.025882133604951158834505067096153, 00226 0.031287306777032798958543119323801, 00227 0.036600169758200798030557240707211, 00228 0.041668873327973686263788305936895, 00229 0.046434821867497674720231880926108, 00230 0.050944573923728691932707670050345, 00231 0.055195105348285994744832372419777, 00232 0.059111400880639572374967220648594, 00233 0.062653237554781168025870122174255, 00234 0.065834597133618422111563556969398, 00235 0.068648672928521619345623411885368, 00236 0.071054423553444068305790361723210, 00237 0.073030690332786667495189417658913, 00238 0.074582875400499188986581418362488, 00239 0.075704497684556674659542775376617, 00240 0.076377867672080736705502835038061, 00241 0.076600711917999656445049901530102 00242 }; 00243 00244 /** \brief Abscissae of the 51-point Kronrod rule */ 00245 static const double qk51_xgk[26] = 00246 { 00247 0.999262104992609834193457486540341, 00248 0.995556969790498097908784946893902, 00249 0.988035794534077247637331014577406, 00250 0.976663921459517511498315386479594, 00251 0.961614986425842512418130033660167, 00252 0.942974571228974339414011169658471, 00253 0.920747115281701561746346084546331, 00254 0.894991997878275368851042006782805, 00255 0.865847065293275595448996969588340, 00256 0.833442628760834001421021108693570, 00257 0.797873797998500059410410904994307, 00258 0.759259263037357630577282865204361, 00259 0.717766406813084388186654079773298, 00260 0.673566368473468364485120633247622, 00261 0.626810099010317412788122681624518, 00262 0.577662930241222967723689841612654, 00263 0.526325284334719182599623778158010, 00264 0.473002731445714960522182115009192, 00265 0.417885382193037748851814394594572, 00266 0.361172305809387837735821730127641, 00267 0.303089538931107830167478909980339, 00268 0.243866883720988432045190362797452, 00269 0.183718939421048892015969888759528, 00270 0.122864692610710396387359818808037, 00271 0.061544483005685078886546392366797, 00272 0.000000000000000000000000000000000 00273 }; 00274 00275 /** \brief Weights of the 25-point Gauss rule */ 00276 static const double qk51_wg[13] = 00277 { 00278 0.011393798501026287947902964113235, 00279 0.026354986615032137261901815295299, 00280 0.040939156701306312655623487711646, 00281 0.054904695975835191925936891540473, 00282 0.068038333812356917207187185656708, 00283 0.080140700335001018013234959669111, 00284 0.091028261982963649811497220702892, 00285 0.100535949067050644202206890392686, 00286 0.108519624474263653116093957050117, 00287 0.114858259145711648339325545869556, 00288 0.119455763535784772228178126512901, 00289 0.122242442990310041688959518945852, 00290 0.123176053726715451203902873079050 00291 }; 00292 00293 /** \brief Weights of the 51-point Kronrod rule */ 00294 static const double qk51_wgk[26] = 00295 { 00296 0.001987383892330315926507851882843, 00297 0.005561932135356713758040236901066, 00298 0.009473973386174151607207710523655, 00299 0.013236229195571674813656405846976, 00300 0.016847817709128298231516667536336, 00301 0.020435371145882835456568292235939, 00302 0.024009945606953216220092489164881, 00303 0.027475317587851737802948455517811, 00304 0.030792300167387488891109020215229, 00305 0.034002130274329337836748795229551, 00306 0.037116271483415543560330625367620, 00307 0.040083825504032382074839284467076, 00308 0.042872845020170049476895792439495, 00309 0.045502913049921788909870584752660, 00310 0.047982537138836713906392255756915, 00311 0.050277679080715671963325259433440, 00312 0.052362885806407475864366712137873, 00313 0.054251129888545490144543370459876, 00314 0.055950811220412317308240686382747, 00315 0.057437116361567832853582693939506, 00316 0.058689680022394207961974175856788, 00317 0.059720340324174059979099291932562, 00318 0.060539455376045862945360267517565, 00319 0.061128509717053048305859030416293, 00320 0.061471189871425316661544131965264, 00321 0.061580818067832935078759824240066 00322 }; 00323 00324 /** \brief Abscissae of the 61-point Kronrod rule */ 00325 static const double qk61_xgk[31] = 00326 { 00327 0.999484410050490637571325895705811, 00328 0.996893484074649540271630050918695, 00329 0.991630996870404594858628366109486, 00330 0.983668123279747209970032581605663, 00331 0.973116322501126268374693868423707, 00332 0.960021864968307512216871025581798, 00333 0.944374444748559979415831324037439, 00334 0.926200047429274325879324277080474, 00335 0.905573307699907798546522558925958, 00336 0.882560535792052681543116462530226, 00337 0.857205233546061098958658510658944, 00338 0.829565762382768397442898119732502, 00339 0.799727835821839083013668942322683, 00340 0.767777432104826194917977340974503, 00341 0.733790062453226804726171131369528, 00342 0.697850494793315796932292388026640, 00343 0.660061064126626961370053668149271, 00344 0.620526182989242861140477556431189, 00345 0.579345235826361691756024932172540, 00346 0.536624148142019899264169793311073, 00347 0.492480467861778574993693061207709, 00348 0.447033769538089176780609900322854, 00349 0.400401254830394392535476211542661, 00350 0.352704725530878113471037207089374, 00351 0.304073202273625077372677107199257, 00352 0.254636926167889846439805129817805, 00353 0.204525116682309891438957671002025, 00354 0.153869913608583546963794672743256, 00355 0.102806937966737030147096751318001, 00356 0.051471842555317695833025213166723, 00357 0.000000000000000000000000000000000 00358 }; 00359 00360 /** \brief Weights of the 30-point Gauss rule */ 00361 static const double qk61_wg[15] = 00362 { 00363 0.007968192496166605615465883474674, 00364 0.018466468311090959142302131912047, 00365 0.028784707883323369349719179611292, 00366 0.038799192569627049596801936446348, 00367 0.048402672830594052902938140422808, 00368 0.057493156217619066481721689402056, 00369 0.065974229882180495128128515115962, 00370 0.073755974737705206268243850022191, 00371 0.080755895229420215354694938460530, 00372 0.086899787201082979802387530715126, 00373 0.092122522237786128717632707087619, 00374 0.096368737174644259639468626351810, 00375 0.099593420586795267062780282103569, 00376 0.101762389748405504596428952168554, 00377 0.102852652893558840341285636705415 00378 }; 00379 00380 /** \brief Weights of the 61-point Kronrod rule */ 00381 static const double qk61_wgk[31] = 00382 { 00383 0.001389013698677007624551591226760, 00384 0.003890461127099884051267201844516, 00385 0.006630703915931292173319826369750, 00386 0.009273279659517763428441146892024, 00387 0.011823015253496341742232898853251, 00388 0.014369729507045804812451432443580, 00389 0.016920889189053272627572289420322, 00390 0.019414141193942381173408951050128, 00391 0.021828035821609192297167485738339, 00392 0.024191162078080601365686370725232, 00393 0.026509954882333101610601709335075, 00394 0.028754048765041292843978785354334, 00395 0.030907257562387762472884252943092, 00396 0.032981447057483726031814191016854, 00397 0.034979338028060024137499670731468, 00398 0.036882364651821229223911065617136, 00399 0.038678945624727592950348651532281, 00400 0.040374538951535959111995279752468, 00401 0.041969810215164246147147541285970, 00402 0.043452539701356069316831728117073, 00403 0.044814800133162663192355551616723, 00404 0.046059238271006988116271735559374, 00405 0.047185546569299153945261478181099, 00406 0.048185861757087129140779492298305, 00407 0.049055434555029778887528165367238, 00408 0.049795683427074206357811569379942, 00409 0.050405921402782346840893085653585, 00410 0.050881795898749606492297473049805, 00411 0.051221547849258772170656282604944, 00412 0.051426128537459025933862879215781, 00413 0.051494729429451567558340433647099 00414 }; 00415 00416 } 00417 00418 #ifndef DOXYGENP 00419 namespace o2scl { 00420 #endif 00421 00422 /** 00423 \brief A GSL integration workspace 00424 00425 \note This class is unused at the moment. This is a future 00426 replacement for the GSL integration workspace. 00427 00428 \comment 00429 For reference: 00430 typedef struct 00431 { 00432 size_t limit; 00433 size_t size; 00434 size_t nrmax; 00435 size_t i; 00436 size_t maximum_level; 00437 double *alist; 00438 double *blist; 00439 double *rlist; 00440 double *elist; 00441 size_t *order; 00442 size_t *level; 00443 } 00444 gsl_integration_workspace; 00445 \endcomment 00446 */ 00447 class gsl_inte_workspace : public gsl_integration_workspace { 00448 00449 public: 00450 00451 gsl_inte_workspace(); 00452 00453 int allocate(size_t sz); 00454 00455 int free(); 00456 00457 /** 00458 \brief Initialize the workspace for an integration with limits 00459 \c a and \c b. 00460 */ 00461 int initialise(double a, double b); 00462 00463 /// Set the result at position zero 00464 int set_initial_result(double result, double error); 00465 00466 /** 00467 \brief Retrieve the ith result from the workspace 00468 00469 The workspace variable \c i is used to specify which 00470 interval is requested. 00471 */ 00472 int retrieve(double *a, double *b, double *r, double *e) const; 00473 00474 /** \brief Sort the workspace 00475 00476 This routine maintains the descending ordering in the list of 00477 the local error estimated resulting from the interval 00478 subdivision process. at each call two error estimates are 00479 inserted using the sequential search method, top-down for the 00480 largest error estimate and bottom-up for the smallest error 00481 estimate. 00482 00483 Originally written in QUADPACK by 00484 \verbatim 00485 piessens, robert, appl. math. & progr. div. - k.u.leuven 00486 de doncker, elise, appl. math. & progr. div. - k.u.leuven 00487 \endverbatim 00488 translated into C for GSL by Brian Gough, and then 00489 rewritten for O2scl. 00490 00491 */ 00492 int qpsrt(); 00493 00494 /// Update workspace with new results and resort 00495 int update(double a1, double b1, double area1, double error1, 00496 double a2, double b2, double area2, double error2); 00497 00498 /// Add up all of the contributions to construct the final result 00499 double sum_results(); 00500 00501 /// Find out if the present subinterval is too small 00502 int subinterval_too_small(double a1, double a2, double b2); 00503 00504 /// Append new results to workspace 00505 int append_interval(double a1, double b1, double area1, double error1); 00506 00507 }; 00508 00509 /** 00510 \brief Base routines for the GSL adaptive integration classes 00511 00512 This class contains several functions for manipulating the GSL 00513 integration workspace. The casual end-user should use the 00514 classes explained in the \ref inte_section section of the User's 00515 guide. 00516 00517 \todo Make the workspace size protected 00518 00519 \future Move gsl_integration_workspace to a separate class 00520 and remove this class, making all children direct descendants 00521 of gsl_inte instead. We'll have to figure out what to 00522 do with the data member \c wkspace though. Some work on this 00523 front is already in gsl_inte_qag_b.h. 00524 00525 QUADPACK workspace documentation: 00526 \verbatim 00527 c parameters (meaning at output) 00528 c limit - integer 00529 c maximum number of error estimates the list 00530 c can contain 00531 c last - integer 00532 c number of error estimates currently in the list 00533 c maxerr - integer 00534 c maxerr points to the nrmax-th largest error 00535 c estimate currently in the list 00536 c ermax - double precision 00537 c nrmax-th largest error estimate 00538 c ermax = elist(maxerr) 00539 c elist - double precision 00540 c vector of dimension last containing 00541 c the error estimates 00542 c iord - integer 00543 c vector of dimension last, the first k elements 00544 c of which contain pointers to the error 00545 c estimates, such that 00546 c elist(iord(1)),..., elist(iord(k)) 00547 c form a decreasing sequence, with 00548 c k = last if last.le.(limit/2+2), and 00549 c k = limit+1-last otherwise 00550 c nrmax - integer 00551 c maxerr = iord(nrmax) 00552 \endverbatim 00553 \verbatim 00554 c alist - real 00555 c vector of dimension at least limit, the first 00556 c last elements of which are the left 00557 c end points of the subintervals in the partition 00558 c of the given integration range (a,b) 00559 c blist - real 00560 c vector of dimension at least limit, the first 00561 c last elements of which are the right 00562 c end points of the subintervals in the partition 00563 c of the given integration range (a,b) 00564 c rlist - real 00565 c vector of dimension at least limit, the first 00566 c last elements of which are the 00567 c integral approximations on the subintervals 00568 c elist - real 00569 c vector of dimension at least limit, the first 00570 c last elements of which are the moduli of the 00571 c absolute error estimates on the subintervals 00572 c iord - integer 00573 c vector of dimension at least limit, the first k 00574 c elements of which are pointers to the 00575 c error estimates over the subintervals, 00576 c such that elist(iord(1)), ..., 00577 c elist(iord(k)) form a decreasing sequence, 00578 c with k = last if last.le.(limit/2+2), and 00579 c k = limit+1-last otherwise 00580 c last - integer 00581 c number of subintervals actually produced in the 00582 c subdivision process 00583 \endverbatim 00584 00585 */ 00586 class gsl_inte_table : public gsl_inte { 00587 00588 public: 00589 00590 gsl_inte_table(); 00591 00592 ~gsl_inte_table(); 00593 00594 /// The integration workspace 00595 gsl_integration_workspace *w; 00596 00597 /// The size of the integration workspace (default 1000) 00598 int wkspace; 00599 00600 /// Set the integration workspace size 00601 int set_wkspace(size_t size); 00602 00603 /** 00604 \brief Initialize the workspace for an integration with limits 00605 \c a and \c b. 00606 */ 00607 void initialise(gsl_integration_workspace *workspace, 00608 double a, double b); 00609 00610 /// Set the result at position zero 00611 void set_initial_result(gsl_integration_workspace *workspace, 00612 double result, double error); 00613 00614 /** 00615 \brief Retrieve the ith result from the workspace 00616 00617 The workspace variable \c i is used to specify which 00618 interval is requested. 00619 */ 00620 void retrieve(const gsl_integration_workspace *workspace, 00621 double *a, double *b, double *r, double *e); 00622 00623 /** \brief Sort the workspace 00624 00625 This routine maintains the descending ordering in the list of 00626 the local error estimated resulting from the interval 00627 subdivision process. at each call two error estimates are 00628 inserted using the sequential search method, top-down for the 00629 largest error estimate and bottom-up for the smallest error 00630 estimate. 00631 00632 Originally written in QUADPACK by 00633 \verbatim 00634 piessens, robert, appl. math. & progr. div. - k.u.leuven 00635 de doncker, elise, appl. math. & progr. div. - k.u.leuven 00636 \endverbatim 00637 translated into C for GSL by Brian Gough, and then 00638 rewritten for O2scl. 00639 00640 */ 00641 void qpsrt(gsl_integration_workspace *workspace); 00642 00643 /// Update workspace with new results and resort 00644 void update(gsl_integration_workspace *workspace, 00645 double a1, double b1, double area1, double error1, 00646 double a2, double b2, double area2, double error2); 00647 00648 /// Add up all of the contributions to construct the final result 00649 double sum_results(const gsl_integration_workspace *workspace); 00650 00651 /// Find out if the present subinterval is too small 00652 int subinterval_too_small(double a1, double a2, double b2); 00653 00654 /// Append new results to workspace 00655 void append_interval(gsl_integration_workspace *workspace, 00656 double a1, double b1, double area1, double error1); 00657 00658 }; 00659 00660 /** 00661 \brief Basic Gauss-Kronrod integration class (GSL) 00662 00663 This class provides the basic Gauss-Kronrod integration function 00664 for some of the GSL-based integration classes. The casual 00665 end-user should use the classes explained in the \ref 00666 inte_section section of the User's guide. 00667 */ 00668 template<class param_t, class func_t> 00669 class gsl_inte_kronrod : public gsl_inte_table, 00670 public inte<param_t,func_t> { 00671 00672 public: 00673 00674 /** 00675 \brief The GSL Gauss-Kronrod integration function 00676 00677 Given abcissas and weights, this performs the integration of 00678 \c func between \c a and \c b, providing a result with 00679 uncertainties. 00680 00681 This function is designed for use with the values given in the 00682 o2scl_inte_qag_coeffs namespace. 00683 00684 This function never calls the error handler. 00685 00686 \future This function, in principle, could be replaced with 00687 a generic integration pointer. 00688 */ 00689 virtual void gsl_integration_qk_o2scl 00690 (func_t &func, const int n, const double xgk[], const double wg[], 00691 const double wgk[], double fv1[], double fv2[], double a, double b, 00692 double *result, double *abserr, double *resabs, double *resasc, 00693 param_t &pa) { 00694 00695 const double center = 0.5 * (a + b); 00696 const double half_length = 0.5 * (b - a); 00697 const double abs_half_length = fabs (half_length); 00698 00699 double f_center; 00700 func(center,f_center,pa); 00701 00702 double result_gauss = 0; 00703 double result_kronrod = f_center * wgk[n - 1]; 00704 00705 double result_abs = fabs (result_kronrod); 00706 double result_asc = 0; 00707 double mean = 0, err = 0; 00708 00709 int j; 00710 00711 if (n % 2 == 0) { 00712 result_gauss = f_center * wg[n / 2 - 1]; 00713 } 00714 00715 // Sum up results from left half of interval 00716 for (j = 0; j < (n - 1) / 2; j++) { 00717 00718 const int jtw = j * 2 + 1; 00719 const double abscissa = half_length * xgk[jtw]; 00720 00721 double fval1, fval2; 00722 func(center-abscissa,fval1,pa); 00723 func(center+abscissa,fval2,pa); 00724 00725 const double fsum = fval1 + fval2; 00726 fv1[jtw] = fval1; 00727 fv2[jtw] = fval2; 00728 result_gauss += wg[j] * fsum; 00729 result_kronrod += wgk[jtw] * fsum; 00730 result_abs += wgk[jtw] * (fabs (fval1) + fabs (fval2)); 00731 } 00732 00733 // Sum up results from right half of interval 00734 for (j = 0; j < n / 2; j++) { 00735 00736 int jtwm1 = j * 2; 00737 const double abscissa = half_length * xgk[jtwm1]; 00738 00739 double fval1, fval2; 00740 func(center-abscissa,fval1,pa); 00741 func(center+abscissa,fval2,pa); 00742 00743 fv1[jtwm1] = fval1; 00744 fv2[jtwm1] = fval2; 00745 result_kronrod += wgk[jtwm1] * (fval1 + fval2); 00746 result_abs += wgk[jtwm1] * (fabs (fval1) + fabs (fval2)); 00747 } 00748 00749 mean = result_kronrod * 0.5; 00750 00751 result_asc = wgk[n - 1] * fabs (f_center - mean); 00752 00753 for (j = 0; j < n - 1; j++) { 00754 result_asc += wgk[j] * (fabs (fv1[j] - mean) + 00755 fabs (fv2[j] - mean)); 00756 } 00757 00758 /* Scale by the width of the integration region */ 00759 00760 err = (result_kronrod - result_gauss) * half_length; 00761 00762 result_kronrod *= half_length; 00763 result_abs *= abs_half_length; 00764 result_asc *= abs_half_length; 00765 00766 *result = result_kronrod; 00767 *resabs = result_abs; 00768 *resasc = result_asc; 00769 *abserr = rescale_error (err, result_abs, result_asc); 00770 00771 return; 00772 } 00773 00774 }; 00775 00776 /** \brief Base class for integrating a function with a 00777 singularity (GSL) 00778 00779 This class contains the extrapolation table mechanics and the 00780 base integration function for singular integrals from GSL. The 00781 casual end-user should use the classes explained in the 00782 \ref inte_section section of the User's guide. 00783 00784 \future Some of the functions inside this class could 00785 be moved out of header files? 00786 */ 00787 template<class param_t, class func_t> class gsl_inte_singular : 00788 public gsl_inte_kronrod<param_t,func_t> { 00789 00790 protected: 00791 00792 /** 00793 \brief A structure for extrapolation for \ref gsl_inte_qags 00794 00795 \future Move this to a new class, with qelg() as a method 00796 */ 00797 struct extrapolation_table 00798 { 00799 /// Index of new element in the first column 00800 size_t n; 00801 /// Lower diagonals of the triangular epsilon table 00802 double rlist2[52]; 00803 /// Number of calls 00804 size_t nres; 00805 /// Three most recent results 00806 double res3la[3]; 00807 }; 00808 00809 /// Initialize the table 00810 void initialise_table(struct extrapolation_table *table) { 00811 table->n = 0; 00812 table->nres = 0; 00813 return; 00814 } 00815 00816 /// Append a result to the table 00817 void append_table(struct extrapolation_table *table, double y) { 00818 size_t n; 00819 n = table->n; 00820 table->rlist2[n] = y; 00821 table->n++; 00822 return; 00823 } 00824 00825 /// Test a result for positivity 00826 inline int test_positivity(double result, double resabs) { 00827 int status = (fabs (result) >= (1 - 50 * GSL_DBL_EPSILON) * resabs); 00828 return status; 00829 } 00830 00831 /** \brief Determines the limit of a given sequence of approximations 00832 00833 This function determines the limit of a given sequence of 00834 approximations, by means of the epsilon algorithm of P. Wynn. 00835 an estimate of the absolute error is also given. the condensed 00836 epsilon table is computed. only those elements needed for the 00837 computation of the next diagonal are preserved. 00838 00839 Quadpack documentation 00840 \verbatim 00841 c 00842 c list of major variables 00843 c ----------------------- 00844 c e0 - the 4 elements on which the computation of a new 00845 c e1 element in the epsilon table is based 00846 c e2 00847 c e3 e0 00848 c e3 e1 new 00849 c e2 00850 c newelm - number of elements to be computed in the new 00851 c diagonal 00852 c error - error = abs(e1-e0)+abs(e2-e1)+abs(new-e2) 00853 c result - the element in the new diagonal with least value 00854 c of error 00855 c 00856 c machine dependent constants 00857 c --------------------------- 00858 c 00859 c epmach is the largest relative spacing. 00860 c oflow is the largest positive magnitude. 00861 c limexp is the maximum number of elements the epsilon 00862 c table can contain. if this number is reached, the upper 00863 c diagonal of the epsilon table is deleted. 00864 c 00865 \endverbatim 00866 */ 00867 void qelg(struct extrapolation_table *table, double *result, 00868 double *abserr) { 00869 00870 double *epstab = table->rlist2; 00871 double *res3la = table->res3la; 00872 const size_t n = table->n - 1; 00873 00874 const double current = epstab[n]; 00875 00876 double absolute = GSL_DBL_MAX; 00877 double relative = 5 * GSL_DBL_EPSILON * fabs (current); 00878 00879 const size_t newelm = n / 2; 00880 const size_t n_orig = n; 00881 size_t n_final = n; 00882 size_t i; 00883 00884 const size_t nres_orig = table->nres; 00885 00886 *result = current; 00887 *abserr = GSL_DBL_MAX; 00888 00889 if (n < 2) { 00890 *result = current; 00891 *abserr = GSL_MAX_DBL (absolute, relative); 00892 return; 00893 } 00894 00895 epstab[n + 2] = epstab[n]; 00896 epstab[n] = GSL_DBL_MAX; 00897 00898 for (i = 0; i < newelm; i++) { 00899 double res = epstab[n - 2 * i + 2]; 00900 double e0 = epstab[n - 2 * i - 2]; 00901 double e1 = epstab[n - 2 * i - 1]; 00902 double e2 = res; 00903 00904 double e1abs = fabs (e1); 00905 double delta2 = e2 - e1; 00906 double err2 = fabs (delta2); 00907 double tol2 = GSL_MAX_DBL (fabs (e2), e1abs) * GSL_DBL_EPSILON; 00908 double delta3 = e1 - e0; 00909 double err3 = fabs (delta3); 00910 double tol3 = GSL_MAX_DBL (e1abs, fabs (e0)) * GSL_DBL_EPSILON; 00911 00912 double e3, delta1, err1, tol1, ss; 00913 00914 if (err2 <= tol2 && err3 <= tol3) 00915 { 00916 /* If e0, e1 and e2 are equal to within machine accuracy, 00917 convergence is assumed. */ 00918 00919 *result = res; 00920 absolute = err2 + err3; 00921 relative = 5 * GSL_DBL_EPSILON * fabs (res); 00922 *abserr = GSL_MAX_DBL (absolute, relative); 00923 return; 00924 } 00925 00926 e3 = epstab[n - 2 * i]; 00927 epstab[n - 2 * i] = e1; 00928 delta1 = e1 - e3; 00929 err1 = fabs (delta1); 00930 tol1 = GSL_MAX_DBL (e1abs, fabs (e3)) * GSL_DBL_EPSILON; 00931 00932 /* If two elements are very close to each other, omit a part of 00933 the table by adjusting the value of n */ 00934 00935 if (err1 <= tol1 || err2 <= tol2 || err3 <= tol3) { 00936 n_final = 2 * i; 00937 break; 00938 } 00939 00940 ss = (1 / delta1 + 1 / delta2) - 1 / delta3; 00941 00942 /* Test to detect irregular behaviour in the table, and 00943 eventually omit a part of the table by adjusting the value of 00944 n. */ 00945 if (fabs (ss * e1) <= 0.0001) { 00946 n_final = 2 * i; 00947 break; 00948 } 00949 /* Compute a new element and eventually adjust the value of 00950 result. */ 00951 00952 res = e1 + 1 / ss; 00953 epstab[n - 2 * i] = res; 00954 00955 { 00956 const double error = err2 + fabs (res - e2) + err3; 00957 00958 if (error <= *abserr) { 00959 *abserr = error; 00960 *result = res; 00961 } 00962 } 00963 } 00964 00965 /* Shift the table */ 00966 00967 { 00968 const size_t limexp = 50 - 1; 00969 00970 if (n_final == limexp) 00971 { 00972 n_final = 2 * (limexp / 2); 00973 } 00974 } 00975 00976 if (n_orig % 2 == 1) { 00977 for (i = 0; i <= newelm; i++) { 00978 epstab[1 + i * 2] = epstab[i * 2 + 3]; 00979 } 00980 } else { 00981 for (i = 0; i <= newelm; i++) { 00982 epstab[i * 2] = epstab[i * 2 + 2]; 00983 } 00984 } 00985 if (n_orig != n_final) { 00986 for (i = 0; i <= n_final; i++) { 00987 epstab[i] = epstab[n_orig - n_final + i]; 00988 } 00989 } 00990 00991 table->n = n_final + 1; 00992 00993 if (nres_orig < 3) { 00994 00995 res3la[nres_orig] = *result; 00996 *abserr = GSL_DBL_MAX; 00997 00998 } else { 00999 /* Compute error estimate */ 01000 *abserr = (fabs (*result - res3la[2]) + fabs (*result - res3la[1]) 01001 + fabs (*result - res3la[0])); 01002 01003 res3la[0] = res3la[1]; 01004 res3la[1] = res3la[2]; 01005 res3la[2] = *result; 01006 } 01007 01008 /* In QUADPACK the variable table->nres is incremented at the top of 01009 qelg, so it increases on every call. This leads to the array 01010 res3la being accessed when its elements are still undefined, so I 01011 have moved the update to this point so that its value more 01012 useful. */ 01013 01014 table->nres = nres_orig + 1; 01015 01016 *abserr = GSL_MAX_DBL (*abserr, 5 * GSL_DBL_EPSILON * fabs (*result)); 01017 01018 return; 01019 } 01020 01021 /// Determine if an interval is large 01022 int large_interval (gsl_integration_workspace * workspace) { 01023 size_t i = workspace->i ; 01024 const size_t * level = workspace->level; 01025 01026 if (level[i] < workspace->maximum_level) { 01027 return 1; 01028 } else { 01029 return 0; 01030 } 01031 } 01032 01033 /// Reset workspace to work on the interval with the largest error 01034 inline void reset_nrmax (gsl_integration_workspace * workspace) { 01035 workspace->nrmax = 0; 01036 workspace->i = workspace->order[0] ; 01037 } 01038 01039 /// Increase workspace 01040 int increase_nrmax (gsl_integration_workspace * workspace) { 01041 int k; 01042 int id = workspace->nrmax; 01043 int jupbnd; 01044 01045 const size_t * level = workspace->level; 01046 const size_t * order = workspace->order; 01047 01048 size_t limit = workspace->limit ; 01049 size_t last = workspace->size - 1 ; 01050 01051 if (last > (1 + limit / 2)) { 01052 jupbnd = limit + 1 - last; 01053 } else { 01054 jupbnd = last; 01055 } 01056 01057 for (k = id; k <= jupbnd; k++) { 01058 size_t i_max = order[workspace->nrmax]; 01059 01060 workspace->i = i_max ; 01061 01062 if (level[i_max] < workspace->maximum_level) { 01063 return 1; 01064 } 01065 01066 workspace->nrmax++; 01067 01068 } 01069 return 0; 01070 } 01071 01072 /** 01073 \brief Integration function 01074 01075 \future Remove goto statements. 01076 */ 01077 int qags(func_t &func, const int qn, const double xgk[], 01078 const double wg[], const double wgk[], double fv1[], 01079 double fv2[], const double a, const double b, 01080 const double l_epsabs, const double l_epsrel, 01081 const size_t limit, double *result, double *abserr, 01082 param_t &pa) { 01083 01084 double area, errsum; 01085 double res_ext, err_ext; 01086 double result0, abserr0, resabs0, resasc0; 01087 double tolerance; 01088 01089 double ertest = 0; 01090 double error_over_large_intervals = 0; 01091 double reseps = 0, abseps = 0, correc = 0; 01092 size_t ktmin = 0; 01093 int roundoff_type1 = 0, roundoff_type2 = 0, roundoff_type3 = 0; 01094 int error_type = 0, error_type2 = 0; 01095 01096 size_t iteration = 0; 01097 01098 int positive_integrand = 0; 01099 int extrapolate = 0; 01100 int disallow_extrapolation = 0; 01101 01102 struct extrapolation_table table; 01103 01104 /* Initialize results */ 01105 01106 initialise(this->w, a, b); 01107 01108 *result = 0; 01109 *abserr = 0; 01110 01111 if (limit > this->w->limit) { 01112 this->last_iter=0; 01113 std::string estr="Iteration limit exceeds workspace "; 01114 estr+="in gsl_inte_kronrod::qags()."; 01115 O2SCL_ERR_RET(estr.c_str(),gsl_einval); 01116 } 01117 01118 /* Test on accuracy */ 01119 if (this->tolx <= 0 && (this->tolf < 50 * GSL_DBL_EPSILON || 01120 this->tolf < 0.5e-28)) { 01121 this->last_iter=0; 01122 01123 std::string estr="Tolerance cannot be achieved with given "; 01124 estr+="value of 'tolx' and 'tolf' in gsl_inte_kronrod::qags()."; 01125 O2SCL_ERR_RET(estr.c_str(),gsl_ebadtol); 01126 } 01127 01128 /* Perform the first integration */ 01129 01130 gsl_integration_qk_o2scl(func,qn,xgk,wg,wgk,fv1,fv2,a,b, 01131 &result0,&abserr0,&resabs0,&resasc0,pa); 01132 01133 set_initial_result (this->w, result0, abserr0); 01134 01135 tolerance = GSL_MAX_DBL (this->tolx, this->tolf * fabs (result0)); 01136 01137 if (abserr0 <= 100 * GSL_DBL_EPSILON * resabs0 && 01138 abserr0 > tolerance) { 01139 01140 *result = result0; 01141 *abserr = abserr0; 01142 01143 this->last_iter=1; 01144 01145 this->last_conv=gsl_eround; 01146 std::string estr="Cannot reach tolerance because of roundoff error "; 01147 estr+="on first attempt in gsl_inte_qag::qag()."; 01148 O2SCL_CONV_RET(estr.c_str(),gsl_eround,this->err_nonconv); 01149 01150 } else if ((abserr0 <= tolerance && 01151 abserr0 != resasc0) || abserr0 == 0.0) { 01152 01153 *result = result0; 01154 *abserr = abserr0; 01155 this->last_iter=1; 01156 return gsl_success; 01157 01158 } else if (limit == 1) { 01159 01160 *result = result0; 01161 *abserr = abserr0; 01162 01163 this->last_iter=1; 01164 this->last_conv=gsl_emaxiter; 01165 O2SCL_CONV2_RET("A maximum of 1 iteration was insufficient ", 01166 "in gsl_inte_qag::qag().", 01167 gsl_emaxiter,this->err_nonconv); 01168 } 01169 01170 /* Initialization */ 01171 01172 initialise_table (&table); 01173 append_table (&table, result0); 01174 01175 area = result0; 01176 errsum = abserr0; 01177 01178 res_ext = result0; 01179 err_ext = GSL_DBL_MAX; 01180 01181 positive_integrand = this->test_positivity (result0, resabs0); 01182 01183 iteration = 1; 01184 01185 do { 01186 size_t current_level; 01187 double a1, b1, a2, b2; 01188 double a_i, b_i, r_i, e_i; 01189 double area1 = 0, area2 = 0, area12 = 0; 01190 double error1 = 0, error2 = 0, error12 = 0; 01191 double resasc1, resasc2; 01192 double resabs1, resabs2; 01193 double last_e_i; 01194 01195 /* Bisect the subinterval with the largest error estimate */ 01196 01197 retrieve (this->w, &a_i, &b_i, &r_i, &e_i); 01198 01199 current_level = this->w->level[this->w->i] + 1; 01200 01201 a1 = a_i; 01202 b1 = 0.5 * (a_i + b_i); 01203 a2 = b1; 01204 b2 = b_i; 01205 01206 iteration++; 01207 01208 gsl_integration_qk_o2scl(func,qn,xgk,wg,wgk,fv1,fv2,a1,b1, 01209 &area1,&error1,&resabs1,&resasc1,pa); 01210 gsl_integration_qk_o2scl(func,qn,xgk,wg,wgk,fv1,fv2,a2,b2, 01211 &area2,&error2,&resabs2,&resasc2,pa); 01212 01213 area12 = area1 + area2; 01214 error12 = error1 + error2; 01215 last_e_i = e_i; 01216 01217 /* Improve previous approximations to the integral and test for 01218 accuracy. 01219 01220 We write these expressions in the same way as the original 01221 QUADPACK code so that the rounding errors are the same, which 01222 makes testing easier. 01223 */ 01224 01225 errsum = errsum + error12 - e_i; 01226 area = area + area12 - r_i; 01227 01228 tolerance = GSL_MAX_DBL (this->tolx, this->tolf * fabs (area)); 01229 if (resasc1 != error1 && resasc2 != error2) { 01230 double delta = r_i - area12; 01231 01232 if (fabs (delta) <= 1.0e-5 * fabs (area12) && 01233 error12 >= 0.99 * e_i) { 01234 if (!extrapolate) { 01235 roundoff_type1++; 01236 } else { 01237 roundoff_type2++; 01238 } 01239 } 01240 if (iteration > 10 && error12 > e_i) { 01241 roundoff_type3++; 01242 } 01243 } 01244 01245 /* Test for roundoff and eventually set error flag */ 01246 01247 if (roundoff_type1 + roundoff_type2 >= 10 || 01248 roundoff_type3 >= 20) { 01249 error_type = 2; /* round off error */ 01250 } 01251 01252 if (roundoff_type2 >= 5) { 01253 error_type2 = 1; 01254 } 01255 01256 /* set error flag in the case of bad integrand behaviour at 01257 a point of the integration range */ 01258 01259 if (this->subinterval_too_small (a1, a2, b2)) { 01260 error_type = 4; 01261 } 01262 01263 /* append the newly-created intervals to the list */ 01264 01265 update (this->w, a1, b1, area1, error1, a2, b2, area2, error2); 01266 01267 if (errsum <= tolerance) { 01268 goto compute_result; 01269 } 01270 01271 if (error_type) { 01272 break; 01273 } 01274 01275 if (iteration >= limit - 1) { 01276 error_type = 1; 01277 break; 01278 } 01279 01280 if (iteration == 2) { 01281 error_over_large_intervals = errsum; 01282 ertest = tolerance; 01283 append_table (&table, area); 01284 continue; 01285 } 01286 01287 if (disallow_extrapolation) { 01288 continue; 01289 } 01290 01291 error_over_large_intervals += -last_e_i; 01292 01293 if (current_level < this->w->maximum_level) { 01294 error_over_large_intervals += error12; 01295 } 01296 01297 if (!extrapolate) { 01298 01299 /* test whether the interval to be bisected next is the 01300 smallest interval. */ 01301 01302 if (large_interval (this->w)) { 01303 continue; 01304 } 01305 01306 extrapolate = 1; 01307 this->w->nrmax = 1; 01308 } 01309 01310 if (!error_type2 && error_over_large_intervals > ertest) { 01311 if (increase_nrmax (this->w)) { 01312 continue; 01313 } 01314 } 01315 01316 /* Perform extrapolation */ 01317 01318 append_table (&table, area); 01319 01320 qelg (&table, &reseps, &abseps); 01321 01322 ktmin++; 01323 01324 if (ktmin > 5 && err_ext < 0.001 * errsum) { 01325 error_type = 5; 01326 } 01327 01328 if (abseps < err_ext) { 01329 ktmin = 0; 01330 err_ext = abseps; 01331 res_ext = reseps; 01332 correc = error_over_large_intervals; 01333 ertest = GSL_MAX_DBL (this->tolx, 01334 this->tolf * fabs (reseps)); 01335 if (err_ext <= ertest) 01336 break; 01337 } 01338 01339 /* Prepare bisection of the smallest interval. */ 01340 01341 if (table.n == 1) { 01342 disallow_extrapolation = 1; 01343 } 01344 01345 if (error_type == 5) { 01346 break; 01347 } 01348 01349 /* work on interval with largest error */ 01350 01351 reset_nrmax (this->w); 01352 extrapolate = 0; 01353 error_over_large_intervals = errsum; 01354 01355 /// Output iteration information 01356 01357 if (this->verbose>0) { 01358 std::cout << this->type(); 01359 std::cout << " Iter: " << iteration; 01360 std::cout.setf(std::ios::showpos); 01361 std::cout << " Res: " << area; 01362 std::cout.unsetf(std::ios::showpos); 01363 std::cout << " Err: " << errsum 01364 << " Tol: " << tolerance << std::endl; 01365 if (this->verbose>1) { 01366 char ch; 01367 std::cout << "Press a key and type enter to continue. " ; 01368 std::cin >> ch; 01369 } 01370 } 01371 01372 } while (iteration < limit); 01373 01374 *result = res_ext; 01375 *abserr = err_ext; 01376 01377 if (err_ext == GSL_DBL_MAX) 01378 goto compute_result; 01379 01380 if (error_type || error_type2) { 01381 if (error_type2) { 01382 err_ext += correc; 01383 } 01384 01385 if (error_type == 0) 01386 error_type = 3; 01387 01388 if (res_ext != 0.0 && area != 0.0) { 01389 if (err_ext / fabs (res_ext) > errsum / fabs (area)) 01390 goto compute_result; 01391 } else if (err_ext > errsum) { 01392 goto compute_result; 01393 } else if (area == 0.0) { 01394 goto return_error; 01395 } 01396 } 01397 01398 /* Test on divergence. */ 01399 01400 { 01401 double max_area = GSL_MAX_DBL (fabs (res_ext), fabs (area)); 01402 01403 if (!positive_integrand && max_area < 0.01 * resabs0) 01404 goto return_error; 01405 } 01406 01407 { 01408 double ratio = res_ext / area; 01409 01410 if (ratio < 0.01 || ratio > 100.0 || errsum > fabs (area)) 01411 error_type = 6; 01412 } 01413 01414 goto return_error; 01415 01416 compute_result: 01417 01418 *result = sum_results (this->w); 01419 *abserr = errsum; 01420 01421 return_error: 01422 01423 if (error_type > 2) { 01424 error_type--; 01425 } 01426 01427 this->last_iter=iteration; 01428 01429 if (error_type == 0) { 01430 return gsl_success; 01431 } else if (error_type == 1) { 01432 this->last_conv=gsl_emaxiter; 01433 std::string estr="Number of iterations was insufficient "; 01434 estr+=" in gsl_inte_kronrod::qags()."; 01435 O2SCL_CONV_RET(estr.c_str(),gsl_emaxiter,this->err_nonconv); 01436 } else if (error_type == 2) { 01437 this->last_conv=gsl_eround; 01438 std::string estr="Roundoff error prevents tolerance "; 01439 estr+="from being achieved in gsl_inte_kronrod::qags()."; 01440 O2SCL_CONV_RET(estr.c_str(),gsl_eround,this->err_nonconv); 01441 } else if (error_type == 3) { 01442 this->last_conv=gsl_esing; 01443 std::string estr="Bad integrand behavior "; 01444 estr+="in gsl_inte_kronrod::qags()."; 01445 O2SCL_CONV_RET(estr.c_str(),gsl_esing,this->err_nonconv); 01446 } else if (error_type == 4) { 01447 this->last_conv=gsl_eround; 01448 std::string estr="Roundoff error detected in extrapolation table "; 01449 estr+="in gsl_inte_kronrod::qags()."; 01450 O2SCL_CONV_RET(estr.c_str(),gsl_eround,this->err_nonconv); 01451 } else if (error_type == 5) { 01452 this->last_conv=gsl_ediverge; 01453 std::string estr="Integral is divergent or slowly convergent "; 01454 estr+="in gsl_inte_kronrod::qags()."; 01455 O2SCL_CONV_RET(estr.c_str(),gsl_ediverge,this->err_nonconv); 01456 } 01457 01458 std::string estr="Could not integrate function in gsl_inte_kronrod"; 01459 estr+="::qags() (it may have returned a non-finite result)."; 01460 O2SCL_ERR(estr.c_str(),gsl_efailed); 01461 01462 return gsl_efailed; 01463 } 01464 01465 }; 01466 01467 /** 01468 \brief Integrate a function with a singularity (GSL) [abstract base] 01469 01470 This class contains the GSL-based integration function for 01471 applying transformations to the user-defined integrand. The 01472 casual end-user should use the classes explained in the 01473 \ref inte_section section of the User's guide. 01474 */ 01475 template<class param_t, class func_t> class gsl_inte_transform : 01476 public gsl_inte_singular<param_t,func_t> { 01477 01478 public: 01479 01480 /// The transformation to apply to the user-supplied function 01481 virtual double transform(func_t &func, double t, param_t &pa)=0; 01482 01483 /** 01484 \brief The basic Gauss-Kronrod integration function 01485 01486 This function is basically just a copy of 01487 gsl_inte_qag::gsl_integration_qk_o2scl() which is rewritten to 01488 call the internal transformed function rather than directly 01489 calling the user-specified function. 01490 01491 This function never calls the error handler. 01492 */ 01493 virtual void gsl_integration_qk_o2scl 01494 (func_t &func, const int n, const double xgk[], const double wg[], 01495 const double wgk[], double fv1[], double fv2[], double a, double b, 01496 double *result, double *abserr, double *resabs, double *resasc, 01497 param_t &pa) { 01498 01499 const double center = 0.5 * (a + b); 01500 const double half_length = 0.5 * (b - a); 01501 const double abs_half_length = fabs (half_length); 01502 01503 const double f_center=transform(func,center,pa); 01504 01505 double result_gauss = 0; 01506 double result_kronrod = f_center * wgk[n - 1]; 01507 01508 double result_abs = fabs (result_kronrod); 01509 double result_asc = 0; 01510 double mean = 0, err = 0; 01511 01512 int j; 01513 01514 if (n % 2 == 0) { 01515 result_gauss = f_center * wg[n / 2 - 1]; 01516 } 01517 01518 // Sum up results from left half of interval 01519 for (j = 0; j < (n - 1) / 2; j++) { 01520 01521 const int jtw = j * 2 + 1; 01522 const double abscissa = half_length * xgk[jtw]; 01523 01524 const double fval1=transform(func,center-abscissa,pa); 01525 const double fval2=transform(func,center+abscissa,pa); 01526 const double fsum = fval1 + fval2; 01527 01528 fv1[jtw] = fval1; 01529 fv2[jtw] = fval2; 01530 result_gauss += wg[j] * fsum; 01531 result_kronrod += wgk[jtw] * fsum; 01532 result_abs += wgk[jtw] * (fabs (fval1) + fabs (fval2)); 01533 } 01534 01535 // Sum up results from right half of interval 01536 for (j = 0; j < n / 2; j++) { 01537 int jtwm1 = j * 2; 01538 const double abscissa = half_length * xgk[jtwm1]; 01539 01540 const double fval1=transform(func,center-abscissa,pa); 01541 const double fval2=transform(func,center+abscissa,pa); 01542 01543 fv1[jtwm1] = fval1; 01544 fv2[jtwm1] = fval2; 01545 result_kronrod += wgk[jtwm1] * (fval1 + fval2); 01546 result_abs += wgk[jtwm1] * (fabs (fval1) + fabs (fval2)); 01547 }; 01548 01549 mean = result_kronrod * 0.5; 01550 01551 result_asc = wgk[n - 1] * fabs (f_center - mean); 01552 01553 for (j = 0; j < n - 1; j++) { 01554 result_asc += wgk[j] * (fabs (fv1[j] - mean) + 01555 fabs (fv2[j] - mean)); 01556 } 01557 01558 // Scale by the width of the integration region 01559 01560 err = (result_kronrod - result_gauss) * half_length; 01561 01562 result_kronrod *= half_length; 01563 result_abs *= abs_half_length; 01564 result_asc *= abs_half_length; 01565 01566 *result = result_kronrod; 01567 *resabs = result_abs; 01568 *resasc = result_asc; 01569 *abserr = this->rescale_error (err, result_abs, result_asc); 01570 01571 } 01572 01573 }; 01574 01575 #ifndef DOXYGENP 01576 } 01577 #endif 01578 01579 #endif
Documentation generated with Doxygen and provided under the GNU Free Documentation License. See License Information for details.
Project hosting provided by
,
O2scl Sourceforge Project Page