3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
9 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
10 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
11 * Copyright (c) 2001-2008, The GROMACS development team,
12 * check out http://www.gromacs.org for more information.
14 * This program is free software; you can redistribute it and/or
15 * modify it under the terms of the GNU General Public License
16 * as published by the Free Software Foundation; either version 2
17 * of the License, or (at your option) any later version.
19 * If you want to redistribute modifications, please consider that
20 * scientific software is very special. Version control is crucial -
21 * bugs must be traceable. We will be happy to consider code for
22 * inclusion in the official distribution, but derived work must not
23 * be called official GROMACS. Details are found in the README & COPYING
24 * files - if they are missing, get the official version at www.gromacs.org.
26 * To help us fund GROMACS development, we humbly ask that you cite
27 * the papers on the package - you can find them in the top README file.
29 * For more info, check our website at http://www.gromacs.org
32 * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
35 #ifndef _GMX_RANDOM_H_
36 #define _GMX_RANDOM_H_
39 #include "types/simple.h"
45 /*! \brief Abstract datatype for a random number generator
47 * This is a handle to the full state of a random number generator.
48 * You can not access anything inside the gmx_rng structure outside this
51 typedef struct gmx_rng
*
55 /*! \brief Returns the size of the RNG integer data structure
57 * Returns the size of the RNG integer data structure.
64 /*! \brief Create a new RNG, seeded from a single integer.
66 * If you dont want to pick a seed, just call it as
67 * rng=gmx_rng_init(gmx_rng_make_seed()) to seed it from
68 * the system time or a random device.
70 * \param seed Random seed, unsigned 32-bit integer.
72 * \return Reference to a random number generator, or NULL if there was an
78 gmx_rng_init(unsigned int seed
);
81 /*! \brief Generate a 'random' RNG seed.
83 * This routine tries to get a seed from /dev/random if present,
84 * and if not it uses time-of-day and process id to generate one.
86 * \return 32-bit unsigned integer random seed.
88 * Tip: If you use this in your code, it is a good idea to write the
89 * returned random seed to a logfile, so you can recreate the exact sequence
90 * of random number if you need to reproduce your run later for one reason
96 gmx_rng_make_seed(void);
99 /*! \brief Initialize a RNG with 624 integers (>32 bits of entropy).
101 * The Mersenne twister RNG used in Gromacs has an extremely long period,
102 * but when you only initialize it with a 32-bit integer there are only
103 * 2^32 different possible sequences of number - much less than the generator
106 * If you really need the full entropy, this routine makes it possible to
107 * initialize the RNG with up to 624 32-bit integers, which will give you
108 * up to 2^19968 bits of entropy.
110 * \param seed Array of unsigned integers to form a seed
111 * \param seed_length Number of integers in the array, up to 624 are used.
113 * \return Reference to a random number generator, or NULL if there was an
119 gmx_rng_init_array(unsigned int seed
[],
123 /*! \brief Release resources of a RNG
125 * This routine destroys a random number generator and releases all
126 * resources allocated by it.
128 * \param rng Handle to random number generator previously returned by
129 * gmx_rng_init() or gmx_rng_init_array().
131 * \threadsafe Function itself is threadsafe, but you should only destroy a
132 * certain RNG once (i.e. from one thread).
135 gmx_rng_destroy(gmx_rng_t rng
);
138 /*! \brief Get the state of a RNG
140 * This routine stores the random state in mt and mti, mt should have
141 * a size of at least 624, mt of 1.
143 * \param rng Handle to random number generator previously returned by
144 * gmx_rng_init() or gmx_rng_init_array().
147 gmx_rng_get_state(gmx_rng_t rng
, unsigned int *mt
,int *mti
);
150 /*! \brief Set the state of a RNG
152 * This routine sets the random state from mt and mti, mt should have
153 * a size of at least 624.
155 * \param rng Handle to random number generator previously returned by
156 * gmx_rng_init() or gmx_rng_init_array().
159 gmx_rng_set_state(gmx_rng_t rng
, unsigned int *mt
,int mti
);
162 /*! \brief Random 32-bit integer from a uniform distribution
164 * This routine returns a random integer from the random number generator
165 * provided, and updates the state of that RNG.
167 * \param rng Handle to random number generator previously returned by
168 * gmx_rng_init() or gmx_rng_init_array().
170 * \return 32-bit unsigned integer from a uniform distribution.
172 * \threadsafe Function yes, input data no. You should not call this function
173 * from two different threads using the same RNG handle at the
174 * same time. For performance reasons we cannot lock the handle
175 * with a mutex every time we need a random number - that would
176 * slow the routine down a factor 2-5. There are two simple
177 * solutions: either use a mutex and lock it before calling
178 * the function, or use a separate RNG handle for each thread.
181 gmx_rng_uniform_uint32(gmx_rng_t rng
);
184 /*! \brief Random gmx_real_t 0<=x<1 from a uniform distribution
186 * This routine returns a random floating-point number from the
187 * random number generator provided, and updates the state of that RNG.
189 * \param rng Handle to random number generator previously returned by
190 * gmx_rng_init() or gmx_rng_init_array().
192 * \return floating-point number 0<=x<1 from a uniform distribution.
194 * \threadsafe Function yes, input data no. You should not call this function
195 * from two different threads using the same RNG handle at the
196 * same time. For performance reasons we cannot lock the handle
197 * with a mutex every time we need a random number - that would
198 * slow the routine down a factor 2-5. There are two simple
199 * solutions: either use a mutex and lock it before calling
200 * the function, or use a separate RNG handle for each thread.
203 gmx_rng_uniform_real(gmx_rng_t rng
);
206 /*! \brief Random gmx_real_t from a gaussian distribution
208 * This routine returns a random floating-point number from the
209 * random number generator provided, and updates the state of that RNG.
211 * The Box-Muller algorithm is used to provide gaussian random numbers. This
212 * is not the fastest known algorithm for gaussian numbers, but in contrast
213 * to the alternatives it is very well studied and you can trust the returned
214 * random numbers to have good properties and no correlations.
216 * \param rng Handle to random number generator previously returned by
217 * gmx_rng_init() or gmx_rng_init_array().
219 * \return Gaussian random floating-point number with average 0.0 and
220 * standard deviation 1.0. You can get any average/mean you want
221 * by first multiplying with the desired average and then adding
222 * the average you want.
224 * \threadsafe Function yes, input data no. You should not call this function
225 * from two different threads using the same RNG handle at the
226 * same time. For performance reasons we cannot lock the handle
227 * with a mutex every time we need a random number - that would
228 * slow the routine down a factor 2-5. There are two simple
229 * solutions: either use a mutex and lock it before calling
230 * the function, or use a separate RNG handle for each thread.
232 * It works perfectly to mix calls to get uniform and gaussian random numbers
233 * from the same generator, but since it will affect the sequence of returned
234 * numbers it is probably better to use separate random number generator
238 gmx_rng_gaussian_real(gmx_rng_t rng
);
242 /* Return a new gaussian random number with expectation value
243 * 0.0 and standard deviation 1.0. This routine uses a table
244 * lookup for maximum speed.
246 * WARNING: The lookup table is 16k by default, which means
247 * the granularity of the random numbers is coarser
248 * than what you get from gmx_rng_gauss_real().
249 * In most cases this is no problem whatsoever,
250 * and it is particularly true for BD/SD integration.
251 * Note that you will NEVER get any really extreme
252 * numbers: the maximum absolute value returned is
258 gmx_rng_gaussian_table(gmx_rng_t rng
);
264 #endif /* _GMX_RANDOM_H_ */