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 #ifdef O2SCL_NEVER_DEFINED 00423 // This is a proposed revision of the class structure 00424 00425 class gsl_inte_workspace : public gsl_integration_workspace { 00426 00427 public: 00428 00429 int set_size(int sz); 00430 00431 int size(); 00432 00433 int resize(int sz); 00434 00435 int initialise(double a, double b); 00436 00437 int set_initial_result(double result, double error); 00438 00439 int retrieve(double *a, double *b, double *r, double *e) const; 00440 00441 int qpsrt(); 00442 00443 int update(double a1, double b1, double area1, double error1, 00444 double a2, double b2, double area2, double error2); 00445 00446 double sum_results(); 00447 00448 int subinterval_too_small(double a1, double a2, double b2); 00449 00450 int append_interval(double a1, double b1, double area1, double error1); 00451 00452 }; 00453 #endif 00454 00455 /** 00456 \brief Base routines for the GSL adaptive integration classes 00457 00458 This class contains several functions for manipulating the GSL 00459 integration workspace. The casual end-user should use the 00460 classes explained in the \ref inte_section section of the User's 00461 guide. 00462 00463 \todo Make the workspace size protected 00464 00465 \future Move gsl_integration_workspace to a separate class 00466 and remove this class, making all children direct descendants 00467 of gsl_inte instead. We'll have to figure out what to 00468 do with the data member \c wkspace though. Some work on this 00469 front is already in gsl_inte_qag_b.h. 00470 00471 QUADPACK workspace documentation: 00472 \verbatim 00473 c parameters (meaning at output) 00474 c limit - integer 00475 c maximum number of error estimates the list 00476 c can contain 00477 c last - integer 00478 c number of error estimates currently in the list 00479 c maxerr - integer 00480 c maxerr points to the nrmax-th largest error 00481 c estimate currently in the list 00482 c ermax - double precision 00483 c nrmax-th largest error estimate 00484 c ermax = elist(maxerr) 00485 c elist - double precision 00486 c vector of dimension last containing 00487 c the error estimates 00488 c iord - integer 00489 c vector of dimension last, the first k elements 00490 c of which contain pointers to the error 00491 c estimates, such that 00492 c elist(iord(1)),..., elist(iord(k)) 00493 c form a decreasing sequence, with 00494 c k = last if last.le.(limit/2+2), and 00495 c k = limit+1-last otherwise 00496 c nrmax - integer 00497 c maxerr = iord(nrmax) 00498 \endverbatim 00499 \verbatim 00500 c alist - real 00501 c vector of dimension at least limit, the first 00502 c last elements of which are the left 00503 c end points of the subintervals in the partition 00504 c of the given integration range (a,b) 00505 c blist - real 00506 c vector of dimension at least limit, the first 00507 c last elements of which are the right 00508 c end points of the subintervals in the partition 00509 c of the given integration range (a,b) 00510 c rlist - real 00511 c vector of dimension at least limit, the first 00512 c last elements of which are the 00513 c integral approximations on the subintervals 00514 c elist - real 00515 c vector of dimension at least limit, the first 00516 c last elements of which are the moduli of the 00517 c absolute error estimates on the subintervals 00518 c iord - integer 00519 c vector of dimension at least limit, the first k 00520 c elements of which are pointers to the 00521 c error estimates over the subintervals, 00522 c such that elist(iord(1)), ..., 00523 c elist(iord(k)) form a decreasing sequence, 00524 c with k = last if last.le.(limit/2+2), and 00525 c k = limit+1-last otherwise 00526 c last - integer 00527 c number of subintervals actually produced in the 00528 c subdivision process 00529 \endverbatim 00530 00531 */ 00532 class gsl_inte_table : public gsl_inte { 00533 00534 public: 00535 00536 gsl_inte_table(); 00537 00538 ~gsl_inte_table(); 00539 00540 /// The integration workspace 00541 gsl_integration_workspace *w; 00542 00543 /// The size of the integration workspace (default 1000) 00544 int wkspace; 00545 00546 /// Set the integration workspace size 00547 int set_wkspace(size_t size); 00548 00549 /** 00550 \brief Initialize the workspace for an integration with limits 00551 \c a and \c b. 00552 */ 00553 void initialise(gsl_integration_workspace *workspace, 00554 double a, double b); 00555 00556 /// Set the result at position zero 00557 void set_initial_result(gsl_integration_workspace *workspace, 00558 double result, double error); 00559 00560 /** 00561 \brief Retrieve the ith result from the workspace 00562 00563 The workspace variable \c i is used to specify which 00564 interval is requested. 00565 */ 00566 void retrieve(const gsl_integration_workspace *workspace, 00567 double *a, double *b, double *r, double *e); 00568 00569 /** \brief Sort the workspace 00570 00571 This routine maintains the descending ordering in the list of 00572 the local error estimated resulting from the interval 00573 subdivision process. at each call two error estimates are 00574 inserted using the sequential search method, top-down for the 00575 largest error estimate and bottom-up for the smallest error 00576 estimate. 00577 00578 Originally written in QUADPACK by 00579 \verbatim 00580 piessens, robert, appl. math. & progr. div. - k.u.leuven 00581 de doncker, elise, appl. math. & progr. div. - k.u.leuven 00582 \endverbatim 00583 translated into C for GSL by Brian Gough, and then 00584 rewritten for O2scl. 00585 00586 */ 00587 void qpsrt(gsl_integration_workspace *workspace); 00588 00589 /// Update workspace with new results and resort 00590 void update(gsl_integration_workspace *workspace, 00591 double a1, double b1, double area1, double error1, 00592 double a2, double b2, double area2, double error2); 00593 00594 /// Add up all of the contributions to construct the final result 00595 double sum_results(const gsl_integration_workspace *workspace); 00596 00597 /// Find out if the present subinterval is too small 00598 int subinterval_too_small(double a1, double a2, double b2); 00599 00600 /// Append new results to workspace 00601 void append_interval(gsl_integration_workspace *workspace, 00602 double a1, double b1, double area1, double error1); 00603 00604 }; 00605 00606 /** 00607 \brief Basic Gauss-Kronrod integration class (GSL) 00608 00609 This class provides the basic Gauss-Kronrod integration function 00610 for some of the GSL-based integration classes. The casual 00611 end-user should use the classes explained in the \ref 00612 inte_section section of the User's guide. 00613 */ 00614 template<class param_t, class func_t> 00615 class gsl_inte_kronrod : public gsl_inte_table, 00616 public inte<param_t,func_t> { 00617 00618 public: 00619 00620 /** 00621 \brief The GSL Gauss-Kronrod integration function 00622 00623 Given abcissas and weights, this performs the integration of 00624 \c func between \c a and \c b, providing a result with 00625 uncertainties. 00626 00627 This function is designed for use with the values given in the 00628 o2scl_inte_qag_coeffs namespace. 00629 00630 This function never calls the error handler. 00631 00632 \future This function, in principle, could be replaced with 00633 a generic integration pointer. 00634 */ 00635 virtual void gsl_integration_qk_o2scl 00636 (func_t &func, const int n, const double xgk[], const double wg[], 00637 const double wgk[], double fv1[], double fv2[], double a, double b, 00638 double *result, double *abserr, double *resabs, double *resasc, 00639 param_t &pa) { 00640 00641 const double center = 0.5 * (a + b); 00642 const double half_length = 0.5 * (b - a); 00643 const double abs_half_length = fabs (half_length); 00644 00645 double f_center; 00646 func(center,f_center,pa); 00647 00648 double result_gauss = 0; 00649 double result_kronrod = f_center * wgk[n - 1]; 00650 00651 double result_abs = fabs (result_kronrod); 00652 double result_asc = 0; 00653 double mean = 0, err = 0; 00654 00655 int j; 00656 00657 if (n % 2 == 0) { 00658 result_gauss = f_center * wg[n / 2 - 1]; 00659 } 00660 00661 // Sum up results from left half of interval 00662 for (j = 0; j < (n - 1) / 2; j++) { 00663 00664 const int jtw = j * 2 + 1; 00665 const double abscissa = half_length * xgk[jtw]; 00666 00667 double fval1, fval2; 00668 func(center-abscissa,fval1,pa); 00669 func(center+abscissa,fval2,pa); 00670 00671 const double fsum = fval1 + fval2; 00672 fv1[jtw] = fval1; 00673 fv2[jtw] = fval2; 00674 result_gauss += wg[j] * fsum; 00675 result_kronrod += wgk[jtw] * fsum; 00676 result_abs += wgk[jtw] * (fabs (fval1) + fabs (fval2)); 00677 } 00678 00679 // Sum up results from right half of interval 00680 for (j = 0; j < n / 2; j++) { 00681 00682 int jtwm1 = j * 2; 00683 const double abscissa = half_length * xgk[jtwm1]; 00684 00685 double fval1, fval2; 00686 func(center-abscissa,fval1,pa); 00687 func(center+abscissa,fval2,pa); 00688 00689 fv1[jtwm1] = fval1; 00690 fv2[jtwm1] = fval2; 00691 result_kronrod += wgk[jtwm1] * (fval1 + fval2); 00692 result_abs += wgk[jtwm1] * (fabs (fval1) + fabs (fval2)); 00693 } 00694 00695 mean = result_kronrod * 0.5; 00696 00697 result_asc = wgk[n - 1] * fabs (f_center - mean); 00698 00699 for (j = 0; j < n - 1; j++) { 00700 result_asc += wgk[j] * (fabs (fv1[j] - mean) + 00701 fabs (fv2[j] - mean)); 00702 } 00703 00704 /* Scale by the width of the integration region */ 00705 00706 err = (result_kronrod - result_gauss) * half_length; 00707 00708 result_kronrod *= half_length; 00709 result_abs *= abs_half_length; 00710 result_asc *= abs_half_length; 00711 00712 *result = result_kronrod; 00713 *resabs = result_abs; 00714 *resasc = result_asc; 00715 *abserr = rescale_error (err, result_abs, result_asc); 00716 00717 return; 00718 } 00719 00720 }; 00721 00722 /** \brief Base class for integrating a function with a 00723 singularity (GSL) 00724 00725 This class contains the extrapolation table mechanics and the 00726 base integration function for singular integrals from GSL. The 00727 casual end-user should use the classes explained in the 00728 \ref inte_section section of the User's guide. 00729 00730 \future Some of the functions inside this class could 00731 be moved out of header files? 00732 */ 00733 template<class param_t, class func_t> class gsl_inte_singular : 00734 public gsl_inte_kronrod<param_t,func_t> { 00735 00736 protected: 00737 00738 /** 00739 \brief A structure for extrapolation for \ref gsl_inte_qags 00740 00741 \future Move this to a new class, with qelg() as a method 00742 */ 00743 struct extrapolation_table 00744 { 00745 /// Index of new element in the first column 00746 size_t n; 00747 /// Lower diagonals of the triangular epsilon table 00748 double rlist2[52]; 00749 /// Number of calls 00750 size_t nres; 00751 /// Three most recent results 00752 double res3la[3]; 00753 }; 00754 00755 /// Initialize the table 00756 void initialise_table(struct extrapolation_table *table) { 00757 table->n = 0; 00758 table->nres = 0; 00759 return; 00760 } 00761 00762 /// Append a result to the table 00763 void append_table(struct extrapolation_table *table, double y) { 00764 size_t n; 00765 n = table->n; 00766 table->rlist2[n] = y; 00767 table->n++; 00768 return; 00769 } 00770 00771 /// Test a result for positivity 00772 inline int test_positivity(double result, double resabs) { 00773 int status = (fabs (result) >= (1 - 50 * GSL_DBL_EPSILON) * resabs); 00774 return status; 00775 } 00776 00777 /** \brief Determines the limit of a given sequence of approximations 00778 00779 This function determines the limit of a given sequence of 00780 approximations, by means of the epsilon algorithm of P. Wynn. 00781 an estimate of the absolute error is also given. the condensed 00782 epsilon table is computed. only those elements needed for the 00783 computation of the next diagonal are preserved. 00784 00785 Quadpack documentation 00786 \verbatim 00787 c 00788 c list of major variables 00789 c ----------------------- 00790 c 00791 c e0 - the 4 elements on which the computation of a new 00792 c e1 element in the epsilon table is based 00793 c e2 00794 c e3 e0 00795 c e3 e1 new 00796 c e2 00797 c newelm - number of elements to be computed in the new 00798 c diagonal 00799 c error - error = abs(e1-e0)+abs(e2-e1)+abs(new-e2) 00800 c result - the element in the new diagonal with least value 00801 c of error 00802 c 00803 c machine dependent constants 00804 c --------------------------- 00805 c 00806 c epmach is the largest relative spacing. 00807 c oflow is the largest positive magnitude. 00808 c limexp is the maximum number of elements the epsilon 00809 c table can contain. if this number is reached, the upper 00810 c diagonal of the epsilon table is deleted. 00811 c 00812 \endverbatim 00813 */ 00814 void qelg(struct extrapolation_table *table, double *result, 00815 double *abserr) { 00816 00817 double *epstab = table->rlist2; 00818 double *res3la = table->res3la; 00819 const size_t n = table->n - 1; 00820 00821 const double current = epstab[n]; 00822 00823 double absolute = GSL_DBL_MAX; 00824 double relative = 5 * GSL_DBL_EPSILON * fabs (current); 00825 00826 const size_t newelm = n / 2; 00827 const size_t n_orig = n; 00828 size_t n_final = n; 00829 size_t i; 00830 00831 const size_t nres_orig = table->nres; 00832 00833 *result = current; 00834 *abserr = GSL_DBL_MAX; 00835 00836 if (n < 2) { 00837 *result = current; 00838 *abserr = GSL_MAX_DBL (absolute, relative); 00839 return; 00840 } 00841 00842 epstab[n + 2] = epstab[n]; 00843 epstab[n] = GSL_DBL_MAX; 00844 00845 for (i = 0; i < newelm; i++) { 00846 double res = epstab[n - 2 * i + 2]; 00847 double e0 = epstab[n - 2 * i - 2]; 00848 double e1 = epstab[n - 2 * i - 1]; 00849 double e2 = res; 00850 00851 double e1abs = fabs (e1); 00852 double delta2 = e2 - e1; 00853 double err2 = fabs (delta2); 00854 double tol2 = GSL_MAX_DBL (fabs (e2), e1abs) * GSL_DBL_EPSILON; 00855 double delta3 = e1 - e0; 00856 double err3 = fabs (delta3); 00857 double tol3 = GSL_MAX_DBL (e1abs, fabs (e0)) * GSL_DBL_EPSILON; 00858 00859 double e3, delta1, err1, tol1, ss; 00860 00861 if (err2 <= tol2 && err3 <= tol3) 00862 { 00863 /* If e0, e1 and e2 are equal to within machine accuracy, 00864 convergence is assumed. */ 00865 00866 *result = res; 00867 absolute = err2 + err3; 00868 relative = 5 * GSL_DBL_EPSILON * fabs (res); 00869 *abserr = GSL_MAX_DBL (absolute, relative); 00870 return; 00871 } 00872 00873 e3 = epstab[n - 2 * i]; 00874 epstab[n - 2 * i] = e1; 00875 delta1 = e1 - e3; 00876 err1 = fabs (delta1); 00877 tol1 = GSL_MAX_DBL (e1abs, fabs (e3)) * GSL_DBL_EPSILON; 00878 00879 /* If two elements are very close to each other, omit a part of 00880 the table by adjusting the value of n */ 00881 00882 if (err1 <= tol1 || err2 <= tol2 || err3 <= tol3) { 00883 n_final = 2 * i; 00884 break; 00885 } 00886 00887 ss = (1 / delta1 + 1 / delta2) - 1 / delta3; 00888 00889 /* Test to detect irregular behaviour in the table, and 00890 eventually omit a part of the table by adjusting the value of 00891 n. */ 00892 if (fabs (ss * e1) <= 0.0001) { 00893 n_final = 2 * i; 00894 break; 00895 } 00896 /* Compute a new element and eventually adjust the value of 00897 result. */ 00898 00899 res = e1 + 1 / ss; 00900 epstab[n - 2 * i] = res; 00901 00902 { 00903 const double error = err2 + fabs (res - e2) + err3; 00904 00905 if (error <= *abserr) { 00906 *abserr = error; 00907 *result = res; 00908 } 00909 } 00910 } 00911 00912 /* Shift the table */ 00913 00914 { 00915 const size_t limexp = 50 - 1; 00916 00917 if (n_final == limexp) 00918 { 00919 n_final = 2 * (limexp / 2); 00920 } 00921 } 00922 00923 if (n_orig % 2 == 1) { 00924 for (i = 0; i <= newelm; i++) { 00925 epstab[1 + i * 2] = epstab[i * 2 + 3]; 00926 } 00927 } else { 00928 for (i = 0; i <= newelm; i++) { 00929 epstab[i * 2] = epstab[i * 2 + 2]; 00930 } 00931 } 00932 if (n_orig != n_final) { 00933 for (i = 0; i <= n_final; i++) { 00934 epstab[i] = epstab[n_orig - n_final + i]; 00935 } 00936 } 00937 00938 table->n = n_final + 1; 00939 00940 if (nres_orig < 3) { 00941 00942 res3la[nres_orig] = *result; 00943 *abserr = GSL_DBL_MAX; 00944 00945 } else { 00946 /* Compute error estimate */ 00947 *abserr = (fabs (*result - res3la[2]) + fabs (*result - res3la[1]) 00948 + fabs (*result - res3la[0])); 00949 00950 res3la[0] = res3la[1]; 00951 res3la[1] = res3la[2]; 00952 res3la[2] = *result; 00953 } 00954 00955 /* In QUADPACK the variable table->nres is incremented at the top of 00956 qelg, so it increases on every call. This leads to the array 00957 res3la being accessed when its elements are still undefined, so I 00958 have moved the update to this point so that its value more 00959 useful. */ 00960 00961 table->nres = nres_orig + 1; 00962 00963 *abserr = GSL_MAX_DBL (*abserr, 5 * GSL_DBL_EPSILON * fabs (*result)); 00964 00965 return; 00966 } 00967 00968 /// Determine if an interval is large 00969 int large_interval (gsl_integration_workspace * workspace) { 00970 size_t i = workspace->i ; 00971 const size_t * level = workspace->level; 00972 00973 if (level[i] < workspace->maximum_level) { 00974 return 1; 00975 } else { 00976 return 0; 00977 } 00978 } 00979 00980 /// Reset workspace to work on the interval with the largest error 00981 inline void reset_nrmax (gsl_integration_workspace * workspace) { 00982 workspace->nrmax = 0; 00983 workspace->i = workspace->order[0] ; 00984 } 00985 00986 /// Increase workspace 00987 int increase_nrmax (gsl_integration_workspace * workspace) { 00988 int k; 00989 int id = workspace->nrmax; 00990 int jupbnd; 00991 00992 const size_t * level = workspace->level; 00993 const size_t * order = workspace->order; 00994 00995 size_t limit = workspace->limit ; 00996 size_t last = workspace->size - 1 ; 00997 00998 if (last > (1 + limit / 2)) { 00999 jupbnd = limit + 1 - last; 01000 } else { 01001 jupbnd = last; 01002 } 01003 01004 for (k = id; k <= jupbnd; k++) { 01005 size_t i_max = order[workspace->nrmax]; 01006 01007 workspace->i = i_max ; 01008 01009 if (level[i_max] < workspace->maximum_level) { 01010 return 1; 01011 } 01012 01013 workspace->nrmax++; 01014 01015 } 01016 return 0; 01017 } 01018 01019 /** 01020 \brief Integration function 01021 01022 \future Remove goto statements? 01023 */ 01024 int qags(func_t &func, const int qn, const double xgk[], 01025 const double wg[], const double wgk[], double fv1[], 01026 double fv2[], const double a, const double b, 01027 const double l_epsabs, const double l_epsrel, 01028 const size_t limit, double *result, double *abserr, param_t &pa) 01029 { 01030 01031 double area, errsum; 01032 double res_ext, err_ext; 01033 double result0, abserr0, resabs0, resasc0; 01034 double tolerance; 01035 01036 double ertest = 0; 01037 double error_over_large_intervals = 0; 01038 double reseps = 0, abseps = 0, correc = 0; 01039 size_t ktmin = 0; 01040 int roundoff_type1 = 0, roundoff_type2 = 0, roundoff_type3 = 0; 01041 int error_type = 0, error_type2 = 0; 01042 01043 size_t iteration = 0; 01044 01045 int positive_integrand = 0; 01046 int extrapolate = 0; 01047 int disallow_extrapolation = 0; 01048 01049 struct extrapolation_table table; 01050 01051 /* Initialize results */ 01052 01053 initialise(this->w, a, b); 01054 01055 *result = 0; 01056 *abserr = 0; 01057 01058 if (limit > this->w->limit) { 01059 this->last_iter=0; 01060 std::string estr="Iteration limit exceeds workspace "; 01061 estr+="in gsl_inte_kronrod::qags()."; 01062 O2SCL_ERR_RET(estr.c_str(),gsl_einval); 01063 } 01064 01065 /* Test on accuracy */ 01066 if (this->tolx <= 0 && (this->tolf < 50 * GSL_DBL_EPSILON || 01067 this->tolf < 0.5e-28)) { 01068 this->last_iter=0; 01069 01070 std::string estr="Tolerance cannot be achieved with given "; 01071 estr+="value of 'tolx' and 'tolf' in gsl_inte_kronrod::qags()."; 01072 O2SCL_ERR_RET(estr.c_str(),gsl_ebadtol); 01073 } 01074 01075 /* Perform the first integration */ 01076 01077 gsl_integration_qk_o2scl(func,qn,xgk,wg,wgk,fv1,fv2,a,b, 01078 &result0,&abserr0,&resabs0,&resasc0,pa); 01079 01080 set_initial_result (this->w, result0, abserr0); 01081 01082 tolerance = GSL_MAX_DBL (this->tolx, this->tolf * fabs (result0)); 01083 01084 if (abserr0 <= 100 * GSL_DBL_EPSILON * resabs0 && 01085 abserr0 > tolerance) { 01086 01087 *result = result0; 01088 *abserr = abserr0; 01089 01090 this->last_iter=1; 01091 std::string estr="Cannot reach tolerance because of roundoff error "; 01092 estr+="on first attempt in gsl_inte_kronrod::qags()."; 01093 O2SCL_ERR_RET(estr.c_str(),gsl_eround); 01094 01095 } else if ((abserr0 <= tolerance && 01096 abserr0 != resasc0) || abserr0 == 0.0) { 01097 01098 *result = result0; 01099 *abserr = abserr0; 01100 this->last_iter=1; 01101 return gsl_success; 01102 01103 } else if (limit == 1) { 01104 01105 *result = result0; 01106 *abserr = abserr0; 01107 01108 this->last_iter=1; 01109 std::string estr="A maximum of 1 iteration was insufficient "; 01110 estr+="in gsl_inte_kronrod::qags()."; 01111 O2SCL_ERR_RET(estr.c_str(),gsl_emaxiter); 01112 } 01113 01114 /* Initialization */ 01115 01116 initialise_table (&table); 01117 append_table (&table, result0); 01118 01119 area = result0; 01120 errsum = abserr0; 01121 01122 res_ext = result0; 01123 err_ext = GSL_DBL_MAX; 01124 01125 positive_integrand = this->test_positivity (result0, resabs0); 01126 01127 iteration = 1; 01128 01129 do { 01130 size_t current_level; 01131 double a1, b1, a2, b2; 01132 double a_i, b_i, r_i, e_i; 01133 double area1 = 0, area2 = 0, area12 = 0; 01134 double error1 = 0, error2 = 0, error12 = 0; 01135 double resasc1, resasc2; 01136 double resabs1, resabs2; 01137 double last_e_i; 01138 01139 /* Bisect the subinterval with the largest error estimate */ 01140 01141 retrieve (this->w, &a_i, &b_i, &r_i, &e_i); 01142 01143 current_level = this->w->level[this->w->i] + 1; 01144 01145 a1 = a_i; 01146 b1 = 0.5 * (a_i + b_i); 01147 a2 = b1; 01148 b2 = b_i; 01149 01150 iteration++; 01151 01152 gsl_integration_qk_o2scl(func,qn,xgk,wg,wgk,fv1,fv2,a1,b1, 01153 &area1,&error1,&resabs1,&resasc1,pa); 01154 gsl_integration_qk_o2scl(func,qn,xgk,wg,wgk,fv1,fv2,a2,b2, 01155 &area2,&error2,&resabs2,&resasc2,pa); 01156 01157 area12 = area1 + area2; 01158 error12 = error1 + error2; 01159 last_e_i = e_i; 01160 01161 /* Improve previous approximations to the integral and test for 01162 accuracy. 01163 01164 We write these expressions in the same way as the original 01165 QUADPACK code so that the rounding errors are the same, which 01166 makes testing easier. 01167 */ 01168 01169 errsum = errsum + error12 - e_i; 01170 area = area + area12 - r_i; 01171 01172 tolerance = GSL_MAX_DBL (this->tolx, this->tolf * fabs (area)); 01173 if (resasc1 != error1 && resasc2 != error2) { 01174 double delta = r_i - area12; 01175 01176 if (fabs (delta) <= 1.0e-5 * fabs (area12) && 01177 error12 >= 0.99 * e_i) { 01178 if (!extrapolate) { 01179 roundoff_type1++; 01180 } else { 01181 roundoff_type2++; 01182 } 01183 } 01184 if (iteration > 10 && error12 > e_i) { 01185 roundoff_type3++; 01186 } 01187 } 01188 01189 /* Test for roundoff and eventually set error flag */ 01190 01191 if (roundoff_type1 + roundoff_type2 >= 10 || 01192 roundoff_type3 >= 20) { 01193 error_type = 2; /* round off error */ 01194 } 01195 01196 if (roundoff_type2 >= 5) { 01197 error_type2 = 1; 01198 } 01199 01200 /* set error flag in the case of bad integrand behaviour at 01201 a point of the integration range */ 01202 01203 if (this->subinterval_too_small (a1, a2, b2)) { 01204 error_type = 4; 01205 } 01206 01207 /* append the newly-created intervals to the list */ 01208 01209 update (this->w, a1, b1, area1, error1, a2, b2, area2, error2); 01210 01211 if (errsum <= tolerance) { 01212 goto compute_result; 01213 } 01214 01215 if (error_type) { 01216 break; 01217 } 01218 01219 if (iteration >= limit - 1) { 01220 error_type = 1; 01221 break; 01222 } 01223 01224 if (iteration == 2) { 01225 error_over_large_intervals = errsum; 01226 ertest = tolerance; 01227 append_table (&table, area); 01228 continue; 01229 } 01230 01231 if (disallow_extrapolation) { 01232 continue; 01233 } 01234 01235 error_over_large_intervals += -last_e_i; 01236 01237 if (current_level < this->w->maximum_level) { 01238 error_over_large_intervals += error12; 01239 } 01240 01241 if (!extrapolate) { 01242 01243 /* test whether the interval to be bisected next is the 01244 smallest interval. */ 01245 01246 if (large_interval (this->w)) { 01247 continue; 01248 } 01249 01250 extrapolate = 1; 01251 this->w->nrmax = 1; 01252 } 01253 01254 if (!error_type2 && error_over_large_intervals > ertest) { 01255 if (increase_nrmax (this->w)) { 01256 continue; 01257 } 01258 } 01259 01260 /* Perform extrapolation */ 01261 01262 append_table (&table, area); 01263 01264 qelg (&table, &reseps, &abseps); 01265 01266 ktmin++; 01267 01268 if (ktmin > 5 && err_ext < 0.001 * errsum) { 01269 error_type = 5; 01270 } 01271 01272 if (abseps < err_ext) { 01273 ktmin = 0; 01274 err_ext = abseps; 01275 res_ext = reseps; 01276 correc = error_over_large_intervals; 01277 ertest = GSL_MAX_DBL (this->tolx, 01278 this->tolf * fabs (reseps)); 01279 if (err_ext <= ertest) 01280 break; 01281 } 01282 01283 /* Prepare bisection of the smallest interval. */ 01284 01285 if (table.n == 1) { 01286 disallow_extrapolation = 1; 01287 } 01288 01289 if (error_type == 5) { 01290 break; 01291 } 01292 01293 /* work on interval with largest error */ 01294 01295 reset_nrmax (this->w); 01296 extrapolate = 0; 01297 error_over_large_intervals = errsum; 01298 01299 /// Output iteration information 01300 01301 if (this->verbose>0) { 01302 std::cout << this->type(); 01303 std::cout << " Iter: " << iteration; 01304 std::cout.setf(std::ios::showpos); 01305 std::cout << " Res: " << area; 01306 std::cout.unsetf(std::ios::showpos); 01307 std::cout << " Err: " << errsum 01308 << " Tol: " << tolerance << std::endl; 01309 if (this->verbose>1) { 01310 char ch; 01311 std::cout << "Press a key and type enter to continue. " ; 01312 std::cin >> ch; 01313 } 01314 } 01315 01316 } while (iteration < limit); 01317 01318 *result = res_ext; 01319 *abserr = err_ext; 01320 01321 if (err_ext == GSL_DBL_MAX) 01322 goto compute_result; 01323 01324 if (error_type || error_type2) { 01325 if (error_type2) { 01326 err_ext += correc; 01327 } 01328 01329 if (error_type == 0) 01330 error_type = 3; 01331 01332 if (res_ext != 0.0 && area != 0.0) { 01333 if (err_ext / fabs (res_ext) > errsum / fabs (area)) 01334 goto compute_result; 01335 } else if (err_ext > errsum) { 01336 goto compute_result; 01337 } else if (area == 0.0) { 01338 goto return_error; 01339 } 01340 } 01341 01342 /* Test on divergence. */ 01343 01344 { 01345 double max_area = GSL_MAX_DBL (fabs (res_ext), fabs (area)); 01346 01347 if (!positive_integrand && max_area < 0.01 * resabs0) 01348 goto return_error; 01349 } 01350 01351 { 01352 double ratio = res_ext / area; 01353 01354 if (ratio < 0.01 || ratio > 100.0 || errsum > fabs (area)) 01355 error_type = 6; 01356 } 01357 01358 goto return_error; 01359 01360 compute_result: 01361 01362 *result = sum_results (this->w); 01363 *abserr = errsum; 01364 01365 return_error: 01366 01367 if (error_type > 2) { 01368 error_type--; 01369 } 01370 01371 this->last_iter=iteration; 01372 01373 if (error_type == 0) { 01374 return gsl_success; 01375 } else if (error_type == 1) { 01376 std::string estr="Number of iterations was insufficient "; 01377 estr+=" in gsl_inte_kronrod::qags()."; 01378 O2SCL_ERR_RET(estr.c_str(),gsl_emaxiter); 01379 } else if (error_type == 2) { 01380 std::string estr="Roundoff error prevents tolerance "; 01381 estr+="from being achieved in gsl_inte_kronrod::qags()."; 01382 O2SCL_ERR_RET(estr.c_str(),gsl_eround); 01383 } else if (error_type == 3) { 01384 std::string estr="Bad integrand behavior "; 01385 estr+="in gsl_inte_kronrod::qags()."; 01386 O2SCL_ERR_RET(estr.c_str(),gsl_esing); 01387 } else if (error_type == 4) { 01388 std::string estr="Roundoff error detected in extrapolation table "; 01389 estr+="in gsl_inte_kronrod::qags()."; 01390 O2SCL_ERR_RET(estr.c_str(),gsl_eround); 01391 } else if (error_type == 5) { 01392 std::string estr="Integral is divergent or slowly convergent "; 01393 estr+="in gsl_inte_kronrod::qags()."; 01394 O2SCL_ERR_RET(estr.c_str(),gsl_ediverge); 01395 } 01396 01397 std::string estr="Could not integrate function in gsl_inte_kronrod"; 01398 estr+="::qags() (it may have returned a non-finite result)."; 01399 O2SCL_ERR(estr.c_str(),gsl_efailed); 01400 01401 return gsl_efailed; 01402 } 01403 01404 }; 01405 01406 /** 01407 \brief Integrate a function with a singularity (GSL) [abstract base] 01408 01409 This class contains the GSL-based integration function for 01410 applying transformations to the user-defined integrand. The 01411 casual end-user should use the classes explained in the 01412 \ref inte_section section of the User's guide. 01413 */ 01414 template<class param_t, class func_t> class gsl_inte_transform : 01415 public gsl_inte_singular<param_t,func_t> { 01416 01417 public: 01418 01419 /// The transformation to apply to the user-supplied function 01420 virtual double transform(func_t &func, double t, param_t &pa)=0; 01421 01422 /** 01423 \brief The basic Gauss-Kronrod integration function 01424 01425 This function is basically just a copy of 01426 gsl_inte_qag::gsl_integration_qk_o2scl() which is rewritten to 01427 call the internal transformed function rather than directly 01428 calling the user-specified function. 01429 01430 This function never calls the error handler. 01431 */ 01432 virtual void gsl_integration_qk_o2scl 01433 (func_t &func, const int n, const double xgk[], const double wg[], 01434 const double wgk[], double fv1[], double fv2[], double a, double b, 01435 double *result, double *abserr, double *resabs, double *resasc, 01436 param_t &pa) { 01437 01438 const double center = 0.5 * (a + b); 01439 const double half_length = 0.5 * (b - a); 01440 const double abs_half_length = fabs (half_length); 01441 01442 const double f_center=transform(func,center,pa); 01443 01444 double result_gauss = 0; 01445 double result_kronrod = f_center * wgk[n - 1]; 01446 01447 double result_abs = fabs (result_kronrod); 01448 double result_asc = 0; 01449 double mean = 0, err = 0; 01450 01451 int j; 01452 01453 if (n % 2 == 0) { 01454 result_gauss = f_center * wg[n / 2 - 1]; 01455 } 01456 01457 // Sum up results from left half of interval 01458 for (j = 0; j < (n - 1) / 2; j++) { 01459 01460 const int jtw = j * 2 + 1; 01461 const double abscissa = half_length * xgk[jtw]; 01462 01463 const double fval1=transform(func,center-abscissa,pa); 01464 const double fval2=transform(func,center+abscissa,pa); 01465 const double fsum = fval1 + fval2; 01466 01467 fv1[jtw] = fval1; 01468 fv2[jtw] = fval2; 01469 result_gauss += wg[j] * fsum; 01470 result_kronrod += wgk[jtw] * fsum; 01471 result_abs += wgk[jtw] * (fabs (fval1) + fabs (fval2)); 01472 } 01473 01474 // Sum up results from right half of interval 01475 for (j = 0; j < n / 2; j++) { 01476 int jtwm1 = j * 2; 01477 const double abscissa = half_length * xgk[jtwm1]; 01478 01479 const double fval1=transform(func,center-abscissa,pa); 01480 const double fval2=transform(func,center+abscissa,pa); 01481 01482 fv1[jtwm1] = fval1; 01483 fv2[jtwm1] = fval2; 01484 result_kronrod += wgk[jtwm1] * (fval1 + fval2); 01485 result_abs += wgk[jtwm1] * (fabs (fval1) + fabs (fval2)); 01486 }; 01487 01488 mean = result_kronrod * 0.5; 01489 01490 result_asc = wgk[n - 1] * fabs (f_center - mean); 01491 01492 for (j = 0; j < n - 1; j++) { 01493 result_asc += wgk[j] * (fabs (fv1[j] - mean) + 01494 fabs (fv2[j] - mean)); 01495 } 01496 01497 // Scale by the width of the integration region 01498 01499 err = (result_kronrod - result_gauss) * half_length; 01500 01501 result_kronrod *= half_length; 01502 result_abs *= abs_half_length; 01503 result_asc *= abs_half_length; 01504 01505 *result = result_kronrod; 01506 *resabs = result_abs; 01507 *resasc = result_asc; 01508 *abserr = this->rescale_error (err, result_abs, result_asc); 01509 01510 } 01511 01512 }; 01513 01514 #ifndef DOXYGENP 01515 } 01516 #endif 01517 01518 #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