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 PetscCall(PetscInitialize(&argc,&argv,NULL,help)); 15 PetscCall(MatCreateDenseCUDA(comm,global_size,global_size,global_size,global_size,NULL,&cuda_matrix)); 16 17 PetscCall(VecCreateSeqCUDA(comm,global_size,&input)); 18 PetscCall(VecDuplicate(input,&output)); 19 PetscCall(VecSet(input,1.)); 20 PetscCall(VecSet(output,2.)); 21 PetscCall(MatMult(cuda_matrix,input,output)); 22 PetscCall(VecNorm(output,NORM_2,&nrm)); 23 PetscCheck(nrm <= PETSC_SMALL,PETSC_COMM_SELF,PETSC_ERR_PLIB,"PETSc generated wrong result. Should be 0, but is %g",(double)nrm); 24 PetscCall(VecDestroy(&input)); 25 PetscCall(VecDestroy(&output)); 26 PetscCall(MatDestroy(&cuda_matrix)); 27 PetscCall(PetscFinalize()); 28 return 0; 29 } 30 31 /*TEST 32 build: 33 requires: cuda 34 35 test: 36 nsize: 1 37 38 TEST*/ 39