xref: /petsc/src/mat/tests/ex151.c (revision ffa8c5705e8ab2cf85ee1d14dbe507a6e2eb5283)
1 static char help[] = "Tests MatPermute() in parallel.\n\n";
2 /* Results:
3    Sequential:
4    - seqaij:   correct permutation
5    - seqbaij:  permutation not supported for this MATTYPE
6    - seqsbaij: permutation not supported for this MATTYPE
7    Parallel:
8    - mpiaij:   correct permutation
9    - mpibaij:  correct permutation
10    - mpisbaij: permutation not supported for this MATTYPE
11  */
12 
13 #include <petscmat.h>
14 
15 int main(int argc,char **argv)
16 {
17   const struct {PetscInt i,j; PetscScalar v;} entries[] = {{0,3,1.},{1,2,2.},{2,1,3.},{2,5,4.},{3,0,5.},{3,6,6.},{4,1,7.},{4,4,8.}};
18   const PetscInt ixrow[5]                               = {4,2,1,0,3},ixcol[7] = {5,3,6,1,2,0,4};
19   Mat            A,B;
20   PetscInt       i,rstart,rend,cstart,cend;
21   IS             isrow,iscol;
22   PetscViewer    viewer;
23   PetscBool      view_sparse;
24 
25   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
26   /* ------- Assemble matrix, --------- */
27   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
28   PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,5,7));
29   PetscCall(MatSetFromOptions(A));
30   PetscCall(MatSetUp(A));
31   PetscCall(MatGetOwnershipRange(A,&rstart,&rend));
32   PetscCall(MatGetOwnershipRangeColumn(A,&cstart,&cend));
33 
34   for (i=0; i<(PetscInt)(sizeof(entries)/sizeof(entries[0])); i++) {
35     PetscCall(MatSetValue(A,entries[i].i,entries[i].j,entries[i].v,INSERT_VALUES));
36   }
37   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
38   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
39 
40   /* ------ Prepare index sets ------ */
41   PetscCall(ISCreateGeneral(PETSC_COMM_WORLD,rend-rstart,ixrow+rstart,PETSC_USE_POINTER,&isrow));
42   PetscCall(ISCreateGeneral(PETSC_COMM_WORLD,cend-cstart,ixcol+cstart,PETSC_USE_POINTER,&iscol));
43   PetscCall(ISSetPermutation(isrow));
44   PetscCall(ISSetPermutation(iscol));
45 
46   PetscCall(PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer));
47   view_sparse = PETSC_FALSE;
48   PetscCall(PetscOptionsGetBool(NULL,NULL, "-view_sparse", &view_sparse, NULL));
49   if (!view_sparse) {
50     PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_DENSE));
51   }
52   PetscCall(PetscViewerASCIIPrintf(viewer,"Original matrix\n"));
53   PetscCall(MatView(A,viewer));
54 
55   PetscCall(MatPermute(A,isrow,iscol,&B));
56   PetscCall(PetscViewerASCIIPrintf(viewer,"Permuted matrix\n"));
57   PetscCall(MatView(B,viewer));
58 
59   if (!view_sparse) {
60     PetscCall(PetscViewerPopFormat(viewer));
61   }
62   PetscCall(PetscViewerASCIIPrintf(viewer,"Row permutation\n"));
63   PetscCall(ISView(isrow,viewer));
64   PetscCall(PetscViewerASCIIPrintf(viewer,"Column permutation\n"));
65   PetscCall(ISView(iscol,viewer));
66 
67   /* Free data structures */
68   PetscCall(ISDestroy(&isrow));
69   PetscCall(ISDestroy(&iscol));
70   PetscCall(MatDestroy(&A));
71   PetscCall(MatDestroy(&B));
72 
73   PetscCall(PetscFinalize());
74   return 0;
75 }
76 
77 /*TEST
78 
79    build:
80       requires: !complex
81 
82    test:
83       args: -view_sparse
84 
85    test:
86       suffix: 2
87       nsize: 2
88       args: -view_sparse
89 
90    test:
91       suffix: 2b
92       nsize: 2
93       args: -mat_type baij -view_sparse
94 
95    test:
96       suffix: 3
97       nsize: 3
98       args: -view_sparse
99 
100    test:
101       suffix: 3b
102       nsize: 3
103       args: -mat_type baij -view_sparse
104 
105    test:
106       suffix: dense
107       args: -mat_type dense
108 
109 TEST*/
110