3 * Copyright (C) 1996, 1997, 1998, 1999, 2000, 2007 Brian Gough
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.
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.
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.
21 #include <gsl/gsl_integration.h>
23 /* Gauss quadrature weights and kronrod quadrature abscissae and
24 weights as evaluated with 80 decimal digit arithmetic by
25 L. W. Fullerton, Bell Labs, Nov. 1981. */
27 static const double xgk[26] = /* abscissae of the 51-point kronrod rule */
29 0.999262104992609834193457486540341,
30 0.995556969790498097908784946893902,
31 0.988035794534077247637331014577406,
32 0.976663921459517511498315386479594,
33 0.961614986425842512418130033660167,
34 0.942974571228974339414011169658471,
35 0.920747115281701561746346084546331,
36 0.894991997878275368851042006782805,
37 0.865847065293275595448996969588340,
38 0.833442628760834001421021108693570,
39 0.797873797998500059410410904994307,
40 0.759259263037357630577282865204361,
41 0.717766406813084388186654079773298,
42 0.673566368473468364485120633247622,
43 0.626810099010317412788122681624518,
44 0.577662930241222967723689841612654,
45 0.526325284334719182599623778158010,
46 0.473002731445714960522182115009192,
47 0.417885382193037748851814394594572,
48 0.361172305809387837735821730127641,
49 0.303089538931107830167478909980339,
50 0.243866883720988432045190362797452,
51 0.183718939421048892015969888759528,
52 0.122864692610710396387359818808037,
53 0.061544483005685078886546392366797,
54 0.000000000000000000000000000000000
57 /* xgk[1], xgk[3], ... abscissae of the 25-point gauss rule.
58 xgk[0], xgk[2], ... abscissae to optimally extend the 25-point gauss rule */
60 static const double wg[13] = /* weights of the 25-point gauss rule */
62 0.011393798501026287947902964113235,
63 0.026354986615032137261901815295299,
64 0.040939156701306312655623487711646,
65 0.054904695975835191925936891540473,
66 0.068038333812356917207187185656708,
67 0.080140700335001018013234959669111,
68 0.091028261982963649811497220702892,
69 0.100535949067050644202206890392686,
70 0.108519624474263653116093957050117,
71 0.114858259145711648339325545869556,
72 0.119455763535784772228178126512901,
73 0.122242442990310041688959518945852,
74 0.123176053726715451203902873079050
77 static const double wgk[26] = /* weights of the 51-point kronrod rule */
79 0.001987383892330315926507851882843,
80 0.005561932135356713758040236901066,
81 0.009473973386174151607207710523655,
82 0.013236229195571674813656405846976,
83 0.016847817709128298231516667536336,
84 0.020435371145882835456568292235939,
85 0.024009945606953216220092489164881,
86 0.027475317587851737802948455517811,
87 0.030792300167387488891109020215229,
88 0.034002130274329337836748795229551,
89 0.037116271483415543560330625367620,
90 0.040083825504032382074839284467076,
91 0.042872845020170049476895792439495,
92 0.045502913049921788909870584752660,
93 0.047982537138836713906392255756915,
94 0.050277679080715671963325259433440,
95 0.052362885806407475864366712137873,
96 0.054251129888545490144543370459876,
97 0.055950811220412317308240686382747,
98 0.057437116361567832853582693939506,
99 0.058689680022394207961974175856788,
100 0.059720340324174059979099291932562,
101 0.060539455376045862945360267517565,
102 0.061128509717053048305859030416293,
103 0.061471189871425316661544131965264,
104 0.061580818067832935078759824240066
107 /* wgk[25] was calculated from the values of wgk[0..24] */
110 gsl_integration_qk51 (const gsl_function * f, double a, double b,
111 double *result, double *abserr,
112 double *resabs, double *resasc)
114 double fv1[26], fv2[26];
115 gsl_integration_qk (26, xgk, wg, wgk, fv1, fv2, f, a, b, result, abserr, resabs, resasc);