SDDSlib
Loading...
Searching...
No Matches
halton.c
Go to the documentation of this file.
1/**
2 * @file halton.c
3 * @brief Implementation of Halton and modified Halton sequences.
4 *
5 * @copyright
6 * - (c) 2002 The University of Chicago, as Operator of Argonne National Laboratory.
7 * - (c) 2002 The Regents of the University of California, as Operator of Los Alamos National Laboratory.
8 *
9 * @license
10 * This file is distributed under the terms of the Software License Agreement
11 * found in the file LICENSE included with this distribution.
12 *
13 * @author M. Borland, R. Soliday, H. Shang, Y. Wang
14 */
15
16#include "mdb.h"
17
18static double *lastPointValue = NULL;
19static long *R = NULL;
20static long sequencesInUse = 0;
21/* These 12 primes work pretty well together.
22 * If more are needed, they are generated on the fly.
23 */
24#define N_SEQ_PREDEFINED 12
25static long Rvalues[N_SEQ_PREDEFINED] = {2, 3, 5, 7, 11, 19, 23, 29, 37, 47, 59, 67};
26
27/**
28 * @brief Initialize and start a new Halton sequence.
29 *
30 * Initializes a new Halton sequence with the given radix and starting value.
31 * If the provided radix is not a prime or is non-positive, a suitable prime is chosen.
32 *
33 * @param radix Pointer to an integer specifying the desired prime radix.
34 * If non-positive, a prime is chosen automatically.
35 * @param value The starting value for the sequence.
36 * @return The sequence ID (1-based) on success, or 0 on failure.
37 */
38int32_t startHaltonSequence(int32_t *radix, double value) {
39 int32_t ID;
40 if ((sequencesInUse == 0 &&
41 (!(lastPointValue = malloc(sizeof(*lastPointValue))) ||
42 !(R = malloc(sizeof(*R))))) ||
43 (!(lastPointValue = realloc(lastPointValue, (sequencesInUse + 1) * sizeof(*lastPointValue))) ||
44 !(R = realloc(R, (sequencesInUse + 1) * sizeof(*R)))))
45 return 0;
46 if (*radix > 0) {
47 if (is_prime(*radix) != 1)
48 return 0;
49 R[sequencesInUse] = *radix;
50 ID = sequencesInUse;
51 } else {
52 /* generate a new, unique prime number for use as the radix */
53 long i, passed;
54 if ((ID = sequencesInUse) < N_SEQ_PREDEFINED)
55 /* try one of the favored values */
56 *radix = Rvalues[ID];
57 else
58 *radix = 2;
59 passed = 0;
60 while (!passed) {
61 passed = 1;
62 for (i = 0; i < sequencesInUse; i++) {
63 if (R[i] == *radix) {
64 passed = 0;
65 (*radix)++;
66 while (is_prime(*radix) != 1)
67 (*radix)++;
68 }
69 }
70 }
71 R[ID] = *radix;
72 }
73 lastPointValue[ID] = value;
74 sequencesInUse++;
75 return ID + 1;
76}
77
78/**
79 * @brief Restart an existing Halton sequence from a new initial value.
80 *
81 * Resets the specified Halton sequence so that subsequent calls will continue from
82 * the given initial value.
83 *
84 * @param ID The sequence ID to restart.
85 * @param value The new starting value.
86 * @return 1 on success, or -1 if the ID is invalid.
87 */
88int32_t restartHaltonSequence(long ID, double value) {
89 ID -= 1;
90
91 if (ID > sequencesInUse || ID < 0)
92 return -1;
93
94 lastPointValue[ID] = value;
95
96 return 1;
97}
98
99/**
100 * @brief Get the next point in a Halton sequence.
101 *
102 * Computes the next value in the specified Halton sequence.
103 *
104 * @param ID The sequence ID.
105 * @return The next point in the sequence, or -1 if the ID is invalid.
106 */
107double nextHaltonSequencePoint(long ID) {
108 double r, f, value;
109
110 ID -= 1;
111
112 if (ID > sequencesInUse || ID < 0)
113 return -1;
114
115 f = 1 - lastPointValue[ID];
116 r = 1. / R[ID];
117 while (f <= r)
118 r = r / R[ID];
119 value = lastPointValue[ID] + (R[ID] + 1) * r - 1;
120 lastPointValue[ID] = value;
121 return value;
122}
123
124/*following code is from outside for improved halton sequence
125 Alogrithm 659, Collected Algorithm from ACM
126 This is the C version of halton sequences
127 Derandom Algorithm is added on 8/12/03
128 by Hongmei CHI (CS/FSU)
129 Modified (9.29.03)
130*/
131#define MAX_D 500
132static int32_t sDim = 12, nextPoint[12];
133static double eError;
134static int32_t prime[] = {
135 2,
136 3,
137 5,
138 7,
139 11,
140 13,
141 17,
142 19,
143 23,
144 29,
145 31,
146 37,
147 41,
148 43,
149 47,
150 53,
151 59,
152 61,
153 67,
154 71,
155 73,
156 79,
157 83,
158 89,
159 97,
160 101,
161 103,
162 107,
163 109,
164 113,
165 127,
166 131,
167 137,
168 139,
169 149,
170 151,
171 157,
172 163,
173 167,
174 173,
175 179,
176 181,
177 191,
178 193,
179 197,
180 199,
181 211,
182 223,
183 227,
184 229,
185 233,
186 239,
187 241,
188 251,
189 257,
190 263,
191 269,
192 271,
193 277,
194 281,
195 283,
196 293,
197 307,
198 311,
199 313,
200 317,
201 331,
202 337,
203 347,
204 349,
205 353,
206 359,
207 367,
208 373,
209 379,
210 383,
211 389,
212 397,
213 401,
214 409,
215 419,
216 421,
217 431,
218 433,
219 439,
220 443,
221 449,
222 457,
223 461,
224 463,
225 467,
226 479,
227 487,
228 491,
229 499,
230};
231static double iprime[] = {
232 2,
233 3,
234 5,
235 7,
236 11,
237 13,
238 17,
239 19,
240 23,
241 29,
242 31,
243 37,
244 41,
245 43,
246 47,
247 53,
248 59,
249 61,
250 67,
251 71,
252 73,
253 79,
254 83,
255 89,
256 97,
257 101,
258 103,
259 107,
260 109,
261 113,
262 127,
263 131,
264 137,
265 139,
266 149,
267 151,
268 157,
269 163,
270 167,
271 173,
272 179,
273 181,
274 191,
275 193,
276 197,
277 199,
278 211,
279 223,
280 227,
281 229,
282 233,
283 239,
284 241,
285 251,
286 257,
287 263,
288 269,
289 271,
290 277,
291 281,
292 283,
293 293,
294 307,
295 311,
296 313,
297 317,
298 331,
299 337,
300 347,
301 349,
302 353,
303 359,
304 367,
305 373,
306 379,
307 383,
308 389,
309 397,
310 401,
311 409,
312 419,
313 421,
314 431,
315 433,
316 439,
317 443,
318 449,
319 457,
320 461,
321 463,
322 467,
323 479,
324 487,
325 491,
326 499,
327};
328
329static int32_t modSequenceInUse = 0;
330static int32_t primroots[][10] = {
331 {1, 2, 3, 3, 8, 11, 12, 14, 7, 18},
332 {12, 13, 17, 18, 29, 14, 18, 43, 41, 44},
333 {40, 30, 47, 65, 71, 28, 40, 60, 79, 89},
334 {56, 50, 52, 61, 108, 56, 66, 63, 60, 66},
335 {104, 76, 111, 142, 71, 154, 118, 84, 127, 142},
336 {84, 105, 186, 178, 188, 152, 165, 159, 103, 205},
337 {166, 173, 188, 181, 91, 233, 210, 217, 153, 212},
338};
339
340static int32_t warnockOpt[] = {
341 1,
342 2,
343 2,
344 5,
345 3,
346 7,
347 3,
348 10,
349 18,
350 11,
351 17,
352 5,
353 17,
354 26,
355 40,
356 14,
357 40,
358 44,
359 12,
360 31,
361 45,
362 70,
363 8,
364 38,
365 82,
366 8,
367 12,
368 38,
369 47,
370 70,
371 29,
372 57,
373 97,
374 110,
375 32,
376 48,
377 84,
378 124,
379 155,
380 26,
381 69,
382 83,
383 157,
384 171,
385 8,
386 22,
387 112,
388 205,
389 15,
390 31,
391 61,
392 105,
393 127,
394 212,
395 12,
396 57,
397 109,
398 133,
399 179,
400 210,
401 231,
402 34,
403 161,
404 199,
405 222,
406 255,
407 59,
408 120,
409 218,
410 237,
411 278,
412 341,
413 54,
414 110,
415 176,
416 218,
417 280,
418 369,
419 17,
420 97,
421 193,
422 221,
423 331,
424 350,
425 419,
426 21,
427 85,
428 173,
429 221,
430 243,
431 288,
432 424,
433 45,
434 78,
435 173,
436 213,
437 288,
438 426,
439 455,
440 138,
441};
442static double *quasi = NULL;
443
444int32_t power(int32_t a, int32_t b, int32_t m) {
445 int32_t i, c = 1;
446 for (i = 0; i < b; i++)
447 c = (c * a) % m;
448 return c;
449}
450
451int32_t primes() {
452 int32_t i, j, a[MAX_D + 1];
453 for (a[1] = 0, i = 2; i <= MAX_D; i++)
454 a[i] = 1;
455 for (i = 2; i <= MAX_D / 2; i++)
456 for (j = 2; j <= MAX_D / i; j++)
457 a[i * j] = 0;
458 for (i = 1, j = 0; i <= MAX_D; i++) {
459 if (a[i]) {
460 prime[j] = i;
461 iprime[j] = (double)i;
462 j++;
463 }
464 }
465 return j;
466}
467
468int inhalt(int dimen, int atmost, double tiny, double *quasi) {
469 double delta;
470 int i /*,m*/;
471
472 /* check dimen */
473 sDim = dimen;
474 if (sDim < 1 || sDim > 1000)
475 return (-1);
476
477 /* compute and check tolerance*/
478
479 eError = 0.9 * (1.0 / (atmost * prime[sDim - 1]) - 10.0 * tiny);
480 delta = 100 * tiny * (double)(atmost + 1) * log10((double)atmost);
481 if (delta >= 0.09 * (eError - 10.0 * tiny))
482 return (-2);
483
484 /* now compute first vector */
485
486 /*m=1;*/
487 for (i = 0; i < sDim; i++) {
488 iprime[i] = 1.0 / prime[i];
489 quasi[i] = iprime[i];
490 /*m=i*prime[i];*/
491 nextPoint[i] = 2;
492 }
493
494 /*printf("largest prime=%d, %d \n",prime[sDim-1], m);*/
495
496 return 0;
497}
498
499/**
500 * @brief Start a modified Halton sequence.
501 *
502 * Initializes a modified Halton sequence with predefined prime bases and internal parameters.
503 *
504 * @param radix Pointer to an integer to store the chosen prime radix.
505 * @param tiny A small tolerance value.
506 * @return The sequence ID (1-based) on success, or -1 on failure.
507 */
508int32_t startModHaltonSequence(int32_t *radix, double tiny) {
509 int32_t modID, dimen = 12, total_points = 100000;
510 tiny = 0;
511 /* check dimen*/
512 if (!modSequenceInUse) {
513 /*generate primes
514 primes(); */
515 /* commented out primes generation and put it as global constants to skip generate it every time when start the sequence*/
516 if (!quasi)
517 quasi = malloc(sizeof(*quasi) * sDim);
518 if (inhalt(dimen, total_points, tiny, quasi) < 0) {
519 fprintf(stderr, "Unable to start modHalton sequence.\n");
520 return -1;
521 }
522 }
523 modID = modSequenceInUse;
524 *radix = prime[modID];
525 modSequenceInUse++;
526 return modID + 1;
527}
528
529/**
530 * @brief Restart an existing modified Halton sequence.
531 *
532 * Resets internal parameters so that the modified Halton sequence restarts.
533 *
534 * @param ID The sequence ID.
535 * @param tiny A small tolerance value.
536 * @return 1 on success, or -1 on failure.
537 */
538int32_t restartModHaltonSequence(long ID, double tiny) {
539 int32_t dimen = 12, total_points = 100000;
540
541 tiny = 0;
542 if (inhalt(dimen, total_points, tiny, quasi) < 0) {
543 fprintf(stderr, "Unable to start modHalton sequence.\n");
544 return -1;
545 }
546
547 return 1;
548}
549
550int32_t generateModHaltSequence(double *quasi, double *dq, double *wq, long ID) {
551 int32_t i = (int32_t)ID, j, k, ytemp[40], xtemp[40], ltemp, mtemp;
552 double t, f, g, h;
553 /* generate quasi one compoment at a time using radix prime[k] for
554 component k */
555
556 t = iprime[i];
557 f = 1.0 - quasi[i];
558 g = 1.0;
559 h = t;
560 while ((f - h) < eError)
561 /* this checks whether q+h>1-eError */
562 {
563 g = h;
564 h *= t;
565 }
566 quasi[i] = g + h - f;
567
568 k = 0;
569 mtemp = nextPoint[i];
570 ltemp = prime[i];
571
572 while (mtemp != 0) {
573 ytemp[k] = mtemp % ltemp;
574 mtemp = mtemp / ltemp;
575 k++;
576 }
577 /*generating Optimal primitive root */
578 for (j = 0; j < k; j++) {
579 /* xtemp[j] = (ytemp[j]*power(primroots[i/10][i%10], nextn%ltemp, ltemp))%ltemp; */
580 xtemp[j] = (warnockOpt[i] * power(primroots[i / 10][i % 10], ytemp[j], ltemp)) % ltemp;
581 xtemp[j] -= ytemp[j];
582 }
583 dq[i] = 0;
584 t = iprime[i];
585 for (j = 0; j < k; j++) {
586 dq[i] += xtemp[j] * t;
587 t *= iprime[i];
588 }
589 dq[i] += quasi[i];
590 /* generating Warnock Optimal sequences */
591 for (j = 0; j < k; j++) {
592 /* if( i%2==0)xtemp[j]= (ytemp[j]*warnockOpt[i])%ltemp;
593 else xtemp[j]= (ytemp[j]*warnockOpt[i]*warnockOpt[i])%ltemp;
594 */
595 xtemp[j] = (ytemp[j] * power(warnockOpt[i], i + 1, ltemp)) % ltemp;
596 /* if(i>=3)
597 xtemp[j]=((prime[i-1]-1)*power(warnockOpt[i],ytemp[j],ltemp))%ltemp;
598 else xtemp[j]= (power(warnockOpt[i],ytemp[j],ltemp))%ltemp;
599 */
600 xtemp[j] -= ytemp[j];
601 }
602
603 wq[i] = 0;
604 t = iprime[i];
605 for (j = 0; j < k; j++) {
606 wq[i] += xtemp[j] * t;
607 t *= iprime[i];
608 }
609 wq[i] += quasi[i];
610
611 nextPoint[i]++;
612
613 return (0);
614}
615
616/**
617 * @brief Retrieve the next point from the modified Halton sequence.
618 *
619 * Generates and returns the next value for a given sequence ID in the modified Halton sequence.
620 *
621 * @param ID The sequence ID.
622 * @return The next point in the modified Halton sequence.
623 */
625 static double *dq = NULL, *wq = NULL;
626
627 ID -= 1;
628 if (dq == NULL)
629 dq = malloc(sizeof(*dq) * sDim);
630 if (wq == NULL)
631 wq = malloc(sizeof(*wq) * sDim);
632 if (!modSequenceInUse) {
633 fprintf(stderr, "ModHalton sequence not started.\n");
634 exit(1);
635 }
636 if (ID < 0 || ID > sDim - 1) {
637 fprintf(stderr, "Invalid ID (%ld) provided\n", ID);
638 exit(1);
639 }
640 generateModHaltSequence(quasi, dq, wq, ID);
641
642 return wq[ID];
643}
int64_t is_prime(int64_t number)
Determine if a number is prime.
Definition factorize.c:26
int32_t startModHaltonSequence(int32_t *radix, double tiny)
Start a modified Halton sequence.
Definition halton.c:508
double nextModHaltonSequencePoint(long ID)
Retrieve the next point from the modified Halton sequence.
Definition halton.c:624
int32_t restartHaltonSequence(long ID, double value)
Restart an existing Halton sequence from a new initial value.
Definition halton.c:88
int32_t startHaltonSequence(int32_t *radix, double value)
Initialize and start a new Halton sequence.
Definition halton.c:38
int32_t restartModHaltonSequence(long ID, double tiny)
Restart an existing modified Halton sequence.
Definition halton.c:538
double nextHaltonSequencePoint(long ID)
Get the next point in a Halton sequence.
Definition halton.c:107