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