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