Actual source code: sprng.c
petsc-3.7.3 2016-08-01
2: #include <../src/sys/classes/random/randomimpl.h>
4: #define USE_MPI
5: #define SIMPLE_SPRNG
6: EXTERN_C_BEGIN
7: #include <sprng.h>
8: EXTERN_C_END
12: PetscErrorCode PetscRandomSeed_Sprng(PetscRandom r)
13: {
15: init_sprng(r->seed,SPRNG_DEFAULT);
16: return(0);
17: }
21: PetscErrorCode PetscRandomGetValue_Sprng(PetscRandom r,PetscScalar *val)
22: {
24: #if defined(PETSC_USE_COMPLEX)
25: if (r->iset) {
26: *val = PetscRealPart(r->width)*sprng() + PetscRealPart(r->low) + (PetscImaginaryPart(r->width)*sprng() + PetscImaginaryPart(r->low)) * PETSC_i;
27: } else {
28: *val = sprng() + sprng()*PETSC_i;
29: }
30: #else
31: if (r->iset) *val = r->width * sprng() + r->low;
32: else *val = sprng();
33: #endif
34: return(0);
35: }
39: PetscErrorCode PetscRandomGetValueReal_Sprng(PetscRandom r,PetscScalar *val)
40: {
42: #if defined(PETSC_USE_COMPLEX)
43: if (r->iset) *val = PetscRealPart(r->width)*sprng() + PetscRealPart(r->low);
44: else *val = sprng();
45: #else
46: if (r->iset) *val = r->width * sprng() + r->low;
47: else *val = sprng();
48: #endif
49: return(0);
50: }
52: static struct _PetscRandomOps PetscRandomOps_Values = {
53: /* 0 */
54: PetscRandomSeed_Sprng,
55: PetscRandomGetValue_Sprng,
56: PetscRandomGetValueReal_Sprng,
57: 0,
58: /* 5 */
59: 0
60: };
62: /*MC
63: PETSCSPRNG- access to the publically available random number generator sprng
65: Options Database Keys:
66: . -random_type <rand,rand48,sprng>
68: Level: beginner
70: PETSc must have been ./configure with the option --download-sprng to use
71: this random number generator.
73: This is NOT currently using a parallel random number generator. Sprng does have
74: an MPI version we should investigate.
76: .seealso: RandomCreate(), RandomSetType(), PETSCRAND, PETSCRAND48
77: M*/
81: PETSC_EXTERN PetscErrorCode PetscRandomCreate_Sprng(PetscRandom r)
82: {
86: PetscMemcpy(r->ops,&PetscRandomOps_Values,sizeof(PetscRandomOps_Values));
87: PetscObjectChangeTypeName((PetscObject)r,PETSCSPRNG);
88: return(0);
89: }