xref: /petsc/src/sys/classes/random/impls/rand48/rand48.c (revision 8cc058d9cd56c1ccb3be12a47760ddfc446aaffc)
15c6c1daeSBarry Smith 
25c6c1daeSBarry Smith #include <../src/sys/classes/random/randomimpl.h>
35c6c1daeSBarry Smith 
45c6c1daeSBarry Smith #undef __FUNCT__
55c6c1daeSBarry Smith #define __FUNCT__ "PetscRandomSeed_Rand48"
65c6c1daeSBarry Smith PetscErrorCode  PetscRandomSeed_Rand48(PetscRandom r)
75c6c1daeSBarry Smith {
85c6c1daeSBarry Smith   PetscFunctionBegin;
95c6c1daeSBarry Smith   srand48(r->seed);
105c6c1daeSBarry Smith   PetscFunctionReturn(0);
115c6c1daeSBarry Smith }
125c6c1daeSBarry Smith 
135c6c1daeSBarry Smith #undef __FUNCT__
145c6c1daeSBarry Smith #define __FUNCT__ "PetscRandomGetValue_Rand48"
155c6c1daeSBarry Smith PetscErrorCode  PetscRandomGetValue_Rand48(PetscRandom r,PetscScalar *val)
165c6c1daeSBarry Smith {
175c6c1daeSBarry Smith   PetscFunctionBegin;
185c6c1daeSBarry Smith #if defined(PETSC_USE_COMPLEX)
195c6c1daeSBarry Smith   if (r->iset) {
20b0043f70SBarry Smith     *val = PetscRealPart(r->width)*(PetscReal)drand48() + PetscRealPart(r->low) + (PetscImaginaryPart(r->width)*(PetscReal)drand48() + PetscImaginaryPart(r->low)) * PETSC_i;
215c6c1daeSBarry Smith   } else {
225c6c1daeSBarry Smith     *val = (PetscReal)drand48() + (PetscReal)drand48()*PETSC_i;
235c6c1daeSBarry Smith   }
245c6c1daeSBarry Smith #else
255c6c1daeSBarry Smith   if (r->iset) *val = r->width * drand48() + r->low;
265c6c1daeSBarry Smith   else         *val = drand48();
275c6c1daeSBarry Smith #endif
285c6c1daeSBarry Smith   PetscFunctionReturn(0);
295c6c1daeSBarry Smith }
305c6c1daeSBarry Smith 
315c6c1daeSBarry Smith #undef __FUNCT__
325c6c1daeSBarry Smith #define __FUNCT__ "PetscRandomGetValueReal_Rand48"
335c6c1daeSBarry Smith PetscErrorCode  PetscRandomGetValueReal_Rand48(PetscRandom r,PetscReal *val)
345c6c1daeSBarry Smith {
355c6c1daeSBarry Smith   PetscFunctionBegin;
365c6c1daeSBarry Smith #if defined(PETSC_USE_COMPLEX)
375c6c1daeSBarry Smith   if (r->iset) *val = PetscRealPart(r->width)*drand48() + PetscRealPart(r->low);
385c6c1daeSBarry Smith   else         *val = drand48();
395c6c1daeSBarry Smith #else
405c6c1daeSBarry Smith   if (r->iset) *val = r->width * drand48() + r->low;
415c6c1daeSBarry Smith   else         *val = drand48();
425c6c1daeSBarry Smith #endif
435c6c1daeSBarry Smith   PetscFunctionReturn(0);
445c6c1daeSBarry Smith }
455c6c1daeSBarry Smith 
465c6c1daeSBarry Smith static struct _PetscRandomOps PetscRandomOps_Values = {
475c6c1daeSBarry Smith   /* 0 */
485c6c1daeSBarry Smith   PetscRandomSeed_Rand48,
495c6c1daeSBarry Smith   PetscRandomGetValue_Rand48,
505c6c1daeSBarry Smith   PetscRandomGetValueReal_Rand48,
515c6c1daeSBarry Smith   0,
525c6c1daeSBarry Smith   /* 5 */
535c6c1daeSBarry Smith   0
545c6c1daeSBarry Smith };
555c6c1daeSBarry Smith 
565c6c1daeSBarry Smith /*MC
575c6c1daeSBarry Smith    PETSCRAND48 - access to the basic Unix drand48() random number generator
585c6c1daeSBarry Smith 
595c6c1daeSBarry Smith    Options Database Keys:
605c6c1daeSBarry Smith . -random_type <rand,rand48,sprng>
615c6c1daeSBarry Smith 
625c6c1daeSBarry Smith   Level: beginner
635c6c1daeSBarry Smith 
645c6c1daeSBarry Smith .seealso: RandomCreate(), RandomSetType(), PETSCRAND, PETSCSPRNG
655c6c1daeSBarry Smith M*/
665c6c1daeSBarry Smith 
675c6c1daeSBarry Smith #undef __FUNCT__
685c6c1daeSBarry Smith #define __FUNCT__ "PetscRandomCreate_Rand48"
69*8cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PetscRandomCreate_Rand48(PetscRandom r)
705c6c1daeSBarry Smith {
715c6c1daeSBarry Smith   PetscErrorCode ierr;
725c6c1daeSBarry Smith 
735c6c1daeSBarry Smith   PetscFunctionBegin;
745c6c1daeSBarry Smith   ierr = PetscMemcpy(r->ops,&PetscRandomOps_Values,sizeof(PetscRandomOps_Values));CHKERRQ(ierr);
755c6c1daeSBarry Smith   /* r->bops->publish   = PetscRandomPublish; */
765c6c1daeSBarry Smith   /*  r->petscnative     = PETSC_TRUE;  */
775c6c1daeSBarry Smith 
785c6c1daeSBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)r,PETSCRAND48);CHKERRQ(ierr);
795c6c1daeSBarry Smith   PetscFunctionReturn(0);
805c6c1daeSBarry Smith }
81