ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
mixmax.cxx
Go to the documentation of this file.
1 // $Id:$
2 // -*- C++ -*-
3 //
4 // -----------------------------------------------------------------------
5 // MixMax Matrix PseudoRandom Number Generator
6 // --- MixMax ---
7 // class header file
8 // -----------------------------------------------------------------------
9 //
10 //
11 // Created by Konstantin Savvidy on Sun Feb 22 2004.
12 // As of version 0.99 and later, the code is being released under
13 // GNU Lesser General Public License v3
14 //
15 // Generator described in
16 // N.Z.Akopov, G.K.Savvidy and N.G.Ter-Arutyunian, Matrix Generator of Pseudorandom Numbers,
17 // J.Comput.Phys. 97, 573 (1991);
18 // Preprint EPI-867(18)-86, Yerevan Jun.1986;
19 //
20 // and
21 //
22 // K.Savvidy
23 // The MIXMAX random number generator
24 // Comp. Phys. Commun. (2015)
25 // http://dx.doi.org/10.1016/j.cpc.2015.06.003
26 //
27 // -----------------------------------------------------------------------
28 
29 #include <stdio.h>
30 #include <stdlib.h>
31 #include <stdint.h>
32 
33 #define __MIXMAX_C // do NOT define it in your own program, just include mixmax.h
34 
35 //#define USE_INLINE_ASM //LM: uncomment if want to use inline asm
36 
37 #include "mixmax.h"
38 
39 
41  X->sumtot = iterate_raw_vec(X->V, X->sumtot);
42  return 0;
43 }
44 
45 #if (SPECIALMUL!=0)
46 inline uint64_t MULWU (uint64_t k){ return (( (k)<<(SPECIALMUL) & M61) | ( (k) >> (BITS-SPECIALMUL)) ) ;}
47 #elif (SPECIALMUL==0)
48 inline uint64_t MULWU (uint64_t k){ (void)k; return 0;}
49 #else
50 #error SPECIALMUL not undefined
51 #endif
52 
54  // operates with a raw vector, uses known sum of elements of Y
55  int i;
56 #ifdef SPECIAL
57  myuint temp2 = Y[1];
58 #endif
59  myuint tempP, tempV;
60  Y[0] = ( tempV = sumtotOld);
61  myuint sumtot = Y[0], ovflow = 0; // will keep a running sum of all new elements (except Y[0])
62  tempP = 0; // will keep a partial sum of all old elements (except Y[0])
63  for (i=1; i<N; i++){
64 #if (SPECIALMUL!=0)
65  myuint tempPO = MULWU(tempP);
66  tempP = modadd(tempP,Y[i]);
67  tempV = MOD_MERSENNE(tempV + tempP + tempPO); // edge cases ?
68 #else
69  tempP = modadd(tempP , Y[i]);
70  tempV = modadd(tempV , tempP);
71 #endif
72  Y[i] = tempV;
73  sumtot += tempV; if (sumtot < tempV) {ovflow++;}
74  }
75 #ifdef SPECIAL
76  temp2 = MOD_MULSPEC(temp2);
77  Y[2] = modadd( Y[2] , temp2 );
78  sumtot += temp2; if (sumtot < temp2) {ovflow++;}
79 #endif
80  return MOD_MERSENNE(MOD_MERSENNE(sumtot) + (ovflow <<3 ));
81 }
82 
84  int i;
85  i=X->counter;
86 
87  if (i<(N) ){
88  X->counter++;
89  return X->V[i];
90  }else{
91  X->sumtot = iterate_raw_vec(X->V, X->sumtot);
92  X->counter=1;
93  return X->V[0];
94  }
95 }
96 
98  /* cast to signed int trick suggested by Andrzej Gòˆrlich */
99  // times for sknuth_Gap test: N = 1, n = 500000000, r = 0, Alpha = 0, Beta = 0.0625
100  //return ((int64_t)get_next(X)) * INV_MERSBASE; // 00:01:17.23
101  //return ( (double)get_next(X)) * INV_MERSBASE; // 00:01:57.89
102  int64_t Z=(int64_t)get_next(X);
103 #if defined(__SSE__) && defined(USE_INLINE_ASM)
104  double F;
105  __asm__ ("pxor %0, %0; "
106  "cvtsi2sdq %1, %0; "
107  :"=x"(F)
108  :"r"(Z)
109  );
110  return F*INV_MERSBASE;
111 #else
112  return Z*INV_MERSBASE;
113 #endif
114 
115 }
116 
117 void fill_array(rng_state_t* X, unsigned int n, double *array)
118 {
119  // Return an array of n random numbers uniformly distributed in (0,1]
120  unsigned int i,j;
121  const int M=N-1;
122  for (i=0; i<(n/M); i++){
123  iterate_and_fill_array(X, array+i*M);
124  }
125  unsigned int rem=(n % M);
126  if (rem) {
127  iterate(X);
128  for (j=0; j< (rem); j++){
129  array[M*i+j] = (int64_t)X->V[j] * (double)(INV_MERSBASE);
130  }
131  X->counter = j; // needed to continue with single fetches from the exact spot, but if you only use fill_array to get numbers then it is not necessary
132  }else{
133  X->counter = N;
134  }
135 }
136 
137 void iterate_and_fill_array(rng_state_t* X, double *array){
138  myuint* Y=X->V;
139  int i;
140  myuint tempP, tempV;
141 #if (SPECIAL != 0)
142  myuint temp2 = Y[1];
143 #endif
144  Y[0] = (tempV = modadd(Y[0] , X->sumtot));
145  //array[0] = (double)tempV * (double)(INV_MERSBASE);
146  myuint sumtot = 0, ovflow = 0; // will keep a running sum of all new elements (except Y[0])
147  tempP = 0; // will keep a partial sum of all old elements (except Y[0])
148  for (i=1; i<N; i++){
149  tempP = modadd(tempP,Y[i]);
150  Y[i] = ( tempV = modadd(tempV,tempP) );
151  sumtot += tempV; if (sumtot < tempV) {ovflow++;}
152  array[i-1] = (int64_t)tempV * (double)(INV_MERSBASE);
153  }
154 #if (SPECIAL != 0)
155  temp2 = MOD_MULSPEC(temp2);
156  Y[2] = modadd( Y[2] , temp2 );
157  sumtot += temp2; if (sumtot < temp2) {ovflow++;}
158 #endif
159  X->sumtot = MOD_MERSENNE(MOD_MERSENNE(sumtot) + (ovflow <<3 ));
160 }
161 
163 #if defined(__x86_64__) && defined(USE_INLINE_ASM)
164  myuint out;
165  /* Assembler trick suggested by Andrzej Gòˆrlich */
166  __asm__ ("addq %2, %0; "
167  "btrq $61, %0; "
168  "adcq $0, %0; "
169  :"=r"(out)
170  :"0"(foo), "r"(bar)
171  );
172  return out;
173 #else
174  return MOD_MERSENNE(foo+bar);
175 #endif
176 }
177 
179 {
180 /* allocate the state */
181  rng_state_t *p = (rng_state_t*)malloc(sizeof(rng_state_t));
182  p->fh=NULL; // by default, set the output file handle to stdout
183  return p;
184 }
185 
186 int rng_free(rng_state_t* X) /* free the memory occupied by the state */
187 {
188  free(X);
189  return 0;
190 }
191 
193 {
194  /* copy the vector stored at Y, and return pointer to the newly allocated and initialized state.
195  It is the user's responsibility to make sure that Y is properly allocated with rng_alloc,
196  then pass Y->V or it can also be an array -- such as myuint Y[N+1] and Y[1]...Y[N] have been set to legal values [0 .. MERSBASE-1]
197  Partial sums on this new state are recalculated, and counter set to zero, so that when get_next is called,
198  it will output the initial vector before any new numbers are produced, call iterate(X) if you want to advance right away */
199  rng_state_t* X = rng_alloc();
200  myuint sumtot=0,ovflow=0;
201  X->counter = 2;
202  int i;
203  for ( i=0; i < N; i++){
204  X->V[i] = Y[i];
205  sumtot += X->V[(i)]; if (sumtot < X->V[(i)]) {ovflow++;}
206 
207  }
208  X->sumtot = MOD_MERSENNE(MOD_MERSENNE(sumtot) + (ovflow <<3 ));
209  return X;
210 }
211 
212 void seed_vielbein(rng_state_t* X, unsigned int index)
213 {
214 int i;
215  if (index<N){
216  for (i=0; i < N; i++){
217  X->V[i] = 0;
218  }
219  X->V[index] = 1;
220  }else{
221  fprintf(stderr, "Out of bounds index, is not ( 0 <= index < N )\n"); exit(ARRAY_INDEX_OUT_OF_BOUNDS);
222  }
223  X->counter = N; // set the counter to N if iteration should happen right away
224  //precalc(X);
225  X->sumtot = 1; //(index ? 1:0);
226  if (X->fh==NULL){X->fh=stdout;}
227 }
228 
230 { // a 64-bit LCG from Knuth line 26, in combination with a bit swap is used to seed
231  const myuint MULT64=6364136223846793005ULL;
232  int i;
233  myuint sumtot=0,ovflow=0;
234  if (seed == 0){
235  fprintf(stderr, " try seeding with nonzero seed next time!\n");
236  exit(SEED_WAS_ZERO);
237  }
238 
239  myuint l = seed;
240 
241  //X->V[0] = l & MERSBASE;
242  if (X->fh==NULL){X->fh=stdout;} // if the filehandle is not yet set, make it stdout
243  for (i=0; i < N; i++){
244  l*=MULT64; l = (l << 32) ^ (l>>32);
245  X->V[i] = l & MERSBASE;
246  sumtot += X->V[(i)]; if (sumtot < X->V[(i)]) {ovflow++;}
247  }
248  X->counter = N; // set the counter to N if iteration should happen right away
249  X->sumtot = MOD_MERSENNE(MOD_MERSENNE(sumtot) + (ovflow <<3 ));
250 }
251 
253  int i;
254  myuint temp;
255  temp = 0;
256  for (i=0; i < N; i++){
257  temp = MOD_MERSENNE(temp + X->V[i]);
258  }
259  X->sumtot = temp;
260  return temp;
261 }
262 
263 
264 int rng_get_N(void){return N;}
265 
266 //#define MASK32 0xFFFFFFFFULL
267 
268 
269 //inline myuint modmulM61(myuint a, myuint b){
270 // // my best modmul so far
271 // __uint128_t temp;
272 // temp = (__uint128_t)a*(__uint128_t)b;
273 // return mod128(temp);
274 //}
275 
276 #if defined(__x86_64__)
277 inline myuint mod128(__uint128_t s){
278  myuint s1;
279  s1 = ( ( ((myuint)s)&MERSBASE ) + ( ((myuint)(s>>64)) * 8 ) + ( ((myuint)s) >>BITS) );
280  return MOD_MERSENNE(s1);
281 }
282 
283 inline myuint fmodmulM61(myuint cum, myuint a, myuint b){
284  __uint128_t temp;
285  temp = (__uint128_t)a*(__uint128_t)b + cum;
286  return mod128(temp);
287 }
288 
289 #else // on all other platforms, including 32-bit linux, PPC and PPC64 and all Windows
290 #define MASK32 0xFFFFFFFFULL
291 
293 {
294  register myuint o,ph,pl,ah,al;
295  o=(s)*a;
296  ph = ((s)>>32);
297  pl = (s) & MASK32;
298  ah = a>>32;
299  al = a & MASK32;
300  o = (o & M61) + ((ph*ah)<<3) + ((ah*pl+al*ph + ((al*pl)>>32))>>29) ;
301  o += cum;
302  o = (o & M61) + ((o>>61));
303  return o;
304 }
305 #endif
306 
308  int j;
309  fprintf(X->fh, "mixmax state, file version 1.0\n" );
310  fprintf(X->fh, "N=%u; V[N]={", rng_get_N() );
311  for (j=0; (j< (rng_get_N()-1) ); j++) {
312  fprintf(X->fh, "%llu, ", X->V[j] );
313  }
314  fprintf(X->fh, "%llu", X->V[rng_get_N()-1] );
315  fprintf(X->fh, "}; " );
316  fprintf(X->fh, "counter=%u; ", X->counter );
317  fprintf(X->fh, "sumtot=%llu;\n", X->sumtot );
318 }
319 
320 void read_state(rng_state_t* X, const char filename[] ){
321  // a function for reading the state from a file, after J. Apostolakis
322  FILE* fin;
323  if( ( fin = fopen(filename, "r") ) ){
324  char l=0;
325  while ( l != '{' ) { // 0x7B = "{"
326  l=fgetc(fin); // proceed until hitting opening bracket
327  }
328  ungetc(' ', fin);
329  }else{
330  fprintf(stderr, "mixmax -> read_state: error reading file %s\n", filename);
332  }
333 
334  myuint vecVal;
335  //printf("mixmax -> read_state: starting to read state from file\n");
336  if (!fscanf(fin, "%llu", &X->V[0]) ) {fprintf(stderr, "mixmax -> read_state: error reading file %s\n", filename); exit(ERROR_READING_STATE_FILE);}
337  //printf("V[%d] = %llu\n",0, X->V[0]);
338  int i;
339  for( i = 1; i < rng_get_N(); i++){
340  if (!fscanf(fin, ", %llu", &vecVal) ) {fprintf(stderr, "mixmax -> read_state: error reading vector component i=%d from file %s\n", i, filename); exit(ERROR_READING_STATE_FILE);}
341  //printf("V[%d] = %llu\n",i, vecVal);
342  if( vecVal <= MERSBASE ){
343  X->V[i] = vecVal;
344  }else{
345  fprintf(stderr, "mixmax -> read_state: Invalid state vector value= %llu"
346  " ( must be less than %llu ) "
347  " obtained from reading file %s\n"
348  , vecVal, MERSBASE, filename);
349 
350  }
351  }
352 
353  unsigned int counter;
354  if (!fscanf( fin, "}; counter=%u; ", &counter)){fprintf(stderr, "mixmax -> read_state: error reading counter from file %s\n", filename); exit(ERROR_READING_STATE_FILE);}
355  if( counter <= N ) {
356  X->counter= counter;
357  }else{
358  fprintf(stderr, "mixmax -> read_state: Invalid counter = %d"
359  " Must be 0 <= counter < %u\n" , counter, N);
360  print_state(X);
362  }
363  precalc(X);
364  myuint sumtot;
365  if (!fscanf( fin, "sumtot=%llu\n", &sumtot)){fprintf(stderr, "mixmax -> read_state: error reading checksum from file %s\n", filename); exit(ERROR_READING_STATE_FILE);}
366 
367  if (X->sumtot != sumtot) {
368  fprintf(stderr, "mixmax -> checksum error while reading state from file %s - corrupted?\n", filename);
370  }
371 // else{fprintf(stderr, "mixmax -> read_state: checksum ok: %llu == %llu\n",X->sumtot, sumtot);}
372  fclose(fin);
373 }
374 
375 
376 #define FUSEDMODMULVEC \
377 { for (i =0; i<N; i++){ \
378 cum[i] = fmodmulM61( cum[i], coeff , Y[i] ) ; \
379 } }
380 
381 
382 #define SKIPISON 1
383 
384 #if ( ( (N==88)||(N==256) ||(N==1000) || (N==3150) ) && BITS==61 && SKIPISON!=0)
385 void seed_uniquestream( rng_state_t* Xin, myID_t clusterID, myID_t machineID, myID_t runID, myID_t streamID ){
386  seed_vielbein(Xin,0);
387  Xin->sumtot = apply_bigskip(Xin->V, Xin->V, clusterID, machineID, runID, streamID );
388  if (Xin->fh==NULL){Xin->fh=stdout;} // if the filehandle is not yet set, make it stdout
389 }
390 
391 void branch_inplace( rng_state_t* Xin, myID_t* IDvec ){
392  Xin->sumtot = apply_bigskip(Xin->V, Xin->V, IDvec[3], IDvec[2], IDvec[1], IDvec[0] );
393 }
394 
395 myuint apply_bigskip(myuint* Vout, myuint* Vin, myID_t clusterID, myID_t machineID, myID_t runID, myID_t streamID ){
396  /*
397  makes a derived state vector, Vout, from the mother state vector Vin
398  by skipping a large number of steps, determined by the given seeding ID's
399 
400  it is mathematically guaranteed that the substreams derived in this way from the SAME (!!!) Vin will not collide provided
401  1) at least one bit of ID is different
402  2) less than 10^100 numbers are drawn from the stream
403  (this is good enough : a single CPU will not exceed this in the lifetime of the universe, 10^19 sec,
404  even if it had a clock cycle of Planch time, 10^44 Hz )
405 
406  Caution: never apply this to a derived vector, just choose some mother vector Vin, for example the unit vector by seed_vielbein(X,0),
407  and use it in all your runs, just change runID to get completely nonoverlapping streams of random numbers on a different day.
408 
409  clusterID and machineID are provided for the benefit of large organizations who wish to ensure that a simulation
410  which is running in parallel on a large number of clusters and machines will have non-colliding source of random numbers.
411 
412  did i repeat it enough times? the non-collision guarantee is absolute, not probabilistic
413 
414  */
415 
416 
417  const myuint skipMat[128][N] =
418 
419 #if (N==88)
420 #include "mixmax_skip_N88.icc" // to make this file, delete all except some chosen 128 rows of the coefficients table
421 #elif (N==256)
422 #include "mixmax_skip_N256.icc"
423 #elif (N==1000)
424 #include "mixmax_skip_N1000.icc"
425 #elif (N==3150)
426 #include "mixmax_skip_N3150.icc"
427 #endif
428  ;
429 
430  myID_t IDvec[4] = {streamID, runID, machineID, clusterID};
431  int r,i,j, IDindex;
432  myID_t id;
433  myuint Y[N], cum[N];
434  myuint coeff;
435  myuint* rowPtr;
436  myuint sumtot=0;
437 
438 
439  for (i=0; i<N; i++) { Y[i] = Vin[i]; sumtot = modadd( sumtot, Vin[i]); } ;
440  for (IDindex=0; IDindex<4; IDindex++) { // go from lower order to higher order ID
441  id=IDvec[IDindex];
442  //printf("now doing ID at level %d, with ID = %d\n", IDindex, id);
443  r = 0;
444  while (id){
445  if (id & 1) {
446  rowPtr = (myuint*)skipMat[r + IDindex*8*sizeof(myID_t)];
447  //printf("free coeff for row %d is %llu\n", r, rowPtr[0]);
448  for (i=0; i<N; i++){ cum[i] = 0; }
449  for (j=0; j<N; j++){ // j is lag, enumerates terms of the poly
450  // for zero lag Y is already given
451  coeff = rowPtr[j]; // same coeff for all i
453  sumtot = iterate_raw_vec(Y, sumtot);
454  }
455  sumtot=0;
456  for (i=0; i<N; i++){ Y[i] = cum[i]; sumtot = modadd( sumtot, cum[i]); } ;
457  }
458  id = (id >> 1); r++; // bring up the r-th bit in the ID
459  }
460  }
461  sumtot=0;
462  for (i=0; i<N; i++){ Vout[i] = Y[i]; sumtot = modadd( sumtot, Y[i]); } ; // returns sumtot, and copy the vector over to Vout
463  return (sumtot) ;
464 }
465 #else
466 #warning For this N, we dont have the skipping coefficients yet, using alternative method to seed
467 
468 void seed_uniquestream( rng_state_t* Xin, myID_t clusterID, myID_t machineID, myID_t runID, myID_t streamID ){
469  Xin->V[0] = (myuint)clusterID;
470  Xin->V[1] = (myuint)machineID;
471  Xin->V[2] = (myuint)runID;
472  Xin->V[3] = (myuint)streamID;
473  Xin->V[4] = (myuint)clusterID << 5;
474  Xin->V[5] = (myuint)machineID << 7;
475  Xin->V[6] = (myuint)runID << 11;
476  Xin->V[7] = (myuint)streamID << 13;
477  precalc(Xin);
478  Xin->sumtot = iterate_raw_vec(Xin->V, Xin->sumtot);
479  Xin->sumtot = iterate_raw_vec(Xin->V, Xin->sumtot);
480 }
481 #endif // SKIPISON
482 
483 
484 
485 
TControlBar * bar
Definition: demos.C:15
rng_state_t * rng_copy(myuint *Y)
Definition: mixmax.cxx:192
myuint V[N]
Definition: mixmax.h:61
#define SPECIALMUL
Definition: mixmax.h:138
#define ERROR_READING_STATE_FILE
Definition: mixmax.h:212
void branch_inplace(rng_state_t *Xin, myID_t *IDvec)
Definition: mixmax.cxx:391
static const char * filename()
#define MOD_MULSPEC(k)
Definition: mixmax.h:140
#define N
int rng_free(rng_state_t *X)
Definition: mixmax.cxx:186
#define SEED_WAS_ZERO
Definition: mixmax.h:211
rng_state_t * rng_alloc()
Definition: mixmax.cxx:178
#define FUSEDMODMULVEC
Definition: mixmax.cxx:376
TSocket * s1
Definition: hserv2.C:36
double get_next_float(rng_state_t *X)
Definition: mixmax.cxx:97
#define MOD_MERSENNE(k)
Definition: mixmax.h:130
#define ARRAY_INDEX_OUT_OF_BOUNDS
Definition: mixmax.h:210
uint64_t myuint
Definition: mixmax.h:52
XFontStruct * id
Definition: TGX11.cxx:108
char * out
Definition: TBase64.cxx:29
#define F(x, y, z)
ROOT::R::TRInterface & r
Definition: Object.C:4
myuint get_next(rng_state_t *X)
Definition: mixmax.cxx:83
void seed_uniquestream(rng_state_t *Xin, myID_t clusterID, myID_t machineID, myID_t runID, myID_t streamID)
Definition: mixmax.cxx:385
FILE * fh
Definition: mixmax.h:64
uint32_t myID_t
Definition: mixmax.h:83
int counter
Definition: mixmax.h:63
void seed_vielbein(rng_state_t *X, unsigned int index)
Definition: mixmax.cxx:212
myuint fmodmulM61(myuint cum, myuint s, myuint a)
Definition: mixmax.cxx:292
TLine * l
Definition: textangle.C:4
tuple pl
Definition: first.py:10
myuint sumtot
Definition: mixmax.h:62
#define M61
Definition: mixmax.h:121
#define BITS
Definition: gifdecode.c:7
#define MERSBASE
Definition: mixmax.h:127
void read_state(rng_state_t *X, const char filename[])
Definition: mixmax.cxx:320
tuple free
Definition: fildir.py:30
#define ERROR_READING_STATE_CHECKSUM
Definition: mixmax.h:214
void seed_spbox(rng_state_t *X, myuint seed)
Definition: mixmax.cxx:229
void print_state(rng_state_t *X)
Definition: mixmax.cxx:307
myuint apply_bigskip(myuint *Vout, myuint *Vin, myID_t clusterID, myID_t machineID, myID_t runID, myID_t streamID)
Definition: mixmax.cxx:395
uint64_t MULWU(uint64_t k)
Definition: mixmax.cxx:48
void iterate_and_fill_array(rng_state_t *X, double *array)
Definition: mixmax.cxx:137
#define INV_MERSBASE
Definition: mixmax.h:132
myuint precalc(rng_state_t *X)
Definition: mixmax.cxx:252
myuint iterate_raw_vec(myuint *Y, myuint sumtotOld)
Definition: mixmax.cxx:53
typedef void((*Func_t)())
int rng_get_N(void)
Definition: mixmax.cxx:264
#define NULL
Definition: Rtypes.h:82
int iterate(rng_state_t *X)
Definition: mixmax.cxx:40
std::complex< float_v > Z
Definition: main.cpp:120
#define ERROR_READING_STATE_COUNTER
Definition: mixmax.h:213
Vc_ALWAYS_INLINE_L T *Vc_ALWAYS_INLINE_R malloc(size_t n)
Allocates memory on the Heap with alignment and padding suitable for vectorized access.
Definition: memory.h:67
const Int_t n
Definition: legend1.C:16
#define MASK32
Definition: mixmax.cxx:290
void fill_array(rng_state_t *X, unsigned int n, double *array)
Definition: mixmax.cxx:117
myuint modadd(myuint foo, myuint bar)
Definition: mixmax.cxx:162