1
  2
  3
  4
  5
  6
  7
  8
  9
 10
 11
 12
 13
 14
 15
 16
 17
 18
 19
 20
 21
 22
 23
 24
 25
 26
 27
 28
 29
 30
 31
 32
 33
 34
 35
 36
 37
 38
 39
 40
 41
 42
 43
 44
 45
 46
 47
 48
 49
 50
 51
 52
 53
 54
 55
 56
 57
 58
 59
 60
 61
 62
 63
 64
 65
 66
 67
 68
 69
 70
 71
 72
 73
 74
 75
 76
 77
 78
 79
 80
 81
 82
 83
 84
 85
 86
 87
 88
 89
 90
 91
 92
 93
 94
 95
 96
 97
 98
 99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
//
// A rust binding for the GSL library by Guillaume Gomez (guillaume1.gomez@gmail.com)
//

#![allow(clippy::upper_case_acronyms)]

use crate::enums;
use crate::Value;
use ffi::FFI;

ffi_wrapper!(IntegrationFixedType, *const sys::gsl_integration_fixed_type);

impl IntegrationFixedType {
    #[doc(alias = "gsl_integration_fixed_legendre")]
    pub fn legendre() -> IntegrationFixedType {
        ffi_wrap!(gsl_integration_fixed_legendre)
    }
    #[doc(alias = "gsl_integration_fixed_chebyshev")]
    pub fn chebyshev() -> IntegrationFixedType {
        ffi_wrap!(gsl_integration_fixed_chebyshev)
    }
    #[doc(alias = "gsl_integration_fixed_chebyshev2")]
    pub fn chebyshev2() -> IntegrationFixedType {
        ffi_wrap!(gsl_integration_fixed_chebyshev2)
    }
    #[doc(alias = "gsl_integration_fixed_gegenbauer")]
    pub fn gegenbauer() -> IntegrationFixedType {
        ffi_wrap!(gsl_integration_fixed_gegenbauer)
    }
    #[doc(alias = "gsl_integration_fixed_jacobi")]
    pub fn jacobi() -> IntegrationFixedType {
        ffi_wrap!(gsl_integration_fixed_jacobi)
    }
    #[doc(alias = "gsl_integration_fixed_laguerre")]
    pub fn laguerre() -> IntegrationFixedType {
        ffi_wrap!(gsl_integration_fixed_laguerre)
    }
    #[doc(alias = "gsl_integration_fixed_hermite")]
    pub fn hermite() -> IntegrationFixedType {
        ffi_wrap!(gsl_integration_fixed_hermite)
    }
    #[doc(alias = "gsl_integration_fixed_exponential")]
    pub fn exponential() -> IntegrationFixedType {
        ffi_wrap!(gsl_integration_fixed_exponential)
    }
    #[doc(alias = "gsl_integration_fixed_rational")]
    pub fn rational() -> IntegrationFixedType {
        ffi_wrap!(gsl_integration_fixed_rational)
    }
}

ffi_wrapper!(
    IntegrationFixedWorkspace,
    *mut sys::gsl_integration_fixed_workspace,
    gsl_integration_fixed_free
);

impl IntegrationFixedWorkspace {
    #[doc(alias = "gsl_integration_fixed_alloc")]
    pub fn new(
        type_: IntegrationFixedType,
        n: usize,
        a: f64,
        b: f64,
        alpha: f64,
        beta: f64,
    ) -> Option<IntegrationFixedWorkspace> {
        let tmp = unsafe {
            sys::gsl_integration_fixed_alloc(type_.unwrap_shared(), n, a, b, alpha, beta)
        };

        if tmp.is_null() {
            None
        } else {
            Some(Self::wrap(tmp))
        }
    }

    #[doc(alias = "gsl_integration_fixed_n")]
    pub fn n(&self) -> usize {
        unsafe { sys::gsl_integration_fixed_n(self.unwrap_shared()) }
    }

    #[doc(alias = "gsl_integration_fixed_nodes")]
    pub fn nodes(&self) -> Option<&[f64]> {
        let tmp = unsafe { sys::gsl_integration_fixed_nodes(self.unwrap_shared()) };
        if tmp.is_null() {
            return None;
        }
        unsafe { Some(::std::slice::from_raw_parts(tmp, self.n())) }
    }

    #[doc(alias = "gsl_integration_fixed_weights")]
    pub fn weights(&self) -> Option<&[f64]> {
        let tmp = unsafe { sys::gsl_integration_fixed_weights(self.unwrap_shared()) };
        if tmp.is_null() {
            return None;
        }
        unsafe { Some(::std::slice::from_raw_parts(tmp, self.n())) }
    }

    #[doc(alias = "gsl_integration_fixed")]
    pub fn fixed<F: Fn(f64) -> f64>(&self, f: F) -> (::Value, f64) {
        let mut result = 0.;
        let function = wrap_callback!(f, F);

        let ret =
            unsafe { sys::gsl_integration_fixed(&function, &mut result, self.unwrap_shared()) };
        (::Value::from(ret), result)
    }
}

ffi_wrapper!(IntegrationWorkspace, *mut sys::gsl_integration_workspace, gsl_integration_workspace_free,
"The QAG algorithm is a simple adaptive integration procedure. The integration region is divided
into subintervals, and on each iteration the subinterval with the largest estimated error is
bisected. This reduces the overall error rapidly, as the subintervals become concentrated
around local difficulties in the integrand. These subintervals are managed by a
gsl_integration_workspace struct, which handles the memory for the subinterval ranges, results
and error estimates.");

impl IntegrationWorkspace {
    /// This function allocates a workspace sufficient to hold n double precision intervals, their
    /// integration results and error estimates. One workspace may be used multiple times as all
    /// necessary reinitialization is performed automatically by the integration routines.
    #[doc(alias = "gsl_integration_workspace_alloc")]
    pub fn new(n: usize) -> Option<IntegrationWorkspace> {
        let tmp = unsafe { sys::gsl_integration_workspace_alloc(n) };

        if tmp.is_null() {
            None
        } else {
            Some(Self::wrap(tmp))
        }
    }

    pub fn limit(&self) -> usize {
        unsafe { (*self.unwrap_shared()).limit }
    }
    pub fn size(&self) -> usize {
        unsafe { (*self.unwrap_shared()).size }
    }
    pub fn nrmax(&self) -> usize {
        unsafe { (*self.unwrap_shared()).nrmax }
    }
    pub fn i(&self) -> usize {
        unsafe { (*self.unwrap_shared()).i }
    }
    pub fn maximum_level(&self) -> usize {
        unsafe { (*self.unwrap_shared()).maximum_level }
    }

    /// This function applies an integration rule adaptively until an estimate of the integral of f
    /// over (a,b) is achieved within the desired absolute and relative error limits, epsabs and
    /// epsrel. The function returns the final approximation, result, and an estimate of the
    /// absolute error, abserr. The integration rule is determined by the value of key, which should
    /// be chosen from the following symbolic names,
    ///
    /// GSL_INTEG_GAUSS15  (key = 1)
    ///
    /// GSL_INTEG_GAUSS21  (key = 2)
    ///
    /// GSL_INTEG_GAUSS31  (key = 3)
    ///
    /// GSL_INTEG_GAUSS41  (key = 4)
    ///
    /// GSL_INTEG_GAUSS51  (key = 5)
    ///
    /// GSL_INTEG_GAUSS61  (key = 6)
    ///
    /// corresponding to the 15f64, 21f64, 31f64, 41f64, 51 and 61 point Gauss-Kronrod rules. The
    /// higher-order rules give better accuracy for smooth functions, while lower-order rules save
    /// time when the function contains local difficulties, such as discontinuities.
    ///
    /// On each iteration the adaptive integration strategy bisects the interval with the largest
    /// error estimate. The subintervals and their results are stored in the memory provided by
    /// workspace. The maximum number of subintervals is given by limit, which may not exceed the
    /// allocated size of the workspace.
    ///
    /// Returns `(result, abs_err)`.
    #[doc(alias = "gsl_integration_qag")]
    pub fn qag<F: Fn(f64) -> f64>(
        &mut self,
        f: F,
        a: f64,
        b: f64,
        epsabs: f64,
        epsrel: f64,
        limit: usize,
        key: enums::GaussKronrodRule,
    ) -> (::Value, f64, f64) {
        let mut result = 0.;
        let mut abs_err = 0.;
        let function = wrap_callback!(f, F);

        let ret = unsafe {
            sys::gsl_integration_qag(
                &function,
                a,
                b,
                epsabs,
                epsrel,
                limit,
                key.into(),
                self.unwrap_unique(),
                &mut result,
                &mut abs_err,
            )
        };
        (::Value::from(ret), result, abs_err)
    }

    /// This function applies the Gauss-Kronrod 21-point integration rule adaptively until an
    /// estimate of the integral of f over (a,b) is achieved within the desired absolute and
    /// relative error limits, epsabs and epsrel. The results are extrapolated using the
    /// epsilon-algorithm, which accelerates the convergence of the integral in the presence of
    /// discontinuities and integrable singularities. The function returns the final approximation
    /// from the extrapolation, result, and an estimate of the absolute error, abserr. The
    /// subintervals and their results are stored in the memory provided by workspace. The maximum
    /// number of subintervals is given by limit, which may not exceed the allocated size of the
    /// workspace.
    ///
    /// Returns `(result, abs_err)`.
    #[doc(alias = "gsl_integration_qags")]
    pub fn qags<F: Fn(f64) -> f64>(
        &mut self,
        f: F,
        a: f64,
        b: f64,
        epsabs: f64,
        epsrel: f64,
        limit: usize,
    ) -> (::Value, f64, f64) {
        let mut result = 0.;
        let mut abs_err = 0.;
        let function = wrap_callback!(f, F);

        let ret = unsafe {
            sys::gsl_integration_qags(
                &function,
                a,
                b,
                epsabs,
                epsrel,
                limit,
                self.unwrap_unique(),
                &mut result,
                &mut abs_err,
            )
        };
        (::Value::from(ret), result, abs_err)
    }

    /// This function applies the adaptive integration algorithm QAGS taking account of the
    /// user-supplied locations of singular points. The array pts of length npts should contain the
    /// endpoints of the integration ranges defined by the integration region and locations of the
    /// singularities.
    ///
    /// For example, to integrate over the region (a,b) with break-points at x_1, x_2, x_3
    /// (where a < x_1 < x_2 < x_3 < b) the following pts array should be used
    ///
    /// ```text
    /// pts[0] = a
    /// pts[1] = x_1
    /// pts[2] = x_2
    /// pts[3] = x_3
    /// pts[4] = b
    /// with npts = 5.
    /// ```
    ///
    /// If you know the locations of the singular points in the integration region then this routine
    /// will be faster than QAGS.
    ///
    /// Returns `(result, abs_err)`.
    #[doc(alias = "gsl_integration_qagp")]
    pub fn qagp<F: Fn(f64) -> f64>(
        &mut self,
        f: F,
        pts: &mut [f64],
        epsabs: f64,
        epsrel: f64,
        limit: usize,
    ) -> (::Value, f64, f64) {
        let mut result = 0.;
        let mut abs_err = 0.;
        let function = wrap_callback!(f, F);

        let ret = unsafe {
            sys::gsl_integration_qagp(
                &function,
                pts.as_mut_ptr(),
                pts.len() as _,
                epsabs,
                epsrel,
                limit,
                self.unwrap_unique(),
                &mut result,
                &mut abs_err,
            )
        };
        (::Value::from(ret), result, abs_err)
    }

    /// This function computes the integral of the function f over the infinite interval
    /// `(-\infty,+\infty)`. The integral is mapped onto the semi-open interval `(0,1]` using the
    /// transformation:
    ///
    /// ```text
    /// x = (1-t)/t,
    ///
    /// \int_{-\infty}^{+\infty} dx f(x) =
    ///      \int_0^1 dt (f((1-t)/t) + f((-1+t)/t))/t^2.
    /// ```
    ///
    /// It is then integrated using the QAGS algorithm. The normal 21-point Gauss-Kronrod rule of
    /// QAGS is replaced by a 15-point rule, because the transformation can generate an integrable
    /// singularity at the origin. In this case a lower-order rule is more efficient.
    ///
    /// Returns `(result, abs_err)`.
    #[doc(alias = "gsl_integration_qagi")]
    pub fn qagi<F: Fn(f64) -> f64>(
        &mut self,
        f: F,
        epsabs: f64,
        epsrel: f64,
        limit: usize,
    ) -> (::Value, f64, f64) {
        let mut result = 0.;
        let mut abs_err = 0.;
        let mut function = wrap_callback!(f, F);

        let ret = unsafe {
            sys::gsl_integration_qagi(
                &mut function,
                epsabs,
                epsrel,
                limit,
                self.unwrap_unique(),
                &mut result,
                &mut abs_err,
            )
        };
        (::Value::from(ret), result, abs_err)
    }

    /// This function computes the integral of the function f over the semi-infinite interval
    /// `(a,+\infty)`. The integral is mapped onto the semi-open interval `(0,1]` using the
    /// transformation:
    ///
    /// ```text
    /// x = a + (1-t)/t,
    ///
    /// \int_{a}^{+\infty} dx f(x) =
    ///      \int_0^1 dt f(a + (1-t)/t)/t^2
    /// ```
    ///
    /// and then integrated using the QAGS algorithm.
    ///
    /// Returns `(result, abs_err)`.
    #[doc(alias = "gsl_integration_qagiu")]
    pub fn qagiu<F: Fn(f64) -> f64>(
        &mut self,
        f: F,
        a: f64,
        epsabs: f64,
        epsrel: f64,
        limit: usize,
    ) -> (::Value, f64, f64) {
        let mut result = 0.;
        let mut abs_err = 0.;
        let mut function = wrap_callback!(f, F);

        let ret = unsafe {
            sys::gsl_integration_qagiu(
                &mut function,
                a,
                epsabs,
                epsrel,
                limit,
                self.unwrap_unique(),
                &mut result,
                &mut abs_err,
            )
        };
        (::Value::from(ret), result, abs_err)
    }

    /// This function computes the integral of the function f over the semi-infinite interval
    /// `(-\infty,b)`. The integral is mapped onto the semi-open interval `(0,1]` using the
    /// transformation:
    ///
    /// ```text
    ///  x = b - (1-t)/t,
    ///
    /// \int_{-\infty}^{b} dx f(x) =
    ///      \int_0^1 dt f(b - (1-t)/t)/t^2
    /// ```
    ///
    /// and then integrated using the QAGS algorithm.
    ///
    /// Returns `(result, abs_err)`.
    #[doc(alias = "gsl_integration_qagil")]
    pub fn qagil<F: Fn(f64) -> f64>(
        &mut self,
        f: F,
        b: f64,
        epsabs: f64,
        epsrel: f64,
        limit: usize,
    ) -> (::Value, f64, f64) {
        let mut result = 0.;
        let mut abs_err = 0.;
        let mut function = wrap_callback!(f, F);

        let ret = unsafe {
            sys::gsl_integration_qagil(
                &mut function,
                b,
                epsabs,
                epsrel,
                limit,
                self.unwrap_unique(),
                &mut result,
                &mut abs_err,
            )
        };
        (::Value::from(ret), result, abs_err)
    }

    /// This function computes the Cauchy principal value of the integral of f over `(a,b)`, with a
    /// singularity at c,
    ///
    /// ```text
    /// I = \int_a^b dx f(x) / (x - c)
    /// ```
    ///
    /// The adaptive bisection algorithm of QAG is used, with modifications to ensure that
    /// subdivisions do not occur at the singular point x = c.
    ///
    /// When a subinterval contains the point x = c or is close to it then a special 25-point
    /// modified Clenshaw-Curtis rule is used to control the singularity. Further away from the
    /// singularity the algorithm uses an ordinary 15-point Gauss-Kronrod integration rule.
    ///
    /// Returns `(result, abs_err)`.
    #[doc(alias = "gsl_integration_qawc")]
    pub fn qawc<F: Fn(f64) -> f64>(
        &mut self,
        f: F,
        a: f64,
        b: f64,
        c: f64,
        epsabs: f64,
        epsrel: f64,
        limit: usize,
    ) -> (::Value, f64, f64) {
        let mut result = 0.;
        let mut abs_err = 0.;
        let mut function = wrap_callback!(f, F);

        let ret = unsafe {
            sys::gsl_integration_qawc(
                &mut function,
                a,
                b,
                c,
                epsabs,
                epsrel,
                limit,
                self.unwrap_unique(),
                &mut result,
                &mut abs_err,
            )
        };
        (::Value::from(ret), result, abs_err)
    }
}

ffi_wrapper!(
    IntegrationQawsTable,
    *mut sys::gsl_integration_qaws_table,
    gsl_integration_qaws_table_free,
    "The QAWS algorithm is designed for integrands with algebraic-logarithmic singularities at the
end-points of an integration region. In order to work efficiently the algorithm requires a
precomputed table of Chebyshev moments."
);

impl IntegrationQawsTable {
    /// This function allocates space for a gsl_integration_qaws_table struct describing a singular
    /// weight function W(x) with the parameters `alpha`, `beta`, `mu` and `nu`,
    ///
    /// ```text
    /// W(x) = (x-a)^alpha (b-x)^beta log^mu (x-a) log^nu (b-x)
    /// ```
    ///
    /// where `alpha > -1f64`, `beta > -1f64`, and `mu = 0, 1`, `nu = 0, 1`. The weight function can
    /// take four different forms depending on the values of `mu` and `nu`,
    ///
    /// ```text
    /// W(x) = (x-a)^alpha (b-x)^beta                   (mu = 0, nu = 0)
    /// W(x) = (x-a)^alpha (b-x)^beta log(x-a)          (mu = 1, nu = 0)
    /// W(x) = (x-a)^alpha (b-x)^beta log(b-x)          (mu = 0, nu = 1)
    /// W(x) = (x-a)^alpha (b-x)^beta log(x-a) log(b-x) (mu = 1, nu = 1)
    /// ```
    ///
    /// The singular points (a,b) do not have to be specified until the integral is computed, where
    /// they are the endpoints of the integration range.
    ///
    /// The function returns a pointer to the newly allocated table gsl_integration_qaws_table if no
    /// errors were detected, and 0 in the case of error.
    #[doc(alias = "gsl_integration_qaws_table_alloc")]
    pub fn new(alpha: f64, beta: f64, mu: i32, nu: i32) -> Option<IntegrationQawsTable> {
        let tmp = unsafe { sys::gsl_integration_qaws_table_alloc(alpha, beta, mu, nu) };

        if tmp.is_null() {
            None
        } else {
            Some(Self::wrap(tmp))
        }
    }

    /// This function modifies the parameters (\alpha, \beta, \mu, \nu)
    #[doc(alias = "gsl_integration_qaws_table_set")]
    pub fn set(&mut self, alpha: f64, beta: f64, mu: i32, nu: i32) -> ::Value {
        ::Value::from(unsafe {
            sys::gsl_integration_qaws_table_set(self.unwrap_unique(), alpha, beta, mu, nu)
        })
    }

    /// This function computes the integral of the function f(x) over the interval (a,b) with the
    /// singular weight function `(x-a)^\alpha (b-x)^\beta \log^\mu (x-a) \log^\nu (b-x)`. The
    /// parameters of the weight function (\alpha, \beta, \mu, \nu) are taken from the table self.
    /// The integral is,
    ///
    /// ```text
    /// I = \int_a^b dx f(x) (x-a)^alpha (b-x)^beta log^mu (x-a) log^nu (b-x).
    /// ```
    ///
    /// The adaptive bisection algorithm of QAG is used. When a subinterval contains one of the
    /// endpoints then a special 25-point modified Clenshaw-Curtis rule is used to control the
    /// singularities. For subintervals which do not include the endpoints an ordinary 15-point
    /// Gauss-Kronrod integration rule is used.
    ///
    /// Returns `(result, abs_err)`
    #[doc(alias = "gsl_integration_qaws")]
    pub fn qaws<F: Fn(f64) -> f64>(
        &mut self,
        f: F,
        a: f64,
        b: f64,
        epsabs: f64,
        epsrel: f64,
        limit: usize,
        workspace: &mut IntegrationWorkspace,
    ) -> (::Value, f64, f64) {
        let mut result = 0.;
        let mut abs_err = 0.;
        let mut function = wrap_callback!(f, F);

        let ret = unsafe {
            sys::gsl_integration_qaws(
                &mut function,
                a,
                b,
                self.unwrap_unique(),
                epsabs,
                epsrel,
                limit,
                workspace.unwrap_unique(),
                &mut result,
                &mut abs_err,
            )
        };
        (::Value::from(ret), result, abs_err)
    }
}

ffi_wrapper!(
    IntegrationQawoTable,
    *mut sys::gsl_integration_qawo_table,
    gsl_integration_qawo_table_free,
    "The QAWO algorithm is designed for integrands with an oscillatory factor, `sin(omega x)` or
`cos(omega x)`. In order to work efficiently the algorithm requires a table of Chebyshev moments
which must be pre-computed with calls to the functions below."
);

impl IntegrationQawoTable {
    /// This function allocates space for a gsl_integration_qawo_table struct and its associated
    /// workspace describing a sine or cosine weight function W(x) with the parameters (\omega, L),
    ///
    /// ```text
    /// W(x) = sin(omega x)
    /// W(x) = cos(omega x)
    /// ```
    ///
    /// The parameter L must be the length of the interval over which the function will be
    /// integrated L = b - a. The choice of sine or cosine is made with the parameter sine which
    /// should be chosen from one of the two following symbolic values:
    ///
    /// ```text
    /// ::Cosine
    /// ::IntegrationQawo::Sine
    /// ```
    ///
    /// The gsl_integration_qawo_table is a table of the trigonometric coefficients required in the
    /// integration process. The parameter n determines the number of levels of coefficients that
    /// are computed. Each level corresponds to one bisection of the interval L, so that n levels
    /// are sufficient for subintervals down to the length L/2^n. The integration routine
    /// gsl_integration_qawo returns the error ::Table if the number of levels is insufficient for
    /// the requested accuracy.
    #[doc(alias = "gsl_integration_qawo_table_alloc")]
    pub fn new(
        omega: f64,
        l: f64,
        sine: ::IntegrationQawo,
        n: usize,
    ) -> Option<IntegrationQawoTable> {
        let tmp = unsafe { sys::gsl_integration_qawo_table_alloc(omega, l, sine.into(), n) };

        if tmp.is_null() {
            None
        } else {
            Some(Self::wrap(tmp))
        }
    }

    /// This function changes the parameters omega, L and sine of the existing self workspace.
    #[doc(alias = "gsl_integration_qawo_table_set")]
    pub fn set(&mut self, omega: f64, l: f64, sine: ::IntegrationQawo) -> ::Value {
        ::Value::from(unsafe {
            sys::gsl_integration_qawo_table_set(self.unwrap_unique(), omega, l, sine.into())
        })
    }

    /// This function allows the length parameter l of the self workspace to be changed.
    #[doc(alias = "gsl_integration_qawo_table_set_length")]
    pub fn set_length(&mut self, l: f64) -> ::Value {
        ::Value::from(unsafe {
            sys::gsl_integration_qawo_table_set_length(self.unwrap_unique(), l)
        })
    }

    /// This function uses an adaptive algorithm to compute the integral of f over (a,b) with the
    /// weight function \sin(\omega x) or \cos(\omega x) defined by the table `wf`,
    ///
    /// I = \int_a^b dx f(x) sin(omega x)
    /// I = \int_a^b dx f(x) cos(omega x)
    ///
    /// The results are extrapolated using the epsilon-algorithm to accelerate the convergence of
    /// the integral. The function returns the final approximation from the extrapolation, result,
    /// and an estimate of the absolute error, abserr. The subintervals and their results are
    /// stored in the memory provided by workspace. The maximum number of subintervals is given by
    /// limit, which may not exceed the allocated size of the workspace.
    ///
    /// Those subintervals with “large” widths d where d\omega > 4 are computed using a 25-point
    /// Clenshaw-Curtis integration rule, which handles the oscillatory behavior. Subintervals with
    /// a "small" widths where d\omega < 4 are computed using a 15-point Gauss-Kronrod integration.
    ///
    /// Returns `(result, abserr)`.
    #[doc(alias = "gsl_integration_qawo")]
    pub fn qawo<F: Fn(f64) -> f64>(
        &mut self,
        f: F,
        a: f64,
        epsabs: f64,
        epsrel: f64,
        limit: usize,
        workspace: &mut IntegrationWorkspace,
    ) -> (::Value, f64, f64) {
        let mut function = wrap_callback!(f, F);
        let mut result = 0.;
        let mut abserr = 0.;

        let ret = unsafe {
            sys::gsl_integration_qawo(
                &mut function,
                a,
                epsabs,
                epsrel,
                limit,
                workspace.unwrap_unique(),
                self.unwrap_unique(),
                &mut result,
                &mut abserr,
            )
        };
        (::Value::from(ret), result, abserr)
    }
}

ffi_wrapper!(CquadWorkspace, *mut sys::gsl_integration_cquad_workspace, gsl_integration_cquad_workspace_free,
"CQUAD is a new doubly-adaptive general-purpose quadrature routine which can handle most types of
singularities, non-numerical function values such as Inf or NaN, as well as some divergent
integrals. It generally requires more function evaluations than the integration routines in
QUADPACK, yet fails less often for difficult integrands.

The underlying algorithm uses a doubly-adaptive scheme in which Clenshaw-Curtis quadrature rules
of increasing degree are used to compute the integral in each interval. The L_2-norm of the
difference between the underlying interpolatory polynomials of two successive rules is used as
an error estimate. The interval is subdivided if the difference between two successive rules is
too large or a rule of maximum degree has been reached.");

impl CquadWorkspace {
    /// This function allocates a workspace sufficient to hold the data for n intervals. The number
    /// n is not the maximum number of intervals that will be evaluated. If the workspace is full,
    /// intervals with smaller error estimates will be discarded. A minimum of 3 intervals
    /// is required and for most functions, a workspace of size 100 is sufficient.
    #[doc(alias = "gsl_integration_cquad_workspace_alloc")]
    pub fn new(n: usize) -> Option<CquadWorkspace> {
        let tmp = unsafe { sys::gsl_integration_cquad_workspace_alloc(n) };

        if tmp.is_null() {
            None
        } else {
            Some(Self::wrap(tmp))
        }
    }

    /// This function computes the integral of f over (a,b) within the desired absolute and relative
    /// error limits, epsabs and epsrel using the CQUAD algorithm. The function returns the final
    /// approximation, result, an estimate of the absolute error, abserr, and the number of function
    /// evaluations required, nevals.
    ///
    /// The CQUAD algorithm divides the integration region into subintervals, and in each iteration,
    /// the subinterval with the largest estimated error is processed. The algorithm uses
    /// Clenshaw-Curits quadrature rules of degree 4, 8, 16 and 32 over 5, 9, 17 and 33 nodes
    /// respectively. Each interval is initialized with the lowest-degree rule. When an interval is
    /// processed, the next-higher degree rule is evaluated and an error estimate is computed based
    /// on the L_2-norm of the difference between the underlying interpolating polynomials of both
    /// rules. If the highest-degree rule has already been used, or the interpolatory polynomials
    /// differ significantly, the interval is bisected.
    ///
    /// The subintervals and their results are stored in the memory provided by workspace. If the
    /// error estimate or the number of function evaluations is not needed, the pointers abserr and
    /// nevals can be set to NULL (not in rgsl).
    ///
    /// Returns `(result, abs_err, n_evals)`.
    #[doc(alias = "gsl_integration_cquad")]
    pub fn cquad<F: Fn(f64) -> f64>(
        &mut self,
        f: F,
        a: f64,
        b: f64,
        epsabs: f64,
        epsrel: f64,
    ) -> (::Value, f64, f64, usize) {
        let function = wrap_callback!(f, F);
        let mut result = 0.;
        let mut abs_err = 0.;
        let mut n_evals = 0;

        let ret = unsafe {
            sys::gsl_integration_cquad(
                &function,
                a,
                b,
                epsabs,
                epsrel,
                self.unwrap_unique(),
                &mut result,
                &mut abs_err,
                &mut n_evals,
            )
        };
        (::Value::from(ret), result, abs_err, n_evals)
    }
}

ffi_wrapper!(GLFixedTable, *mut sys::gsl_integration_glfixed_table, gsl_integration_glfixed_table_free,
"The fixed-order Gauss-Legendre integration routines are provided for fast integration of smooth
functions with known polynomial order. The n-point Gauss-Legendre rule is exact for polynomials
of order 2*n-1 or less. For example, these rules are useful when integrating basis functions to
form mass matrices for the Galerkin method. Unlike other numerical integration routines within
the library, these routines do not accept absolute or relative error bounds.");

impl GLFixedTable {
    /// This function determines the Gauss-Legendre abscissae and weights necessary for an n-point
    /// fixed order integration scheme. If possible, high precision precomputed coefficients are
    /// used. If precomputed weights are not available, lower precision coefficients are computed
    /// on the fly.
    #[doc(alias = "gsl_integration_glfixed_table_alloc")]
    pub fn new(n: usize) -> Option<GLFixedTable> {
        let tmp = unsafe { sys::gsl_integration_glfixed_table_alloc(n) };

        if tmp.is_null() {
            None
        } else {
            Some(Self::wrap(tmp))
        }
    }

    /// For i in [0, …, t->n - 1], this function obtains the i-th Gauss-Legendre point xi and weight
    /// wi on the interval [a,b]. The points and weights are ordered by increasing point value. A
    /// function f may be integrated on [a,b] by summing wi * f(xi) over i.
    ///
    /// Returns `(Value, xi, wi)`.
    #[doc(alias = "gsl_integration_glfixed_point")]
    pub fn point(&self, a: f64, b: f64, i: usize) -> (::Value, f64, f64) {
        let mut xi = 0.;
        let mut wi = 0.;
        let ret = unsafe {
            sys::gsl_integration_glfixed_point(a, b, i, &mut xi, &mut wi, self.unwrap_shared())
        };
        (::Value::from(ret), xi, wi)
    }

    /// This function applies the Gauss-Legendre integration rule contained in table self and
    /// returns the result.
    #[doc(alias = "gsl_integration_glfixed")]
    pub fn glfixed<F: Fn(f64) -> f64>(&self, f: F, a: f64, b: f64) -> f64 {
        let function = wrap_callback!(f, F);
        unsafe { sys::gsl_integration_glfixed(&function, a, b, self.unwrap_shared()) }
    }

    #[doc(alias = "gsl_integration_glfixed_point")]
    pub fn glfixed_point(&self, a: f64, b: f64, xi: &mut [f64], wi: &mut [f64]) -> Value {
        assert!(xi.len() == wi.len());

        Value::from(unsafe {
            sys::gsl_integration_glfixed_point(
                a,
                b,
                xi.len() as _,
                xi.as_mut_ptr(),
                wi.as_mut_ptr(),
                self.unwrap_shared(),
            )
        })
    }
}