Actual source code: ex59.c

  2: static char help[] = "Tests MatGetSubmatrix() in parallel.";

 4:  #include petscmat.h

  8: int main(int argc,char **args)
  9: {
 10:   Mat            C,A;
 11:   PetscInt       i,j,m = 3,n = 2,rstart,rend;
 12:   PetscMPIInt    size,rank;
 14:   PetscScalar    v;
 15:   IS             isrow,iscol;
 16:   PetscTruth     flg;
 17:   char           type[256];

 19:   PetscInitialize(&argc,&args,(char *)0,help);
 20:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
 21:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
 22:   n = 2*size;

 24:   PetscStrcpy(type,MATSAME);
 25:   PetscOptionsGetString(PETSC_NULL,"-mat_type",type,256,PETSC_NULL);

 27:   PetscStrcmp(type,MATMPIDENSE,&flg);
 28:   if (flg) {
 29:     MatCreateMPIDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,
 30:            m*n,m*n,PETSC_NULL,&C);
 31:   } else {
 32:     MatCreateMPIAIJ(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,
 33:            m*n,m*n,5,PETSC_NULL,5,PETSC_NULL,&C);
 34:   }

 36:   /*
 37:         This is JUST to generate a nice test matrix, all processors fill up
 38:     the entire matrix. This is not something one would ever do in practice.
 39:   */
 40:   for (i=0; i<m*n; i++) {
 41:     for (j=0; j<m*n; j++) {
 42:       v = i + j + 1;
 43:       MatSetValues(C,1,&i,1,&j,&v,INSERT_VALUES);
 44:     }
 45:   }
 46:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
 47:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
 48:   PetscViewerSetFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_COMMON);
 49:   MatView(C,PETSC_VIEWER_STDOUT_WORLD);

 51:   /* 
 52:      Generate a new matrix consisting of every second row and column of
 53:    the original matrix
 54:   */
 55:   MatGetOwnershipRange(C,&rstart,&rend);
 56:   /* list the rows we want on THIS processor */
 57:   ISCreateStride(PETSC_COMM_WORLD,(rend-rstart)/2,rstart,2,&isrow);
 58:   /* list ALL the columns we want */
 59:   ISCreateStride(PETSC_COMM_WORLD,(m*n)/2,0,2,&iscol);
 60:   MatGetSubMatrix(C,isrow,iscol,PETSC_DECIDE,MAT_INITIAL_MATRIX,&A);
 61:   MatView(A,PETSC_VIEWER_STDOUT_WORLD);

 63:   MatGetSubMatrix(C,isrow,iscol,PETSC_DECIDE,MAT_REUSE_MATRIX,&A);
 64:   MatView(A,PETSC_VIEWER_STDOUT_WORLD);

 66:   ISDestroy(isrow);
 67:   ISDestroy(iscol);
 68:   MatDestroy(A);
 69:   MatDestroy(C);
 70:   PetscFinalize();
 71:   return 0;
 72: }