xref: /petsc/src/sys/classes/random/impls/curand/curand.c (revision 1179163e0bf5c4dd309079707fd3c0dfe8d44eee)
1a4af0ceeSJacob Faibussowitsch #include <petsc/private/deviceimpl.h>
2d6cc7855SJacob Faibussowitsch #include <petsc/private/randomimpl.h>
3808ba619SStefano Zampini #include <curand.h>
4808ba619SStefano Zampini 
5808ba619SStefano Zampini typedef struct {
6808ba619SStefano Zampini   curandGenerator_t gen;
7808ba619SStefano Zampini } PetscRandom_CURAND;
8808ba619SStefano Zampini 
9808ba619SStefano Zampini PetscErrorCode PetscRandomSeed_CURAND(PetscRandom r)
10808ba619SStefano Zampini {
11808ba619SStefano Zampini   PetscRandom_CURAND *curand = (PetscRandom_CURAND*)r->data;
12808ba619SStefano Zampini 
13808ba619SStefano Zampini   PetscFunctionBegin;
149566063dSJacob Faibussowitsch   PetscCallCURAND(curandSetPseudoRandomGeneratorSeed(curand->gen,r->seed));
15808ba619SStefano Zampini   PetscFunctionReturn(0);
16808ba619SStefano Zampini }
17808ba619SStefano Zampini 
18808ba619SStefano Zampini PETSC_INTERN PetscErrorCode PetscRandomCurandScale_Private(PetscRandom,size_t,PetscReal*,PetscBool);
19808ba619SStefano Zampini 
20808ba619SStefano Zampini PetscErrorCode  PetscRandomGetValuesReal_CURAND(PetscRandom r, PetscInt n, PetscReal *val)
21808ba619SStefano Zampini {
22808ba619SStefano Zampini   PetscRandom_CURAND *curand = (PetscRandom_CURAND*)r->data;
235162e2cfSBarry Smith   size_t             nn = n < 0 ? (size_t)(-2*n) : n; /* handle complex case */
24808ba619SStefano Zampini 
25808ba619SStefano Zampini   PetscFunctionBegin;
26808ba619SStefano Zampini #if defined(PETSC_USE_REAL_SINGLE)
279566063dSJacob Faibussowitsch   PetscCallCURAND(curandGenerateUniform(curand->gen,val,nn));
28808ba619SStefano Zampini #else
299566063dSJacob Faibussowitsch   PetscCallCURAND(curandGenerateUniformDouble(curand->gen,val,nn));
30808ba619SStefano Zampini #endif
31808ba619SStefano Zampini   if (r->iset) {
329566063dSJacob Faibussowitsch     PetscCall(PetscRandomCurandScale_Private(r,nn,val,(PetscBool)(n<0)));
33808ba619SStefano Zampini   }
34808ba619SStefano Zampini   PetscFunctionReturn(0);
35808ba619SStefano Zampini }
36808ba619SStefano Zampini 
37808ba619SStefano Zampini PetscErrorCode PetscRandomGetValues_CURAND(PetscRandom r, PetscInt n, PetscScalar *val)
38808ba619SStefano Zampini {
39808ba619SStefano Zampini   PetscFunctionBegin;
40808ba619SStefano Zampini #if defined(PETSC_USE_COMPLEX)
41808ba619SStefano Zampini   /* pass negative size to flag complex scaling (if needed) */
429566063dSJacob Faibussowitsch   PetscCall(PetscRandomGetValuesReal_CURAND(r,-n,(PetscReal*)val));
43808ba619SStefano Zampini #else
449566063dSJacob Faibussowitsch   PetscCall(PetscRandomGetValuesReal_CURAND(r,n,val));
45808ba619SStefano Zampini #endif
46808ba619SStefano Zampini   PetscFunctionReturn(0);
47808ba619SStefano Zampini }
48808ba619SStefano Zampini 
49808ba619SStefano Zampini PetscErrorCode PetscRandomDestroy_CURAND(PetscRandom r)
50808ba619SStefano Zampini {
51808ba619SStefano Zampini   PetscRandom_CURAND *curand = (PetscRandom_CURAND*)r->data;
52808ba619SStefano Zampini 
53808ba619SStefano Zampini   PetscFunctionBegin;
549566063dSJacob Faibussowitsch   PetscCallCURAND(curandDestroyGenerator(curand->gen));
559566063dSJacob Faibussowitsch   PetscCall(PetscFree(r->data));
56808ba619SStefano Zampini   PetscFunctionReturn(0);
57808ba619SStefano Zampini }
58808ba619SStefano Zampini 
59808ba619SStefano Zampini static struct _PetscRandomOps PetscRandomOps_Values = {
60267267bdSJacob Faibussowitsch   PetscDesignatedInitializer(seed,PetscRandomSeed_CURAND),
61267267bdSJacob Faibussowitsch   PetscDesignatedInitializer(getvalue,NULL),
62267267bdSJacob Faibussowitsch   PetscDesignatedInitializer(getvaluereal,NULL),
63267267bdSJacob Faibussowitsch   PetscDesignatedInitializer(getvalues,PetscRandomGetValues_CURAND),
64267267bdSJacob Faibussowitsch   PetscDesignatedInitializer(getvaluesreal,PetscRandomGetValuesReal_CURAND),
65267267bdSJacob Faibussowitsch   PetscDesignatedInitializer(destroy,PetscRandomDestroy_CURAND),
66808ba619SStefano Zampini };
67808ba619SStefano Zampini 
68808ba619SStefano Zampini /*MC
69808ba619SStefano Zampini    PETSCCURAND - access to the CUDA random number generator
70808ba619SStefano Zampini 
71*1179163eSBarry Smith   PETSc must be ./configure with the option --with-cuda to use this random number generator.
72*1179163eSBarry Smith 
73808ba619SStefano Zampini   Level: beginner
74808ba619SStefano Zampini 
75db781477SPatrick Sanan .seealso: `PetscRandomCreate()`, `PetscRandomSetType()`
76808ba619SStefano Zampini M*/
77808ba619SStefano Zampini 
78808ba619SStefano Zampini PETSC_EXTERN PetscErrorCode PetscRandomCreate_CURAND(PetscRandom r)
79808ba619SStefano Zampini {
80808ba619SStefano Zampini   PetscRandom_CURAND *curand;
81808ba619SStefano Zampini 
82808ba619SStefano Zampini   PetscFunctionBegin;
839566063dSJacob Faibussowitsch   PetscCall(PetscDeviceInitialize(PETSC_DEVICE_CUDA));
849566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(r,&curand));
859566063dSJacob Faibussowitsch   PetscCallCURAND(curandCreateGenerator(&curand->gen,CURAND_RNG_PSEUDO_DEFAULT));
86808ba619SStefano Zampini   /* https://docs.nvidia.com/cuda/curand/host-api-overview.html#performance-notes2 */
879566063dSJacob Faibussowitsch   PetscCallCURAND(curandSetGeneratorOrdering(curand->gen,CURAND_ORDERING_PSEUDO_SEEDED));
889566063dSJacob Faibussowitsch   PetscCall(PetscMemcpy(r->ops,&PetscRandomOps_Values,sizeof(PetscRandomOps_Values)));
899566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)r,PETSCCURAND));
90808ba619SStefano Zampini   r->data = curand;
91808ba619SStefano Zampini   r->seed = 1234ULL; /* taken from example */
929566063dSJacob Faibussowitsch   PetscCall(PetscRandomSeed_CURAND(r));
93808ba619SStefano Zampini   PetscFunctionReturn(0);
94808ba619SStefano Zampini }
95