rgsl/types/qrng.rs
1//
2// A rust binding for the GSL library by Guillaume Gomez (guillaume1.gomez@gmail.com)
3//
4
5/*!
6# Quasi-Random Sequences
7
8This chapter describes functions for generating quasi-random sequences in arbitrary dimensions. A quasi-random sequence progressively
9covers a d-dimensional space with a set of points that are uniformly distributed. Quasi-random sequences are also known as low-discrepancy
10sequences. The quasi-random sequence generators use an interface that is similar to the interface for random number generators, except
11that seeding is not required—each generator produces a single sequence.
12
13## References
14
15The implementations of the quasi-random sequence routines are based on the algorithms described in the following paper,
16
17P. Bratley and B.L. Fox and H. Niederreiter, “Algorithm 738: Programs to Generate Niederreiter’s Low-discrepancy Sequences”, ACM
18Transactions on Mathematical Software, Vol. 20, No. 4, December, 1994, p. 494–495.
19!*/
20
21use crate::Value;
22use ffi::FFI;
23
24ffi_wrapper!(QRng, *mut sys::gsl_qrng, gsl_qrng_free);
25
26impl QRng {
27 /// This function returns a pointer to a newly-created instance of a quasi-random sequence
28 /// generator of type T and dimension d. If there is insufficient memory to create the generator
29 /// then the function returns a null pointer and the error handler is invoked with an error code
30 /// of [`Value::NoMemory`].
31 #[doc(alias = "gsl_qrng_alloc")]
32 pub fn new(t: QRngType, d: u32) -> Option<Self> {
33 let tmp = unsafe { sys::gsl_qrng_alloc(t.unwrap_shared(), d) };
34
35 if tmp.is_null() {
36 None
37 } else {
38 Some(Self::wrap(tmp))
39 }
40 }
41
42 /// This function reinitializes the generator self to its starting point. Note that quasi-random
43 /// sequences do not use a seed and always produce the same set of values.
44 #[doc(alias = "gsl_qrng_init")]
45 pub fn init(&mut self) {
46 unsafe { sys::gsl_qrng_init(self.unwrap_unique()) }
47 }
48
49 /// This function stores the next point from the sequence generator self in the array x. The
50 /// space available for x must match the dimension of the generator. The point x will lie in the
51 /// range 0 < x_i < 1 for each x_i.
52 #[doc(alias = "gsl_qrng_get")]
53 pub fn get(&self, x: &mut [f64]) -> Result<(), Value> {
54 let ret = unsafe { sys::gsl_qrng_get(self.unwrap_shared(), x.as_mut_ptr()) };
55 result_handler!(ret, ())
56 }
57
58 /// This function returns a pointer to the name of the generator.
59 #[doc(alias = "gsl_qrng_name")]
60 pub fn name(&self) -> Option<String> {
61 let tmp = unsafe { sys::gsl_qrng_name(self.unwrap_shared()) };
62
63 if tmp.is_null() {
64 None
65 } else {
66 unsafe {
67 Some(
68 String::from_utf8_lossy(::std::ffi::CStr::from_ptr(tmp).to_bytes()).to_string(),
69 )
70 }
71 }
72 }
73
74 /// These functions return a pointer to the state of generator r and its size.
75 #[doc(alias = "gsl_qrng_size")]
76 pub fn size(&self) -> usize {
77 unsafe { sys::gsl_qrng_size(self.unwrap_shared()) }
78 }
79
80 /// This function returns a pointer to the state of generator `self`.
81 #[doc(alias = "gsl_qrng_state")]
82 pub fn state(&mut self) -> Option<&[i8]> {
83 let tmp = unsafe { sys::gsl_qrng_state(self.unwrap_shared()) };
84
85 if tmp.is_null() {
86 None
87 } else {
88 Some(unsafe { std::slice::from_raw_parts(tmp as _, self.size()) })
89 }
90 }
91
92 /// This function returns a pointer to the state of generator `self`.
93 // checker:ignore
94 #[doc(alias = "gsl_qrng_state")]
95 pub fn state_mut(&mut self) -> Option<&mut [i8]> {
96 let tmp = unsafe { sys::gsl_qrng_state(self.unwrap_shared()) };
97
98 if tmp.is_null() {
99 None
100 } else {
101 Some(unsafe { std::slice::from_raw_parts_mut(tmp as _, self.size()) })
102 }
103 }
104
105 /// This function copies the quasi-random sequence generator src into the pre-existing generator
106 /// `dest`, making dest into an exact copy of `self`. The two generators must be of the same
107 /// type.
108 // checker:ignore
109 #[doc(alias = "gsl_qrng_memcpy")]
110 pub fn copy(&self, dest: &mut QRng) -> Result<(), Value> {
111 let ret = unsafe { sys::gsl_qrng_memcpy(dest.unwrap_unique(), self.unwrap_shared()) };
112 result_handler!(ret, ())
113 }
114}
115
116impl Clone for QRng {
117 /// This function returns a pointer to a newly created generator which is an exact copy of the
118 /// generator `self`.
119 #[doc(alias = "gsl_qrng_clone")]
120 fn clone(&self) -> Self {
121 unsafe { Self::wrap(sys::gsl_qrng_clone(self.unwrap_shared())) }
122 }
123}
124
125ffi_wrapper!(QRngType, *const sys::gsl_qrng_type);
126
127impl QRngType {
128 /// This generator uses the algorithm described in Bratley, Fox, Niederreiter, ACM Trans. Model.
129 /// Comp. Sim. 2, 195 (1992). It is valid up to 12 dimensions.
130 #[doc(alias = "gsl_qrng_niederreiter_2")]
131 pub fn niederreiter_2() -> QRngType {
132 ffi_wrap!(gsl_qrng_niederreiter_2)
133 }
134
135 /// This generator uses the Sobol sequence described in Antonov, Saleev, USSR Comput. Maths.
136 /// Math. Phys. 19, 252 (1980). It is valid up to 40 dimensions.
137 #[doc(alias = "gsl_qrng_sobol")]
138 pub fn sobol() -> QRngType {
139 ffi_wrap!(gsl_qrng_sobol)
140 }
141
142 /// These generators use the Halton and reverse Halton sequences described in J.H. Halton,
143 /// Numerische Mathematik 2, 84-90 (1960) and B. Vandewoestyne and R. Cools Computational and
144 /// Applied Mathematics 189, 1&2, 341-361 (2006). They are valid up to 1229 dimensions.
145 #[doc(alias = "gsl_qrng_halton")]
146 pub fn halton() -> QRngType {
147 ffi_wrap!(gsl_qrng_halton)
148 }
149
150 #[doc(alias = "gsl_qrng_reversehalton")]
151 pub fn reversehalton() -> QRngType {
152 ffi_wrap!(gsl_qrng_reversehalton)
153 }
154}