xref: /petsc/src/mat/tests/ex239.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1 static char help[] = "Test device/host memory allocation in MatDenseSeqCUDA()\n\n";
2 
3 /* Contributed by: Victor Eijkhout <eijkhout@tacc.utexas.edu> */
4 
5 #include <petscmat.h>
6 int main(int argc, char** argv)
7 {
8   PetscInt  global_size = 100;
9   Mat       cuda_matrix;
10   Vec       input,output;
11   MPI_Comm  comm        = PETSC_COMM_SELF;
12   PetscReal nrm         = 1;
13 
14   CHKERRQ(PetscInitialize(&argc,&argv,NULL,help));
15   CHKERRQ(MatCreateDenseCUDA(comm,global_size,global_size,global_size,global_size,NULL,&cuda_matrix));
16   CHKERRQ(MatAssemblyBegin(cuda_matrix,MAT_FINAL_ASSEMBLY));
17   CHKERRQ(MatAssemblyEnd(cuda_matrix,MAT_FINAL_ASSEMBLY));
18 
19   CHKERRQ(VecCreateSeqCUDA(comm,global_size,&input));
20   CHKERRQ(VecDuplicate(input,&output));
21   CHKERRQ(VecSet(input,1.));
22   CHKERRQ(VecSet(output,2.));
23   CHKERRQ(MatMult(cuda_matrix,input,output));
24   CHKERRQ(VecNorm(output,NORM_2,&nrm));
25   PetscCheckFalse(nrm > PETSC_SMALL,PETSC_COMM_SELF,PETSC_ERR_PLIB,"PETSc generated wrong result. Should be 0, but is %g",(double)nrm);
26   CHKERRQ(VecDestroy(&input));
27   CHKERRQ(VecDestroy(&output));
28   CHKERRQ(MatDestroy(&cuda_matrix));
29   CHKERRQ(PetscFinalize());
30   return 0;
31 }
32 
33 /*TEST
34    build:
35      requires: cuda
36 
37    test:
38     nsize: 1
39 
40 TEST*/
41