/*$Id: ex42.c,v 1.25 2001/08/07 03:03:07 balay Exp $*/ static char help[] = "Tests MatIncreaseOverlap() and MatGetSubmatrices() 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 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 "petscsles.h" #undef __FUNCT__ #define __FUNCT__ "main" int main(int argc,char **args) { int ierr,nd = 2,ov=1,i,j,size,m,n,rank,*idx; PetscTruth flg; Mat A,B,*submatA,*submatB; char file[128]; 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,127,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,PETSC_BINARY_RDONLY,&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,PETSC_BINARY_RDONLY,&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,RANDOM_DEFAULT,&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(int),&idx);CHKERRQ(ierr); /* Create the random Index Sets */ for (i=0; i