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