Actual source code: ex87.c

petsc-3.6.4 2016-04-12
Report Typos and Errors
  2: static char help[] = "Tests MatGetSubMatrices() for SBAIJ matrices\n\n";

  4: #include <petscmat.h>

  8: int main(int argc,char **args)
  9: {
 10:   Mat            BAIJ,SBAIJ,*subBAIJ,*subSBAIJ;
 11:   PetscViewer    viewer;
 12:   char           file[PETSC_MAX_PATH_LEN];
 13:   PetscBool      flg;
 15:   PetscInt       n = 2,issize,M,N;
 16:   PetscMPIInt    rank;
 17:   IS             isrow,iscol,irow[n],icol[n];

 19:   PetscInitialize(&argc,&args,(char*)0,help);
 20:   PetscOptionsGetString(NULL,"-f",file,PETSC_MAX_PATH_LEN,&flg);
 21:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&viewer);

 23:   MatCreate(PETSC_COMM_WORLD,&BAIJ);
 24:   MatSetType(BAIJ,MATMPIBAIJ);
 25:   MatLoad(BAIJ,viewer);
 26:   PetscViewerDestroy(&viewer);

 28:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&viewer);
 29:   MatCreate(PETSC_COMM_WORLD,&SBAIJ);
 30:   MatSetType(SBAIJ,MATMPISBAIJ);
 31:   MatLoad(SBAIJ,viewer);
 32:   PetscViewerDestroy(&viewer);

 34:   MatGetSize(BAIJ,&M,&N);
 35:   issize = M/4;
 36:   ISCreateStride(PETSC_COMM_SELF,issize,0,1,&isrow);
 37:   irow[0] = irow[1] = isrow;
 38:   issize = N/2;
 39:   ISCreateStride(PETSC_COMM_SELF,issize,0,1,&iscol);
 40:   icol[0] = icol[1] = iscol;
 41:   MatGetSubMatrices(BAIJ,n,irow,icol,MAT_INITIAL_MATRIX,&subBAIJ);

 43:   /* irow and icol must be same for SBAIJ matrices! */
 44:   icol[0] = icol[1] = isrow;
 45:   MatGetSubMatrices(SBAIJ,n,irow,icol,MAT_INITIAL_MATRIX,&subSBAIJ);

 47:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
 48:   if (!rank) {
 49:     MatView(subBAIJ[0],PETSC_VIEWER_STDOUT_SELF);
 50:     MatView(subSBAIJ[0],PETSC_VIEWER_STDOUT_SELF);
 51:   }

 53:   /* Free data structures */
 54:   ISDestroy(&isrow);
 55:   ISDestroy(&iscol);
 56:   MatDestroyMatrices(n,&subBAIJ);
 57:   MatDestroyMatrices(n,&subSBAIJ);
 58:   MatDestroy(&BAIJ);
 59:   MatDestroy(&SBAIJ);

 61:   PetscFinalize();
 62:   return 0;
 63: }