xref: /petsc/src/mat/tests/ex239.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
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   PetscErrorCode ierr;
9a78de62eSJunchao Zhang   PetscInt       global_size=100;
10a78de62eSJunchao Zhang   Mat            cuda_matrix;
11a78de62eSJunchao Zhang   Vec            input,output;
12a78de62eSJunchao Zhang   MPI_Comm       comm = PETSC_COMM_SELF;
13a78de62eSJunchao Zhang   PetscReal      nrm = 1;
14a78de62eSJunchao Zhang 
15a78de62eSJunchao Zhang   ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr;
16*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateDenseCUDA(comm,global_size,global_size,global_size,global_size,NULL,&cuda_matrix));
17*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(cuda_matrix,MAT_FINAL_ASSEMBLY));
18*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(cuda_matrix,MAT_FINAL_ASSEMBLY));
19a78de62eSJunchao Zhang 
20*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateSeqCUDA(comm,global_size,&input));
21*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(input,&output));
22*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(input,1.));
23*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(output,2.));
24*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMult(cuda_matrix,input,output));
25*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecNorm(output,NORM_2,&nrm));
262c71b3e2SJacob Faibussowitsch   PetscCheckFalse(nrm > PETSC_SMALL,PETSC_COMM_SELF,PETSC_ERR_PLIB,"PETSc generated wrong result. Should be 0, but is %g",(double)nrm);
27*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&input));
28*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&output));
29*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&cuda_matrix));
30a78de62eSJunchao Zhang   ierr = PetscFinalize();
31a78de62eSJunchao Zhang   return ierr;
32a78de62eSJunchao Zhang }
33a78de62eSJunchao Zhang 
34a78de62eSJunchao Zhang /*TEST
35a78de62eSJunchao Zhang    build:
36a78de62eSJunchao Zhang      requires: cuda
37a78de62eSJunchao Zhang 
38a78de62eSJunchao Zhang    test:
39a78de62eSJunchao Zhang     nsize: 1
40a78de62eSJunchao Zhang 
41a78de62eSJunchao Zhang TEST*/
42