gsl_inte_qawo.h

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

Documentation generated with Doxygen and provided under the GNU Free Documentation License. See License Information for details.

Project hosting provided by SourceForge.net Logo, O2scl Sourceforge Project Page