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