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_QAWO_H 00024 #define O2SCL_GSL_INTE_QAWO_H 00025 00026 #include <o2scl/inte.h> 00027 #include <o2scl/gsl_inte_qawc.h> 00028 00029 #ifndef DOXYGENP 00030 namespace o2scl { 00031 #endif 00032 00033 /** 00034 \brief Adaptive integration for oscillatory integrals (GSL) 00035 00036 The number of subdivisions of the original interval which 00037 this class is allowed to make is dictated by the workspace 00038 size for the integration class, which can be set using 00039 \ref gsl_inte_table::set_wkspace() . 00040 00041 \todo Improve documentation 00042 */ 00043 template<class param_t, class func_t> class gsl_inte_qawo_sin : 00044 public gsl_inte_cheb<param_t,func_t> { 00045 00046 public: 00047 00048 gsl_inte_qawo_sin() { 00049 tab_size=4; 00050 omega=1.0; 00051 } 00052 00053 virtual ~gsl_inte_qawo_sin() {} 00054 00055 /** \brief Integrate function \c func from \c a to \c b. 00056 */ 00057 virtual double integ(func_t &func, double a, double b, param_t &pa) { 00058 double res, err; 00059 integ_err(func,a,b,pa,res,err); 00060 this->interror=err; 00061 return res; 00062 } 00063 00064 /// Desc 00065 double omega; 00066 00067 /// Desc 00068 size_t tab_size; 00069 00070 /** \brief Integrate function \c func from \c a to \c b and place 00071 the result in \c res and the error in \c err 00072 */ 00073 virtual int integ_err(func_t &func, double a, double b, 00074 param_t &pa, double &res, double &err2) { 00075 00076 otable=gsl_integration_qawo_table_alloc 00077 (omega,b-a,GSL_INTEG_SINE,tab_size); 00078 00079 int status=qawo(func,a,this->tolx,this->tolf,this->wkspace, 00080 this->w,this->otable, 00081 &res,&err2,pa); 00082 00083 gsl_integration_qawo_table_free(otable); 00084 00085 return status; 00086 00087 } 00088 00089 #ifndef DOXYGEN_INTERNAL 00090 00091 protected: 00092 00093 /// The integration workspace 00094 gsl_integration_qawo_table *otable; 00095 00096 /** 00097 \brief The full GSL integration routine called by integ_err() 00098 */ 00099 int qawo(func_t &func, const double a, 00100 const double epsabs, const double epsrel, const size_t limit, 00101 gsl_integration_workspace *loc_w, 00102 gsl_integration_qawo_table *wf, 00103 double *result, double *abserr, param_t &pa) { 00104 00105 double area, errsum; 00106 double res_ext, err_ext; 00107 double result0, abserr0, resabs0, resasc0; 00108 double tolerance; 00109 00110 double ertest = 0; 00111 double error_over_large_intervals = 0; 00112 double reseps = 0, abseps = 0, correc = 0; 00113 size_t ktmin = 0; 00114 int roundoff_type1 = 0, roundoff_type2 = 0, roundoff_type3 = 0; 00115 int error_type = 0, error_type2 = 0; 00116 00117 size_t iteration = 0; 00118 00119 int positive_integrand = 0; 00120 int extrapolate = 0; 00121 int extall = 0; 00122 int disallow_extrapolation = 0; 00123 00124 struct extrapolation_table table; 00125 00126 double b = a + wf->L ; 00127 double abs_omega = fabs (wf->omega) ; 00128 00129 /* Initialize results */ 00130 00131 this->initialise (loc_w, a, b); 00132 00133 *result = 0; 00134 *abserr = 0; 00135 00136 if (limit > loc_w->limit) 00137 { 00138 std::string estr="Iteration limit exceeds workspace "; 00139 estr+="in gsl_inte_qawf::qawf()."; 00140 O2SCL_ERR_RET(estr.c_str(),gsl_einval); 00141 } 00142 00143 /* Test on accuracy */ 00144 00145 if (epsabs <= 0 && (epsrel < 50 * GSL_DBL_EPSILON || epsrel < 0.5e-28)) 00146 { 00147 std::string estr="Tolerance cannot be achieved with given "; 00148 estr+="value of 'tolx' and 'tolf' in gsl_inte_qawo::qawo()."; 00149 O2SCL_ERR_RET(estr.c_str(),gsl_ebadtol); 00150 } 00151 00152 /* Perform the first integration */ 00153 00154 qc25f (func, a, b, wf, 0, 00155 &result0, &abserr0, &resabs0, &resasc0, pa); 00156 00157 this->set_initial_result (loc_w, result0, abserr0); 00158 00159 tolerance = GSL_MAX_DBL (epsabs, epsrel * fabs (result0)); 00160 00161 if (abserr0 <= 100 * GSL_DBL_EPSILON * resabs0 && 00162 abserr0 > tolerance) { 00163 *result = result0; 00164 *abserr = abserr0; 00165 00166 std::string estr="Cannot reach tolerance because of roundoff error "; 00167 estr+="on first attempt in gsl_inte_qawo::qawo()."; 00168 O2SCL_ERR_RET(estr.c_str(),gsl_eround); 00169 } else if ((abserr0 <= tolerance && abserr0 != resasc0) || 00170 abserr0 == 0.0) { 00171 *result = result0; 00172 *abserr = abserr0; 00173 00174 return GSL_SUCCESS; 00175 } else if (limit == 1) { 00176 *result = result0; 00177 *abserr = abserr0; 00178 00179 std::string estr="A maximum of 1 iteration was insufficient "; 00180 estr+="in gsl_inte_qawo::qawo()."; 00181 O2SCL_ERR_RET(estr.c_str(),gsl_emaxiter); 00182 } 00183 00184 /* Initialization */ 00185 00186 initialise_table (&table); 00187 00188 if (0.5 * abs_omega * fabs(b - a) <= 2) 00189 { 00190 append_table (&table, result0); 00191 extall = 1; 00192 } 00193 00194 area = result0; 00195 errsum = abserr0; 00196 00197 res_ext = result0; 00198 err_ext = GSL_DBL_MAX; 00199 00200 positive_integrand = this->test_positivity (result0, resabs0); 00201 00202 iteration = 1; 00203 00204 do 00205 { 00206 size_t current_level; 00207 double a1, b1, a2, b2; 00208 double a_i, b_i, r_i, e_i; 00209 double area1 = 0, area2 = 0, area12 = 0; 00210 double error1 = 0, error2 = 0, error12 = 0; 00211 double resasc1, resasc2; 00212 double resabs1, resabs2; 00213 double last_e_i; 00214 00215 /* Bisect the subinterval with the largest error estimate */ 00216 00217 this->retrieve (loc_w, &a_i, &b_i, &r_i, &e_i); 00218 00219 current_level = loc_w->level[loc_w->i] + 1; 00220 00221 if (current_level >= wf->n) 00222 { 00223 error_type = -1 ; /* exceeded limit of table */ 00224 break ; 00225 } 00226 00227 a1 = a_i; 00228 b1 = 0.5 * (a_i + b_i); 00229 a2 = b1; 00230 b2 = b_i; 00231 00232 iteration++; 00233 00234 qc25f (func, a1, b1, wf, current_level, 00235 &area1, &error1, &resabs1, &resasc1, pa); 00236 qc25f (func, a2, b2, wf, current_level, 00237 &area2, &error2, &resabs2, &resasc2, pa); 00238 00239 area12 = area1 + area2; 00240 error12 = error1 + error2; 00241 last_e_i = e_i; 00242 00243 /* Improve previous approximations to the integral and test for 00244 accuracy. 00245 00246 We write these expressions in the same way as the original 00247 QUADPACK code so that the rounding errors are the same, which 00248 makes testing easier. */ 00249 00250 errsum = errsum + error12 - e_i; 00251 area = area + area12 - r_i; 00252 00253 tolerance = GSL_MAX_DBL (epsabs, epsrel * fabs (area)); 00254 00255 if (resasc1 != error1 && resasc2 != error2) 00256 { 00257 double delta = r_i - area12; 00258 00259 if (fabs (delta) <= 1.0e-5 * fabs (area12) && 00260 error12 >= 0.99 * e_i) 00261 { 00262 if (!extrapolate) 00263 { 00264 roundoff_type1++; 00265 } 00266 else 00267 { 00268 roundoff_type2++; 00269 } 00270 } 00271 if (iteration > 10 && error12 > e_i) 00272 { 00273 roundoff_type3++; 00274 } 00275 } 00276 00277 /* Test for roundoff and eventually set error flag */ 00278 00279 if (roundoff_type1 + roundoff_type2 >= 10 || roundoff_type3 >= 20) 00280 { 00281 error_type = 2; /* round off error */ 00282 } 00283 00284 if (roundoff_type2 >= 5) 00285 { 00286 error_type2 = 1; 00287 } 00288 00289 /* set error flag in the case of bad integrand behaviour at 00290 a point of the integration range */ 00291 00292 if (this->subinterval_too_small (a1, a2, b2)) 00293 { 00294 error_type = 4; 00295 } 00296 00297 /* append the newly-created intervals to the list */ 00298 00299 this->update (loc_w, a1, b1, area1, error1, a2, b2, area2, error2); 00300 00301 if (errsum <= tolerance) 00302 { 00303 goto compute_result; 00304 } 00305 00306 if (error_type) 00307 { 00308 break; 00309 } 00310 00311 if (iteration >= limit - 1) 00312 { 00313 error_type = 1; 00314 break; 00315 } 00316 00317 /* set up variables on first iteration */ 00318 00319 if (iteration == 2 && extall) 00320 { 00321 error_over_large_intervals = errsum; 00322 ertest = tolerance; 00323 append_table (&table, area); 00324 continue; 00325 } 00326 00327 if (disallow_extrapolation) 00328 { 00329 continue; 00330 } 00331 00332 if (extall) 00333 { 00334 error_over_large_intervals += -last_e_i; 00335 00336 if (current_level < loc_w->maximum_level) 00337 { 00338 error_over_large_intervals += error12; 00339 } 00340 00341 if (extrapolate) 00342 goto label70; 00343 } 00344 00345 if (this->large_interval(loc_w)) 00346 { 00347 continue; 00348 } 00349 00350 if (extall) 00351 { 00352 extrapolate = 1; 00353 loc_w->nrmax = 1; 00354 } 00355 else 00356 { 00357 /* test whether the interval to be bisected next is the 00358 smallest interval. */ 00359 size_t i = loc_w->i; 00360 double width = loc_w->blist[i] - loc_w->alist[i]; 00361 00362 if (0.25 * fabs(width) * abs_omega > 2) 00363 continue; 00364 00365 extall = 1; 00366 error_over_large_intervals = errsum; 00367 ertest = tolerance; 00368 continue; 00369 } 00370 00371 label70: 00372 if (!error_type2 && error_over_large_intervals > ertest) 00373 { 00374 if (this->increase_nrmax (loc_w)) 00375 continue; 00376 } 00377 00378 /* Perform extrapolation */ 00379 00380 append_table (&table, area); 00381 00382 if (table.n < 3) 00383 { 00384 this->reset_nrmax(loc_w); 00385 extrapolate = 0; 00386 error_over_large_intervals = errsum; 00387 continue; 00388 } 00389 00390 qelg (&table, &reseps, &abseps); 00391 00392 ktmin++; 00393 00394 if (ktmin > 5 && err_ext < 0.001 * errsum) 00395 { 00396 error_type = 5; 00397 } 00398 00399 if (abseps < err_ext) 00400 { 00401 ktmin = 0; 00402 err_ext = abseps; 00403 res_ext = reseps; 00404 correc = error_over_large_intervals; 00405 ertest = GSL_MAX_DBL (epsabs, epsrel * fabs (reseps)); 00406 if (err_ext <= ertest) 00407 break; 00408 } 00409 00410 /* Prepare bisection of the smallest interval. */ 00411 00412 if (table.n == 1) 00413 { 00414 disallow_extrapolation = 1; 00415 } 00416 00417 if (error_type == 5) 00418 { 00419 break; 00420 } 00421 00422 /* work on interval with largest error */ 00423 00424 this->reset_nrmax (loc_w); 00425 extrapolate = 0; 00426 error_over_large_intervals = errsum; 00427 00428 } 00429 while (iteration < limit); 00430 00431 *result = res_ext; 00432 *abserr = err_ext; 00433 00434 if (err_ext == GSL_DBL_MAX) 00435 goto compute_result; 00436 00437 if (error_type || error_type2) 00438 { 00439 if (error_type2) 00440 { 00441 err_ext += correc; 00442 } 00443 00444 if (error_type == 0) 00445 error_type = 3; 00446 00447 if (result != 0 && area != 0) 00448 { 00449 if (err_ext / fabs (res_ext) > errsum / fabs (area)) 00450 goto compute_result; 00451 } 00452 else if (err_ext > errsum) 00453 { 00454 goto compute_result; 00455 } 00456 else if (area == 0.0) 00457 { 00458 goto return_error; 00459 } 00460 } 00461 00462 /* Test on divergence. */ 00463 00464 { 00465 double max_area = GSL_MAX_DBL (fabs (res_ext), fabs (area)); 00466 00467 if (!positive_integrand && max_area < 0.01 * resabs0) 00468 goto return_error; 00469 } 00470 00471 { 00472 double ratio = res_ext / area; 00473 00474 if (ratio < 0.01 || ratio > 100 || errsum > fabs (area)) 00475 error_type = 6; 00476 } 00477 00478 goto return_error; 00479 00480 compute_result: 00481 00482 *result = this->sum_results (loc_w); 00483 *abserr = errsum; 00484 00485 return_error: 00486 00487 if (error_type > 2) 00488 error_type--; 00489 00490 if (error_type == 0) 00491 { 00492 return GSL_SUCCESS; 00493 } 00494 else if (error_type == 1) 00495 { 00496 std::string estr="Number of iterations was insufficient "; 00497 estr+=" in gsl_inte_qawo::qawo()."; 00498 O2SCL_ERR_RET(estr.c_str(),gsl_emaxiter); 00499 } 00500 else if (error_type == 2) 00501 { 00502 std::string estr="Roundoff error prevents tolerance "; 00503 estr+="from being achieved in gsl_inte_qawo::qawo()."; 00504 O2SCL_ERR_RET(estr.c_str(),gsl_eround); 00505 } 00506 else if (error_type == 3) 00507 { 00508 std::string estr="Bad integrand behavior "; 00509 estr+=" in gsl_inte_qawo::qawo()."; 00510 O2SCL_ERR_RET(estr.c_str(),gsl_esing); 00511 } 00512 else if (error_type == 4) 00513 { 00514 std::string estr="Roundoff error detected in extrapolation table "; 00515 estr+="in gsl_inte_qawo::qawo()."; 00516 O2SCL_ERR_RET(estr.c_str(),gsl_eround); 00517 } 00518 else if (error_type == 5) 00519 { 00520 std::string estr="Integral is divergent or slowly convergent "; 00521 estr+="in gsl_inte_qawo::qawo()."; 00522 O2SCL_ERR_RET(estr.c_str(),gsl_ediverge); 00523 } 00524 else if (error_type == -1) 00525 { 00526 std::string estr="Exceeded limit of trigonometric table "; 00527 estr+="gsl_inte_qawo::qawo()"; 00528 O2SCL_ERR_RET(estr.c_str(),gsl_etable); 00529 } 00530 else 00531 { 00532 std::string estr="Could not integrate function in gsl_inte_qawo"; 00533 estr+="::qawo() (it may have returned a non-finite result)."; 00534 O2SCL_ERR_RET(estr.c_str(),gsl_efailed); 00535 } 00536 00537 } 00538 00539 /// 25-point quadrature for oscillating functions 00540 void qc25f(func_t &func, double a, double b, 00541 gsl_integration_qawo_table *wf, size_t level, 00542 double *result, double *abserr, double *resabs, 00543 double *resasc, param_t &pa) { 00544 00545 const double center = 0.5 * (a + b); 00546 const double half_length = 0.5 * (b - a); 00547 //const double omega = wf->omega ; 00548 00549 const double par = omega * half_length; 00550 00551 if (fabs (par) < 2) 00552 { 00553 00554 double fv1[8], fv2[8]; 00555 gsl_integration_qk_o2scl(func,8,o2scl_inte_qag_coeffs::qk15_xgk, 00556 o2scl_inte_qag_coeffs::qk15_wg, 00557 o2scl_inte_qag_coeffs::qk15_wgk, 00558 fv1,fv2,a,b,result,abserr,resabs,resasc, 00559 pa); 00560 00561 return; 00562 } 00563 else 00564 { 00565 double *moment; 00566 double cheb12[13], cheb24[25]; 00567 double result_abs, res12_cos, res12_sin, res24_cos, res24_sin; 00568 double est_cos, est_sin; 00569 double c, s; 00570 size_t i; 00571 00572 this->gsl_integration_qcheb (func, a, b, cheb12, cheb24, pa); 00573 00574 if (level >= wf->n) { 00575 /* table overflow should not happen, check before calling */ 00576 O2SCL_ERR("Table overflow in gsl_inte_qawo::qc25f().", 00577 gsl_esanity); 00578 return; 00579 } 00580 00581 /* obtain moments from the table */ 00582 00583 moment = wf->chebmo + 25 * level; 00584 00585 res12_cos = cheb12[12] * moment[12]; 00586 res12_sin = 0 ; 00587 00588 for (i = 0; i < 6; i++) 00589 { 00590 size_t k = 10 - 2 * i; 00591 res12_cos += cheb12[k] * moment[k]; 00592 res12_sin += cheb12[k + 1] * moment[k + 1]; 00593 } 00594 00595 res24_cos = cheb24[24] * moment[24]; 00596 res24_sin = 0 ; 00597 00598 result_abs = fabs(cheb24[24]) ; 00599 00600 for (i = 0; i < 12; i++) 00601 { 00602 size_t k = 22 - 2 * i; 00603 res24_cos += cheb24[k] * moment[k]; 00604 res24_sin += cheb24[k + 1] * moment[k + 1]; 00605 result_abs += fabs(cheb24[k]) + fabs(cheb24[k+1]); 00606 } 00607 00608 est_cos = fabs(res24_cos - res12_cos); 00609 est_sin = fabs(res24_sin - res12_sin); 00610 00611 c = half_length * cos(center * omega); 00612 s = half_length * sin(center * omega); 00613 00614 if (wf->sine == GSL_INTEG_SINE) 00615 { 00616 *result = c * res24_sin + s * res24_cos; 00617 *abserr = fabs(c * est_sin) + fabs(s * est_cos); 00618 } 00619 else 00620 { 00621 *result = c * res24_cos - s * res24_sin; 00622 *abserr = fabs(c * est_cos) + fabs(s * est_sin); 00623 } 00624 00625 *resabs = result_abs * half_length; 00626 *resasc = GSL_DBL_MAX; 00627 00628 return; 00629 } 00630 } 00631 00632 /// Add the oscillating part to the integrand 00633 virtual double transform(func_t &func, double x, param_t &pa) { 00634 double wx = omega * x; 00635 double sinwx = sin(wx) ; 00636 double y; 00637 func(x,y,pa); 00638 return y*sinwx; 00639 } 00640 00641 #endif 00642 00643 /// Return string denoting type ("gsl_inte_qawo_sin") 00644 const char *type() { return "gsl_inte_qawo_sin"; } 00645 00646 }; 00647 00648 /** \brief Adaptive integration a function with finite limits of 00649 integration (GSL) 00650 00651 The number of subdivisions of the original interval which 00652 this class is allowed to make is dictated by the workspace 00653 size for the integration class, which can be set using 00654 \ref gsl_inte_table::set_wkspace() . 00655 00656 \todo Verbose output has been setup for this class, but this 00657 needs to be done for the other GSL-like integrators 00658 */ 00659 template<class param_t, class func_t> class gsl_inte_qawo_cos : 00660 public gsl_inte_qawo_sin<param_t,func_t> { 00661 00662 public: 00663 00664 gsl_inte_qawo_cos() { 00665 } 00666 00667 virtual ~gsl_inte_qawo_cos() {} 00668 00669 /** \brief Integrate function \c func from \c a to \c b and place 00670 the result in \c res and the error in \c err 00671 */ 00672 virtual int integ_err(func_t &func, double a, double b, 00673 param_t &pa, double &res, double &err2) { 00674 00675 this->otable=gsl_integration_qawo_table_alloc 00676 (this->omega,b-a,GSL_INTEG_COSINE,this->tab_size); 00677 00678 int status=qawo(func,a,this->tolx,this->tolf,this->wkspace, 00679 this->w,this->otable,&res,&err2,pa); 00680 00681 gsl_integration_qawo_table_free(this->otable); 00682 00683 return status; 00684 00685 } 00686 00687 #ifndef DOXYGEN_INTERNAL 00688 00689 protected: 00690 00691 /// Add the oscillating part to the integrand 00692 virtual double transform(func_t &func, double x, param_t &pa) { 00693 double wx = this->omega * x; 00694 double coswx = cos(wx) ; 00695 double y; 00696 func(x,y,pa); 00697 return y*coswx; 00698 } 00699 00700 #endif 00701 00702 /// Return string denoting type ("gsl_inte_qawo_cos") 00703 const char *type() { return "gsl_inte_qawo_cos"; } 00704 00705 }; 00706 00707 #ifndef DOXYGENP 00708 } 00709 #endif 00710 00711 #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