Actual source code: ex151.c

petsc-3.8.4 2018-03-24
Report Typos and Errors
  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:  #include <petscmat.h>

 14: int main(int argc,char **argv)
 15: {
 16:   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.}};
 17:   const PetscInt ixrow[5]                               = {4,2,1,0,3},ixcol[7] = {5,3,6,1,2,0,4};
 18:   Mat            A,B;
 20:   PetscInt       i,rstart,rend,cstart,cend;
 21:   IS             isrow,iscol;
 22:   PetscViewer    viewer;
 23:   PetscBool      view_sparse;

 25:   PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
 26:   /* ------- Assemble matrix, --------- */
 27:   MatCreate(PETSC_COMM_WORLD,&A);
 28:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,5,7);
 29:   MatSetFromOptions(A);
 30:   MatSetUp(A);
 31:   MatGetOwnershipRange(A,&rstart,&rend);
 32:   MatGetOwnershipRangeColumn(A,&cstart,&cend);

 34:   for (i=0; i<(PetscInt)(sizeof(entries)/sizeof(entries[0])); i++) {
 35:     MatSetValue(A,entries[i].i,entries[i].j,entries[i].v,INSERT_VALUES);
 36:   }
 37:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 38:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

 40:   /* ------ Prepare index sets ------ */
 41:   ISCreateGeneral(PETSC_COMM_WORLD,rend-rstart,ixrow+rstart,PETSC_USE_POINTER,&isrow);
 42:   ISCreateGeneral(PETSC_COMM_WORLD,cend-cstart,ixcol+cstart,PETSC_USE_POINTER,&iscol);
 43:   ISSetPermutation(isrow);
 44:   ISSetPermutation(iscol);

 46:   PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
 47:   view_sparse = PETSC_FALSE;
 48:   PetscOptionsGetBool(NULL,NULL, "-view_sparse", &view_sparse, NULL);
 49:   if (!view_sparse) {
 50:     PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_DENSE);
 51:   }
 52:   PetscViewerASCIIPrintf(viewer,"Original matrix\n");
 53:   MatView(A,viewer);

 55:   MatPermute(A,isrow,iscol,&B);
 56:   PetscViewerASCIIPrintf(viewer,"Permuted matrix\n");
 57:   MatView(B,viewer);

 59:   if (!view_sparse) {
 60:     PetscViewerPopFormat(viewer);
 61:   }
 62:   PetscViewerASCIIPrintf(viewer,"Row permutation\n");
 63:   ISView(isrow,viewer);
 64:   PetscViewerASCIIPrintf(viewer,"Column permutation\n");
 65:   ISView(iscol,viewer);

 67:   /* Free data structures */
 68:   ISDestroy(&isrow);
 69:   ISDestroy(&iscol);
 70:   MatDestroy(&A);
 71:   MatDestroy(&B);

 73:   PetscFinalize();
 74:   return ierr;
 75: }