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
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
//
// A rust binding for the GSL library by Guillaume Gomez (guillaume1.gomez@gmail.com)
//

/*!
#Monte Carlo Integration

This chapter describes routines for multidimensional Monte Carlo integration. These include the traditional Monte Carlo method and adaptive 
algorithms such as VEGAS and MISER which use importance sampling and stratified sampling techniques. Each algorithm computes an estimate of 
a multidimensional definite integral of the form,

I = \int_xl^xu dx \int_yl^yu  dy ...  f(x, y, ...)

over a hypercubic region ((x_l,x_u), (y_l,y_u), ...) using a fixed number of function calls. The routines also provide a statistical estimate 
of the error on the result. This error estimate should be taken as a guide rather than as a strict error bound—random sampling of the region 
may not uncover all the important features of the function, resulting in an underestimate of the error.

##Interface

All of the Monte Carlo integration routines use the same general form of interface. There is an allocator to allocate memory for control 
variables and workspace, a routine to initialize those control variables, the integrator itself, and a function to free the space when done.

Each integration function requires a random number generator to be supplied, and returns an estimate of the integral and its standard deviation. 
The accuracy of the result is determined by the number of function calls specified by the user. If a known level of accuracy is required this 
can be achieved by calling the integrator several times and averaging the individual results until the desired accuracy is obtained.

Random sample points used within the Monte Carlo routines are always chosen strictly within the integration region, so that endpoint singularities 
are automatically avoided.

##VEGAS

The VEGAS algorithm of Lepage is based on importance sampling. It samples points from the probability distribution described by the function 
|f|, so that the points are concentrated in the regions that make the largest contribution to the integral.

In general, if the Monte Carlo integral of f is sampled with points distributed according to a probability distribution described by the function 
g, we obtain an estimate E_g(f; N),

E_g(f; N) = E(f/g; N)
with a corresponding variance,

\Var_g(f; N) = \Var(f/g; N).
If the probability distribution is chosen as g = |f|/I(|f|) then it can be shown that the variance V_g(f; N) vanishes, and the error in the 
estimate will be zero. In practice it is not possible to sample from the exact distribution g for an arbitrary function, so importance sampling 
algorithms aim to produce efficient approximations to the desired distribution.

The VEGAS algorithm approximates the exact distribution by making a number of passes over the integration region while histogramming the 
function f. Each histogram is used to define a sampling distribution for the next pass. Asymptotically this procedure converges to the desired 
distribution. In order to avoid the number of histogram bins growing like K^d the probability distribution is approximated by a separable 
function: g(x_1, x_2, ...) = g_1(x_1) g_2(x_2) ... so that the number of bins required is only Kd. This is equivalent to locating the 
peaks of the function from the projections of the integrand onto the coordinate axes. The efficiency of VEGAS depends on the validity of 
this assumption. It is most efficient when the peaks of the integrand are well-localized. If an integrand can be rewritten in a form which 
is approximately separable this will increase the efficiency of integration with VEGAS.

VEGAS incorporates a number of additional features, and combines both stratified sampling and importance sampling. The integration region 
is divided into a number of “boxes”, with each box getting a fixed number of points (the goal is 2). Each box can then have a fractional 
number of bins, but if the ratio of bins-per-box is less than two, Vegas switches to a kind variance reduction (rather than importance 
sampling).

The VEGAS algorithm computes a number of independent estimates of the integral internally, according to the iterations parameter described 
below, and returns their weighted average. Random sampling of the integrand can occasionally produce an estimate where the error is zero, 
particularly if the function is constant in some regions. An estimate with zero error causes the weighted average to break down and must 
be handled separately. In the original Fortran implementations of VEGAS the error estimate is made non-zero by substituting a small value 
(typically 1e-30). The implementation in GSL differs from this and avoids the use of an arbitrary constant—it either assigns the value a 
weight which is the average weight of the preceding estimates or discards it according to the following procedure,

current estimate has zero error, weighted average has finite error
The current estimate is assigned a weight which is the average weight of the preceding estimates.

current estimate has finite error, previous estimates had zero error
The previous estimates are discarded and the weighted averaging procedure begins with the current estimate.

current estimate has zero error, previous estimates had zero error
The estimates are averaged using the arithmetic mean, but no error is computed.
!*/

use ffi;
use libc::{c_void, c_double, size_t};
use std::slice;
use std::mem::transmute;
use std::marker::PhantomData;

/// The plain Monte Carlo algorithm samples points randomly from the integration region to estimate the integral and its error. Using this algorithm
/// the estimate of the integral E(f; N) for N randomly distributed points x_i is given by,
///
/// E(f; N) = =  V <f> = (V / N) \sum_i^N f(x_i)
/// where V is the volume of the integration region. The error on this estimate \sigma(E;N) is calculated from the estimated variance of the mean,
///
/// \sigma^2 (E; N) = (V^2 / N^2) \sum_i^N (f(x_i) -  <f>)^2.
///
/// For large N this variance decreases asymptotically as \Var(f)/N, where \Var(f) is the true variance of the function over the integration region.
/// The error estimate itself should decrease as \sigma(f)/\sqrt{N}. The familiar law of errors decreasing as 1/\sqrt{N} applies—to reduce the
/// error by a factor of 10 requires a 100-fold increase in the number of sample points.
pub struct PlainMonteCarlo {
    s: *mut ffi::gsl_monte_plain_state,
}

impl PlainMonteCarlo {
    /// This function allocates and initializes a workspace for Monte Carlo integration in dim dimensions.
    pub fn new(dim: usize) -> Option<PlainMonteCarlo> {
        let tmp = unsafe { ffi::gsl_monte_plain_alloc(dim) };

        if tmp.is_null() {
            None
        } else {
            Some(PlainMonteCarlo { s: tmp })
        }
    }

    /// This function initializes a previously allocated integration state. This allows an existing workspace to be reused for different
    /// integrations.
    pub fn init(&mut self) -> ::Value {
        unsafe { ffi::gsl_monte_plain_init(self.s) }
    }

    /// This routines uses the plain Monte Carlo algorithm to integrate the function f over the dim-dimensional hypercubic region defined
    /// by the lower and upper limits in the arrays xl and xu, each of the same size. The integration uses a fixed number of function calls
    /// calls, and obtains random sampling points using the random number generator r. A previously allocated workspace s must be supplied.
    /// The result of the integration is returned in result, with an estimated absolute error abserr.
    ///
    /// In C, the function takes a `gsl_monte_function` as first argument. In here, you have to
    /// pass the `dim` argument and the function pointer (which became a closure) directly to the
    /// function.
    ///
    /// It returns either Ok((result, abserr)) or Err(enums::Value).
    pub fn integrate<F: FnMut(&[f64]) -> f64>(&mut self,
                                              dim: usize,
                                              f: F,
                                              xl: &[f64],
                                              xu: &[f64],
                                              t_calls: usize,
                                              r: &mut ::Rng)
                                              -> Result<(f64, f64), ::Value> {
        unsafe {
            assert!(xl.len() == xu.len());
            let mut result = 0f64;
            let mut abserr = 0f64;
            let f: Box<Box<FnMut(&[f64]) -> f64>> = Box::new(Box::new(f));
            let mut func = ffi::gsl_monte_function {
                f: transmute(monte_trampoline as usize),
                dim: dim,
                params: Box::into_raw(f) as *mut _,
            };
            let ret = ffi::gsl_monte_plain_integrate(&mut func as *mut _ as *mut c_void,
                                                     xl.as_ptr(),
                                                     xu.as_ptr(),
                                                     xl.len(),
                                                     t_calls,
                                                     ffi::FFI::unwrap_unique(r),
                                                     self.s,
                                                     (&mut result) as *mut c_double,
                                                     (&mut abserr) as *mut c_double);

            if ret == ::Value::Success {
                Ok((result, abserr))
            } else {
                Err(ret)
            }
        }
    }
}

impl Drop for PlainMonteCarlo {
    fn drop(&mut self) {
        unsafe { ffi::gsl_monte_plain_free(self.s) };
        self.s = ::std::ptr::null_mut();
    }
}

impl ffi::FFI<ffi::gsl_monte_plain_state> for PlainMonteCarlo {
    fn wrap(s: *mut ffi::gsl_monte_plain_state) -> PlainMonteCarlo {
        PlainMonteCarlo { s: s }
    }

    fn soft_wrap(s: *mut ffi::gsl_monte_plain_state) -> PlainMonteCarlo {
        Self::wrap(s)
    }

    fn unwrap_shared(s: &PlainMonteCarlo) -> *const ffi::gsl_monte_plain_state {
        s.s as *const _
    }

    fn unwrap_unique(s: &mut PlainMonteCarlo) -> *mut ffi::gsl_monte_plain_state {
        s.s
    }
}

/// The MISER algorithm of Press and Farrar is based on recursive stratified sampling. This technique aims to reduce the overall integration error
/// by concentrating integration points in the regions of highest variance.
///
/// The idea of stratified sampling begins with the observation that for two disjoint regions a and b with Monte Carlo estimates of the integral
/// E_a(f) and E_b(f) and variances \sigma_a^2(f) and \sigma_b^2(f), the variance \Var(f) of the combined estimate E(f) = (1/2) (E_a(f) + E_b(f))
/// is given by,
///
/// \Var(f) = (\sigma_a^2(f) / 4 N_a) + (\sigma_b^2(f) / 4 N_b).
///
/// It can be shown that this variance is minimized by distributing the points such that,
///
/// N_a / (N_a + N_b) = \sigma_a / (\sigma_a + \sigma_b).
///
/// Hence the smallest error estimate is obtained by allocating sample points in proportion to the standard deviation of the function in each
/// sub-region.
///
/// The MISER algorithm proceeds by bisecting the integration region along one coordinate axis to give two sub-regions at each step. The direction
/// is chosen by examining all d possible bisections and selecting the one which will minimize the combined variance of the two sub-regions. The
/// variance in the sub-regions is estimated by sampling with a fraction of the total number of points available to the current step. The same
/// procedure is then repeated recursively for each of the two half-spaces from the best bisection. The remaining sample points are allocated to
/// the sub-regions using the formula for N_a and N_b. This recursive allocation of integration points continues down to a user-specified depth
/// where each sub-region is integrated using a plain Monte Carlo estimate. These individual values and their error estimates are then combined
/// upwards to give an overall result and an estimate of its error.
pub struct MiserMonteCarlo {
    s: *mut ffi::gsl_monte_miser_state,
}

impl MiserMonteCarlo {
    /// This function allocates and initializes a workspace for Monte Carlo integration in dim dimensions. The workspace is used to maintain
    /// the state of the integration.
    pub fn new(dim: usize) -> Option<MiserMonteCarlo> {
        let tmp_pointer = unsafe { ffi::gsl_monte_miser_alloc(dim) };

        if tmp_pointer.is_null() {
            None
        } else {
            Some(MiserMonteCarlo { s: tmp_pointer })
        }
    }

    /// This function initializes a previously allocated integration state. This allows an existing workspace to be reused for different integrations.
    pub fn init(&mut self) -> ::Value {
        unsafe { ffi::gsl_monte_miser_init(self.s) }
    }

    /// This routines uses the MISER Monte Carlo algorithm to integrate the function f over the dim-dimensional hypercubic region defined by
    /// the lower and upper limits in the arrays xl and xu, each of size dim. The integration uses a fixed number of function calls calls,
    /// and obtains random sampling points using the random number generator r. A previously allocated workspace s must be supplied. The result
    /// of the integration is returned in result, with an estimated absolute error abserr.
    ///
    /// In C, the function takes a `gsl_monte_function` as first argument. In here, you have to
    /// pass the `dim` argument and the function pointer (which became a closure) directly to the
    /// function.
    ///
    /// It returns either Ok((result, abserr)) or Err(enums::Value).
    pub fn integrate<F: FnMut(&[f64]) -> f64>(&mut self,
                                              dim: usize,
                                              f: F,
                                              xl: &[f64],
                                              xu: &[f64],
                                              t_calls: usize,
                                              r: &mut ::Rng)
                                              -> Result<(f64, f64), ::Value> {
        unsafe {
            assert!(xl.len() == xu.len());
            let mut result = 0f64;
            let mut abserr = 0f64;
            let f: Box<Box<FnMut(&[f64]) -> f64>> = Box::new(Box::new(f));
            let mut func = ffi::gsl_monte_function {
                f: transmute(monte_trampoline as usize),
                dim: dim,
                params: Box::into_raw(f) as *mut _,
            };
            let ret = ffi::gsl_monte_miser_integrate(&mut func as *mut _ as *mut c_void,
                                                     xl.as_ptr(),
                                                     xu.as_ptr(),
                                                     xl.len(),
                                                     t_calls,
                                                     ffi::FFI::unwrap_unique(r),
                                                     self.s,
                                                     (&mut result) as *mut c_double,
                                                     (&mut abserr) as *mut c_double);

            if ret == ::Value::Success {
                Ok((result, abserr))
            } else {
                Err(ret)
            }
        }
    }

    /// This function copies the parameters of the integrator state into the user-supplied params structure.
    pub fn get_params(&self) -> MiserParams {
        let mut m = MiserParams {
            estimate_frac: 0f64,
            min_calls: 0,
            min_calls_per_bisection: 0,
            alpha: 0f64,
            dither: 0f64,
        };

        unsafe {
            ffi::gsl_monte_miser_params_get(self.s, &mut m as *mut MiserParams);
        }
        m
    }

    /// This function sets the integrator parameters based on values provided in the params structure.
    pub fn set_params(&mut self, params: &MiserParams) {
        unsafe {
            ffi::gsl_monte_miser_params_set(self.s, params as *const MiserParams);
        }
    }
}

impl Drop for MiserMonteCarlo {
    fn drop(&mut self) {
        unsafe { ffi::gsl_monte_miser_free(self.s) };
        self.s = ::std::ptr::null_mut();
    }
}

impl ffi::FFI<ffi::gsl_monte_miser_state> for MiserMonteCarlo {
    fn wrap(s: *mut ffi::gsl_monte_miser_state) -> MiserMonteCarlo {
        MiserMonteCarlo { s: s }
    }

    fn soft_wrap(s: *mut ffi::gsl_monte_miser_state) -> MiserMonteCarlo {
        Self::wrap(s)
    }

    fn unwrap_shared(s: &MiserMonteCarlo) -> *const ffi::gsl_monte_miser_state {
        s.s as *const _
    }

    fn unwrap_unique(s: &mut MiserMonteCarlo) -> *mut ffi::gsl_monte_miser_state {
        s.s
    }
}

#[repr(C)]
pub struct MiserParams {
    /// This parameter specifies the fraction of the currently available number of function calls which
    /// are allocated to estimating the variance at each recursive step. The default value is 0.1.
    pub estimate_frac: f64,
    /// This parameter specifies the minimum number of function calls required for each estimate of the
    /// variance. If the number of function calls allocated to the estimate using estimate_frac falls
    /// below min_calls then min_calls are used instead. This ensures that each estimate maintains a
    /// reasonable level of accuracy. The default value of min_calls is 16 * dim.
    pub min_calls: usize,
    /// This parameter specifies the minimum number of function calls required to proceed with a bisection
    /// step. When a recursive step has fewer calls available than min_calls_per_bisection it performs
    /// a plain Monte Carlo estimate of the current sub-region and terminates its branch of the recursion.
    /// The default value of this parameter is 32 * min_calls.
    pub min_calls_per_bisection: usize,
    /// This parameter controls how the estimated variances for the two sub-regions of a bisection are
    /// combined when allocating points. With recursive sampling the overall variance should scale better
    /// than 1/N, since the values from the sub-regions will be obtained using a procedure which explicitly
    /// minimizes their variance. To accommodate this behavior the MISER algorithm allows the total variance
    /// to depend on a scaling parameter \alpha,
    ///
    /// \Var(f) = {\sigma_a \over N_a^\alpha} + {\sigma_b \over N_b^\alpha}.
    ///
    /// The authors of the original paper describing MISER recommend the value \alpha = 2 as a good choice,
    /// obtained from numerical experiments, and this is used as the default value in this implementation.
    pub alpha: f64,
    /// This parameter introduces a random fractional variation of size dither into each bisection, which
    /// can be used to break the symmetry of integrands which are concentrated near the exact center of
    /// the hypercubic integration region. The default value of dither is zero, so no variation is introduced.
    /// If needed, a typical value of dither is 0.1.
    pub dither: f64,
}

/// The VEGAS algorithm of Lepage is based on importance sampling. It samples points from the probability
/// distribution described by the function |f|, so that the points are concentrated in the regions that
/// make the largest contribution to the integral.
///
/// In general, if the Monte Carlo integral of f is sampled with points distributed according to a
/// probability distribution described by the function g, we obtain an estimate E_g(f; N),
///
/// E_g(f; N) = E(f/g; N)
///
/// with a corresponding variance,
///
/// \Var_g(f; N) = \Var(f/g; N).
///
/// If the probability distribution is chosen as g = |f|/I(|f|) then it can be shown that the variance
/// V_g(f; N) vanishes, and the error in the estimate will be zero. In practice it is not possible to
/// sample from the exact distribution g for an arbitrary function, so importance sampling algorithms
/// aim to produce efficient approximations to the desired distribution.
///
/// The VEGAS algorithm approximates the exact distribution by making a number of passes over the
/// integration region while histogramming the function f. Each histogram is used to define a sampling
/// distribution for the next pass. Asymptotically this procedure converges to the desired distribution.
/// In order to avoid the number of histogram bins growing like K^d the probability distribution is
/// approximated by a separable function: g(x_1, x_2, ...) = g_1(x_1) g_2(x_2) ... so that the number
/// of bins required is only Kd. This is equivalent to locating the peaks of the function from the
/// projections of the integrand onto the coordinate axes. The efficiency of VEGAS depends on the
/// validity of this assumption. It is most efficient when the peaks of the integrand are well-localized.
/// If an integrand can be rewritten in a form which is approximately separable this will increase
/// the efficiency of integration with VEGAS.
///
/// VEGAS incorporates a number of additional features, and combines both stratified sampling and
/// importance sampling. The integration region is divided into a number of “boxes”, with each box
/// getting a fixed number of points (the goal is 2). Each box can then have a fractional number of
/// bins, but if the ratio of bins-per-box is less than two, Vegas switches to a kind variance reduction
/// (rather than importance sampling).
///
/// The VEGAS algorithm computes a number of independent estimates of the integral internally, according
/// to the iterations parameter described below, and returns their weighted average. Random sampling of
/// the integrand can occasionally produce an estimate where the error is zero, particularly if the function
/// is constant in some regions. An estimate with zero error causes the weighted average to break down and
/// must be handled separately. In the original Fortran implementations of VEGAS the error estimate is made
/// non-zero by substituting a small value (typically 1e-30). The implementation in GSL differs from this
/// and avoids the use of an arbitrary constant—it either assigns the value a weight which is the average
/// weight of the preceding estimates or discards it according to the following procedure,
///
/// current estimate has zero error, weighted average has finite error
///
/// * The current estimate is assigned a weight which is the average weight of the preceding estimates.
/// current estimate has finite error, previous estimates had zero error
///
/// * The previous estimates are discarded and the weighted averaging procedure begins with the current estimate.
/// current estimate has zero error, previous estimates had zero error
///
/// * The estimates are averaged using the arithmetic mean, but no error is computed.
pub struct VegasMonteCarlo {
    s: *mut ffi::gsl_monte_vegas_state,
}

impl VegasMonteCarlo {
    /// This function allocates and initializes a workspace for Monte Carlo integration in dim dimensions.
    /// The workspace is used to maintain the state of the integration.
    pub fn new(dim: usize) -> Option<VegasMonteCarlo> {
        let tmp_pointer = unsafe { ffi::gsl_monte_vegas_alloc(dim) };

        if tmp_pointer.is_null() {
            None
        } else {
            Some(VegasMonteCarlo { s: tmp_pointer })
        }
    }

    /// This function initializes a previously allocated integration state. This allows an existing workspace
    /// to be reused for different integrations.
    pub fn init(&mut self) -> ::Value {
        unsafe { ffi::gsl_monte_vegas_init(self.s) }
    }

    /// This routines uses the VEGAS Monte Carlo algorithm to integrate the function f over the dim-dimensional
    /// hypercubic region defined by the lower and upper limits in the arrays xl and xu, each of size dim.
    /// The integration uses a fixed number of function calls calls, and obtains random sampling points using
    /// the random number generator r. A previously allocated workspace s must be supplied. The result of the
    /// integration is returned in result, with an estimated absolute error abserr. The result and its error
    /// estimate are based on a weighted average of independent samples. The chi-squared per degree of freedom
    /// for the weighted average is returned via the state struct component, s->chisq, and must be consistent
    /// with 1 for the weighted average to be reliable.
    ///
    /// In C, the function takes a `gsl_monte_function` as first argument. In here, you have to
    /// pass the `dim` argument and the function pointer (which became a closure) directly to the
    /// function.
    ///
    /// It returns either Ok((result, abserr)) or Err(enums::Value).
    pub fn integrate<F: FnMut(&[f64]) -> f64>(&mut self,
                                              dim: usize,
                                              f: F,
                                              xl: &[f64],
                                              xu: &[f64],
                                              t_calls: usize,
                                              r: &mut ::Rng)
                                              -> Result<(f64, f64), ::Value> {
        unsafe {
            assert!(xl.len() == xu.len());
            let mut result = 0f64;
            let mut abserr = 0f64;
            let f: Box<Box<FnMut(&[f64]) -> f64>> = Box::new(Box::new(f));
            let mut func = ffi::gsl_monte_function {
                f: transmute(monte_trampoline as usize),
                dim: dim,
                params: Box::into_raw(f) as *mut _,
            };
            let ret = ffi::gsl_monte_vegas_integrate(&mut func as *mut _ as *mut c_void,
                                                     xl.as_ptr(),
                                                     xu.as_ptr(),
                                                     xl.len(),
                                                     t_calls,
                                                     ffi::FFI::unwrap_unique(r),
                                                     self.s,
                                                     (&mut result) as *mut c_double,
                                                     (&mut abserr) as *mut c_double);

            if ret == ::Value::Success {
                Ok((result, abserr))
            } else {
                Err(ret)
            }
        }
    }

    /// This function returns the chi-squared per degree of freedom for the weighted estimate of the integral.
    /// The returned value should be close to 1. A value which differs significantly from 1 indicates that
    /// the values from different iterations are inconsistent. In this case the weighted error will be
    /// under-estimated, and further iterations of the algorithm are needed to obtain reliable results.
    pub fn chisq(&mut self) -> f64 {
        unsafe { ffi::gsl_monte_vegas_chisq(self.s) }
    }

    /// This function returns the raw (unaveraged) values of the integral result and its error sigma from
    /// the most recent iteration of the algorithm.
    pub fn runval(&mut self, result: &mut f64, sigma: &mut f64) {
        unsafe {
            ffi::gsl_monte_vegas_runval(self.s, result as *mut c_double, sigma as *mut c_double)
        }
    }

    pub fn get_params(&self) -> VegasParams {
        let mut params = VegasParams::default();
        unsafe {
            ffi::gsl_monte_vegas_params_get(self.s, &mut params.inner as *mut _);
        }
        params
    }

    pub fn set_params(&mut self, params: &VegasParams) {
        unsafe {
            ffi::gsl_monte_vegas_params_set(self.s, &params.inner as *const _);
        }
    }
}

pub struct VegasParams<'a> {
    inner: ffi::gsl_monte_vegas_params,
    lt: PhantomData<&'a ()>,
}

impl<'a> VegasParams<'a> {
    /// alpha: The parameter alpha controls the stiffness of the rebinning algorithm. It is typically
    /// set between one and two. A value of zero prevents rebinning of the grid. The default
    /// value is 1.5.
    ///
    /// iterations: The number of iterations to perform for each call to the routine. The default value
    /// is 5 iterations.
    ///
    /// stage: Setting this determines the stage of the calculation. Normally, stage = 0 which begins
    /// with a new uniform grid and empty weighted average. Calling vegas with stage =
    /// 1 retains the grid from the previous run but discards the weighted average, so that
    /// one can “tune” the grid using a relatively small number of points and then do a large
    /// run with stage = 1 on the optimized grid. Setting stage = 2 keeps the grid and the
    /// weighted average from the previous run, but may increase (or decrease) the number
    /// of histogram bins in the grid depending on the number of calls available. Choosing
    /// stage = 3 enters at the main loop, so that nothing is changed, and is equivalent to
    /// performing additional iterations in a previous call.
    ///
    /// mode: The possible choices are GSL_VEGAS_MODE_IMPORTANCE, GSL_VEGAS_MODE_
    /// STRATIFIED, GSL_VEGAS_MODE_IMPORTANCE_ONLY. This determines whether vegas
    /// will use importance sampling or stratified sampling, or whether it can pick on
    /// its own. In low dimensions vegas uses strict stratified sampling (more precisely,
    /// stratified sampling is chosen if there are fewer than 2 bins per box).
    ///
    /// verbosity + stream: These parameters set the level of information printed by vegas.
    pub fn new(alpha: f64,
               iterations: usize,
               stage: i32,
               mode: ::VegasMode,
               verbosity: VegasVerbosity,
               stream: Option<&'a mut ::IOStream>)
               -> Result<VegasParams, String> {
        if !verbosity.is_off() && stream.is_none() {
            return Err("rust-GSL: need to provide an input stream for Vegas Monte Carlo \
                        integration if verbosity is not 'Off'"
                .to_string());
        } else if verbosity.is_off() && stream.is_some() {
            return Err("rust-GSL: need to provide the verbosity flag for Vegas Monta Carlo \
                        integration, currently set to 'Off'"
                .to_string());
        }

        let stream = if let Some(stream) = stream {
            if !stream.write_mode() {
                return Err("rust-GSL: input stream not flagged as 'write' mode".to_string());
            }
            stream.as_raw()
        } else {
            ::std::ptr::null_mut()
        };
        Ok(VegasParams {
            inner: ffi::gsl_monte_vegas_params {
                alpha: alpha,
                iterations: iterations,
                stage: stage,
                mode: mode,
                verbose: verbosity.to_int(),
                ostream: stream,
            },
            lt: PhantomData,
        })
    }
}

impl<'a> ::std::default::Default for VegasParams<'a> {
    fn default() -> VegasParams<'a> {
        VegasParams {
            inner: ffi::gsl_monte_vegas_params {
                alpha: 1.5,
                iterations: 5,
                stage: 0,
                mode: ::VegasMode::ImportanceOnly,
                verbose: -1,
                ostream: ::std::ptr::null_mut(),
            },
            lt: PhantomData,
        }
    }
}

/// The default setting of verbose is `Off`, which turns off all output.
/// A verbose value of `Summary` prints summary information about the weighted average
/// and final result, while a value of `Grid` also displays the grid coordinates.
/// A value of 'Rebinning' prints information from the rebinning procedure for each iteration.
#[derive(Clone, Copy)]
pub enum VegasVerbosity {
    Off, // -1
    Summary, // 0
    Grid, // 1
    Rebinning, // 2
}

impl VegasVerbosity {
    fn to_int(&self) -> i32 {
        match *self {
            VegasVerbosity::Off => -1,
            VegasVerbosity::Summary => 0,
            VegasVerbosity::Grid => 1,
            VegasVerbosity::Rebinning => 2,
        }
    }

    fn is_off(&self) -> bool {
        match *self {
            VegasVerbosity::Off => true,
            _ => false,
        }
    }
}

impl Drop for VegasMonteCarlo {
    fn drop(&mut self) {
        unsafe { ffi::gsl_monte_vegas_free(self.s) };
        self.s = ::std::ptr::null_mut();
    }
}

impl ffi::FFI<ffi::gsl_monte_vegas_state> for VegasMonteCarlo {
    fn wrap(s: *mut ffi::gsl_monte_vegas_state) -> VegasMonteCarlo {
        VegasMonteCarlo { s: s }
    }

    fn soft_wrap(s: *mut ffi::gsl_monte_vegas_state) -> VegasMonteCarlo {
        Self::wrap(s)
    }

    fn unwrap_shared(s: &VegasMonteCarlo) -> *const ffi::gsl_monte_vegas_state {
        s.s as *const _
    }

    fn unwrap_unique(s: &mut VegasMonteCarlo) -> *mut ffi::gsl_monte_vegas_state {
        s.s
    }
}

unsafe extern "C" fn monte_trampoline(x: *mut c_double,
                                      dim: size_t,
                                      param: *mut c_void)
                                      -> c_double {
    let f: &mut Box<FnMut(&[f64]) -> f64> = transmute(param);
    f(slice::from_raw_parts(x, dim as usize))
}

// The following tests have been made and tested against the following C code:
//
// ```ignore
// double
// g (double *k, size_t dim, void *params)
// {
//   (void)(dim);
//   (void)(params);
//   double A = 1.0 / (M_PI * M_PI * M_PI);
//   return A / (1.0 - cos (k[0]) * cos (k[1]) * cos (k[2]));
// }
//
// void
// display_results (char *title, double result, double error)
// {
//   printf ("%s ==================\n", title);
//   printf ("result = % .6f\n", result);
//   printf ("sigma  = % .6f\n", error);
//   printf ("exact  = % .6f\n", exact);
//   printf ("error  = % .6f = %.2g sigma\n", result - exact,
//       fabs (result - exact) / error);
// }
//
// int main(void) {
//   double res, err;
//
//   double xl[3] = { 0, 0, 0 };
//   double xu[3] = { M_PI, M_PI, M_PI };
//
//   const gsl_rng_type *T;
//   gsl_rng *r;
//
//   gsl_monte_function G = { &g, 3, 0 };
//
//   size_t calls = 500000;
//
//   gsl_rng_env_setup ();
//
//   T = gsl_rng_default;
//   r = gsl_rng_alloc (T);
//
//   {
//     gsl_monte_vegas_state *s = gsl_monte_vegas_alloc (3);
//
//     gsl_monte_vegas_integrate (&G, xl, xu, 3, 10000, r, s,
//                    &res, &err);
//     display_results ("vegas warm-up", res, err);
//
//     printf ("converging...\n");
//
//         do
//       {
//         gsl_monte_vegas_integrate (&G, xl, xu, 3, calls/5, r, s,
//                        &res, &err);
//         printf ("result = % .6f sigma = % .6f "
//             "chisq/dof = %.1f\n", res, err, gsl_monte_vegas_chisq (s));
//       }
//     while (fabs (gsl_monte_vegas_chisq (s) - 1.0) > 0.5);
//
//     display_results ("vegas final", res, err);
//
//     gsl_monte_vegas_free (s);
//   }
//
//   gsl_rng_free (r);
//
//   return 0;
// }
// ```
#[test]
fn plain() {
    use std::f64::consts::PI;
    fn g(k: &[f64]) -> f64 {
        let a = 1f64 / (PI * PI * PI);

        a / (1.0 - k[0].cos() * k[1].cos() * k[2].cos())
    }

    let xl: [f64; 3] = [0f64; 3];
    let xu: [f64; 3] = [PI, PI, PI];

    let calls = 500000;

    ::RngType::env_setup();
    let t: ::RngType = ::rng::default();
    let mut r = ::Rng::new(&t).unwrap();

    {
        let mut s = PlainMonteCarlo::new(3).unwrap();

        let (res, err) = s.integrate(3, g, &xl, &xu, calls, &mut r).unwrap();
        assert_eq!(&format!("{:.6}", res), "1.412209");
        assert_eq!(&format!("{:.6}", err), "0.013436");
    }
}

#[test]
fn miser() {
    use std::f64::consts::PI;
    fn g(k: &[f64]) -> f64 {
        let a = 1f64 / (PI * PI * PI);

        a / (1.0 - k[0].cos() * k[1].cos() * k[2].cos())
    }

    let xl: [f64; 3] = [0f64; 3];
    let xu: [f64; 3] = [PI, PI, PI];

    let calls = 500000;

    ::RngType::env_setup();
    let t: ::RngType = ::rng::default();
    let mut r = ::Rng::new(&t).unwrap();

    {
        let mut s = MiserMonteCarlo::new(3).unwrap();

        let (res, err) = s.integrate(3, g, &xl, &xu, calls, &mut r).unwrap();
        assert_eq!(&format!("{:.6}", res), "1.389530");
        assert_eq!(&format!("{:.6}", err), "0.005011");
    }
}

#[test]
fn miser_closure() {
    use std::f64::consts::PI;

    let xl: [f64; 3] = [0f64; 3];
    let xu: [f64; 3] = [PI, PI, PI];

    let calls = 500000;

    ::RngType::env_setup();
    let t: ::RngType = ::rng::default();
    let mut r = ::Rng::new(&t).unwrap();

    {
        let mut s = MiserMonteCarlo::new(3).unwrap();

        let (res, err) = s.integrate(3,
                       |k| {
                           let a = 1f64 / (PI * PI * PI);

                           a / (1.0 - k[0].cos() * k[1].cos() * k[2].cos())
                       },
                       &xl,
                       &xu,
                       calls,
                       &mut r)
            .unwrap();
        assert_eq!(&format!("{:.6}", res), "1.389530");
        assert_eq!(&format!("{:.6}", err), "0.005011");
    }
}

#[test]
fn vegas_warm_up() {
    use std::f64::consts::PI;
    fn g(k: &[f64]) -> f64 {
        let a = 1f64 / (PI * PI * PI);

        a / (1.0 - k[0].cos() * k[1].cos() * k[2].cos())
    }

    let xl: [f64; 3] = [0f64; 3];
    let xu: [f64; 3] = [PI, PI, PI];

    ::RngType::env_setup();
    let t: ::RngType = ::rng::default();
    let mut r = ::Rng::new(&t).unwrap();

    {
        let mut s = VegasMonteCarlo::new(3).unwrap();

        let (res, err) = s.integrate(3, g, &xl, &xu, 10000, &mut r).unwrap();
        assert_eq!(&format!("{:.6}", res), "1.385603");
        assert_eq!(&format!("{:.6}", err), "0.002212");
    }
}

#[test]
fn vegas() {
    use std::f64::consts::PI;
    fn g(k: &[f64]) -> f64 {
        let a = 1f64 / (PI * PI * PI);

        a / (1.0 - k[0].cos() * k[1].cos() * k[2].cos())
    }

    let calls = 500000;

    let xl: [f64; 3] = [0f64; 3];
    let xu: [f64; 3] = [PI, PI, PI];

    ::RngType::env_setup();
    let t: ::RngType = ::rng::default();
    let mut r = ::Rng::new(&t).unwrap();

    {
        let mut s = VegasMonteCarlo::new(3).unwrap();

        s.integrate(3, g, &xl, &xu, 10000, &mut r).unwrap();
        let mut res;
        let mut err;
        loop {
            let (_res, _err) = s.integrate(3, g, &xl, &xu, calls / 5, &mut r).unwrap();
            res = _res;
            err = _err;
            println!("result = {:.6} sigma = {:.6} chisq/dof = {:.1}",
                     res,
                     err,
                     s.chisq());
            if (s.chisq() - 1f64).abs() <= 0.5f64 {
                break;
            }
        }
        assert_eq!(&format!("{:.6}", res), "1.393307");
        assert_eq!(&format!("{:.6}", err), "0.000335");
    }
}