static char help[] = "Tests MatIncreaseOverlap() - the parallel case. This example\n\ is similar to ex40.c; here the index sets used are random. 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,m,n,*idx,lsize; PetscErrorCode ierr; PetscMPIInt rank; PetscBool flg; Mat A,B; char file[PETSC_MAX_PATH_LEN]; PetscViewer fd; IS *is1,*is2; PetscRandom r; PetscScalar rand; ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr); ierr = PetscOptionsGetString(NULL,NULL,"-f",file,sizeof(file),NULL);CHKERRQ(ierr); ierr = PetscOptionsGetInt(NULL,NULL,"-nd",&nd,NULL);CHKERRQ(ierr); ierr = PetscOptionsGetInt(NULL,NULL,"-ov",&ov,NULL);CHKERRQ(ierr); /* Read matrix and RHS */ ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd);CHKERRQ(ierr); ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); ierr = MatSetType(A,MATMPIAIJ);CHKERRQ(ierr); ierr = MatLoad(A,fd);CHKERRQ(ierr); ierr = PetscViewerDestroy(&fd);CHKERRQ(ierr); /* Read the matrix again as a seq matrix */ ierr = PetscViewerBinaryOpen(PETSC_COMM_SELF,file,FILE_MODE_READ,&fd);CHKERRQ(ierr); ierr = MatCreate(PETSC_COMM_SELF,&B);CHKERRQ(ierr); ierr = MatSetType(B,MATSEQAIJ);CHKERRQ(ierr); ierr = MatLoad(B,fd);CHKERRQ(ierr); ierr = PetscViewerDestroy(&fd);CHKERRQ(ierr); /* Create the Random no generator */ ierr = MatGetSize(A,&m,&n);CHKERRQ(ierr); ierr = PetscRandomCreate(PETSC_COMM_SELF,&r);CHKERRQ(ierr); ierr = PetscRandomSetFromOptions(r);CHKERRQ(ierr); /* Create the IS corresponding to subdomains */ ierr = PetscMalloc1(nd,&is1);CHKERRQ(ierr); ierr = PetscMalloc1(nd,&is2);CHKERRQ(ierr); ierr = PetscMalloc1(m ,&idx);CHKERRQ(ierr); /* Create the random Index Sets */ for (i=0; i