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