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 a 5X5 example of the 5-pt. stencil,\n\ use the file petsc/src/mat/examples/matbinary.ex\n\ -nd : > 0 no of domains per processor \n\ -ov : >=0 amount of overlap between domains\n\n"; #include "petscksp.h" #undef __FUNCT__ #define __FUNCT__ "main" int main(int argc,char **args) { PetscInt nd = 2,ov=1,i,j,m,n,*idx,lsize; PetscErrorCode ierr; PetscMPIInt rank; PetscTruth flg; Mat A,B; char file[PETSC_MAX_PATH_LEN]; PetscViewer fd; IS *is1,*is2; PetscRandom r; PetscScalar rand; PetscInitialize(&argc,&args,(char *)0,help); #if defined(PETSC_USE_COMPLEX) SETERRQ(1,"This example does not work with complex numbers"); #else ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr); ierr = PetscOptionsGetString(PETSC_NULL,"-f",file,PETSC_MAX_PATH_LEN-1,PETSC_NULL);CHKERRQ(ierr); ierr = PetscOptionsGetInt(PETSC_NULL,"-nd",&nd,PETSC_NULL);CHKERRQ(ierr); ierr = PetscOptionsGetInt(PETSC_NULL,"-ov",&ov,PETSC_NULL);CHKERRQ(ierr); /* Read matrix and RHS */ ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd);CHKERRQ(ierr); ierr = MatLoad(fd,MATMPIAIJ,&A);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 = MatLoad(fd,MATSEQAIJ,&B);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 = PetscMalloc(nd*sizeof(IS **),&is1);CHKERRQ(ierr); ierr = PetscMalloc(nd*sizeof(IS **),&is2);CHKERRQ(ierr); ierr = PetscMalloc(m *sizeof(PetscInt),&idx);CHKERRQ(ierr); /* Create the random Index Sets */ for (i=0; i