gsl_integration__qk61.c 4.6 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127
  1. /* integration/qk61.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[31] = /* abscissae of the 61-point kronrod rule */
  25. {
  26. 0.999484410050490637571325895705811,
  27. 0.996893484074649540271630050918695,
  28. 0.991630996870404594858628366109486,
  29. 0.983668123279747209970032581605663,
  30. 0.973116322501126268374693868423707,
  31. 0.960021864968307512216871025581798,
  32. 0.944374444748559979415831324037439,
  33. 0.926200047429274325879324277080474,
  34. 0.905573307699907798546522558925958,
  35. 0.882560535792052681543116462530226,
  36. 0.857205233546061098958658510658944,
  37. 0.829565762382768397442898119732502,
  38. 0.799727835821839083013668942322683,
  39. 0.767777432104826194917977340974503,
  40. 0.733790062453226804726171131369528,
  41. 0.697850494793315796932292388026640,
  42. 0.660061064126626961370053668149271,
  43. 0.620526182989242861140477556431189,
  44. 0.579345235826361691756024932172540,
  45. 0.536624148142019899264169793311073,
  46. 0.492480467861778574993693061207709,
  47. 0.447033769538089176780609900322854,
  48. 0.400401254830394392535476211542661,
  49. 0.352704725530878113471037207089374,
  50. 0.304073202273625077372677107199257,
  51. 0.254636926167889846439805129817805,
  52. 0.204525116682309891438957671002025,
  53. 0.153869913608583546963794672743256,
  54. 0.102806937966737030147096751318001,
  55. 0.051471842555317695833025213166723,
  56. 0.000000000000000000000000000000000
  57. };
  58. /* xgk[1], xgk[3], ... abscissae of the 30-point gauss rule.
  59. xgk[0], xgk[2], ... abscissae to optimally extend the 30-point gauss rule */
  60. static const double wg[15] = /* weights of the 30-point gauss rule */
  61. {
  62. 0.007968192496166605615465883474674,
  63. 0.018466468311090959142302131912047,
  64. 0.028784707883323369349719179611292,
  65. 0.038799192569627049596801936446348,
  66. 0.048402672830594052902938140422808,
  67. 0.057493156217619066481721689402056,
  68. 0.065974229882180495128128515115962,
  69. 0.073755974737705206268243850022191,
  70. 0.080755895229420215354694938460530,
  71. 0.086899787201082979802387530715126,
  72. 0.092122522237786128717632707087619,
  73. 0.096368737174644259639468626351810,
  74. 0.099593420586795267062780282103569,
  75. 0.101762389748405504596428952168554,
  76. 0.102852652893558840341285636705415
  77. };
  78. static const double wgk[31] = /* weights of the 61-point kronrod rule */
  79. {
  80. 0.001389013698677007624551591226760,
  81. 0.003890461127099884051267201844516,
  82. 0.006630703915931292173319826369750,
  83. 0.009273279659517763428441146892024,
  84. 0.011823015253496341742232898853251,
  85. 0.014369729507045804812451432443580,
  86. 0.016920889189053272627572289420322,
  87. 0.019414141193942381173408951050128,
  88. 0.021828035821609192297167485738339,
  89. 0.024191162078080601365686370725232,
  90. 0.026509954882333101610601709335075,
  91. 0.028754048765041292843978785354334,
  92. 0.030907257562387762472884252943092,
  93. 0.032981447057483726031814191016854,
  94. 0.034979338028060024137499670731468,
  95. 0.036882364651821229223911065617136,
  96. 0.038678945624727592950348651532281,
  97. 0.040374538951535959111995279752468,
  98. 0.041969810215164246147147541285970,
  99. 0.043452539701356069316831728117073,
  100. 0.044814800133162663192355551616723,
  101. 0.046059238271006988116271735559374,
  102. 0.047185546569299153945261478181099,
  103. 0.048185861757087129140779492298305,
  104. 0.049055434555029778887528165367238,
  105. 0.049795683427074206357811569379942,
  106. 0.050405921402782346840893085653585,
  107. 0.050881795898749606492297473049805,
  108. 0.051221547849258772170656282604944,
  109. 0.051426128537459025933862879215781,
  110. 0.051494729429451567558340433647099
  111. };
  112. void
  113. gsl_integration_qk61 (const gsl_function * f, double a, double b,
  114. double *result, double *abserr,
  115. double *resabs, double *resasc)
  116. {
  117. double fv1[31], fv2[31];
  118. gsl_integration_qk (31, xgk, wg, wgk, fv1, fv2, f, a, b, result, abserr, resabs, resasc);
  119. }