SDDSlib
Loading...
Searching...
No Matches
gaussQuad.c File Reference

Implements recursive Gaussian quadrature for numerical integration. More...

#include <math.h>
#include "mdb.h"

Go to the source code of this file.

Macros

#define _CGQ   .57735026918962576451
 
#define gauss_quad(fn, z, a, w)
 
#define MAXSTACK   16384
 

Functions

long gaussianQuadrature (double(*fn)(), double a, double b, long n, double err, double *result)
 Integrate a function using recursive Gaussian quadrature.
 

Detailed Description

Implements recursive Gaussian quadrature for numerical integration.

This file contains the implementation of the recursive Gaussian quadrature technique for integrating functions numerically. The primary function gaussianQuadrature divides the integration interval into panels and recursively refines the estimate to achieve the desired error tolerance.

License
This file is distributed under the terms of the Software License Agreement found in the file LICENSE included with this distribution.
Author
M. Borland, R. Soliday

Definition in file gaussQuad.c.

Macro Definition Documentation

◆ _CGQ

#define _CGQ   .57735026918962576451

Definition at line 24 of file gaussQuad.c.

◆ gauss_quad

#define gauss_quad ( fn,
z,
a,
w )
Value:
(((*fn)(z + a) + (*fn)(a - z)) * w)

Definition at line 25 of file gaussQuad.c.

◆ MAXSTACK

#define MAXSTACK   16384

Definition at line 27 of file gaussQuad.c.

Function Documentation

◆ gaussianQuadrature()

long gaussianQuadrature ( double(* fn )(),
double a,
double b,
long n,
double err,
double * result )

Integrate a function using recursive Gaussian quadrature.

Performs numerical integration of the specified function over the interval [a, b] using a recursive Gaussian quadrature technique. The interval is initially divided into n panels, and each panel is recursively subdivided until the desired fractional error err is achieved.

Parameters
fnPointer to the function to be integrated.
aLower limit of integration.
bUpper limit of integration.
nInitial number of panels to divide the interval into.
errFractional error permissible on any quadrature.
resultPointer to a double where the result of the integration will be stored.
Returns
The total number of function evaluations performed, or -1 if an error occurs.

Definition at line 45 of file gaussQuad.c.

50 {
51 register long s; /* stack pointer for panel limits */
52 double ave, i0, i1, i2, i3, idiv;
53 long totalEvals;
54 double d_ab, i_tot, d_ab2;
55 double ab, z, width;
56 /* stacks of upper,lower limits, and integral values from two
57 * point quadratures on each interval */
58 double a_stack[MAXSTACK];
59 double b_stack[MAXSTACK];
60 double i_stack[MAXSTACK];
61
62 if (n > MAXSTACK)
63 return -1;
64
65 d_ab = (b - a) / n;
66 d_ab2 = d_ab / 2;
67 totalEvals = 0;
68 for (s = 0; s < n; s++) {
69 b = b_stack[s] = (a_stack[s] = a) + d_ab;
70 ave = a + d_ab2;
71 z = d_ab2 * _CGQ;
72 i_stack[s] = gauss_quad(fn, z, ave, d_ab2);
73 totalEvals += 2;
74 a = b;
75 }
76
77 s = n - 1;
78 i_tot = 0;
79 while (s >= 0) {
80 if (s == MAXSTACK)
81 return -1;
82 a = a_stack[s];
83 b = b_stack[s];
84 ab = (a + b) / 2;
85 i0 = i_stack[s];
86
87 ave = (a + ab) / 2;
88 z = (width = (ab - a) / 2) * _CGQ;
89 i1 = gauss_quad(fn, z, ave, width);
90 totalEvals += 2;
91
92 ave = (ab + b) / 2;
93 z = (width = (b - ab) / 2) * _CGQ;
94 i2 = gauss_quad(fn, z, ave, width);
95 totalEvals += 2;
96
97 i3 = idiv = i1 + i2;
98 if (i3 == 0)
99 idiv = i0;
100 if (idiv != 0) {
101 if (fabs((i3 - i0) / idiv) > err) {
102 b_stack[s] = ab;
103 i_stack[s] = i1;
104 a_stack[++s] = ab;
105 b_stack[s] = b;
106 i_stack[s] = i2;
107 } else {
108 i_tot += i3;
109 s--;
110 }
111 } else {
112 i_tot += i3;
113 s--;
114 }
115 }
116 *result = i_tot;
117 return totalEvals;
118}