Actual source code: rand48.c
1: #define PETSC_DESIRE_FEATURE_TEST_MACROS /* for drand48() */
2: #include <petsc/private/randomimpl.h>
4: static PetscErrorCode PetscRandomSeed_Rand48(PetscRandom r)
5: {
6: PetscFunctionBegin;
7: srand48(r->seed);
8: PetscFunctionReturn(PETSC_SUCCESS);
9: }
11: static PetscErrorCode PetscRandomGetValue_Rand48(PetscRandom r, PetscScalar *val)
12: {
13: PetscFunctionBegin;
14: #if defined(PETSC_USE_COMPLEX)
15: if (r->iset) {
16: *val = PetscRealPart(r->width) * (PetscReal)drand48() + PetscRealPart(r->low) + (PetscImaginaryPart(r->width) * (PetscReal)drand48() + PetscImaginaryPart(r->low)) * PETSC_i;
17: } else {
18: *val = (PetscReal)drand48() + (PetscReal)drand48() * PETSC_i;
19: }
20: #else
21: if (r->iset) *val = r->width * drand48() + r->low;
22: else *val = drand48();
23: #endif
24: PetscFunctionReturn(PETSC_SUCCESS);
25: }
27: static PetscErrorCode PetscRandomGetValueReal_Rand48(PetscRandom r, PetscReal *val)
28: {
29: PetscFunctionBegin;
30: #if defined(PETSC_USE_COMPLEX)
31: if (r->iset) *val = PetscRealPart(r->width) * drand48() + PetscRealPart(r->low);
32: else *val = drand48();
33: #else
34: if (r->iset) *val = r->width * drand48() + r->low;
35: else *val = drand48();
36: #endif
37: PetscFunctionReturn(PETSC_SUCCESS);
38: }
40: static struct _PetscRandomOps PetscRandomOps_Values = {
41: PetscDesignatedInitializer(seed, PetscRandomSeed_Rand48),
42: PetscDesignatedInitializer(getvalue, PetscRandomGetValue_Rand48),
43: PetscDesignatedInitializer(getvaluereal, PetscRandomGetValueReal_Rand48),
44: };
46: /*MC
47: PETSCRAND48 - access to the basic Unix `drand48()` random number generator
49: Options Database Key:
50: . -random_type <rand,rand48,sprng> - select the random number generator at runtime
52: Level: beginner
54: Note:
55: Not recommended because it may produce different results on different systems.
57: .seealso: `PetscRandomCreate()`, `PetscRandomSetType()`, `PETSCRAND`, `PETSCSPRNG`, `PetscRandomSetFromOptions()`
58: M*/
60: PETSC_EXTERN PetscErrorCode PetscRandomCreate_Rand48(PetscRandom r)
61: {
62: PetscFunctionBegin;
63: r->ops[0] = PetscRandomOps_Values;
64: PetscCall(PetscObjectChangeTypeName((PetscObject)r, PETSCRAND48));
65: PetscFunctionReturn(PETSC_SUCCESS);
66: }