Actual source code: ex40.c

  1: /*$Id: ex40.c,v 1.25 2001/08/07 03:03:07 balay Exp $*/

  3: static char help[] = "Tests the parallel case for MatIncreaseOverlap(). Input arguments are:\n\
  4:   -f <input_file> : file to load.  For a 5X5 example of the 5-pt. stencil,\n\
  5:                        use the file petsc/src/mat/examples/matbinary.ex\n\
  6:   -nd <size>      : > 0  number of domains per processor \n\
  7:   -ov <overlap>   : >=0  amount of overlap between domains\n\n";

 9:  #include petscsles.h

 13: int main(int argc,char **args)
 14: {
 15:   int         ierr,nd = 2,ov=1,i,size,start,m,n,end,rank;
 16:   PetscTruth  flg;
 17:   Mat         A,B;
 18:   char        file[128];
 19:   PetscViewer fd;
 20:   IS          *is1,*is2;
 21:   PetscRandom r;
 22:   PetscScalar rand;

 24:   PetscInitialize(&argc,&args,(char *)0,help);
 25: #if defined(PETSC_USE_COMPLEX)
 26:   SETERRQ(1,"This example does not work with complex numbers");
 27: #else
 28: 
 29:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
 30:   PetscOptionsGetString(PETSC_NULL,"-f",file,127,PETSC_NULL);
 31:   PetscOptionsGetInt(PETSC_NULL,"-nd",&nd,PETSC_NULL);
 32:   PetscOptionsGetInt(PETSC_NULL,"-ov",&ov,PETSC_NULL);

 34:   /* Read matrix and RHS */
 35:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,PETSC_BINARY_RDONLY,&fd);
 36:   MatLoad(fd,MATMPIAIJ,&A);
 37:   PetscViewerDestroy(fd);

 39:   /* Read the matrix again as a sequential matrix */
 40:   PetscViewerBinaryOpen(PETSC_COMM_SELF,file,PETSC_BINARY_RDONLY,&fd);
 41:   MatLoad(fd,MATSEQAIJ,&B);
 42:   PetscViewerDestroy(fd);
 43: 
 44:   /* Create the IS corresponding to subdomains */
 45:   PetscMalloc(nd*sizeof(IS **),&is1);
 46:   PetscMalloc(nd*sizeof(IS **),&is2);

 48:   /* Create the random Index Sets */
 49:   MatGetSize(A,&m,&n);
 50:   PetscRandomCreate(PETSC_COMM_SELF,RANDOM_DEFAULT,&r);
 51:   for (i=0; i<nd; i++) {
 52:     PetscRandomGetValue(r,&rand);
 53:     start = (int)(rand*m);
 54:     PetscRandomGetValue(r,&rand);
 55:     end  = (int)(rand*m);
 56:     size =  end - start;
 57:     if (start > end) { start = end; size = -size ;}
 58:     ISCreateStride(PETSC_COMM_SELF,size,start,1,is1+i);
 59:     ISCreateStride(PETSC_COMM_SELF,size,start,1,is2+i);
 60:   }
 61:   MatIncreaseOverlap(A,nd,is1,ov);
 62:   MatIncreaseOverlap(B,nd,is2,ov);

 64:   /* Now see if the serial and parallel case have the same answers */
 65:   for (i=0; i<nd; ++i) {
 66:     ISEqual(is1[i],is2[i],&flg);
 67:     PetscPrintf(PETSC_COMM_SELF,"proc:[%d], i=%d, flg =%d\n",rank,i,flg);
 68:   }

 70:   /* Free allocated memory */
 71:   for (i=0; i<nd; ++i) {
 72:     ISDestroy(is1[i]);
 73:     ISDestroy(is2[i]);
 74:   }
 75:   PetscFree(is1);
 76:   PetscFree(is2);
 77:   PetscRandomDestroy(r);
 78:   MatDestroy(A);
 79:   MatDestroy(B);

 81:   PetscFinalize();
 82: #endif
 83:   return 0;
 84: }