rgsl/integration.rs
1//
2// A rust binding for the GSL library by Guillaume Gomez (guillaume1.gomez@gmail.com)
3//
4
5/*!
6## Introduction
7
8Each algorithm computes an approximation to a definite integral of the form,
9
10I = \int_a^b f(x) w(x) dx
11where w(x) is a weight function (for general integrands w(x)=1). The user provides absolute and relative error bounds (epsabs, epsrel) which
12specify the following accuracy requirement,
13
14|RESULT - I| <= max(epsabs, epsrel |I|)
15
16where RESULT is the numerical approximation obtained by the algorithm. The algorithms attempt to estimate the absolute error ABSERR = |RESULT
17- I| in such a way that the following inequality holds,
18
19|RESULT - I| <= ABSERR <= max(epsabs, epsrel |I|)
20
21In short, the routines return the first approximation which has an absolute error smaller than epsabs or a relative error smaller than epsrel.
22
23Note that this is an either-or constraint, not simultaneous. To compute to a specified absolute error, set epsrel to zero. To compute to a
24specified relative error, set epsabs to zero. The routines will fail to converge if the error bounds are too stringent, but always return the
25best approximation obtained up to that stage.
26
27The algorithms in QUADPACK use a naming convention based on the following letters,
28
29Q - quadrature routine
30
31N - non-adaptive integrator
32A - adaptive integrator
33
34G - general integrand (user-defined)
35W - weight function with integrand
36
37S - singularities can be more readily integrated
38P - points of special difficulty can be supplied
39I - infinite range of integration
40O - oscillatory weight function, cos or sin
41F - Fourier integral
42C - Cauchy principal value
43The algorithms are built on pairs of quadrature rules, a higher order rule and a lower order rule. The higher order rule is used to compute the
44best approximation to an integral over a small range. The difference between the results of the higher order rule and the lower order rule gives
45an estimate of the error in the approximation.
46
47 * [Integrands without weight functions](http://www.gnu.org/software/gsl/manual/html_node/Integrands-without-weight-functions.html#Integrands-without-weight-functions)
48 * [Integrands with weight functions](http://www.gnu.org/software/gsl/manual/html_node/Integrands-with-weight-functions.html#Integrands-with-weight-functions)
49 * [Integrands with singular weight functions](http://www.gnu.org/software/gsl/manual/html_node/Integrands-with-singular-weight-functions.html#Integrands-with-singular-weight-functions)
50
51## QNG non-adaptive Gauss-Kronrod integration
52
53The QNG algorithm is a non-adaptive procedure which uses fixed Gauss-Kronrod-Patterson abscissae to sample the integrand at a maximum of 87
54points. It is provided for fast integration of smooth functions.
55
56## QAG adaptive integration
57
58The QAG algorithm is a simple adaptive integration procedure. The integration region is divided into subintervals, and on each iteration the
59subinterval with the largest estimated error is bisected. This reduces the overall error rapidly, as the subintervals become concentrated
60around local difficulties in the integrand. These subintervals are managed by a gsl_integration_workspace struct, which handles the memory
61for the subinterval ranges, results and error estimates.
62
63## QAGS adaptive integration with singularities
64
65The presence of an integrable singularity in the integration region causes an adaptive routine to concentrate new subintervals around the
66singularity. As the subintervals decrease in size the successive approximations to the integral converge in a limiting fashion. This
67approach to the limit can be accelerated using an extrapolation procedure. The QAGS algorithm combines adaptive bisection with the Wynn
68epsilon-algorithm to speed up the integration of many types of integrable singularities.
69
70## References and Further Reading
71
72The following book is the definitive reference for QUADPACK, and was written by the original authors. It provides descriptions of the
73algorithms, program listings, test programs and examples. It also includes useful advice on numerical integration and many references
74to the numerical integration literature used in developing QUADPACK.
75
76R. Piessens, E. de Doncker-Kapenga, C.W. Ueberhuber, D.K. Kahaner. QUADPACK A subroutine package for automatic integration Springer Verlag, 1983.
77The CQUAD integration algorithm is described in the following paper:
78
79P. Gonnet, “Increasing the Reliability of Adaptive Quadrature Using Explicit Interpolants”, ACM Transactions on Mathematical Software, Volume 37
80(2010), Issue 3, Article 26.
81!*/
82
83use crate::Value;
84use ffi::FFI;
85
86/// This function applies the Gauss-Kronrod 10-point, 21-point, 43-point and 87-point integration
87/// rules in succession until an estimate of the integral of f over (a,b) is achieved within the
88/// desired absolute and relative error limits, eps_abs and eps_rel. The function returns the final
89/// approximation, result, an estimate of the absolute error, abserr and the number of function
90/// evaluations used, neval. The Gauss-Kronrod rules are designed in such a way that each rule uses
91/// all the results of its predecessors, in order to minimize the total number of function
92/// evaluations.
93///
94/// Returns `(result, abs_err, n_eval)`.
95#[doc(alias = "gsl_integration_qng")]
96pub fn qng<F: Fn(f64) -> f64>(
97 f: F,
98 a: f64,
99 b: f64,
100 eps_abs: f64,
101 eps_rel: f64,
102) -> Result<(f64, f64, usize), Value> {
103 let function = wrap_callback!(f, F);
104 let mut result = 0.;
105 let mut abs_err = 0.;
106 let mut n_eval = 0;
107
108 let ret = unsafe {
109 sys::gsl_integration_qng(
110 &function,
111 a,
112 b,
113 eps_abs,
114 eps_rel,
115 &mut result,
116 &mut abs_err,
117 &mut n_eval,
118 )
119 };
120 result_handler!(ret, (result, abs_err, n_eval))
121}
122
123/// Gauss quadrature weights and kronrod quadrature abscissae and weights as evaluated with 80
124/// decimal digit arithmetic by L. W.
125///
126/// Fullerton, Bell Labs, Nov. 1981.
127///
128/// Returns `(result, abs_err, resabs, resasc)`.
129#[doc(alias = "gsl_integration_qk15")]
130pub fn qk15<F: Fn(f64) -> f64>(f: F, a: f64, b: f64) -> (f64, f64, f64, f64) {
131 let function = wrap_callback!(f, F);
132 let mut result = 0.;
133 let mut abs_err = 0.;
134 let mut resabs = 0.;
135 let mut resasc = 0.;
136
137 unsafe {
138 sys::gsl_integration_qk15(
139 &function,
140 a,
141 b,
142 &mut result,
143 &mut abs_err,
144 &mut resabs,
145 &mut resasc,
146 )
147 };
148 (result, abs_err, resabs, resasc)
149}
150
151/// Returns `(result, abs_err, resabs, resasc)`.
152#[doc(alias = "gsl_integration_qk21")]
153pub fn qk21<F: Fn(f64) -> f64>(f: F, a: f64, b: f64) -> (f64, f64, f64, f64) {
154 let function = wrap_callback!(f, F);
155 let mut result = 0.;
156 let mut abs_err = 0.;
157 let mut resabs = 0.;
158 let mut resasc = 0.;
159
160 unsafe {
161 sys::gsl_integration_qk21(
162 &function,
163 a,
164 b,
165 &mut result,
166 &mut abs_err,
167 &mut resabs,
168 &mut resasc,
169 )
170 };
171 (result, abs_err, resabs, resasc)
172}
173
174/// Returns `(result, abs_err, resabs, resasc)`.
175#[doc(alias = "gsl_integration_qk31")]
176pub fn qk31<F: Fn(f64) -> f64>(f: F, a: f64, b: f64) -> (f64, f64, f64, f64) {
177 let function = wrap_callback!(f, F);
178 let mut result = 0.;
179 let mut abs_err = 0.;
180 let mut resabs = 0.;
181 let mut resasc = 0.;
182
183 unsafe {
184 sys::gsl_integration_qk31(
185 &function,
186 a,
187 b,
188 &mut result,
189 &mut abs_err,
190 &mut resabs,
191 &mut resasc,
192 )
193 };
194 (result, abs_err, resabs, resasc)
195}
196
197/// Returns `(result, abs_err, resabs, resasc)`.
198#[doc(alias = "gsl_integration_qk41")]
199pub fn qk41<F: Fn(f64) -> f64>(f: F, a: f64, b: f64) -> (f64, f64, f64, f64) {
200 let function = wrap_callback!(f, F);
201 let mut result = 0.;
202 let mut abs_err = 0.;
203 let mut resabs = 0.;
204 let mut resasc = 0.;
205
206 unsafe {
207 sys::gsl_integration_qk41(
208 &function,
209 a,
210 b,
211 &mut result,
212 &mut abs_err,
213 &mut resabs,
214 &mut resasc,
215 )
216 };
217 (result, abs_err, resabs, resasc)
218}
219
220/// Returns `(result, abs_err, resabs, resasc)`.
221#[doc(alias = "gsl_integration_qk51")]
222pub fn qk51<F: Fn(f64) -> f64>(f: F, a: f64, b: f64) -> (f64, f64, f64, f64) {
223 let function = wrap_callback!(f, F);
224 let mut result = 0.;
225 let mut abs_err = 0.;
226 let mut resabs = 0.;
227 let mut resasc = 0.;
228
229 unsafe {
230 sys::gsl_integration_qk51(
231 &function,
232 a,
233 b,
234 &mut result,
235 &mut abs_err,
236 &mut resabs,
237 &mut resasc,
238 )
239 };
240 (result, abs_err, resabs, resasc)
241}
242
243/// Returns `(result, abs_err, resabs, resasc)`.
244#[doc(alias = "gsl_integration_qk61")]
245pub fn qk61<F: Fn(f64) -> f64>(f: F, a: f64, b: f64) -> (f64, f64, f64, f64) {
246 let function = wrap_callback!(f, F);
247 let mut result = 0.;
248 let mut abs_err = 0.;
249 let mut resabs = 0.;
250 let mut resasc = 0.;
251
252 unsafe {
253 sys::gsl_integration_qk61(
254 &function,
255 a,
256 b,
257 &mut result,
258 &mut abs_err,
259 &mut resabs,
260 &mut resasc,
261 )
262 };
263 (result, abs_err, resabs, resasc)
264}
265
266/// Returns `(result, abs_err, resabs, resasc)`.
267#[doc(alias = "gsl_integration_qk")]
268pub fn qk<F: Fn(f64) -> f64>(
269 xgk: &[f64],
270 wg: &[f64],
271 wgk: &[f64],
272 fv1: &mut [f64],
273 fv2: &mut [f64],
274 f: F,
275 a: f64,
276 b: f64,
277) -> (f64, f64, f64, f64) {
278 assert!(xgk.len() == wg.len());
279 assert!(xgk.len() == wgk.len());
280 assert!(xgk.len() == fv1.len());
281 assert!(xgk.len() == fv2.len());
282
283 let function = wrap_callback!(f, F);
284 let mut result = 0.;
285 let mut abs_err = 0.;
286 let mut resabs = 0.;
287 let mut resasc = 0.;
288
289 unsafe {
290 sys::gsl_integration_qk(
291 xgk.len() as _,
292 xgk.as_ptr(),
293 wg.as_ptr(),
294 wgk.as_ptr(),
295 fv1.as_mut_ptr(),
296 fv2.as_mut_ptr(),
297 &function,
298 a,
299 b,
300 &mut result,
301 &mut abs_err,
302 &mut resabs,
303 &mut resasc,
304 );
305 }
306 (result, abs_err, resabs, resasc)
307}
308
309/// This function attempts to compute a Fourier integral of the function f over the semi-infinite
310/// interval `[a,+\infty)`.
311///
312/// ```text
313/// I = \int_a^{+\infty} dx f(x) sin(omega x)
314/// I = \int_a^{+\infty} dx f(x) cos(omega x)
315/// ```
316///
317/// The parameter \omega and choice of \sin or \cos is taken from the table wf (the length L can
318/// take any value, since it is overridden by this function to a value appropriate for the Fourier
319/// integration). The integral is computed using the QAWO algorithm over each of the subintervals,
320///
321/// ```text
322/// C_1 = [a, a + c]
323/// C_2 = [a + c, a + 2 c]
324/// ... = ...
325/// C_k = [a + (k-1) c, a + k c]
326/// ```
327///
328/// where c = (2 floor(|\omega|) + 1) \pi/|\omega|. The width c is chosen to cover an odd number of
329/// periods so that the contributions from the intervals alternate in sign and are monotonically
330/// decreasing when f is positive and monotonically decreasing. The sum of this sequence of
331/// contributions is accelerated using the epsilon-algorithm.
332///
333/// This function works to an overall absolute tolerance of abserr. The following strategy is used:
334/// on each interval C_k the algorithm tries to achieve the tolerance
335///
336/// ```text
337/// TOL_k = u_k abserr
338/// ```
339///
340/// where u_k = (1 - p)p^{k-1} and p = 9/10. The sum of the geometric series of contributions from
341/// each interval gives an overall tolerance of abserr.
342///
343/// If the integration of a subinterval leads to difficulties then the accuracy requirement for
344/// subsequent intervals is relaxed,
345///
346/// ```text
347/// TOL_k = u_k max(abserr, max_{i<k}{E_i})
348/// ```
349///
350/// where E_k is the estimated error on the interval C_k.
351///
352/// The subintervals and their results are stored in the memory provided by workspace. The maximum
353/// number of subintervals is given by limit, which may not exceed the allocated size of the
354/// workspace. The integration over each subinterval uses the memory provided by cycle_workspace
355/// as workspace for the QAWO algorithm.
356///
357/// Returns `(result, abs_err)`.
358#[doc(alias = "gsl_integration_qawf")]
359pub fn qawf<F: Fn(f64) -> f64>(
360 f: F,
361 a: f64,
362 epsabs: f64,
363 limit: usize,
364 workspace: &mut crate::IntegrationWorkspace,
365 cycle_workspace: &mut crate::IntegrationWorkspace,
366 wf: &mut crate::IntegrationQawoTable,
367) -> Result<(f64, f64), Value> {
368 let mut result = 0.;
369 let mut abs_err = 0.;
370
371 let mut function = wrap_callback!(f, F);
372 let ret = unsafe {
373 sys::gsl_integration_qawf(
374 &mut function,
375 a,
376 epsabs,
377 limit,
378 FFI::unwrap_unique(workspace),
379 FFI::unwrap_unique(cycle_workspace),
380 FFI::unwrap_unique(wf),
381 &mut result,
382 &mut abs_err,
383 )
384 };
385 result_handler!(ret, (result, abs_err))
386}