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