gsl_integration__qk31.c 3.2 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990
  1. /* integration/qk31.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[16] = /* abscissae of the 31-point kronrod rule */
  25. {
  26. 0.998002298693397060285172840152271,
  27. 0.987992518020485428489565718586613,
  28. 0.967739075679139134257347978784337,
  29. 0.937273392400705904307758947710209,
  30. 0.897264532344081900882509656454496,
  31. 0.848206583410427216200648320774217,
  32. 0.790418501442465932967649294817947,
  33. 0.724417731360170047416186054613938,
  34. 0.650996741297416970533735895313275,
  35. 0.570972172608538847537226737253911,
  36. 0.485081863640239680693655740232351,
  37. 0.394151347077563369897207370981045,
  38. 0.299180007153168812166780024266389,
  39. 0.201194093997434522300628303394596,
  40. 0.101142066918717499027074231447392,
  41. 0.000000000000000000000000000000000
  42. };
  43. /* xgk[1], xgk[3], ... abscissae of the 15-point gauss rule.
  44. xgk[0], xgk[2], ... abscissae to optimally extend the 15-point gauss rule */
  45. static const double wg[8] = /* weights of the 15-point gauss rule */
  46. {
  47. 0.030753241996117268354628393577204,
  48. 0.070366047488108124709267416450667,
  49. 0.107159220467171935011869546685869,
  50. 0.139570677926154314447804794511028,
  51. 0.166269205816993933553200860481209,
  52. 0.186161000015562211026800561866423,
  53. 0.198431485327111576456118326443839,
  54. 0.202578241925561272880620199967519
  55. };
  56. static const double wgk[16] = /* weights of the 31-point kronrod rule */
  57. {
  58. 0.005377479872923348987792051430128,
  59. 0.015007947329316122538374763075807,
  60. 0.025460847326715320186874001019653,
  61. 0.035346360791375846222037948478360,
  62. 0.044589751324764876608227299373280,
  63. 0.053481524690928087265343147239430,
  64. 0.062009567800670640285139230960803,
  65. 0.069854121318728258709520077099147,
  66. 0.076849680757720378894432777482659,
  67. 0.083080502823133021038289247286104,
  68. 0.088564443056211770647275443693774,
  69. 0.093126598170825321225486872747346,
  70. 0.096642726983623678505179907627589,
  71. 0.099173598721791959332393173484603,
  72. 0.100769845523875595044946662617570,
  73. 0.101330007014791549017374792767493
  74. };
  75. void
  76. gsl_integration_qk31 (const gsl_function * f, double a, double b,
  77. double *result, double *abserr,
  78. double *resabs, double *resasc)
  79. {
  80. double fv1[16], fv2[16];
  81. gsl_integration_qk (16, xgk, wg, wgk, fv1, fv2, f, a, b, result, abserr, resabs, resasc);
  82. }