xref: /petsc/src/mat/tests/ex166.c (revision ebead697dbf761eb322f829370bbe90b3bd93fa3)
1 static char help[] = "Tests MatPermute() for a square matrix in parallel.\n\n";
2 #include <petscmat.h>
3 /* Results:
4    Sequential:
5    - seqaij:   correct permutation
6    - seqbaij:  permutation not supported for this MATTYPE
7    - seqsbaij: permutation not supported for this MATTYPE
8    Parallel:
9    - mpiaij:   incorrect permutation (disable this op for now)
10    - seqbaij:  correct permutation, but all manner of memory leaks
11                (disable this op for now, because it appears broken for nonsquare matrices; see ex151)
12    - seqsbaij: permutation not supported for this MATTYPE
13 
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,4,4.},{3,0,5.},{3,3,6.},{4,1,7.},{4,4,8.}};
18   const PetscInt ixrow[5]                               = {4,2,1,3,0},ixcol[5] = {3,2,1,4,0};
19   Mat            A,B;
20   PetscInt       i,rstart,rend,cstart,cend;
21   IS             isrow,iscol;
22   PetscViewer    viewer,sviewer;
23   PetscBool      view_sparse;
24 
25   PetscFunctionBeginUser;
26   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
27   /* ------- Assemble matrix, --------- */
28   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
29   PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,5,5));
30   PetscCall(MatSetFromOptions(A));
31   PetscCall(MatSetUp(A));
32   PetscCall(MatGetOwnershipRange(A,&rstart,&rend));
33   PetscCall(MatGetOwnershipRangeColumn(A,&cstart,&cend));
34 
35   for (i=0; i<(PetscInt)PETSC_STATIC_ARRAY_LENGTH(entries); i++) {
36     PetscCall(MatSetValue(A,entries[i].i,entries[i].j,entries[i].v,INSERT_VALUES));
37   }
38   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
39   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
40 
41   /* ------ Prepare index sets ------ */
42   PetscCall(ISCreateGeneral(PETSC_COMM_WORLD,rend-rstart,ixrow+rstart,PETSC_USE_POINTER,&isrow));
43   PetscCall(ISCreateGeneral(PETSC_COMM_SELF,5,ixcol,PETSC_USE_POINTER,&iscol));
44   PetscCall(ISSetPermutation(isrow));
45   PetscCall(ISSetPermutation(iscol));
46 
47   PetscCall(PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer));
48   view_sparse = PETSC_FALSE;
49   PetscCall(PetscOptionsGetBool(NULL,NULL, "-view_sparse", &view_sparse, NULL));
50   if (!view_sparse) {
51     PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_DENSE));
52   }
53   PetscCall(PetscViewerASCIIPrintf(viewer,"Original matrix\n"));
54   PetscCall(MatView(A,viewer));
55 
56   PetscCall(MatPermute(A,isrow,iscol,&B));
57   PetscCall(PetscViewerASCIIPrintf(viewer,"Permuted matrix\n"));
58   PetscCall(MatView(B,viewer));
59   if (!view_sparse) {
60     PetscCall(PetscViewerPopFormat(viewer));
61   }
62 
63   PetscCall(PetscViewerASCIIPrintf(viewer,"Row permutation\n"));
64   PetscCall(ISView(isrow,viewer));
65   PetscCall(PetscViewerASCIIPrintf(viewer,"Column permutation\n"));
66   PetscCall(PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer));
67   PetscCall(ISView(iscol,sviewer));
68   PetscCall(PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer));
69 
70   /* Free data structures */
71   PetscCall(ISDestroy(&isrow));
72   PetscCall(ISDestroy(&iscol));
73   PetscCall(MatDestroy(&A));
74   PetscCall(MatDestroy(&B));
75 
76   PetscCall(PetscFinalize());
77   return 0;
78 }
79