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