Blame integration/qk41.c

Packit 67cb25
/* integration/qk41.c
Packit 67cb25
 * 
Packit 67cb25
 * Copyright (C) 1996, 1997, 1998, 1999, 2000, 2007 Brian Gough
Packit 67cb25
 * 
Packit 67cb25
 * This program is free software; you can redistribute it and/or modify
Packit 67cb25
 * it under the terms of the GNU General Public License as published by
Packit 67cb25
 * the Free Software Foundation; either version 3 of the License, or (at
Packit 67cb25
 * your option) any later version.
Packit 67cb25
 * 
Packit 67cb25
 * This program is distributed in the hope that it will be useful, but
Packit 67cb25
 * WITHOUT ANY WARRANTY; without even the implied warranty of
Packit 67cb25
 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
Packit 67cb25
 * General Public License for more details.
Packit 67cb25
 * 
Packit 67cb25
 * You should have received a copy of the GNU General Public License
Packit 67cb25
 * along with this program; if not, write to the Free Software
Packit 67cb25
 * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
Packit 67cb25
 */
Packit 67cb25
Packit 67cb25
#include <config.h>
Packit 67cb25
#include <gsl/gsl_integration.h>
Packit 67cb25
Packit 67cb25
/* Gauss quadrature weights and kronrod quadrature abscissae and
Packit 67cb25
   weights as evaluated with 80 decimal digit arithmetic by
Packit 67cb25
   L. W. Fullerton, Bell Labs, Nov. 1981. */
Packit 67cb25
Packit 67cb25
static const double xgk[21] =   /* abscissae of the 41-point kronrod rule */
Packit 67cb25
{
Packit 67cb25
  0.998859031588277663838315576545863,
Packit 67cb25
  0.993128599185094924786122388471320,
Packit 67cb25
  0.981507877450250259193342994720217,
Packit 67cb25
  0.963971927277913791267666131197277,
Packit 67cb25
  0.940822633831754753519982722212443,
Packit 67cb25
  0.912234428251325905867752441203298,
Packit 67cb25
  0.878276811252281976077442995113078,
Packit 67cb25
  0.839116971822218823394529061701521,
Packit 67cb25
  0.795041428837551198350638833272788,
Packit 67cb25
  0.746331906460150792614305070355642,
Packit 67cb25
  0.693237656334751384805490711845932,
Packit 67cb25
  0.636053680726515025452836696226286,
Packit 67cb25
  0.575140446819710315342946036586425,
Packit 67cb25
  0.510867001950827098004364050955251,
Packit 67cb25
  0.443593175238725103199992213492640,
Packit 67cb25
  0.373706088715419560672548177024927,
Packit 67cb25
  0.301627868114913004320555356858592,
Packit 67cb25
  0.227785851141645078080496195368575,
Packit 67cb25
  0.152605465240922675505220241022678,
Packit 67cb25
  0.076526521133497333754640409398838,
Packit 67cb25
  0.000000000000000000000000000000000
Packit 67cb25
};
Packit 67cb25
Packit 67cb25
/* xgk[1], xgk[3], ... abscissae of the 20-point gauss rule. 
Packit 67cb25
   xgk[0], xgk[2], ... abscissae to optimally extend the 20-point gauss rule */
Packit 67cb25
Packit 67cb25
static const double wg[10] =    /* weights of the 20-point gauss rule */
Packit 67cb25
{
Packit 67cb25
  0.017614007139152118311861962351853,
Packit 67cb25
  0.040601429800386941331039952274932,
Packit 67cb25
  0.062672048334109063569506535187042,
Packit 67cb25
  0.083276741576704748724758143222046,
Packit 67cb25
  0.101930119817240435036750135480350,
Packit 67cb25
  0.118194531961518417312377377711382,
Packit 67cb25
  0.131688638449176626898494499748163,
Packit 67cb25
  0.142096109318382051329298325067165,
Packit 67cb25
  0.149172986472603746787828737001969,
Packit 67cb25
  0.152753387130725850698084331955098
Packit 67cb25
};
Packit 67cb25
Packit 67cb25
static const double wgk[21] =   /* weights of the 41-point kronrod rule */
Packit 67cb25
{
Packit 67cb25
  0.003073583718520531501218293246031,
Packit 67cb25
  0.008600269855642942198661787950102,
Packit 67cb25
  0.014626169256971252983787960308868,
Packit 67cb25
  0.020388373461266523598010231432755,
Packit 67cb25
  0.025882133604951158834505067096153,
Packit 67cb25
  0.031287306777032798958543119323801,
Packit 67cb25
  0.036600169758200798030557240707211,
Packit 67cb25
  0.041668873327973686263788305936895,
Packit 67cb25
  0.046434821867497674720231880926108,
Packit 67cb25
  0.050944573923728691932707670050345,
Packit 67cb25
  0.055195105348285994744832372419777,
Packit 67cb25
  0.059111400880639572374967220648594,
Packit 67cb25
  0.062653237554781168025870122174255,
Packit 67cb25
  0.065834597133618422111563556969398,
Packit 67cb25
  0.068648672928521619345623411885368,
Packit 67cb25
  0.071054423553444068305790361723210,
Packit 67cb25
  0.073030690332786667495189417658913,
Packit 67cb25
  0.074582875400499188986581418362488,
Packit 67cb25
  0.075704497684556674659542775376617,
Packit 67cb25
  0.076377867672080736705502835038061,
Packit 67cb25
  0.076600711917999656445049901530102
Packit 67cb25
};
Packit 67cb25
Packit 67cb25
void
Packit 67cb25
gsl_integration_qk41 (const gsl_function * f, double a, double b,
Packit 67cb25
                      double *result, double *abserr,
Packit 67cb25
                      double *resabs, double *resasc)
Packit 67cb25
{
Packit 67cb25
  double fv1[21], fv2[21];
Packit 67cb25
  gsl_integration_qk (21, xgk, wg, wgk, fv1, fv2, f, a, b, result, abserr, resabs, resasc);
Packit 67cb25
}
Packit 67cb25