gsl_integration__qk51.c 4.2 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118
  1. /* integration/qk51.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[26] = /* abscissae of the 51-point kronrod rule */
  25. {
  26. 0.999262104992609834193457486540341,
  27. 0.995556969790498097908784946893902,
  28. 0.988035794534077247637331014577406,
  29. 0.976663921459517511498315386479594,
  30. 0.961614986425842512418130033660167,
  31. 0.942974571228974339414011169658471,
  32. 0.920747115281701561746346084546331,
  33. 0.894991997878275368851042006782805,
  34. 0.865847065293275595448996969588340,
  35. 0.833442628760834001421021108693570,
  36. 0.797873797998500059410410904994307,
  37. 0.759259263037357630577282865204361,
  38. 0.717766406813084388186654079773298,
  39. 0.673566368473468364485120633247622,
  40. 0.626810099010317412788122681624518,
  41. 0.577662930241222967723689841612654,
  42. 0.526325284334719182599623778158010,
  43. 0.473002731445714960522182115009192,
  44. 0.417885382193037748851814394594572,
  45. 0.361172305809387837735821730127641,
  46. 0.303089538931107830167478909980339,
  47. 0.243866883720988432045190362797452,
  48. 0.183718939421048892015969888759528,
  49. 0.122864692610710396387359818808037,
  50. 0.061544483005685078886546392366797,
  51. 0.000000000000000000000000000000000
  52. };
  53. /* xgk[1], xgk[3], ... abscissae of the 25-point gauss rule.
  54. xgk[0], xgk[2], ... abscissae to optimally extend the 25-point gauss rule */
  55. static const double wg[13] = /* weights of the 25-point gauss rule */
  56. {
  57. 0.011393798501026287947902964113235,
  58. 0.026354986615032137261901815295299,
  59. 0.040939156701306312655623487711646,
  60. 0.054904695975835191925936891540473,
  61. 0.068038333812356917207187185656708,
  62. 0.080140700335001018013234959669111,
  63. 0.091028261982963649811497220702892,
  64. 0.100535949067050644202206890392686,
  65. 0.108519624474263653116093957050117,
  66. 0.114858259145711648339325545869556,
  67. 0.119455763535784772228178126512901,
  68. 0.122242442990310041688959518945852,
  69. 0.123176053726715451203902873079050
  70. };
  71. static const double wgk[26] = /* weights of the 51-point kronrod rule */
  72. {
  73. 0.001987383892330315926507851882843,
  74. 0.005561932135356713758040236901066,
  75. 0.009473973386174151607207710523655,
  76. 0.013236229195571674813656405846976,
  77. 0.016847817709128298231516667536336,
  78. 0.020435371145882835456568292235939,
  79. 0.024009945606953216220092489164881,
  80. 0.027475317587851737802948455517811,
  81. 0.030792300167387488891109020215229,
  82. 0.034002130274329337836748795229551,
  83. 0.037116271483415543560330625367620,
  84. 0.040083825504032382074839284467076,
  85. 0.042872845020170049476895792439495,
  86. 0.045502913049921788909870584752660,
  87. 0.047982537138836713906392255756915,
  88. 0.050277679080715671963325259433440,
  89. 0.052362885806407475864366712137873,
  90. 0.054251129888545490144543370459876,
  91. 0.055950811220412317308240686382747,
  92. 0.057437116361567832853582693939506,
  93. 0.058689680022394207961974175856788,
  94. 0.059720340324174059979099291932562,
  95. 0.060539455376045862945360267517565,
  96. 0.061128509717053048305859030416293,
  97. 0.061471189871425316661544131965264,
  98. 0.061580818067832935078759824240066
  99. };
  100. /* wgk[25] was calculated from the values of wgk[0..24] */
  101. void
  102. gsl_integration_qk51 (const gsl_function * f, double a, double b,
  103. double *result, double *abserr,
  104. double *resabs, double *resasc)
  105. {
  106. double fv1[26], fv2[26];
  107. gsl_integration_qk (26, xgk, wg, wgk, fv1, fv2, f, a, b, result, abserr, resabs, resasc);
  108. }