xref: /petsc/src/mat/tutorials/ex18cu.cu (revision 3307d110e72ee4e6d2468971620073eb5ff93529)
1 #include <petscdevice.h>
2 #include "ex18.h"
3 
4 __global__ void FillValues(PetscInt n, PetscScalar *v)
5 {
6   PetscInt     i = blockIdx.x * blockDim.x + threadIdx.x;
7   PetscScalar  *s;
8   if (i < n) {
9     s = &v[3*3*i];
10     for (PetscInt vi=0; vi<3; vi++) {
11       for (PetscInt vj=0; vj<3; vj++) {
12         s[vi*3+vj] = vi+2*vj;
13       }
14     }
15   }
16 }
17 
18 PetscErrorCode FillMatrixCUDACOO(FEStruct *fe,Mat A)
19 {
20   PetscScalar *v;
21 
22   PetscFunctionBeginUser;
23   PetscCallCUDA(cudaMalloc((void**)&v,3*3*fe->Ne*sizeof(PetscScalar)));
24   FillValues<<<(fe->Ne+255)/256,256>>>(fe->Ne,v);
25   PetscCall(MatSetValuesCOO(A,v,INSERT_VALUES));
26   PetscCallCUDA(cudaFree(v));
27   PetscFunctionReturn(0);
28 }
29