24 #ifndef O2SCL_GSL_INTE_SINGULAR_H
25 #define O2SCL_GSL_INTE_SINGULAR_H
33 #include <gsl/gsl_integration.h>
35 #include <o2scl/string_conv.h>
36 #include <o2scl/inte.h>
37 #include <o2scl/inte_kronrod_gsl.h>
39 #ifndef DOXYGEN_NO_O2NS
95 double dbl_eps=std::numeric_limits<double>::epsilon();
96 int status=(fabs(result) >= (1-50*dbl_eps)*resabs);
144 double *epstab = table->
rlist2;
145 double *res3la = table->
res3la;
146 const size_t n = table->
n - 1;
148 const double current = epstab[n];
150 double absolute = GSL_DBL_MAX;
151 double relative = 5 * GSL_DBL_EPSILON * fabs (current);
153 const size_t newelm = n / 2;
154 const size_t n_orig = n;
158 const size_t nres_orig = table->
nres;
161 *abserr = GSL_DBL_MAX;
165 *abserr = GSL_MAX_DBL (absolute, relative);
169 epstab[n + 2] = epstab[n];
170 epstab[n] = GSL_DBL_MAX;
172 for (i = 0; i < newelm; i++) {
173 double res = epstab[n - 2 * i + 2];
174 double e0 = epstab[n - 2 * i - 2];
175 double e1 = epstab[n - 2 * i - 1];
178 double e1abs = fabs (e1);
179 double delta2 = e2 - e1;
180 double err2 = fabs (delta2);
181 double tol2 = GSL_MAX_DBL (fabs (e2), e1abs) * GSL_DBL_EPSILON;
182 double delta3 = e1 - e0;
183 double err3 = fabs (delta3);
184 double tol3 = GSL_MAX_DBL (e1abs, fabs (e0)) * GSL_DBL_EPSILON;
186 double e3, delta1, err1, tol1, ss;
188 if (err2 <= tol2 && err3 <= tol3) {
193 absolute = err2 + err3;
194 relative = 5 * GSL_DBL_EPSILON * fabs (res);
195 *abserr = GSL_MAX_DBL (absolute, relative);
199 e3 = epstab[n - 2 * i];
200 epstab[n - 2 * i] = e1;
202 err1 = fabs (delta1);
203 tol1 = GSL_MAX_DBL (e1abs, fabs (e3)) * GSL_DBL_EPSILON;
208 if (err1 <= tol1 || err2 <= tol2 || err3 <= tol3) {
213 ss = (1 / delta1 + 1 / delta2) - 1 / delta3;
218 if (fabs (ss * e1) <= 0.0001) {
226 epstab[n - 2 * i] = res;
229 const double error = err2 + fabs (res - e2) + err3;
231 if (error <= *abserr) {
241 const size_t limexp = 50 - 1;
243 if (n_final == limexp) {
244 n_final = 2 * (limexp / 2);
248 if (n_orig % 2 == 1) {
249 for (i = 0; i <= newelm; i++) {
250 epstab[1 + i * 2] = epstab[i * 2 + 3];
253 for (i = 0; i <= newelm; i++) {
254 epstab[i * 2] = epstab[i * 2 + 2];
257 if (n_orig != n_final) {
258 for (i = 0; i <= n_final; i++) {
259 epstab[i] = epstab[n_orig - n_final + i];
263 table->
n = n_final + 1;
267 res3la[nres_orig] = *result;
268 *abserr = GSL_DBL_MAX;
272 *abserr = (fabs (*result - res3la[2]) + fabs (*result - res3la[1])
273 + fabs (*result - res3la[0]));
275 res3la[0] = res3la[1];
276 res3la[1] = res3la[2];
286 table->
nres = nres_orig + 1;
288 *abserr = GSL_MAX_DBL (*abserr, 5 * GSL_DBL_EPSILON * fabs (*result));
295 size_t i = workspace->
i ;
296 const size_t * level = workspace->
level;
307 workspace->
nrmax = 0;
308 workspace->
i = workspace->
order[0] ;
314 int id = workspace->
nrmax;
317 const size_t * level = workspace->
level;
318 const size_t * order = workspace->
order;
320 size_t limit = workspace->
limit ;
321 size_t last = workspace->
size - 1 ;
323 if (last > (1 + limit / 2)) {
324 jupbnd = limit + 1 - last;
329 for (k =
id; k <= jupbnd; k++) {
330 size_t i_max = order[workspace->
nrmax];
332 workspace->
i = i_max ;
349 int qags(func_t &func,
const double a,
const double b,
350 const double l_epsabs,
const double l_epsrel,
351 double *result,
double *abserr) {
354 double res_ext, err_ext;
355 double result0, abserr0, resabs0, resasc0;
359 double error_over_large_intervals = 0;
360 double reseps = 0, abseps = 0, correc = 0;
362 int roundoff_type1 = 0, roundoff_type2 = 0, roundoff_type3 = 0;
363 int error_type = 0, error_type2 = 0;
365 size_t iteration = 0;
367 int positive_integrand = 0;
369 int disallow_extrapolation = 0;
380 size_t limit=this->
w->
limit;
384 if (this->
tol_abs <= 0 && (this->
tol_rel < 50 * GSL_DBL_EPSILON ||
388 std::string estr=
"Tolerance cannot be achieved with given ";
389 estr+=
"value of tol_abs, "+
dtos(l_epsabs)+
", and tol_rel, "+
390 dtos(l_epsrel)+
", in inte_singular_gsl::qags().";
396 this->
gauss_kronrod(func,a,b,&result0,&abserr0,&resabs0,&resasc0);
400 tolerance = GSL_MAX_DBL (this->
tol_abs, this->
tol_rel * fabs (result0));
402 if (abserr0 <= 100 * GSL_DBL_EPSILON * resabs0 &&
403 abserr0 > tolerance) {
410 std::string estr=
"Cannot reach tolerance because of roundoff error ";
411 estr+=
"on first attempt in inte_singular_gsl::qags().";
414 }
else if ((abserr0 <= tolerance &&
415 abserr0 != resasc0) || abserr0 == 0.0) {
422 }
else if (limit == 1) {
429 "in inte_singular_gsl::qags().",
442 err_ext = GSL_DBL_MAX;
452 std::cout << this->
type();
453 std::cout <<
" Iter: " << iteration;
454 std::cout.setf(std::ios::showpos);
455 std::cout <<
" Res: " << area;
456 std::cout.unsetf(std::ios::showpos);
457 std::cout <<
" Err: " << errsum
458 <<
" Tol: " << tolerance << std::endl;
461 std::cout <<
"Press a key and type enter to continue. " ;
466 size_t current_level;
467 double a1, b1, a2, b2;
468 double a_i, b_i, r_i, e_i;
469 double area1 = 0, area2 = 0, area12 = 0;
470 double error1 = 0, error2 = 0, error12 = 0;
471 double resasc1, resasc2;
472 double resabs1, resabs2;
477 this->
w->
retrieve (&a_i, &b_i, &r_i, &e_i);
479 current_level = this->
w->
level[this->
w->
i] + 1;
482 b1 = 0.5 * (a_i + b_i);
488 this->
gauss_kronrod(func,a1,b1,&area1,&error1,&resabs1,&resasc1);
489 this->
gauss_kronrod(func,a2,b2,&area2,&error2,&resabs2,&resasc2);
491 area12 = area1 + area2;
492 error12 = error1 + error2;
503 errsum = errsum + error12 - e_i;
504 area = area + area12 - r_i;
506 tolerance = GSL_MAX_DBL (this->
tol_abs, this->
tol_rel * fabs (area));
507 if (resasc1 != error1 && resasc2 != error2) {
508 double delta = r_i - area12;
510 if (fabs (delta) <= 1.0e-5 * fabs (area12) &&
511 error12 >= 0.99 * e_i) {
518 if (iteration > 10 && error12 > e_i) {
525 if (roundoff_type1 + roundoff_type2 >= 10 ||
526 roundoff_type3 >= 20) {
531 if (roundoff_type2 >= 5) {
544 this->
w->
update(a1,b1,area1,error1,a2,b2,area2,error2);
546 if (errsum <= tolerance) {
550 std::cout << this->
type();
551 std::cout <<
" Iter: " << iteration;
552 std::cout.setf(std::ios::showpos);
553 std::cout <<
" Res: " << area;
554 std::cout.unsetf(std::ios::showpos);
555 std::cout <<
" Err: " << errsum
556 <<
" Tol: " << tolerance << std::endl;
559 std::cout <<
"Press a key and type enter to continue. " ;
571 if (iteration >= limit - 1) {
576 if (iteration == 2) {
577 error_over_large_intervals = errsum;
583 if (disallow_extrapolation) {
587 error_over_large_intervals += -last_e_i;
590 error_over_large_intervals += error12;
606 if (!error_type2 && error_over_large_intervals > ertest) {
616 qelg(&table,&reseps,&abseps);
620 if (ktmin > 5 && err_ext < 0.001 * errsum) {
624 if (abseps < err_ext) {
628 correc = error_over_large_intervals;
629 ertest = GSL_MAX_DBL (this->
tol_abs,
630 this->
tol_rel * fabs (reseps));
631 if (err_ext <= ertest) {
639 disallow_extrapolation = 1;
642 if (error_type == 5) {
650 error_over_large_intervals = errsum;
652 }
while (iteration < limit);
657 if (err_ext == GSL_DBL_MAX)
660 if (error_type || error_type2) {
668 if (res_ext != 0.0 && area != 0.0) {
669 if (err_ext / fabs (res_ext) > errsum / fabs (area))
671 }
else if (err_ext > errsum) {
673 }
else if (area == 0.0) {
681 double max_area = GSL_MAX_DBL (fabs (res_ext),fabs (area));
683 if (!positive_integrand && max_area < 0.01 * resabs0)
688 double ratio = res_ext / area;
690 if (ratio < 0.01 || ratio > 100.0 || errsum > fabs (area))
703 if (error_type > 2) {
709 if (error_type == 0) {
711 }
else if (error_type == 1) {
712 std::string estr=
"Maximum number of subdivisions ("+
itos(iteration);
713 estr+=
") reached in inte_singular_gsl::qags().";
715 }
else if (error_type == 2) {
716 std::string estr=
"Roundoff error prevents tolerance ";
717 estr+=
"from being achieved in inte_singular_gsl::qags().";
719 }
else if (error_type == 3) {
720 std::string estr=
"Bad integrand behavior ";
721 estr+=
"in inte_singular_gsl::qags().";
723 }
else if (error_type == 4) {
724 std::string estr=
"Roundoff error detected in extrapolation table ";
725 estr+=
"in inte_singular_gsl::qags().";
727 }
else if (error_type == 5) {
728 std::string estr=
"Integral is divergent or slowly convergent ";
729 estr+=
"in inte_singular_gsl::qags().";
733 std::string estr=
"Could not integrate function in inte_kronrod_gsl";
734 estr+=
"::qags() (it may have returned a non-finite result).";
756 virtual double transform(
double t, func_t &func)=0;
762 (func_t &func,
double a,
double b,
763 double *result,
double *abserr,
double *resabs,
double *resasc) {
765 funct11 fmp=std::bind(std::mem_fn<
double(
double,func_t &)>
767 this,std::placeholders::_1,func);
770 (fmp,a,b,result,abserr,resabs,resasc);
775 #ifndef DOXYGEN_NO_O2NS
Integration workspace for the GSL integrators.
void initialise_table(struct extrapolation_table *table)
Initialize the table.
#define O2SCL_CONV_RET(d, n, b)
Set a "convergence" error and return the error value.
int qags(func_t &func, const double a, const double b, const double l_epsabs, const double l_epsrel, double *result, double *abserr)
Integration function.
std::function< double(double)> funct11
One-dimensional function typedef.
size_t * level
Numbers of subdivisions made.
void append_table(struct extrapolation_table *table, double y)
Append a result to the table.
struct o2scl::inte_singular_gsl::extrapolation_table extrap_table
A structure for extrapolation for inte_qags_gsl.
apparent singularity detected
exceeded max number of iterations
Base class for integrating a function with a singularity (GSL)
int large_interval(inte_workspace_gsl *workspace)
Determine if an interval is large.
#define O2SCL_CONV2_RET(d, d2, n, b)
Set an error and return the error value, two-string version.
size_t * order
Linear ordering vector for sort routine.
size_t last_iter
The most recent number of iterations taken.
bool err_nonconv
If true, call the error handler if the routine does not converge or reach the desired tolerance (defa...
int update(double a1, double b1, double area1, double error1, double a2, double b2, double area2, double error2)
Determine which new subinterval to add to the workspace stack and perform update. ...
size_t i
Index of current subinterval.
virtual void gauss_kronrod(func_t &func, double a, double b, double *result, double *abserr, double *resabs, double *resasc)
Integration wrapper for user-specified function type.
int subinterval_too_small(double a1, double a2, double b2)
Test whether the proposed subdivision falls before floating-point precision.
double tol_abs
The maximum absolute uncertainty in the value of the integral (default )
int test_positivity(double result, double resabs)
Test if the integrand satisfies .
int initialise(double a, double b)
Initialize the workspace for an integration with limits a and b.
int set_initial_result(double result, double error)
Update the workspace with the result and error from the first integration.
std::string dtos(double x, int prec=6, bool auto_prec=false)
Convert a double to a string.
#define O2SCL_ERR(d, n)
Set an error with message d and code n.
size_t size
Current number of subintervals being used.
int retrieve(double *a, double *b, double *r, double *e) const
Retrieve the ith result from the workspace stack.
Basic Gauss-Kronrod integration class (GSL)
size_t nrmax
Counter for extrapolation routine.
void qelg(struct extrapolation_table *table, double *result, double *abserr)
Determines the limit of a given sequence of approximations.
inte_workspace_gsl * w
The integration workspace.
void gauss_kronrod_base(func2_t &func, double a, double b, double *result, double *abserr, double *resabs, double *resasc)
The base Gauss-Kronrod integration function template.
double sum_results()
Add up all of the contributions to construct the final result.
double tol_rel
The maximum relative uncertainty in the value of the integral (default )
void reset_nrmax(inte_workspace_gsl *workspace)
Reset workspace to work on the interval with the largest error.
int increase_nrmax(inte_workspace_gsl *workspace)
Increase workspace.
size_t limit
Maximum number of subintervals allocated.
size_t maximum_level
Depth of subdivisions reached.
std::string itos(int x)
Convert an integer to a string.
user specified an invalid tolerance
failed because of roundoff error
integral or series is divergent
virtual const char * type()
Return string denoting type ("inte")