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