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_
44 #include <types/simple.h>
46 /*! \brief Abstract datatype for a random number generator
48 * This is a handle to the full state of a random number generator.
49 * You can not access anything inside the gmx_rng structure outside this
52 typedef struct gmx_rng
*
56 /*! \brief Returns the size of the RNG integer data structure
58 * Returns the size of the RNG integer data structure.
65 /*! \brief Create a new RNG, seeded from a single integer.
67 * If you dont want to pick a seed, just call it as
68 * rng=gmx_rng_init(gmx_rng_make_seed()) to seed it from
69 * the system time or a random device.
71 * \param seed Random seed, unsigned 32-bit integer.
73 * \return Reference to a random number generator, or NULL if there was an
79 gmx_rng_init(unsigned int seed
);
82 /*! \brief Generate a 'random' RNG seed.
84 * This routine tries to get a seed from /dev/random if present,
85 * and if not it uses time-of-day and process id to generate one.
87 * \return 32-bit unsigned integer random seed.
89 * Tip: If you use this in your code, it is a good idea to write the
90 * returned random seed to a logfile, so you can recreate the exact sequence
91 * of random number if you need to reproduce your run later for one reason
97 gmx_rng_make_seed(void);
100 /*! \brief Initialize a RNG with 624 integers (>32 bits of entropy).
102 * The Mersenne twister RNG used in Gromacs has an extremely long period,
103 * but when you only initialize it with a 32-bit integer there are only
104 * 2^32 different possible sequences of number - much less than the generator
107 * If you really need the full entropy, this routine makes it possible to
108 * initialize the RNG with up to 624 32-bit integers, which will give you
109 * up to 2^19968 bits of entropy.
111 * \param seed Array of unsigned integers to form a seed
112 * \param seed_length Number of integers in the array, up to 624 are used.
114 * \return Reference to a random number generator, or NULL if there was an
120 gmx_rng_init_array(unsigned int seed
[],
124 /*! \brief Release resources of a RNG
126 * This routine destroys a random number generator and releases all
127 * resources allocated by it.
129 * \param rng Handle to random number generator previously returned by
130 * gmx_rng_init() or gmx_rng_init_array().
132 * \threadsafe Function itself is threadsafe, but you should only destroy a
133 * certain RNG once (i.e. from one thread).
136 gmx_rng_destroy(gmx_rng_t rng
);
139 /*! \brief Get the state of a RNG
141 * This routine stores the random state in mt and mti, mt should have
142 * a size of at least 624, mt of 1.
144 * \param rng Handle to random number generator previously returned by
145 * gmx_rng_init() or gmx_rng_init_array().
148 gmx_rng_get_state(gmx_rng_t rng
, unsigned int *mt
,int *mti
);
151 /*! \brief Set the state of a RNG
153 * This routine sets the random state from mt and mti, mt should have
154 * a size of at least 624.
156 * \param rng Handle to random number generator previously returned by
157 * gmx_rng_init() or gmx_rng_init_array().
160 gmx_rng_set_state(gmx_rng_t rng
, unsigned int *mt
,int mti
);
163 /*! \brief Random 32-bit integer from a uniform distribution
165 * This routine returns a random integer from the random number generator
166 * provided, and updates the state of that RNG.
168 * \param rng Handle to random number generator previously returned by
169 * gmx_rng_init() or gmx_rng_init_array().
171 * \return 32-bit unsigned integer from a uniform distribution.
173 * \threadsafe Function yes, input data no. You should not call this function
174 * from two different threads using the same RNG handle at the
175 * same time. For performance reasons we cannot lock the handle
176 * with a mutex every time we need a random number - that would
177 * slow the routine down a factor 2-5. There are two simple
178 * solutions: either use a mutex and lock it before calling
179 * the function, or use a separate RNG handle for each thread.
182 gmx_rng_uniform_uint32(gmx_rng_t rng
);
185 /*! \brief Random gmx_real_t 0<=x<1 from a uniform distribution
187 * This routine returns a random floating-point number from the
188 * random number generator provided, and updates the state of that RNG.
190 * \param rng Handle to random number generator previously returned by
191 * gmx_rng_init() or gmx_rng_init_array().
193 * \return floating-point number 0<=x<1 from a uniform distribution.
195 * \threadsafe Function yes, input data no. You should not call this function
196 * from two different threads using the same RNG handle at the
197 * same time. For performance reasons we cannot lock the handle
198 * with a mutex every time we need a random number - that would
199 * slow the routine down a factor 2-5. There are two simple
200 * solutions: either use a mutex and lock it before calling
201 * the function, or use a separate RNG handle for each thread.
204 gmx_rng_uniform_real(gmx_rng_t rng
);
207 /*! \brief Random gmx_real_t from a gaussian distribution
209 * This routine returns a random floating-point number from the
210 * random number generator provided, and updates the state of that RNG.
212 * The Box-Muller algorithm is used to provide gaussian random numbers. This
213 * is not the fastest known algorithm for gaussian numbers, but in contrast
214 * to the alternatives it is very well studied and you can trust the returned
215 * random numbers to have good properties and no correlations.
217 * \param rng Handle to random number generator previously returned by
218 * gmx_rng_init() or gmx_rng_init_array().
220 * \return Gaussian random floating-point number with average 0.0 and
221 * standard deviation 1.0. You can get any average/mean you want
222 * by first multiplying with the desired average and then adding
223 * the average you want.
225 * \threadsafe Function yes, input data no. You should not call this function
226 * from two different threads using the same RNG handle at the
227 * same time. For performance reasons we cannot lock the handle
228 * with a mutex every time we need a random number - that would
229 * slow the routine down a factor 2-5. There are two simple
230 * solutions: either use a mutex and lock it before calling
231 * the function, or use a separate RNG handle for each thread.
233 * It works perfectly to mix calls to get uniform and gaussian random numbers
234 * from the same generator, but since it will affect the sequence of returned
235 * numbers it is probably better to use separate random number generator
239 gmx_rng_gaussian_real(gmx_rng_t rng
);
243 /* Return a new gaussian random number with expectation value
244 * 0.0 and standard deviation 1.0. This routine uses a table
245 * lookup for maximum speed.
247 * WARNING: The lookup table is 16k by default, which means
248 * the granularity of the random numbers is coarser
249 * than what you get from gmx_rng_gauss_real().
250 * In most cases this is no problem whatsoever,
251 * and it is particularly true for BD/SD integration.
252 * Note that you will NEVER get any really extreme
253 * numbers: the maximum absolute value returned is
259 gmx_rng_gaussian_table(gmx_rng_t rng
);
262 #endif /* _GMX_RANDOM_H_ */