Actual source code: ex87.c

petsc-3.4.5 2014-06-29
  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;
 16:   PetscMPIInt    rank;
 17:   IS             is,iss[2];

 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);
 22:   MatCreate(PETSC_COMM_WORLD,&BAIJ);
 23:   MatSetType(BAIJ,MATMPIBAIJ);
 24:   MatLoad(BAIJ,viewer);
 25:   PetscViewerDestroy(&viewer);
 26:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&viewer);
 27:   MatCreate(PETSC_COMM_WORLD,&SBAIJ);
 28:   MatSetType(SBAIJ,MATMPISBAIJ);
 29:   MatLoad(SBAIJ,viewer);
 30:   PetscViewerDestroy(&viewer);

 32:   MatGetSize(BAIJ,&issize,0);
 33:   issize = 9;
 34:   ISCreateStride(PETSC_COMM_SELF,issize,0,1,&is);
 35:   iss[0] = is;iss[1] = is;
 36:   MatGetSubMatrices(BAIJ,n,iss,iss,MAT_INITIAL_MATRIX,&subBAIJ);
 37:   MatGetSubMatrices(SBAIJ,n,iss,iss,MAT_INITIAL_MATRIX,&subSBAIJ);

 39:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
 40: #if defined(PETSC_USE_SOCKET_VIEWER)
 41:   if (!rank) {
 42:     MatView(subBAIJ[0],PETSC_VIEWER_SOCKET_SELF);
 43:     MatView(subSBAIJ[0],PETSC_VIEWER_SOCKET_SELF);
 44:   }
 45: #endif

 47:   /* Free data structures */
 48:   ISDestroy(&is);
 49:   MatDestroyMatrices(n,&subBAIJ);
 50:   MatDestroyMatrices(n,&subSBAIJ);
 51:   MatDestroy(&BAIJ);
 52:   MatDestroy(&SBAIJ);

 54:   PetscFinalize();
 55:   return 0;
 56: }