1 static char help[] = "Test MatSetValuesCOO for MPIAIJ and its subclasses \n\n"; 2 3 #include <petscmat.h> 4 int main(int argc,char **args) 5 { 6 PetscErrorCode ierr; 7 Mat A,B; 8 PetscInt k; 9 const PetscInt M = 18,N = 18; 10 PetscMPIInt rank,size; 11 PetscBool equal; 12 PetscScalar *vals; 13 PetscBool flg = PETSC_FALSE; 14 15 /* Construct 18 x 18 matrices, which are big enough to have complex communication patterns but still small enough for debugging */ 16 PetscInt i0[] = {7, 7, 8, 8, 9, 16, 17, 9, 10, 1, 1, -2, 2, 3, 3, 14, 4, 5, 10, 13, 9, 9, 10, 1, 0, 0, 5, 5, 6, 6, 13, 13, 14, -14, 4, 4, 5, 11, 11, 12, 15, 15, 16}; 17 PetscInt j0[] = {1, 6, 2, 4, 10, 15, 13, 16, 11, 2, 7, 3, 8, 4, 9, 13, 5, 2, 15, 14, 10, 16, 11, 2, 0, 1, 6, -11, 0, 7, 15, 17, 11, 13, 5, 8, 2, 12, 17, 13, 3, 16, 9}; 18 19 PetscInt i1[] = {8, 5, 15, 16, 6, 13, 4, 17, 8, 9, 9, 10, -6, 12, 7, 3, -4, 1, 1, 2, 5, 5, 6, 14, 17, 8, 9, 9, 10, 4, 5, 10, 11, 1, 2}; 20 PetscInt j1[] = {2, 3, 16, 9, 5, 17, 1, 13, 4, 10, 16, 11, -5, 12, 1, 7, -1, 2, 7, 3, 6, 11, 0, 11, 13, 4, 10, 16, 11, 8, -2, 15, 12, 7, 3}; 21 22 PetscInt i2[] = {3, 4, 1, 10, 0, 1, 1, 2, 1, 1, 2, 2, 3, 3, 4, 4, 1, 2, 5, 5, 6, 4, 17, 0, 1, 1, 8, 5, 5, 6, 4, 7, 8, 5}; 23 PetscInt j2[] = {7, 1, 2, 11, 5, 2, 7, 3, 2, 7, 3, 8, 4, 9, 3, 5, 7, 3, 6, 11, 0, 1, 13, 5, 2, 7, 4, 6, 11, 0, 1, 3, 4, 2}; 24 25 struct { 26 PetscInt *i,*j,n; 27 } coo[3] = {{i0,j0,sizeof(i0)/sizeof(PetscInt)}, {i1,j1,sizeof(i1)/sizeof(PetscInt)}, {i2,j2,sizeof(i2)/sizeof(PetscInt)}}; 28 29 ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 30 ierr = PetscOptionsGetBool(NULL,NULL,"-ignore_remote",&flg,NULL);CHKERRQ(ierr); 31 ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr); 32 ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr); 33 34 PetscCheckFalse(size > 3,PETSC_COMM_WORLD,PETSC_ERR_ARG_WRONG,"This test requires at most 3 processes"); 35 36 ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 37 ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,M,N); 38 ierr = MatSetType(A,MATAIJ);CHKERRQ(ierr); 39 ierr = MatSeqAIJSetPreallocation(A,2,NULL); 40 ierr = MatMPIAIJSetPreallocation(A,2,NULL,2,NULL);CHKERRQ(ierr); 41 ierr = MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr); 42 ierr = MatSetOption(A,MAT_IGNORE_OFF_PROC_ENTRIES,flg);CHKERRQ(ierr); 43 44 for (k=0; k<coo[rank].n; k++) { 45 PetscScalar val = coo[rank].j[k]; 46 ierr = MatSetValue(A,coo[rank].i[k],coo[rank].j[k],val,ADD_VALUES);CHKERRQ(ierr); 47 } 48 ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 49 ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 50 51 ierr = MatCreate(PETSC_COMM_WORLD,&B);CHKERRQ(ierr); 52 ierr = MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,M,N); 53 ierr = MatSetFromOptions(B);CHKERRQ(ierr); 54 ierr = MatSetOption(B,MAT_IGNORE_OFF_PROC_ENTRIES,flg);CHKERRQ(ierr); 55 ierr = MatSetPreallocationCOO(B,coo[rank].n,coo[rank].i,coo[rank].j);CHKERRQ(ierr); 56 57 ierr = PetscMalloc1(coo[rank].n,&vals);CHKERRQ(ierr); 58 for (k=0; k<coo[rank].n; k++) vals[k] = coo[rank].j[k]; 59 ierr = MatSetValuesCOO(B,vals,ADD_VALUES);CHKERRQ(ierr); 60 61 ierr = MatEqual(A,B,&equal);CHKERRQ(ierr); 62 63 if (!equal) { 64 ierr = MatView(A,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 65 ierr = MatView(B,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 66 SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"MatSetValuesCOO() failed"); 67 } 68 69 ierr = PetscFree(vals);CHKERRQ(ierr); 70 ierr = MatDestroy(&A);CHKERRQ(ierr); 71 ierr = MatDestroy(&B);CHKERRQ(ierr); 72 73 ierr = PetscFinalize(); 74 return ierr; 75 } 76 77 /*TEST 78 79 testset: 80 output_file: output/ex254_1.out 81 nsize: {{1 2 3}} 82 args: -ignore_remote {{0 1}} 83 84 test: 85 suffix: kokkos 86 requires: kokkos_kernels 87 args: -mat_type aijkokkos 88 89 test: 90 suffix: cuda 91 requires: cuda 92 args: -mat_type aijcusparse 93 94 test: 95 suffix: aij 96 args: -mat_type aij 97 98 TEST*/ 99 100