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.exn
  6:   -nd <size>      : > 0  number of domains per processor n
  7:   -ov <overlap>   : >=0  amount of overlap between domainsnn";

 9:  #include petscsles.h

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

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

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

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

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

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

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

 79:   PetscFinalize();
 80: #endif
 81:   return 0;
 82: }