xref: /petsc/src/mat/tests/ex239.c (revision 08401ef684002a709c6d3db98a0c9f54a8bcf1ec)
1a78de62eSJunchao Zhang static char help[] = "Test device/host memory allocation in MatDenseSeqCUDA()\n\n";
2a78de62eSJunchao Zhang 
3a78de62eSJunchao Zhang /* Contributed by: Victor Eijkhout <eijkhout@tacc.utexas.edu> */
4a78de62eSJunchao Zhang 
5a78de62eSJunchao Zhang #include <petscmat.h>
6a78de62eSJunchao Zhang int main(int argc, char** argv)
7a78de62eSJunchao Zhang {
8a78de62eSJunchao Zhang   PetscInt  global_size = 100;
9a78de62eSJunchao Zhang   Mat       cuda_matrix;
10a78de62eSJunchao Zhang   Vec       input,output;
11a78de62eSJunchao Zhang   MPI_Comm  comm        = PETSC_COMM_SELF;
12a78de62eSJunchao Zhang   PetscReal nrm         = 1;
13a78de62eSJunchao Zhang 
149566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,NULL,help));
159566063dSJacob Faibussowitsch   PetscCall(MatCreateDenseCUDA(comm,global_size,global_size,global_size,global_size,NULL,&cuda_matrix));
169566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(cuda_matrix,MAT_FINAL_ASSEMBLY));
179566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(cuda_matrix,MAT_FINAL_ASSEMBLY));
18a78de62eSJunchao Zhang 
199566063dSJacob Faibussowitsch   PetscCall(VecCreateSeqCUDA(comm,global_size,&input));
209566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(input,&output));
219566063dSJacob Faibussowitsch   PetscCall(VecSet(input,1.));
229566063dSJacob Faibussowitsch   PetscCall(VecSet(output,2.));
239566063dSJacob Faibussowitsch   PetscCall(MatMult(cuda_matrix,input,output));
249566063dSJacob Faibussowitsch   PetscCall(VecNorm(output,NORM_2,&nrm));
25*08401ef6SPierre Jolivet   PetscCheck(nrm <= PETSC_SMALL,PETSC_COMM_SELF,PETSC_ERR_PLIB,"PETSc generated wrong result. Should be 0, but is %g",(double)nrm);
269566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&input));
279566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&output));
289566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&cuda_matrix));
299566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
30b122ec5aSJacob Faibussowitsch   return 0;
31a78de62eSJunchao Zhang }
32a78de62eSJunchao Zhang 
33a78de62eSJunchao Zhang /*TEST
34a78de62eSJunchao Zhang    build:
35a78de62eSJunchao Zhang      requires: cuda
36a78de62eSJunchao Zhang 
37a78de62eSJunchao Zhang    test:
38a78de62eSJunchao Zhang     nsize: 1
39a78de62eSJunchao Zhang 
40a78de62eSJunchao Zhang TEST*/
41