Actual source code: ex42.c
2: static char help[] = "Tests MatIncreaseOverlap() and MatGetSubmatrices() for the parallel case.\n\
3: This example is similar to ex40.c; here the index sets used are random.\n\
4: 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.ex\n\
7: -nd <size> : > 0 no of domains per processor \n\
8: -ov <overlap> : >=0 amount of overlap between domains\n\n";
10: #include petscksp.h
14: int main(int argc,char **args)
15: {
17: PetscInt nd = 2,ov=1,i,j,lsize,m,n,*idx;
18: PetscMPIInt rank;
19: PetscTruth flg;
20: Mat A,B,*submatA,*submatB;
21: char file[PETSC_MAX_PATH_LEN];
22: PetscViewer fd;
23: IS *is1,*is2;
24: PetscRandom r;
25: PetscScalar rand;
27: PetscInitialize(&argc,&args,(char *)0,help);
28: #if defined(PETSC_USE_COMPLEX)
29: SETERRQ(1,"This example does not work with complex numbers");
30: #else
31:
32: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
33: PetscOptionsGetString(PETSC_NULL,"-f",file,PETSC_MAX_PATH_LEN-1,PETSC_NULL);
34: PetscOptionsGetInt(PETSC_NULL,"-nd",&nd,PETSC_NULL);
35: PetscOptionsGetInt(PETSC_NULL,"-ov",&ov,PETSC_NULL);
37: /* Read matrix and RHS */
38: PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd);
39: MatLoad(fd,MATMPIAIJ,&A);
40: PetscViewerDestroy(fd);
42: /* Read the matrix again as a seq matrix */
43: PetscViewerBinaryOpen(PETSC_COMM_SELF,file,FILE_MODE_READ,&fd);
44: MatLoad(fd,MATSEQAIJ,&B);
45: PetscViewerDestroy(fd);
46:
47: /* Create the Random no generator */
48: MatGetSize(A,&m,&n);
49: PetscRandomCreate(PETSC_COMM_SELF,&r);
50: PetscRandomSetFromOptions(r);
52: /* Create the IS corresponding to subdomains */
53: PetscMalloc(nd*sizeof(IS **),&is1);
54: PetscMalloc(nd*sizeof(IS **),&is2);
55: PetscMalloc(m *sizeof(PetscInt),&idx);
56:
57: /* Create the random Index Sets */
58: for (i=0; i<nd; i++) {
59: /* Skip a few,so that the IS on different procs are diffeent*/
60: for (j=0; j<rank; j++) {
61: PetscRandomGetValue(r,&rand);
62: }
63: PetscRandomGetValue(r,&rand);
64: lsize = (PetscInt)(rand*m);
65: for (j=0; j<lsize; j++) {
66: PetscRandomGetValue(r,&rand);
67: idx[j] = (PetscInt)(rand*m);
68: }
69: PetscSortInt(lsize,idx);
70: ISCreateGeneral(PETSC_COMM_SELF,lsize,idx,is1+i);
71: ISCreateGeneral(PETSC_COMM_SELF,lsize,idx,is2+i);
72: }
74: MatIncreaseOverlap(A,nd,is1,ov);
75: MatIncreaseOverlap(B,nd,is2,ov);
77: for (i=0; i<nd; ++i) {
78: ISSort(is1[i]);
79: ISSort(is2[i]);
80: }
81:
82: MatGetSubMatrices(A,nd,is1,is1,MAT_INITIAL_MATRIX,&submatA);
83: MatGetSubMatrices(B,nd,is2,is2,MAT_INITIAL_MATRIX,&submatB);
84:
85: /* Now see if the serial and parallel case have the same answers */
86: for (i=0; i<nd; ++i) {
87: MatEqual(submatA[i],submatB[i],&flg);
88: PetscPrintf(PETSC_COMM_SELF,"proc:[%d], i=%D, flg =%d\n",rank,i,(int)flg);
89: }
91: /* Free Allocated Memory */
92: for (i=0; i<nd; ++i) {
93: ISDestroy(is1[i]);
94: ISDestroy(is2[i]);
95: MatDestroy(submatA[i]);
96: MatDestroy(submatB[i]);
97: }
98: PetscFree(submatA);
99: PetscFree(submatB);
100: PetscRandomDestroy(r);
101: PetscFree(is1);
102: PetscFree(is2);
103: MatDestroy(A);
104: MatDestroy(B);
105: PetscFree(idx);
107: PetscFinalize();
108: #endif
109: return 0;
110: }