Actual source code: sprng.c
2: #include <petsc/private/randomimpl.h>
4: #define USE_MPI
5: #define SIMPLE_SPRNG
6: EXTERN_C_BEGIN
7: #include <sprng.h>
8: EXTERN_C_END
10: PetscErrorCode PetscRandomSeed_Sprng(PetscRandom r)
11: {
12: init_sprng(r->seed,SPRNG_DEFAULT);
13: return 0;
14: }
16: PetscErrorCode PetscRandomGetValue_Sprng(PetscRandom r,PetscScalar *val)
17: {
18: #if defined(PETSC_USE_COMPLEX)
19: if (r->iset) {
20: *val = PetscRealPart(r->width)*sprng() + PetscRealPart(r->low) + (PetscImaginaryPart(r->width)*sprng() + PetscImaginaryPart(r->low)) * PETSC_i;
21: } else {
22: *val = sprng() + sprng()*PETSC_i;
23: }
24: #else
25: if (r->iset) *val = r->width * sprng() + r->low;
26: else *val = sprng();
27: #endif
28: return 0;
29: }
31: PetscErrorCode PetscRandomGetValueReal_Sprng(PetscRandom r,PetscReal *val)
32: {
33: #if defined(PETSC_USE_COMPLEX)
34: if (r->iset) *val = PetscRealPart(r->width)*sprng() + PetscRealPart(r->low);
35: else *val = sprng();
36: #else
37: if (r->iset) *val = r->width * sprng() + r->low;
38: else *val = sprng();
39: #endif
40: return 0;
41: }
43: static struct _PetscRandomOps PetscRandomOps_Values = {
44: PetscDesignatedInitializer(seed,PetscRandomSeed_Sprng),
45: PetscDesignatedInitializer(getvalue,PetscRandomGetValue_Sprng),
46: PetscDesignatedInitializer(getvaluereal,PetscRandomGetValueReal_Sprng),
47: };
49: /*MC
50: PETSCSPRNG - access to the publically available random number generator sprng
52: Options Database Keys:
53: . -random_type <rand,rand48,sprng>
55: Level: beginner
57: PETSc must have been ./configure with the option --download-sprng to use
58: this random number generator.
60: This is NOT currently using a parallel random number generator. Sprng does have
61: an MPI version we should investigate.
63: .seealso: RandomCreate(), RandomSetType(), PETSCRAND, PETSCRAND48
64: M*/
66: PETSC_EXTERN PetscErrorCode PetscRandomCreate_Sprng(PetscRandom r)
67: {
68: PetscMemcpy(r->ops,&PetscRandomOps_Values,sizeof(PetscRandomOps_Values));
69: PetscObjectChangeTypeName((PetscObject)r,PETSCSPRNG);
70: return 0;
71: }