Actual source code: sprng.c
petsc-3.14.6 2021-03-30
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
10: PetscErrorCode PetscRandomSeed_Sprng(PetscRandom r)
11: {
13: init_sprng(r->seed,SPRNG_DEFAULT);
14: return(0);
15: }
17: PetscErrorCode PetscRandomGetValue_Sprng(PetscRandom r,PetscScalar *val)
18: {
20: #if defined(PETSC_USE_COMPLEX)
21: if (r->iset) {
22: *val = PetscRealPart(r->width)*sprng() + PetscRealPart(r->low) + (PetscImaginaryPart(r->width)*sprng() + PetscImaginaryPart(r->low)) * PETSC_i;
23: } else {
24: *val = sprng() + sprng()*PETSC_i;
25: }
26: #else
27: if (r->iset) *val = r->width * sprng() + r->low;
28: else *val = sprng();
29: #endif
30: return(0);
31: }
33: PetscErrorCode PetscRandomGetValueReal_Sprng(PetscRandom r,PetscReal *val)
34: {
36: #if defined(PETSC_USE_COMPLEX)
37: if (r->iset) *val = PetscRealPart(r->width)*sprng() + PetscRealPart(r->low);
38: else *val = sprng();
39: #else
40: if (r->iset) *val = r->width * sprng() + r->low;
41: else *val = sprng();
42: #endif
43: return(0);
44: }
46: static struct _PetscRandomOps PetscRandomOps_Values = {
47: /* 0 */
48: PetscRandomSeed_Sprng,
49: PetscRandomGetValue_Sprng,
50: PetscRandomGetValueReal_Sprng,
51: 0,
52: /* 5 */
53: 0
54: };
56: /*MC
57: PETSCSPRNG- access to the publically available random number generator sprng
59: Options Database Keys:
60: . -random_type <rand,rand48,sprng>
62: Level: beginner
64: PETSc must have been ./configure with the option --download-sprng to use
65: this random number generator.
67: This is NOT currently using a parallel random number generator. Sprng does have
68: an MPI version we should investigate.
70: .seealso: RandomCreate(), RandomSetType(), PETSCRAND, PETSCRAND48
71: M*/
73: PETSC_EXTERN PetscErrorCode PetscRandomCreate_Sprng(PetscRandom r)
74: {
78: PetscMemcpy(r->ops,&PetscRandomOps_Values,sizeof(PetscRandomOps_Values));
79: PetscObjectChangeTypeName((PetscObject)r,PETSCSPRNG);
80: return(0);
81: }