gsl_integration__qk21.c 2.7 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778
  1. /* integration/qk21.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 "gsl_integration.h"
  21. /* Gauss quadrature weights and kronrod quadrature abscissae and
  22. weights as evaluated with 80 decimal digit arithmetic by
  23. L. W. Fullerton, Bell Labs, Nov. 1981. */
  24. static const double xgk[11] = /* abscissae of the 21-point kronrod rule */
  25. {
  26. 0.995657163025808080735527280689003,
  27. 0.973906528517171720077964012084452,
  28. 0.930157491355708226001207180059508,
  29. 0.865063366688984510732096688423493,
  30. 0.780817726586416897063717578345042,
  31. 0.679409568299024406234327365114874,
  32. 0.562757134668604683339000099272694,
  33. 0.433395394129247190799265943165784,
  34. 0.294392862701460198131126603103866,
  35. 0.148874338981631210884826001129720,
  36. 0.000000000000000000000000000000000
  37. };
  38. /* xgk[1], xgk[3], ... abscissae of the 10-point gauss rule.
  39. xgk[0], xgk[2], ... abscissae to optimally extend the 10-point gauss rule */
  40. static const double wg[5] = /* weights of the 10-point gauss rule */
  41. {
  42. 0.066671344308688137593568809893332,
  43. 0.149451349150580593145776339657697,
  44. 0.219086362515982043995534934228163,
  45. 0.269266719309996355091226921569469,
  46. 0.295524224714752870173892994651338
  47. };
  48. static const double wgk[11] = /* weights of the 21-point kronrod rule */
  49. {
  50. 0.011694638867371874278064396062192,
  51. 0.032558162307964727478818972459390,
  52. 0.054755896574351996031381300244580,
  53. 0.075039674810919952767043140916190,
  54. 0.093125454583697605535065465083366,
  55. 0.109387158802297641899210590325805,
  56. 0.123491976262065851077958109831074,
  57. 0.134709217311473325928054001771707,
  58. 0.142775938577060080797094273138717,
  59. 0.147739104901338491374841515972068,
  60. 0.149445554002916905664936468389821
  61. };
  62. void
  63. gsl_integration_qk21 (const gsl_function * f, double a, double b,
  64. double *result, double *abserr,
  65. double *resabs, double *resasc)
  66. {
  67. double fv1[11], fv2[11];
  68. gsl_integration_qk (11, xgk, wg, wgk, fv1, fv2, f, a, b, result, abserr, resabs, resasc);
  69. }