xref: /petsc/src/sys/classes/random/impls/rand48/rand48.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1 #define PETSC_DESIRE_FEATURE_TEST_MACROS /* for drand48() */
2 #include <../src/sys/classes/random/randomimpl.h>
3 
4 PetscErrorCode  PetscRandomSeed_Rand48(PetscRandom r)
5 {
6   PetscFunctionBegin;
7   srand48(r->seed);
8   PetscFunctionReturn(0);
9 }
10 
11 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(0);
25 }
26 
27 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(0);
38 }
39 
40 static struct _PetscRandomOps PetscRandomOps_Values = {
41   /* 0 */
42   PetscRandomSeed_Rand48,
43   PetscRandomGetValue_Rand48,
44   PetscRandomGetValueReal_Rand48,
45   NULL,
46   /* 5 */
47   NULL
48 };
49 
50 /*MC
51    PETSCRAND48 - access to the basic Unix drand48() random number generator
52 
53    Options Database Keys:
54 . -random_type <rand,rand48,sprng>
55 
56   Level: beginner
57 
58 .seealso: RandomCreate(), RandomSetType(), PETSCRAND, PETSCSPRNG
59 M*/
60 
61 PETSC_EXTERN PetscErrorCode PetscRandomCreate_Rand48(PetscRandom r)
62 {
63   PetscErrorCode ierr;
64 
65   PetscFunctionBegin;
66   ierr = PetscMemcpy(r->ops,&PetscRandomOps_Values,sizeof(PetscRandomOps_Values));CHKERRQ(ierr);
67   ierr = PetscObjectChangeTypeName((PetscObject)r,PETSCRAND48);CHKERRQ(ierr);
68   PetscFunctionReturn(0);
69 }
70