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