gsl_integration__qag.c 6.8 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266
  1. /* integration/qag.c
  2. *
  3. * Copyright (C) 1996, 1997, 1998, 1999, 2000, 2007 Brian Gough
  4. *
  5. * This program is free software; you can redistribute it and/or modify
  6. * it under the terms of the GNU General Public License as published by
  7. * the Free Software Foundation; either version 3 of the License, or (at
  8. * your option) any later version.
  9. *
  10. * This program is distributed in the hope that it will be useful, but
  11. * WITHOUT ANY WARRANTY; without even the implied warranty of
  12. * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
  13. * General Public License for more details.
  14. *
  15. * You should have received a copy of the GNU General Public License
  16. * along with this program; if not, write to the Free Software
  17. * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
  18. */
  19. #include "gsl__config.h"
  20. #include <stdlib.h>
  21. #include "gsl_math.h"
  22. #include "gsl_errno.h"
  23. #include "gsl_integration.h"
  24. #include "gsl_integration__initialise.c"
  25. #include "gsl_integration__set_initial.c"
  26. #include "gsl_integration__qpsrt.c"
  27. #include "gsl_integration__util.c"
  28. static int
  29. qag (const gsl_function *f,
  30. const double a, const double b,
  31. const double epsabs, const double epsrel,
  32. const size_t limit,
  33. gsl_integration_workspace * workspace,
  34. double * result, double * abserr,
  35. gsl_integration_rule * q) ;
  36. int
  37. gsl_integration_qag (const gsl_function *f,
  38. double a, double b,
  39. double epsabs, double epsrel, size_t limit,
  40. int key,
  41. gsl_integration_workspace * workspace,
  42. double * result, double * abserr)
  43. {
  44. int status ;
  45. gsl_integration_rule * integration_rule = gsl_integration_qk15 ;
  46. if (key < GSL_INTEG_GAUSS15)
  47. {
  48. key = GSL_INTEG_GAUSS15 ;
  49. }
  50. else if (key > GSL_INTEG_GAUSS61)
  51. {
  52. key = GSL_INTEG_GAUSS61 ;
  53. }
  54. switch (key)
  55. {
  56. case GSL_INTEG_GAUSS15:
  57. integration_rule = gsl_integration_qk15 ;
  58. break ;
  59. case GSL_INTEG_GAUSS21:
  60. integration_rule = gsl_integration_qk21 ;
  61. break ;
  62. case GSL_INTEG_GAUSS31:
  63. integration_rule = gsl_integration_qk31 ;
  64. break ;
  65. case GSL_INTEG_GAUSS41:
  66. integration_rule = gsl_integration_qk41 ;
  67. break ;
  68. case GSL_INTEG_GAUSS51:
  69. integration_rule = gsl_integration_qk51 ;
  70. break ;
  71. case GSL_INTEG_GAUSS61:
  72. integration_rule = gsl_integration_qk61 ;
  73. break ;
  74. default:
  75. GSL_ERROR("value of key does specify a known integration rule",
  76. GSL_EINVAL) ;
  77. }
  78. status = qag (f, a, b, epsabs, epsrel, limit,
  79. workspace,
  80. result, abserr,
  81. integration_rule) ;
  82. return status ;
  83. }
  84. static int
  85. qag (const gsl_function * f,
  86. const double a, const double b,
  87. const double epsabs, const double epsrel,
  88. const size_t limit,
  89. gsl_integration_workspace * workspace,
  90. double *result, double *abserr,
  91. gsl_integration_rule * q)
  92. {
  93. double area, errsum;
  94. double result0, abserr0, resabs0, resasc0;
  95. double tolerance;
  96. size_t iteration = 0;
  97. int roundoff_type1 = 0, roundoff_type2 = 0, error_type = 0;
  98. double round_off;
  99. /* Initialize results */
  100. initialise (workspace, a, b);
  101. *result = 0;
  102. *abserr = 0;
  103. if (limit > workspace->limit)
  104. {
  105. GSL_ERROR ("iteration limit exceeds available workspace", GSL_EINVAL) ;
  106. }
  107. if (epsabs <= 0 && (epsrel < 50 * GSL_DBL_EPSILON || epsrel < 0.5e-28))
  108. {
  109. GSL_ERROR ("tolerance cannot be acheived with given epsabs and epsrel",
  110. GSL_EBADTOL);
  111. }
  112. /* perform the first integration */
  113. q (f, a, b, &result0, &abserr0, &resabs0, &resasc0);
  114. set_initial_result (workspace, result0, abserr0);
  115. /* Test on accuracy */
  116. tolerance = GSL_MAX_DBL (epsabs, epsrel * fabs (result0));
  117. /* need IEEE rounding here to match original quadpack behavior */
  118. round_off = GSL_COERCE_DBL (50 * GSL_DBL_EPSILON * resabs0);
  119. if (abserr0 <= round_off && abserr0 > tolerance)
  120. {
  121. *result = result0;
  122. *abserr = abserr0;
  123. GSL_ERROR ("cannot reach tolerance because of roundoff error "
  124. "on first attempt", GSL_EROUND);
  125. }
  126. else if ((abserr0 <= tolerance && abserr0 != resasc0) || abserr0 == 0.0)
  127. {
  128. *result = result0;
  129. *abserr = abserr0;
  130. return GSL_SUCCESS;
  131. }
  132. else if (limit == 1)
  133. {
  134. *result = result0;
  135. *abserr = abserr0;
  136. GSL_ERROR ("a maximum of one iteration was insufficient", GSL_EMAXITER);
  137. }
  138. area = result0;
  139. errsum = abserr0;
  140. iteration = 1;
  141. do
  142. {
  143. double a1, b1, a2, b2;
  144. double a_i, b_i, r_i, e_i;
  145. double area1 = 0, area2 = 0, area12 = 0;
  146. double error1 = 0, error2 = 0, error12 = 0;
  147. double resasc1, resasc2;
  148. double resabs1, resabs2;
  149. /* Bisect the subinterval with the largest error estimate */
  150. retrieve (workspace, &a_i, &b_i, &r_i, &e_i);
  151. a1 = a_i;
  152. b1 = 0.5 * (a_i + b_i);
  153. a2 = b1;
  154. b2 = b_i;
  155. q (f, a1, b1, &area1, &error1, &resabs1, &resasc1);
  156. q (f, a2, b2, &area2, &error2, &resabs2, &resasc2);
  157. area12 = area1 + area2;
  158. error12 = error1 + error2;
  159. errsum += (error12 - e_i);
  160. area += area12 - r_i;
  161. if (resasc1 != error1 && resasc2 != error2)
  162. {
  163. double delta = r_i - area12;
  164. if (fabs (delta) <= 1.0e-5 * fabs (area12) && error12 >= 0.99 * e_i)
  165. {
  166. roundoff_type1++;
  167. }
  168. if (iteration >= 10 && error12 > e_i)
  169. {
  170. roundoff_type2++;
  171. }
  172. }
  173. tolerance = GSL_MAX_DBL (epsabs, epsrel * fabs (area));
  174. if (errsum > tolerance)
  175. {
  176. if (roundoff_type1 >= 6 || roundoff_type2 >= 20)
  177. {
  178. error_type = 2; /* round off error */
  179. }
  180. /* set error flag in the case of bad integrand behaviour at
  181. a point of the integration range */
  182. if (subinterval_too_small (a1, a2, b2))
  183. {
  184. error_type = 3;
  185. }
  186. }
  187. update (workspace, a1, b1, area1, error1, a2, b2, area2, error2);
  188. retrieve (workspace, &a_i, &b_i, &r_i, &e_i);
  189. iteration++;
  190. }
  191. while (iteration < limit && !error_type && errsum > tolerance);
  192. *result = sum_results (workspace);
  193. *abserr = errsum;
  194. if (errsum <= tolerance)
  195. {
  196. return GSL_SUCCESS;
  197. }
  198. else if (error_type == 2)
  199. {
  200. GSL_ERROR ("roundoff error prevents tolerance from being achieved",
  201. GSL_EROUND);
  202. }
  203. else if (error_type == 3)
  204. {
  205. GSL_ERROR ("bad integrand behavior found in the integration interval",
  206. GSL_ESING);
  207. }
  208. else if (iteration == limit)
  209. {
  210. GSL_ERROR ("maximum number of subdivisions reached", GSL_EMAXITER);
  211. }
  212. else
  213. {
  214. GSL_ERROR ("could not integrate function", GSL_EFAILED);
  215. }
  216. }