15c6c1daeSBarry Smith 2b6bfaf9aSBarry Smith #include <../src/sys/classes/random/randomimpl.h> 35c6c1daeSBarry Smith 45c6c1daeSBarry Smith #define USE_MPI 55c6c1daeSBarry Smith #define SIMPLE_SPRNG 65c6c1daeSBarry Smith EXTERN_C_BEGIN 75c6c1daeSBarry Smith #include <sprng.h> 85c6c1daeSBarry Smith EXTERN_C_END 95c6c1daeSBarry Smith 105c6c1daeSBarry Smith #undef __FUNCT__ 115c6c1daeSBarry Smith #define __FUNCT__ "PetscRandomSeed_Sprng" 125c6c1daeSBarry Smith PetscErrorCode PetscRandomSeed_Sprng(PetscRandom r) 135c6c1daeSBarry Smith { 145c6c1daeSBarry Smith PetscFunctionBegin; 155c6c1daeSBarry Smith init_sprng(r->seed,SPRNG_DEFAULT); 165c6c1daeSBarry Smith PetscFunctionReturn(0); 175c6c1daeSBarry Smith } 185c6c1daeSBarry Smith 195c6c1daeSBarry Smith #undef __FUNCT__ 205c6c1daeSBarry Smith #define __FUNCT__ "PetscRandomGetValue_Sprng" 215c6c1daeSBarry Smith PetscErrorCode PetscRandomGetValue_Sprng(PetscRandom r,PetscScalar *val) 225c6c1daeSBarry Smith { 235c6c1daeSBarry Smith PetscFunctionBegin; 245c6c1daeSBarry Smith #if defined(PETSC_USE_COMPLEX) 255c6c1daeSBarry Smith if (r->iset) { 26b0043f70SBarry Smith *val = PetscRealPart(r->width)*sprng() + PetscRealPart(r->low) + (PetscImaginaryPart(r->width)*sprng() + PetscImaginaryPart(r->low)) * PETSC_i; 275c6c1daeSBarry Smith } else { 285c6c1daeSBarry Smith *val = sprng() + sprng()*PETSC_i; 295c6c1daeSBarry Smith } 305c6c1daeSBarry Smith #else 315c6c1daeSBarry Smith if (r->iset) *val = r->width * sprng() + r->low; 325c6c1daeSBarry Smith else *val = sprng(); 335c6c1daeSBarry Smith #endif 345c6c1daeSBarry Smith PetscFunctionReturn(0); 355c6c1daeSBarry Smith } 365c6c1daeSBarry Smith 375c6c1daeSBarry Smith #undef __FUNCT__ 385c6c1daeSBarry Smith #define __FUNCT__ "PetscRandomGetValue_Sprng" 395c6c1daeSBarry Smith PetscErrorCode PetscRandomGetValueReal_Sprng(PetscRandom r,PetscScalar *val) 405c6c1daeSBarry Smith { 415c6c1daeSBarry Smith PetscFunctionBegin; 425c6c1daeSBarry Smith #if defined(PETSC_USE_COMPLEX) 435c6c1daeSBarry Smith if (r->iset) *val = PetscRealPart(r->width)*sprng() + PetscRealPart(r->low); 445c6c1daeSBarry Smith else *val = sprng(); 455c6c1daeSBarry Smith #else 465c6c1daeSBarry Smith if (r->iset) *val = r->width * sprng() + r->low; 475c6c1daeSBarry Smith else *val = sprng(); 485c6c1daeSBarry Smith #endif 495c6c1daeSBarry Smith PetscFunctionReturn(0); 505c6c1daeSBarry Smith } 515c6c1daeSBarry Smith 525c6c1daeSBarry Smith static struct _PetscRandomOps PetscRandomOps_Values = { 535c6c1daeSBarry Smith /* 0 */ 545c6c1daeSBarry Smith PetscRandomSeed_Sprng, 555c6c1daeSBarry Smith PetscRandomGetValue_Sprng, 565c6c1daeSBarry Smith PetscRandomGetValueReal_Sprng, 575c6c1daeSBarry Smith 0, 585c6c1daeSBarry Smith /* 5 */ 595c6c1daeSBarry Smith 0 605c6c1daeSBarry Smith }; 615c6c1daeSBarry Smith 625c6c1daeSBarry Smith /*MC 635c6c1daeSBarry Smith PETSCSPRNG- access to the publically available random number generator sprng 645c6c1daeSBarry Smith 655c6c1daeSBarry Smith Options Database Keys: 665c6c1daeSBarry Smith . -random_type <rand,rand48,sprng> 675c6c1daeSBarry Smith 685c6c1daeSBarry Smith Level: beginner 695c6c1daeSBarry Smith 705c6c1daeSBarry Smith PETSc must have been ./configure with the option --download-sprng to use 715c6c1daeSBarry Smith this random number generator. 725c6c1daeSBarry Smith 735c6c1daeSBarry Smith This is NOT currently using a parallel random number generator. Sprng does have 745c6c1daeSBarry Smith an MPI version we should investigate. 755c6c1daeSBarry Smith 765c6c1daeSBarry Smith .seealso: RandomCreate(), RandomSetType(), PETSCRAND, PETSCRAND48 775c6c1daeSBarry Smith M*/ 785c6c1daeSBarry Smith 795c6c1daeSBarry Smith #undef __FUNCT__ 805c6c1daeSBarry Smith #define __FUNCT__ "PetscRandomCreate_Sprng" 81*8cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PetscRandomCreate_Sprng(PetscRandom r) 825c6c1daeSBarry Smith { 835c6c1daeSBarry Smith PetscErrorCode ierr; 845c6c1daeSBarry Smith 855c6c1daeSBarry Smith PetscFunctionBegin; 865c6c1daeSBarry Smith ierr = PetscMemcpy(r->ops,&PetscRandomOps_Values,sizeof(PetscRandomOps_Values));CHKERRQ(ierr); 875c6c1daeSBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)r,PETSCSPRNG);CHKERRQ(ierr); 885c6c1daeSBarry Smith PetscFunctionReturn(0); 895c6c1daeSBarry Smith } 90