static char help[] = "Tests MatIncreaseOverlap() and MatCreateSubmatrices() for the parallel case.\n\ This example is similar to ex40.c; here the index sets used are random.\n\ Input arguments are:\n\ -f : file to load. For example see $PETSC_DIR/share/petsc/datafiles/matrices\n\ -nd : > 0 no of domains per processor \n\ -ov : >=0 amount of overlap between domains\n\n"; #include int main(int argc,char **args) { PetscInt nd = 2,ov=1,i,j,lsize,m,n,*idx,bs; PetscMPIInt rank, size; PetscBool flg; Mat A,B,*submatA,*submatB; char file[PETSC_MAX_PATH_LEN]; PetscViewer fd; IS *is1,*is2; PetscRandom r; PetscBool test_unsorted = PETSC_FALSE; PetscScalar rand; CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help)); CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); CHKERRQ(PetscOptionsGetString(NULL,NULL,"-f",file,sizeof(file),NULL)); CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-nd",&nd,NULL)); CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-ov",&ov,NULL)); CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-test_unsorted",&test_unsorted,NULL)); /* Read matrix A and RHS */ CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd)); CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A)); CHKERRQ(MatSetType(A,MATAIJ)); CHKERRQ(MatSetFromOptions(A)); CHKERRQ(MatLoad(A,fd)); CHKERRQ(PetscViewerDestroy(&fd)); /* Read the same matrix as a seq matrix B */ CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_SELF,file,FILE_MODE_READ,&fd)); CHKERRQ(MatCreate(PETSC_COMM_SELF,&B)); CHKERRQ(MatSetType(B,MATSEQAIJ)); CHKERRQ(MatSetFromOptions(B)); CHKERRQ(MatLoad(B,fd)); CHKERRQ(PetscViewerDestroy(&fd)); CHKERRQ(MatGetBlockSize(A,&bs)); /* Create the Random no generator */ CHKERRQ(MatGetSize(A,&m,&n)); CHKERRQ(PetscRandomCreate(PETSC_COMM_SELF,&r)); CHKERRQ(PetscRandomSetFromOptions(r)); /* Create the IS corresponding to subdomains */ CHKERRQ(PetscMalloc1(nd,&is1)); CHKERRQ(PetscMalloc1(nd,&is2)); CHKERRQ(PetscMalloc1(m ,&idx)); for (i = 0; i < m; i++) {idx[i] = i;} /* Create the random Index Sets */ for (i=0; i