Actual source code: ex41.c
1: /*$Id: ex41.c,v 1.25 2001/08/07 03:03:07 balay Exp $*/
3: static char help[] = "Tests MatIncreaseOverlap() - the parallel case. This examplen
4: is similar to ex40.c; here the index sets used are random. Input arguments are:n
5: -f <input_file> : file to load. For a 5X5 example of the 5-pt. stencil,n
6: use the file petsc/src/mat/examples/matbinary.exn
7: -nd <size> : > 0 no of domains per processor n
8: -ov <overlap> : >=0 amount of overlap between domainsnn";
10: #include petscsles.h
12: int main(int argc,char **args)
13: {
14: int ierr,nd = 2,ov=1,i,j,size,m,n,rank,*idx;
15: PetscTruth flg;
16: Mat A,B;
17: char file[128];
18: PetscViewer fd;
19: IS *is1,*is2;
20: PetscRandom r;
21: PetscScalar rand;
23: PetscInitialize(&argc,&args,(char *)0,help);
24: #if defined(PETSC_USE_COMPLEX)
25: SETERRQ(1,"This example does not work with complex numbers");
26: #else
27:
28: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
29: PetscOptionsGetString(PETSC_NULL,"-f",file,127,PETSC_NULL);
30: PetscOptionsGetInt(PETSC_NULL,"-nd",&nd,PETSC_NULL);
31: PetscOptionsGetInt(PETSC_NULL,"-ov",&ov,PETSC_NULL);
33: /* Read matrix and RHS */
34: PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,PETSC_BINARY_RDONLY,&fd);
35: MatLoad(fd,MATMPIAIJ,&A);
36: PetscViewerDestroy(fd);
38: /* Read the matrix again as a seq matrix */
39: PetscViewerBinaryOpen(PETSC_COMM_SELF,file,PETSC_BINARY_RDONLY,&fd);
40: MatLoad(fd,MATSEQAIJ,&B);
41: PetscViewerDestroy(fd);
42:
43: /* Create the Random no generator */
44: MatGetSize(A,&m,&n);
45: PetscRandomCreate(PETSC_COMM_SELF,RANDOM_DEFAULT,&r);
46:
47: /* Create the IS corresponding to subdomains */
48: PetscMalloc(nd*sizeof(IS **),&is1);
49: PetscMalloc(nd*sizeof(IS **),&is2);
50: PetscMalloc(m *sizeof(int),&idx);
52: /* Create the random Index Sets */
53: for (i=0; i<nd; i++) {
54: for (j=0; j<rank; j++) {
55: ierr = PetscRandomGetValue(r,&rand);
56: }
57: ierr = PetscRandomGetValue(r,&rand);
58: size = (int)(rand*m);
59: for (j=0; j<size; j++) {
60: ierr = PetscRandomGetValue(r,&rand);
61: idx[j] = (int)(rand*m);
62: }
63: ISCreateGeneral(PETSC_COMM_SELF,size,idx,is1+i);
64: ISCreateGeneral(PETSC_COMM_SELF,size,idx,is2+i);
65: }
66:
67: MatIncreaseOverlap(A,nd,is1,ov);
68: MatIncreaseOverlap(B,nd,is2,ov);
69:
70: /* Now see if the serial and parallel case have the same answers */
71: for (i=0; i<nd; ++i) {
72: int sz1,sz2;
73: ISEqual(is1[i],is2[i],&flg);
74: ISGetSize(is1[i],&sz1);
75: ISGetSize(is2[i],&sz2);
76: PetscPrintf(PETSC_COMM_SELF,"proc:[%d], i=%d, flg =%d sz1 = %d sz2 = %dn",rank,i,flg,sz1,sz2);
77: /* ISView(is1[i],PETSC_VIEWER_STDOUT_SELF);
78: ISView(is2[i],PETSC_VIEWER_STDOUT_SELF); */
79: }
81: /* Free Allocated Memory */
82: for (i=0; i<nd; ++i) {
83: ISDestroy(is1[i]);
84: ISDestroy(is2[i]);
85: }
86: PetscRandomDestroy(r);
87: PetscFree(is1);
88: PetscFree(is2);
89: MatDestroy(A);
90: MatDestroy(B);
91: PetscFree(idx);
93: PetscFinalize();
94: #endif
95: return 0;
96: }