Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-02-21 10:03:48

0001 /* integration/gsl_integration.h
0002  * 
0003  * Copyright (C) 1996, 1997, 1998, 1999, 2000, 2007 Brian Gough
0004  * 
0005  * This program is free software; you can redistribute it and/or modify
0006  * it under the terms of the GNU General Public License as published by
0007  * the Free Software Foundation; either version 3 of the License, or (at
0008  * your option) any later version.
0009  * 
0010  * This program is distributed in the hope that it will be useful, but
0011  * WITHOUT ANY WARRANTY; without even the implied warranty of
0012  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
0013  * General Public License for more details.
0014  * 
0015  * You should have received a copy of the GNU General Public License
0016  * along with this program; if not, write to the Free Software
0017  * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
0018  */
0019 
0020 #ifndef __GSL_INTEGRATION_H__
0021 #define __GSL_INTEGRATION_H__
0022 #include <stdlib.h>
0023 #include <gsl/gsl_math.h>
0024 
0025 #undef __BEGIN_DECLS
0026 #undef __END_DECLS
0027 #ifdef __cplusplus
0028 # define __BEGIN_DECLS extern "C" {
0029 # define __END_DECLS }
0030 #else
0031 # define __BEGIN_DECLS /* empty */
0032 # define __END_DECLS /* empty */
0033 #endif
0034 
0035 __BEGIN_DECLS
0036 
0037 /* Workspace for adaptive integrators */
0038 
0039 typedef struct
0040   {
0041     size_t limit;
0042     size_t size;
0043     size_t nrmax;
0044     size_t i;
0045     size_t maximum_level;
0046     double *alist;
0047     double *blist;
0048     double *rlist;
0049     double *elist;
0050     size_t *order;
0051     size_t *level;
0052   }
0053 gsl_integration_workspace;
0054 
0055 gsl_integration_workspace *
0056   gsl_integration_workspace_alloc (const size_t n);
0057 
0058 void
0059   gsl_integration_workspace_free (gsl_integration_workspace * w);
0060 
0061 
0062 /* Workspace for QAWS integrator */
0063 
0064 typedef struct
0065 {
0066   double alpha;
0067   double beta;
0068   int mu;
0069   int nu;
0070   double ri[25];
0071   double rj[25];
0072   double rg[25];
0073   double rh[25];
0074 }
0075 gsl_integration_qaws_table;
0076 
0077 gsl_integration_qaws_table * 
0078 gsl_integration_qaws_table_alloc (double alpha, double beta, int mu, int nu);
0079 
0080 int
0081 gsl_integration_qaws_table_set (gsl_integration_qaws_table * t,
0082                                 double alpha, double beta, int mu, int nu);
0083 
0084 void
0085 gsl_integration_qaws_table_free (gsl_integration_qaws_table * t);
0086 
0087 /* Workspace for QAWO integrator */
0088 
0089 enum gsl_integration_qawo_enum { GSL_INTEG_COSINE, GSL_INTEG_SINE };
0090 
0091 typedef struct
0092 {
0093   size_t n;
0094   double omega;
0095   double L;
0096   double par;
0097   enum gsl_integration_qawo_enum sine;
0098   double *chebmo;
0099 }
0100 gsl_integration_qawo_table;
0101 
0102 gsl_integration_qawo_table * 
0103 gsl_integration_qawo_table_alloc (double omega, double L, 
0104                                   enum gsl_integration_qawo_enum sine,
0105                                   size_t n);
0106 
0107 int
0108 gsl_integration_qawo_table_set (gsl_integration_qawo_table * t,
0109                                 double omega, double L,
0110                                 enum gsl_integration_qawo_enum sine);
0111 
0112 int
0113 gsl_integration_qawo_table_set_length (gsl_integration_qawo_table * t,
0114                                        double L);
0115 
0116 void
0117 gsl_integration_qawo_table_free (gsl_integration_qawo_table * t);
0118 
0119 
0120 /* Definition of an integration rule */
0121 
0122 typedef void gsl_integration_rule (const gsl_function * f,
0123                                    double a, double b,
0124                                    double *result, double *abserr,
0125                                    double *defabs, double *resabs);
0126 
0127 void gsl_integration_qk15 (const gsl_function * f, double a, double b,
0128                            double *result, double *abserr,
0129                            double *resabs, double *resasc);
0130 
0131 void gsl_integration_qk21 (const gsl_function * f, double a, double b,
0132                            double *result, double *abserr,
0133                            double *resabs, double *resasc);
0134 
0135 void gsl_integration_qk31 (const gsl_function * f, double a, double b,
0136                            double *result, double *abserr,
0137                            double *resabs, double *resasc);
0138 
0139 void gsl_integration_qk41 (const gsl_function * f, double a, double b,
0140                            double *result, double *abserr,
0141                            double *resabs, double *resasc);
0142 
0143 void gsl_integration_qk51 (const gsl_function * f, double a, double b,
0144                            double *result, double *abserr,
0145                            double *resabs, double *resasc);
0146 
0147 void gsl_integration_qk61 (const gsl_function * f, double a, double b,
0148                            double *result, double *abserr,
0149                            double *resabs, double *resasc);
0150 
0151 void gsl_integration_qcheb (gsl_function * f, double a, double b, 
0152                             double *cheb12, double *cheb24);
0153 
0154 /* The low-level integration rules in QUADPACK are identified by small
0155    integers (1-6). We'll use symbolic constants to refer to them.  */
0156 
0157 enum
0158   {
0159     GSL_INTEG_GAUSS15 = 1,      /* 15 point Gauss-Kronrod rule */
0160     GSL_INTEG_GAUSS21 = 2,      /* 21 point Gauss-Kronrod rule */
0161     GSL_INTEG_GAUSS31 = 3,      /* 31 point Gauss-Kronrod rule */
0162     GSL_INTEG_GAUSS41 = 4,      /* 41 point Gauss-Kronrod rule */
0163     GSL_INTEG_GAUSS51 = 5,      /* 51 point Gauss-Kronrod rule */
0164     GSL_INTEG_GAUSS61 = 6       /* 61 point Gauss-Kronrod rule */
0165   };
0166 
0167 void 
0168 gsl_integration_qk (const int n, const double xgk[], 
0169                     const double wg[], const double wgk[],
0170                     double fv1[], double fv2[],
0171                     const gsl_function *f, double a, double b,
0172                     double * result, double * abserr, 
0173                     double * resabs, double * resasc);
0174 
0175 
0176 int gsl_integration_qng (const gsl_function * f,
0177                          double a, double b,
0178                          double epsabs, double epsrel,
0179                          double *result, double *abserr,
0180                          size_t * neval);
0181 
0182 int gsl_integration_qag (const gsl_function * f,
0183                          double a, double b,
0184                          double epsabs, double epsrel, size_t limit,
0185                          int key,
0186                          gsl_integration_workspace * workspace,
0187                          double *result, double *abserr);
0188 
0189 int gsl_integration_qagi (gsl_function * f,
0190                           double epsabs, double epsrel, size_t limit,
0191                           gsl_integration_workspace * workspace,
0192                           double *result, double *abserr);
0193 
0194 int gsl_integration_qagiu (gsl_function * f,
0195                            double a,
0196                            double epsabs, double epsrel, size_t limit,
0197                            gsl_integration_workspace * workspace,
0198                            double *result, double *abserr);
0199 
0200 int gsl_integration_qagil (gsl_function * f,
0201                            double b,
0202                            double epsabs, double epsrel, size_t limit,
0203                            gsl_integration_workspace * workspace,
0204                            double *result, double *abserr);
0205 
0206 
0207 int gsl_integration_qags (const gsl_function * f,
0208                           double a, double b,
0209                           double epsabs, double epsrel, size_t limit,
0210                           gsl_integration_workspace * workspace,
0211                           double *result, double *abserr);
0212 
0213 int gsl_integration_qagp (const gsl_function * f,
0214                           double *pts, size_t npts,
0215                           double epsabs, double epsrel, size_t limit,
0216                           gsl_integration_workspace * workspace,
0217                           double *result, double *abserr);
0218 
0219 int gsl_integration_qawc (gsl_function *f,
0220                           const double a, const double b, const double c,
0221                           const double epsabs, const double epsrel, const size_t limit,
0222                           gsl_integration_workspace * workspace,
0223                           double * result, double * abserr);
0224 
0225 int gsl_integration_qaws (gsl_function * f,
0226                           const double a, const double b,
0227                           gsl_integration_qaws_table * t,
0228                           const double epsabs, const double epsrel,
0229                           const size_t limit,
0230                           gsl_integration_workspace * workspace,
0231                           double *result, double *abserr);
0232 
0233 int gsl_integration_qawo (gsl_function * f,
0234                           const double a,
0235                           const double epsabs, const double epsrel,
0236                           const size_t limit,
0237                           gsl_integration_workspace * workspace,
0238                           gsl_integration_qawo_table * wf,
0239                           double *result, double *abserr);
0240 
0241 int gsl_integration_qawf (gsl_function * f,
0242                           const double a,
0243                           const double epsabs,
0244                           const size_t limit,
0245                           gsl_integration_workspace * workspace,
0246                           gsl_integration_workspace * cycle_workspace,
0247                           gsl_integration_qawo_table * wf,
0248                           double *result, double *abserr);
0249 
0250 /* Workspace for fixed-order Gauss-Legendre integration */
0251 
0252 typedef struct
0253   {
0254     size_t n;         /* number of points */
0255     double *x;        /* Gauss abscissae/points */
0256     double *w;        /* Gauss weights for each abscissae */
0257     int precomputed;  /* high precision abscissae/weights precomputed? */
0258   }
0259 gsl_integration_glfixed_table;
0260 
0261 
0262 gsl_integration_glfixed_table * gsl_integration_glfixed_table_alloc (size_t n);
0263 
0264 void gsl_integration_glfixed_table_free (gsl_integration_glfixed_table * t);
0265 
0266 /* Routine for fixed-order Gauss-Legendre integration */
0267 
0268 double gsl_integration_glfixed (const gsl_function *f,
0269                                 double a,
0270                                 double b,
0271                                 const gsl_integration_glfixed_table * t);
0272 
0273 /* Routine to retrieve the i-th Gauss-Legendre point and weight from t */
0274 
0275 int gsl_integration_glfixed_point (double a,
0276                                    double b,
0277                                    size_t i,
0278                                    double *xi,
0279                                    double *wi,
0280                                    const gsl_integration_glfixed_table * t);
0281 
0282 
0283 /* Cquad integration - Pedro Gonnet */
0284 
0285 /* Data of a single interval */
0286 typedef struct
0287 {
0288   double a, b;
0289   double c[64];
0290   double fx[33];
0291   double igral, err;
0292   int depth, rdepth, ndiv;
0293 } gsl_integration_cquad_ival;
0294 
0295 
0296 /* The workspace is just a collection of intervals */
0297 typedef struct
0298 {
0299   size_t size;
0300   gsl_integration_cquad_ival *ivals;
0301   size_t *heap;
0302 } gsl_integration_cquad_workspace;
0303 
0304 gsl_integration_cquad_workspace *
0305 gsl_integration_cquad_workspace_alloc (const size_t n);
0306 
0307 void
0308 gsl_integration_cquad_workspace_free (gsl_integration_cquad_workspace * w);
0309 
0310 int
0311 gsl_integration_cquad (const gsl_function * f, double a, double b,
0312                            double epsabs, double epsrel,
0313                            gsl_integration_cquad_workspace * ws,
0314                            double *result, double *abserr, size_t * nevals);
0315 
0316 /* Romberg integration workspace and routines */
0317 
0318 typedef struct
0319 {
0320   size_t n;       /* maximum number of steps */
0321   double *work1;  /* workspace for a row of R matrix, size n */
0322   double *work2;  /* workspace for a row of R matrix, size n */
0323 } gsl_integration_romberg_workspace;
0324 
0325 gsl_integration_romberg_workspace *gsl_integration_romberg_alloc(const size_t n);
0326 void gsl_integration_romberg_free(gsl_integration_romberg_workspace * w);
0327 int gsl_integration_romberg(const gsl_function * f, const double a, const double b,
0328                             const double epsabs, const double epsrel, double * result,
0329                             size_t * neval, gsl_integration_romberg_workspace * w);
0330 
0331 /* IQPACK related structures and routines */
0332 
0333 typedef struct
0334 {
0335   double alpha;
0336   double beta;
0337   double a;
0338   double b;
0339   double zemu;
0340   double shft;
0341   double slp;
0342   double al;
0343   double be;
0344 } gsl_integration_fixed_params;
0345 
0346 typedef struct
0347 {
0348   int (*check)(const size_t n, const gsl_integration_fixed_params * params);
0349   int (*init)(const size_t n, double * diag, double * subdiag, gsl_integration_fixed_params * params);
0350 } gsl_integration_fixed_type;
0351 
0352 typedef struct
0353 {
0354   size_t n;        /* number of nodes/weights */
0355   double *weights; /* quadrature weights */
0356   double *x;       /* quadrature nodes */
0357   double *diag;    /* diagonal of Jacobi matrix */
0358   double *subdiag; /* subdiagonal of Jacobi matrix */
0359   const gsl_integration_fixed_type * type;
0360 } gsl_integration_fixed_workspace;
0361 
0362 /* IQPACK integral types */
0363 GSL_VAR const gsl_integration_fixed_type * gsl_integration_fixed_legendre;
0364 GSL_VAR const gsl_integration_fixed_type * gsl_integration_fixed_chebyshev;
0365 GSL_VAR const gsl_integration_fixed_type * gsl_integration_fixed_gegenbauer;
0366 GSL_VAR const gsl_integration_fixed_type * gsl_integration_fixed_jacobi;
0367 GSL_VAR const gsl_integration_fixed_type * gsl_integration_fixed_laguerre;
0368 GSL_VAR const gsl_integration_fixed_type * gsl_integration_fixed_hermite;
0369 GSL_VAR const gsl_integration_fixed_type * gsl_integration_fixed_exponential;
0370 GSL_VAR const gsl_integration_fixed_type * gsl_integration_fixed_rational;
0371 GSL_VAR const gsl_integration_fixed_type * gsl_integration_fixed_chebyshev2;
0372 
0373 gsl_integration_fixed_workspace *
0374 gsl_integration_fixed_alloc(const gsl_integration_fixed_type * type, const size_t n,
0375                             const double a, const double b, const double alpha, const double beta);
0376 
0377 void gsl_integration_fixed_free(gsl_integration_fixed_workspace * w);
0378 
0379 size_t gsl_integration_fixed_n(const gsl_integration_fixed_workspace * w);
0380 
0381 double *gsl_integration_fixed_nodes(const gsl_integration_fixed_workspace * w);
0382 
0383 double *gsl_integration_fixed_weights(const gsl_integration_fixed_workspace * w);
0384 
0385 int gsl_integration_fixed(const gsl_function * func, double * result,
0386                           const gsl_integration_fixed_workspace * w);
0387 
0388 /* Lebedev quadrature */
0389 
0390 typedef struct
0391 {
0392   size_t n;        /* number of nodes/weights */
0393   double *weights; /* quadrature weights */
0394   double *x;       /* x quadrature nodes */
0395   double *y;       /* y quadrature nodes */
0396   double *z;       /* z quadrature nodes */
0397   double *theta;   /* theta quadrature nodes */
0398   double *phi;     /* phi quadrature nodes */
0399 } gsl_integration_lebedev_workspace;
0400 
0401 gsl_integration_lebedev_workspace * gsl_integration_lebedev_alloc(const size_t n);
0402 
0403 void gsl_integration_lebedev_free(gsl_integration_lebedev_workspace * w);
0404 
0405 size_t gsl_integration_lebedev_n(const gsl_integration_lebedev_workspace * w);
0406 
0407 __END_DECLS
0408 
0409 #endif /* __GSL_INTEGRATION_H__ */