xref: /petsc/src/mat/tests/ex166.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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   CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help));
26   /* ------- Assemble matrix, --------- */
27   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A));
28   CHKERRQ(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,5,5));
29   CHKERRQ(MatSetFromOptions(A));
30   CHKERRQ(MatSetUp(A));
31   CHKERRQ(MatGetOwnershipRange(A,&rstart,&rend));
32   CHKERRQ(MatGetOwnershipRangeColumn(A,&cstart,&cend));
33 
34   for (i=0; i<(PetscInt)(sizeof(entries)/sizeof(entries[0])); i++) {
35     CHKERRQ(MatSetValue(A,entries[i].i,entries[i].j,entries[i].v,INSERT_VALUES));
36   }
37   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
38   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
39 
40   /* ------ Prepare index sets ------ */
41   CHKERRQ(ISCreateGeneral(PETSC_COMM_WORLD,rend-rstart,ixrow+rstart,PETSC_USE_POINTER,&isrow));
42   CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF,5,ixcol,PETSC_USE_POINTER,&iscol));
43   CHKERRQ(ISSetPermutation(isrow));
44   CHKERRQ(ISSetPermutation(iscol));
45 
46   CHKERRQ(PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer));
47   view_sparse = PETSC_FALSE;
48   CHKERRQ(PetscOptionsGetBool(NULL,NULL, "-view_sparse", &view_sparse, NULL));
49   if (!view_sparse) {
50     CHKERRQ(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_DENSE));
51   }
52   CHKERRQ(PetscViewerASCIIPrintf(viewer,"Original matrix\n"));
53   CHKERRQ(MatView(A,viewer));
54 
55   CHKERRQ(MatPermute(A,isrow,iscol,&B));
56   CHKERRQ(PetscViewerASCIIPrintf(viewer,"Permuted matrix\n"));
57   CHKERRQ(MatView(B,viewer));
58   if (!view_sparse) {
59     CHKERRQ(PetscViewerPopFormat(viewer));
60   }
61 
62   CHKERRQ(PetscViewerASCIIPrintf(viewer,"Row permutation\n"));
63   CHKERRQ(ISView(isrow,viewer));
64   CHKERRQ(PetscViewerASCIIPrintf(viewer,"Column permutation\n"));
65   CHKERRQ(PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer));
66   CHKERRQ(ISView(iscol,sviewer));
67   CHKERRQ(PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer));
68 
69   /* Free data structures */
70   CHKERRQ(ISDestroy(&isrow));
71   CHKERRQ(ISDestroy(&iscol));
72   CHKERRQ(MatDestroy(&A));
73   CHKERRQ(MatDestroy(&B));
74 
75   CHKERRQ(PetscFinalize());
76   return 0;
77 }
78