xref: /petsc/src/sys/classes/random/impls/rand48/rand48.c (revision be37439ebbbdb2f81c3420c175a94aa72e59929c)
10039db0dSBarry Smith #define PETSC_DESIRE_FEATURE_TEST_MACROS /* for drand48() */
2d6cc7855SJacob Faibussowitsch #include <petsc/private/randomimpl.h>
35c6c1daeSBarry Smith 
PetscRandomSeed_Rand48(PetscRandom r)466976f2fSJacob Faibussowitsch static PetscErrorCode PetscRandomSeed_Rand48(PetscRandom r)
5d71ae5a4SJacob Faibussowitsch {
65c6c1daeSBarry Smith   PetscFunctionBegin;
75c6c1daeSBarry Smith   srand48(r->seed);
83ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
95c6c1daeSBarry Smith }
105c6c1daeSBarry Smith 
PetscRandomGetValue_Rand48(PetscRandom r,PetscScalar * val)1166976f2fSJacob Faibussowitsch static PetscErrorCode PetscRandomGetValue_Rand48(PetscRandom r, PetscScalar *val)
12d71ae5a4SJacob Faibussowitsch {
135c6c1daeSBarry Smith   PetscFunctionBegin;
145c6c1daeSBarry Smith #if defined(PETSC_USE_COMPLEX)
155c6c1daeSBarry Smith   if (r->iset) {
16b0043f70SBarry Smith     *val = PetscRealPart(r->width) * (PetscReal)drand48() + PetscRealPart(r->low) + (PetscImaginaryPart(r->width) * (PetscReal)drand48() + PetscImaginaryPart(r->low)) * PETSC_i;
175c6c1daeSBarry Smith   } else {
185c6c1daeSBarry Smith     *val = (PetscReal)drand48() + (PetscReal)drand48() * PETSC_i;
195c6c1daeSBarry Smith   }
205c6c1daeSBarry Smith #else
21*835f2295SStefano Zampini   if (r->iset) *val = r->width * drand48() + r->low;
225c6c1daeSBarry Smith   else *val = drand48();
235c6c1daeSBarry Smith #endif
243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
255c6c1daeSBarry Smith }
265c6c1daeSBarry Smith 
PetscRandomGetValueReal_Rand48(PetscRandom r,PetscReal * val)2766976f2fSJacob Faibussowitsch static PetscErrorCode PetscRandomGetValueReal_Rand48(PetscRandom r, PetscReal *val)
28d71ae5a4SJacob Faibussowitsch {
295c6c1daeSBarry Smith   PetscFunctionBegin;
305c6c1daeSBarry Smith #if defined(PETSC_USE_COMPLEX)
315c6c1daeSBarry Smith   if (r->iset) *val = PetscRealPart(r->width) * drand48() + PetscRealPart(r->low);
325c6c1daeSBarry Smith   else *val = drand48();
335c6c1daeSBarry Smith #else
34*835f2295SStefano Zampini   if (r->iset) *val = r->width * (PetscReal)drand48() + r->low;
356497c311SBarry Smith   else *val = (PetscReal)drand48();
365c6c1daeSBarry Smith #endif
373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
385c6c1daeSBarry Smith }
395c6c1daeSBarry Smith 
405c6c1daeSBarry Smith static struct _PetscRandomOps PetscRandomOps_Values = {
41267267bdSJacob Faibussowitsch   PetscDesignatedInitializer(seed, PetscRandomSeed_Rand48),
42267267bdSJacob Faibussowitsch   PetscDesignatedInitializer(getvalue, PetscRandomGetValue_Rand48),
43267267bdSJacob Faibussowitsch   PetscDesignatedInitializer(getvaluereal, PetscRandomGetValueReal_Rand48),
448a403008SStefano Zampini   PetscDesignatedInitializer(getvalues, NULL),
458a403008SStefano Zampini   PetscDesignatedInitializer(getvaluesreal, NULL),
468a403008SStefano Zampini   PetscDesignatedInitializer(destroy, NULL),
478a403008SStefano Zampini   PetscDesignatedInitializer(setfromoptions, NULL),
485c6c1daeSBarry Smith };
495c6c1daeSBarry Smith 
505c6c1daeSBarry Smith /*MC
51c31d2375SBarry Smith    PETSCRAND48 - access to the basic Unix `drand48()` random number generator
525c6c1daeSBarry Smith 
53c31d2375SBarry Smith    Options Database Key:
541179163eSBarry Smith . -random_type <rand,rand48,sprng> - select the random number generator at runtime
555c6c1daeSBarry Smith 
565c6c1daeSBarry Smith   Level: beginner
575c6c1daeSBarry Smith 
58811af0c4SBarry Smith   Note:
591179163eSBarry Smith   Not recommended because it may produce different results on different systems.
601179163eSBarry Smith 
611179163eSBarry Smith .seealso: `PetscRandomCreate()`, `PetscRandomSetType()`, `PETSCRAND`, `PETSCSPRNG`, `PetscRandomSetFromOptions()`
625c6c1daeSBarry Smith M*/
635c6c1daeSBarry Smith 
PetscRandomCreate_Rand48(PetscRandom r)64d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscRandomCreate_Rand48(PetscRandom r)
65d71ae5a4SJacob Faibussowitsch {
665c6c1daeSBarry Smith   PetscFunctionBegin;
67aea10558SJacob Faibussowitsch   r->ops[0] = PetscRandomOps_Values;
689566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)r, PETSCRAND48));
693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
705c6c1daeSBarry Smith }
71