xref: /petsc/src/sys/classes/random/impls/curand/curand.c (revision 4dfa11a44d5adf2389f1d3acbc8f3c1116dc6c3a)
1a4af0ceeSJacob Faibussowitsch #include <petsc/private/deviceimpl.h>
2d6cc7855SJacob Faibussowitsch #include <petsc/private/randomimpl.h>
30e6b6b59SJacob Faibussowitsch #include <petscdevice_cuda.h>
4808ba619SStefano Zampini #include <curand.h>
5808ba619SStefano Zampini 
6808ba619SStefano Zampini typedef struct {
7808ba619SStefano Zampini   curandGenerator_t gen;
8808ba619SStefano Zampini } PetscRandom_CURAND;
9808ba619SStefano Zampini 
109371c9d4SSatish Balay PetscErrorCode PetscRandomSeed_CURAND(PetscRandom r) {
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 
209371c9d4SSatish Balay PetscErrorCode PetscRandomGetValuesReal_CURAND(PetscRandom r, PetscInt n, PetscReal *val) {
21808ba619SStefano Zampini   PetscRandom_CURAND *curand = (PetscRandom_CURAND *)r->data;
225162e2cfSBarry Smith   size_t              nn     = n < 0 ? (size_t)(-2 * n) : n; /* handle complex case */
23808ba619SStefano Zampini 
24808ba619SStefano Zampini   PetscFunctionBegin;
25808ba619SStefano Zampini #if defined(PETSC_USE_REAL_SINGLE)
269566063dSJacob Faibussowitsch   PetscCallCURAND(curandGenerateUniform(curand->gen, val, nn));
27808ba619SStefano Zampini #else
289566063dSJacob Faibussowitsch   PetscCallCURAND(curandGenerateUniformDouble(curand->gen, val, nn));
29808ba619SStefano Zampini #endif
3048a46eb9SPierre Jolivet   if (r->iset) PetscCall(PetscRandomCurandScale_Private(r, nn, val, (PetscBool)(n < 0)));
31808ba619SStefano Zampini   PetscFunctionReturn(0);
32808ba619SStefano Zampini }
33808ba619SStefano Zampini 
349371c9d4SSatish Balay PetscErrorCode PetscRandomGetValues_CURAND(PetscRandom r, PetscInt n, PetscScalar *val) {
35808ba619SStefano Zampini   PetscFunctionBegin;
36808ba619SStefano Zampini #if defined(PETSC_USE_COMPLEX)
37808ba619SStefano Zampini   /* pass negative size to flag complex scaling (if needed) */
389566063dSJacob Faibussowitsch   PetscCall(PetscRandomGetValuesReal_CURAND(r, -n, (PetscReal *)val));
39808ba619SStefano Zampini #else
409566063dSJacob Faibussowitsch   PetscCall(PetscRandomGetValuesReal_CURAND(r, n, val));
41808ba619SStefano Zampini #endif
42808ba619SStefano Zampini   PetscFunctionReturn(0);
43808ba619SStefano Zampini }
44808ba619SStefano Zampini 
459371c9d4SSatish Balay PetscErrorCode PetscRandomDestroy_CURAND(PetscRandom r) {
46808ba619SStefano Zampini   PetscRandom_CURAND *curand = (PetscRandom_CURAND *)r->data;
47808ba619SStefano Zampini 
48808ba619SStefano Zampini   PetscFunctionBegin;
499566063dSJacob Faibussowitsch   PetscCallCURAND(curandDestroyGenerator(curand->gen));
509566063dSJacob Faibussowitsch   PetscCall(PetscFree(r->data));
51808ba619SStefano Zampini   PetscFunctionReturn(0);
52808ba619SStefano Zampini }
53808ba619SStefano Zampini 
54808ba619SStefano Zampini static struct _PetscRandomOps PetscRandomOps_Values = {
55267267bdSJacob Faibussowitsch   PetscDesignatedInitializer(seed, PetscRandomSeed_CURAND),
56267267bdSJacob Faibussowitsch   PetscDesignatedInitializer(getvalue, NULL),
57267267bdSJacob Faibussowitsch   PetscDesignatedInitializer(getvaluereal, NULL),
58267267bdSJacob Faibussowitsch   PetscDesignatedInitializer(getvalues, PetscRandomGetValues_CURAND),
59267267bdSJacob Faibussowitsch   PetscDesignatedInitializer(getvaluesreal, PetscRandomGetValuesReal_CURAND),
60267267bdSJacob Faibussowitsch   PetscDesignatedInitializer(destroy, PetscRandomDestroy_CURAND),
61808ba619SStefano Zampini };
62808ba619SStefano Zampini 
63808ba619SStefano Zampini /*MC
64811af0c4SBarry Smith    PETSCCURAND - access to the CUDA random number generator from a `PetscRandom` object
65808ba619SStefano Zampini 
661179163eSBarry Smith   PETSc must be ./configure with the option --with-cuda to use this random number generator.
671179163eSBarry Smith 
68808ba619SStefano Zampini   Level: beginner
69808ba619SStefano Zampini 
70811af0c4SBarry Smith .seealso: `PetscRandomCreate()`, `PetscRandomSetType()`, `PetscRandomType`
71808ba619SStefano Zampini M*/
72808ba619SStefano Zampini 
739371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode PetscRandomCreate_CURAND(PetscRandom r) {
74808ba619SStefano Zampini   PetscRandom_CURAND *curand;
75808ba619SStefano Zampini 
76808ba619SStefano Zampini   PetscFunctionBegin;
779566063dSJacob Faibussowitsch   PetscCall(PetscDeviceInitialize(PETSC_DEVICE_CUDA));
78*4dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&curand));
799566063dSJacob Faibussowitsch   PetscCallCURAND(curandCreateGenerator(&curand->gen, CURAND_RNG_PSEUDO_DEFAULT));
80808ba619SStefano Zampini   /* https://docs.nvidia.com/cuda/curand/host-api-overview.html#performance-notes2 */
819566063dSJacob Faibussowitsch   PetscCallCURAND(curandSetGeneratorOrdering(curand->gen, CURAND_ORDERING_PSEUDO_SEEDED));
829566063dSJacob Faibussowitsch   PetscCall(PetscMemcpy(r->ops, &PetscRandomOps_Values, sizeof(PetscRandomOps_Values)));
839566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)r, PETSCCURAND));
84808ba619SStefano Zampini   r->data = curand;
85808ba619SStefano Zampini   r->seed = 1234ULL; /* taken from example */
869566063dSJacob Faibussowitsch   PetscCall(PetscRandomSeed_CURAND(r));
87808ba619SStefano Zampini   PetscFunctionReturn(0);
88808ba619SStefano Zampini }
89