gsl_integration__qawf.c 6.2 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282
  1. /* integration/qawf.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 <math.h>
  21. #include <float.h>
  22. #include "gsl_math.h"
  23. #include "gsl_errno.h"
  24. #include "gsl_integration.h"
  25. #include "gsl_integration__initialise.c"
  26. #include "gsl_integration__append.c"
  27. #include "gsl_integration__qelg.c"
  28. int
  29. gsl_integration_qawf (gsl_function * f,
  30. const double a,
  31. const double epsabs,
  32. const size_t limit,
  33. gsl_integration_workspace * workspace,
  34. gsl_integration_workspace * cycle_workspace,
  35. gsl_integration_qawo_table * wf,
  36. double *result, double *abserr)
  37. {
  38. double area, errsum;
  39. double res_ext, err_ext;
  40. double correc, total_error = 0.0, truncation_error;
  41. size_t ktmin = 0;
  42. size_t iteration = 0;
  43. struct extrapolation_table table;
  44. double cycle;
  45. double omega = wf->omega;
  46. const double p = 0.9;
  47. double factor = 1;
  48. double initial_eps, eps;
  49. int error_type = 0;
  50. /* Initialize results */
  51. initialise (workspace, a, a);
  52. *result = 0;
  53. *abserr = 0;
  54. if (limit > workspace->limit)
  55. {
  56. GSL_ERROR ("iteration limit exceeds available workspace", GSL_EINVAL) ;
  57. }
  58. /* Test on accuracy */
  59. if (epsabs <= 0)
  60. {
  61. GSL_ERROR ("absolute tolerance epsabs must be positive", GSL_EBADTOL) ;
  62. }
  63. if (omega == 0.0)
  64. {
  65. if (wf->sine == GSL_INTEG_SINE)
  66. {
  67. /* The function sin(w x) f(x) is always zero for w = 0 */
  68. *result = 0;
  69. *abserr = 0;
  70. return GSL_SUCCESS;
  71. }
  72. else
  73. {
  74. /* The function cos(w x) f(x) is always f(x) for w = 0 */
  75. int status = gsl_integration_qagiu (f, a, epsabs, 0.0,
  76. cycle_workspace->limit,
  77. cycle_workspace,
  78. result, abserr);
  79. return status;
  80. }
  81. }
  82. if (epsabs > GSL_DBL_MIN / (1 - p))
  83. {
  84. eps = epsabs * (1 - p);
  85. }
  86. else
  87. {
  88. eps = epsabs;
  89. }
  90. initial_eps = eps;
  91. area = 0;
  92. errsum = 0;
  93. res_ext = 0;
  94. err_ext = GSL_DBL_MAX;
  95. correc = 0;
  96. cycle = (2 * floor (fabs (omega)) + 1) * M_PI / fabs (omega);
  97. gsl_integration_qawo_table_set_length (wf, cycle);
  98. initialise_table (&table);
  99. for (iteration = 0; iteration < limit; iteration++)
  100. {
  101. double area1, error1, reseps, erreps;
  102. double a1 = a + iteration * cycle;
  103. double b1 = a1 + cycle;
  104. double epsabs1 = eps * factor;
  105. int status = gsl_integration_qawo (f, a1, epsabs1, 0.0, limit,
  106. cycle_workspace, wf,
  107. &area1, &error1);
  108. append_interval (workspace, a1, b1, area1, error1);
  109. factor *= p;
  110. area = area + area1;
  111. errsum = errsum + error1;
  112. /* estimate the truncation error as 50 times the final term */
  113. truncation_error = 50 * fabs (area1);
  114. total_error = errsum + truncation_error;
  115. if (total_error < epsabs && iteration > 4)
  116. {
  117. goto compute_result;
  118. }
  119. if (error1 > correc)
  120. {
  121. correc = error1;
  122. }
  123. if (status)
  124. {
  125. eps = GSL_MAX_DBL (initial_eps, correc * (1.0 - p));
  126. }
  127. if (status && total_error < 10 * correc && iteration > 3)
  128. {
  129. goto compute_result;
  130. }
  131. append_table (&table, area);
  132. if (table.n < 2)
  133. {
  134. continue;
  135. }
  136. qelg (&table, &reseps, &erreps);
  137. ktmin++;
  138. if (ktmin >= 15 && err_ext < 0.001 * total_error)
  139. {
  140. error_type = 4;
  141. }
  142. if (erreps < err_ext)
  143. {
  144. ktmin = 0;
  145. err_ext = erreps;
  146. res_ext = reseps;
  147. if (err_ext + 10 * correc <= epsabs)
  148. break;
  149. if (err_ext <= epsabs && 10 * correc >= epsabs)
  150. break;
  151. }
  152. }
  153. if (iteration == limit)
  154. error_type = 1;
  155. if (err_ext == GSL_DBL_MAX)
  156. goto compute_result;
  157. err_ext = err_ext + 10 * correc;
  158. *result = res_ext;
  159. *abserr = err_ext;
  160. if (error_type == 0)
  161. {
  162. return GSL_SUCCESS ;
  163. }
  164. if (res_ext != 0.0 && area != 0.0)
  165. {
  166. if (err_ext / fabs (res_ext) > errsum / fabs (area))
  167. goto compute_result;
  168. }
  169. else if (err_ext > errsum)
  170. {
  171. goto compute_result;
  172. }
  173. else if (area == 0.0)
  174. {
  175. goto return_error;
  176. }
  177. if (error_type == 4)
  178. {
  179. err_ext = err_ext + truncation_error;
  180. }
  181. goto return_error;
  182. compute_result:
  183. *result = area;
  184. *abserr = total_error;
  185. return_error:
  186. if (error_type > 2)
  187. error_type--;
  188. if (error_type == 0)
  189. {
  190. return GSL_SUCCESS;
  191. }
  192. else if (error_type == 1)
  193. {
  194. GSL_ERROR ("number of iterations was insufficient", GSL_EMAXITER);
  195. }
  196. else if (error_type == 2)
  197. {
  198. GSL_ERROR ("cannot reach tolerance because of roundoff error",
  199. GSL_EROUND);
  200. }
  201. else if (error_type == 3)
  202. {
  203. GSL_ERROR ("bad integrand behavior found in the integration interval",
  204. GSL_ESING);
  205. }
  206. else if (error_type == 4)
  207. {
  208. GSL_ERROR ("roundoff error detected in the extrapolation table",
  209. GSL_EROUND);
  210. }
  211. else if (error_type == 5)
  212. {
  213. GSL_ERROR ("integral is divergent, or slowly convergent",
  214. GSL_EDIVERGE);
  215. }
  216. else
  217. {
  218. GSL_ERROR ("could not integrate function", GSL_EFAILED);
  219. }
  220. }