rtnorm_rust/lib.rs
1// Pseudorandom numbers from a truncated Gaussian distribution.
2//
3// This is a port of c++ code.
4//
5// This implements an extension of Chopin's algorithm detailed in
6// N. Chopin, "Fast simulation of truncated Gaussian distributions",
7// Stat Comput (2011) 21:275-288
8//
9// Original copyright holder:
10// Copyright (C) 2012 Guillaume Dollé, Vincent Mazet
11// (LSIIT, CNRS/Université de Strasbourg)
12// Version 2012-07-04, Contact: vincent.mazet@unistra.fr
13//
14// Port copyright holder:
15// Copyright (C) 2021 Nozomu Shimaoka
16//
17// 07/10/2021:
18// - Ported from c++ version
19// 06/07/2012:
20// - first launch of rtnorm.cpp
21//
22// Licence: GNU General Public License Version 2
23// This program is free software; you can redistribute it and/or modify it
24// under the terms of the GNU General Public License as published by the
25// Free Software Foundation; either version 2 of the License, or (at your
26// option) any later version. This program is distributed in the hope that
27// it will be useful, but WITHOUT ANY WARRANTY; without even the implied
28// warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
29// GNU General Public License for more details. You should have received a
30// copy of the GNU General Public License along with this program; if not,
31// see http://www.gnu.org/licenses/old-licenses/gpl-2.0.txt
32//
33// Depends: LibGSL
34// OS: Unix based system
35mod table;
36
37use num::ToPrimitive;
38use rgsl::{self, error::erf, Rng};
39
40use crate::table::{NCELL, X, YU};
41
42/// Index of the right tail
43const N: usize = 4001;
44
45#[cfg(test)]
46mod test {
47 use rgsl::{Rng, RngType};
48
49 use crate::rtnorm;
50
51 #[test]
52 fn generate() -> Result<(), ()> {
53 let a = 1.0; // Left bound
54 let b = 9.0; // Right bound
55 let mu = 2.0; // Mean
56 let sigma = 3.0; // Standard deviation
57 const K: i64 = 100000; // Number of random variables to generate
58
59 //--- GSL random init ---
60 // Read variable environnement
61 RngType::env_setup();
62 // Rand generator allocation
63 // Default algorithm 'twister'
64 let mut gen =
65 Rng::new(RngType::default()).expect("could not initialize a random number generator");
66
67 //--- generate the random numbers ---
68 // println!("# x p(x)");
69 // println!("{} {}", a, b);
70 // println!("{} {}", mu, sigma);
71 for _ in 0..K {
72 let (x, _p) = rtnorm(&mut gen, a, b, mu, sigma);
73 assert!(a <= x && x <= b);
74 // println!("{} {}", s.0, s.1);
75 }
76
77 Ok(())
78 // GSL rand generator will be deallocated
79 }
80}
81
82/// Compute y_l from y_k
83#[inline]
84const fn yl(k: usize) -> f64 {
85 // y_l of the leftmost rectangle
86 const YL0: f64 = 0.053513975472;
87 // y_l of the rightmost rectangle
88 const YLN: f64 = 0.000914116389555;
89
90 if k == 0 {
91 YL0
92 } else if k == N - 1 {
93 YLN
94 } else if k <= 1953 {
95 // TODO: always k > 0 ?
96 YU[k - 1]
97 } else {
98 YU[k + 1]
99 }
100}
101
102/// Rejection algorithm with a truncated exponential proposal
103#[inline]
104fn rtexp(gen: &mut Rng, a: f64, b: f64) -> f64 {
105 let twoasq = 2.0 * a.powi(2);
106 let expab = (-a * (b - a)).exp() - 1.0;
107 let mut z;
108 let mut e;
109
110 loop {
111 z = (1.0 + gen.uniform() * expab).ln();
112 e = -(gen.uniform()).ln();
113 if twoasq * e > z.powi(2) {
114 break;
115 }
116 }
117 return a - z / a;
118}
119
120/// Pseudorandom numbers from a truncated Gaussian distribution.
121/// The Gaussian has parameters mu and sigma
122/// and is truncated on the interval \[a,b\].
123/// Returns the random variable x and its probability p(x).
124#[inline]
125pub fn rtnorm(gen: &mut Rng, mut a: f64, mut b: f64, mu: f64, sigma: f64) -> (f64, f64) {
126 // Design variables
127 const XMIN: f64 = -2.00443204036; // Left bound
128 const XMAX: f64 = 3.48672170399; // Right bound
129 const KMIN: i64 = 5; // if kb-ka < kmin then use a rejection algorithm
130 const INVH: f64 = 1631.73284006; // = 1/h, h being the minimal interval range
131 const I0: i64 = 3271; // = - floor(x(0)/h)
132 const ALPHA: f64 = 1.837877066409345; // = log(2*pi)
133 const SQ2: f64 = 7.071067811865475e-1; // = 1/sqrt(2)
134 const SQPI: f64 = 1.772453850905516; // = sqrt(pi)
135
136 // Scaling
137 if mu != 0.0 || sigma != 1.0 {
138 a = (a - mu) / sigma;
139 b = (b - mu) / sigma;
140 }
141
142 // Check if a < b
143 assert!(a < b, "B must be greater than A");
144 // Check if |a| < |b|
145 let r = if a.abs() > b.abs() {
146 // Tuple (r,p)
147 -rtnorm(gen, -b, -a, 0.0, 1.0).0
148 } else if a > XMAX {
149 // If a in the right tail (a > xmax), use rejection algorithm with a truncated exponential proposal
150 rtexp(gen, a, b)
151 } else if a < XMIN {
152 // If a in the left tail (a < xmin), use rejection algorithm with a Gaussian proposal
153 let mut r;
154 loop {
155 r = gen.gaussian(1.0);
156 if (r >= a) && (r <= b) {
157 break;
158 }
159 }
160 r
161 } else {
162 // In other cases (xmin < a < xmax), use Chopin's algorithm
163 // Compute ka
164 let ka = NCELL[(I0 + (a * INVH).floor() as i64).to_usize().unwrap()];
165
166 // Compute kb
167 let kb = if b >= XMAX {
168 N as i64
169 } else {
170 NCELL[(I0 + (b * INVH).floor() as i64).to_usize().unwrap()]
171 };
172
173 // If |b-a| is small, use rejection algorithm with a truncated exponential proposal
174 if (kb - ka).abs() < KMIN {
175 rtexp(gen, a, b)
176 } else {
177 loop {
178 // Sample integer between ka and kb
179 let k = (gen.uniform() * (kb - ka + 1) as f64).floor() as i64 + ka;
180 let k = k.to_usize().unwrap();
181
182 if k == N {
183 // Right tail
184 let lbound = X[X.len() - 1];
185 let z = -gen.uniform().ln() / lbound;
186 let e = -gen.uniform().ln();
187
188 if (z.powi(2) <= 2.0 * e) && (z < b - lbound) {
189 // Accept this proposition, otherwise reject
190 break lbound + z;
191 }
192 } else if (k as i64 <= ka + 1) || (k as i64 >= kb - 1 && b < XMAX) {
193 // Two leftmost and rightmost regions
194 let sim = X[k] + (X[k + 1] - X[k]) * gen.uniform();
195
196 if (sim >= a) && (sim <= b) {
197 // Accept this proposition, otherwise reject
198 let simy = YU[k] * gen.uniform();
199 if (simy < yl(k)) || (sim * sim + 2.0 * simy.ln() + ALPHA) < 0.0 {
200 break sim;
201 }
202 }
203 } else {
204 // All the other boxes
205
206 let u = gen.uniform();
207 let simy = YU[k] * u;
208 let d = X[k + 1] - X[k];
209 if simy < yl(k)
210 // That's what happens most of the time
211 {
212 break X[k] + u * d * YU[k] / yl(k);
213 } else {
214 let sim = X[k] + d * gen.uniform();
215
216 // Otherwise, check you're below the pdf curve
217 if (sim * sim + 2.0 * simy.ln() + ALPHA) < 0.0 {
218 break sim;
219 }
220 }
221 }
222 }
223 }
224 };
225
226 let r = if mu != 0.0 || sigma != 1.0 {
227 // Scaling
228 r * sigma + mu
229 } else {
230 r
231 };
232
233 // Compute the probability
234 let large_z = SQPI * SQ2 * sigma * (erf(b * SQ2) - erf(a * SQ2));
235 let p = (-((r - mu) / sigma).powi(2) / 2.0).exp() / large_z;
236
237 return (r, p);
238}