1 static char help[] = "Tests MatPermute() in parallel.\n\n"; 2 /* Results: 3 Sequential: 4 - seqaij: correct permutation 5 - seqbaij: permutation not supported for this MATTYPE 6 - seqsbaij: permutation not supported for this MATTYPE 7 Parallel: 8 - mpiaij: correct permutation 9 - mpibaij: correct permutation 10 - mpisbaij: permutation not supported for this MATTYPE 11 */ 12 13 #include <petscmat.h> 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,5,4.},{3,0,5.},{3,6,6.},{4,1,7.},{4,4,8.}}; 18 const PetscInt ixrow[5] = {4,2,1,0,3},ixcol[7] = {5,3,6,1,2,0,4}; 19 Mat A,B; 20 PetscInt i,rstart,rend,cstart,cend; 21 IS isrow,iscol; 22 PetscViewer viewer; 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,7)); 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_WORLD,cend-cstart,ixcol+cstart,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 60 if (!view_sparse) { 61 PetscCall(PetscViewerPopFormat(viewer)); 62 } 63 PetscCall(PetscViewerASCIIPrintf(viewer,"Row permutation\n")); 64 PetscCall(ISView(isrow,viewer)); 65 PetscCall(PetscViewerASCIIPrintf(viewer,"Column permutation\n")); 66 PetscCall(ISView(iscol,viewer)); 67 68 /* Free data structures */ 69 PetscCall(ISDestroy(&isrow)); 70 PetscCall(ISDestroy(&iscol)); 71 PetscCall(MatDestroy(&A)); 72 PetscCall(MatDestroy(&B)); 73 74 PetscCall(PetscFinalize()); 75 return 0; 76 } 77 78 /*TEST 79 80 build: 81 requires: !complex 82 83 test: 84 args: -view_sparse 85 86 test: 87 suffix: 2 88 nsize: 2 89 args: -view_sparse 90 91 test: 92 suffix: 2b 93 nsize: 2 94 args: -mat_type baij -view_sparse 95 96 test: 97 suffix: 3 98 nsize: 3 99 args: -view_sparse 100 101 test: 102 suffix: 3b 103 nsize: 3 104 args: -mat_type baij -view_sparse 105 106 test: 107 suffix: dense 108 args: -mat_type dense 109 110 TEST*/ 111