xref: /petsc/src/mat/tests/ex92.c (revision 28b400f66ebc7ae0049166a2294dfcd3df27e64b)
1 
2 static char help[] = "Tests MatIncreaseOverlap(), MatCreateSubMatrices() for parallel MatSBAIJ format.\n";
3 /* Example of usage:
4       mpiexec -n 2 ./ex92 -nd 2 -ov 3 -mat_block_size 2 -view_id 0 -test_overlap -test_submat
5 */
6 #include <petscmat.h>
7 
8 int main(int argc,char **args)
9 {
10   Mat            A,Atrans,sA,*submatA,*submatsA;
11   PetscErrorCode ierr;
12   PetscMPIInt    size,rank;
13   PetscInt       bs=1,mbs=10,ov=1,i,j,k,*rows,*cols,nd=2,*idx,rstart,rend,sz,M,N,Mbs;
14   PetscScalar    *vals,rval,one=1.0;
15   IS             *is1,*is2;
16   PetscRandom    rand;
17   PetscBool      flg,TestOverlap,TestSubMat,TestAllcols,test_sorted=PETSC_FALSE;
18   PetscInt       vid = -1;
19 #if defined(PETSC_USE_LOG)
20   PetscLogStage  stages[2];
21 #endif
22 
23   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
24   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
25   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
26 
27   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-mat_block_size",&bs,NULL));
28   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-mat_mbs",&mbs,NULL));
29   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-ov",&ov,NULL));
30   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-nd",&nd,NULL));
31   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-view_id",&vid,NULL));
32   CHKERRQ(PetscOptionsHasName(NULL,NULL, "-test_overlap", &TestOverlap));
33   CHKERRQ(PetscOptionsHasName(NULL,NULL, "-test_submat", &TestSubMat));
34   CHKERRQ(PetscOptionsHasName(NULL,NULL, "-test_allcols", &TestAllcols));
35   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-test_sorted",&test_sorted,NULL));
36 
37   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A));
38   CHKERRQ(MatSetSizes(A,mbs*bs,mbs*bs,PETSC_DECIDE,PETSC_DECIDE));
39   CHKERRQ(MatSetType(A,MATBAIJ));
40   CHKERRQ(MatSeqBAIJSetPreallocation(A,bs,PETSC_DEFAULT,NULL));
41   CHKERRQ(MatMPIBAIJSetPreallocation(A,bs,PETSC_DEFAULT,NULL,PETSC_DEFAULT,NULL));
42 
43   CHKERRQ(PetscRandomCreate(PETSC_COMM_WORLD,&rand));
44   CHKERRQ(PetscRandomSetFromOptions(rand));
45 
46   CHKERRQ(MatGetOwnershipRange(A,&rstart,&rend));
47   CHKERRQ(MatGetSize(A,&M,&N));
48   Mbs  = M/bs;
49 
50   CHKERRQ(PetscMalloc1(bs,&rows));
51   CHKERRQ(PetscMalloc1(bs,&cols));
52   CHKERRQ(PetscMalloc1(bs*bs,&vals));
53   CHKERRQ(PetscMalloc1(M,&idx));
54 
55   /* Now set blocks of values */
56   for (j=0; j<bs*bs; j++) vals[j] = 0.0;
57   for (i=0; i<Mbs; i++) {
58     cols[0] = i*bs; rows[0] = i*bs;
59     for (j=1; j<bs; j++) {
60       rows[j] = rows[j-1]+1;
61       cols[j] = cols[j-1]+1;
62     }
63     CHKERRQ(MatSetValues(A,bs,rows,bs,cols,vals,ADD_VALUES));
64   }
65   /* second, add random blocks */
66   for (i=0; i<20*bs; i++) {
67     CHKERRQ(PetscRandomGetValue(rand,&rval));
68     cols[0] = bs*(PetscInt)(PetscRealPart(rval)*Mbs);
69     CHKERRQ(PetscRandomGetValue(rand,&rval));
70     rows[0] = rstart + bs*(PetscInt)(PetscRealPart(rval)*mbs);
71     for (j=1; j<bs; j++) {
72       rows[j] = rows[j-1]+1;
73       cols[j] = cols[j-1]+1;
74     }
75 
76     for (j=0; j<bs*bs; j++) {
77       CHKERRQ(PetscRandomGetValue(rand,&rval));
78       vals[j] = rval;
79     }
80     CHKERRQ(MatSetValues(A,bs,rows,bs,cols,vals,ADD_VALUES));
81   }
82 
83   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
84   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
85 
86   /* make A a symmetric matrix: A <- A^T + A */
87   CHKERRQ(MatTranspose(A,MAT_INITIAL_MATRIX, &Atrans));
88   CHKERRQ(MatAXPY(A,one,Atrans,DIFFERENT_NONZERO_PATTERN));
89   CHKERRQ(MatDestroy(&Atrans));
90   CHKERRQ(MatTranspose(A,MAT_INITIAL_MATRIX, &Atrans));
91   CHKERRQ(MatEqual(A, Atrans, &flg));
92   if (flg) {
93     CHKERRQ(MatSetOption(A,MAT_SYMMETRIC,PETSC_TRUE));
94   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"A+A^T is non-symmetric");
95   CHKERRQ(MatDestroy(&Atrans));
96 
97   /* create a SeqSBAIJ matrix sA (= A) */
98   CHKERRQ(MatConvert(A,MATSBAIJ,MAT_INITIAL_MATRIX,&sA));
99   if (vid >= 0 && vid < size) {
100     CHKERRQ(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"A:\n"));
101     CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
102     CHKERRQ(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"sA:\n"));
103     CHKERRQ(MatView(sA,PETSC_VIEWER_STDOUT_WORLD));
104   }
105 
106   /* Test sA==A through MatMult() */
107   CHKERRQ(MatMultEqual(A,sA,10,&flg));
108   PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Error in MatConvert(): A != sA");
109 
110   /* Test MatIncreaseOverlap() */
111   CHKERRQ(PetscMalloc1(nd,&is1));
112   CHKERRQ(PetscMalloc1(nd,&is2));
113 
114   for (i=0; i<nd; i++) {
115     if (!TestAllcols) {
116       CHKERRQ(PetscRandomGetValue(rand,&rval));
117       sz   = (PetscInt)((0.5+0.2*PetscRealPart(rval))*mbs); /* 0.5*mbs < sz < 0.7*mbs */
118 
119       for (j=0; j<sz; j++) {
120         CHKERRQ(PetscRandomGetValue(rand,&rval));
121         idx[j*bs] = bs*(PetscInt)(PetscRealPart(rval)*Mbs);
122         for (k=1; k<bs; k++) idx[j*bs+k] = idx[j*bs]+k;
123       }
124       CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF,sz*bs,idx,PETSC_COPY_VALUES,is1+i));
125       CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF,sz*bs,idx,PETSC_COPY_VALUES,is2+i));
126       if (rank == vid) {
127         CHKERRQ(PetscPrintf(PETSC_COMM_SELF," [%d] IS sz[%" PetscInt_FMT "]: %" PetscInt_FMT "\n",rank,i,sz));
128         CHKERRQ(ISView(is2[i],PETSC_VIEWER_STDOUT_SELF));
129       }
130     } else { /* Test all rows and columns */
131       sz   = M;
132       CHKERRQ(ISCreateStride(PETSC_COMM_SELF,sz,0,1,is1+i));
133       CHKERRQ(ISCreateStride(PETSC_COMM_SELF,sz,0,1,is2+i));
134 
135       if (rank == vid) {
136         PetscBool colflag;
137         CHKERRQ(ISIdentity(is2[i],&colflag));
138         CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"[%d] is2[%" PetscInt_FMT "], colflag %d\n",rank,i,colflag));
139         CHKERRQ(ISView(is2[i],PETSC_VIEWER_STDOUT_SELF));
140       }
141     }
142   }
143 
144   CHKERRQ(PetscLogStageRegister("MatOv_SBAIJ",&stages[0]));
145   CHKERRQ(PetscLogStageRegister("MatOv_BAIJ",&stages[1]));
146 
147   /* Test MatIncreaseOverlap */
148   if (TestOverlap) {
149     CHKERRQ(PetscLogStagePush(stages[0]));
150     CHKERRQ(MatIncreaseOverlap(sA,nd,is2,ov));
151     CHKERRQ(PetscLogStagePop());
152 
153     CHKERRQ(PetscLogStagePush(stages[1]));
154     CHKERRQ(MatIncreaseOverlap(A,nd,is1,ov));
155     CHKERRQ(PetscLogStagePop());
156 
157     if (rank == vid) {
158       CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"\n[%d] IS from BAIJ:\n",rank));
159       CHKERRQ(ISView(is1[0],PETSC_VIEWER_STDOUT_SELF));
160       CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"\n[%d] IS from SBAIJ:\n",rank));
161       CHKERRQ(ISView(is2[0],PETSC_VIEWER_STDOUT_SELF));
162     }
163 
164     for (i=0; i<nd; ++i) {
165       CHKERRQ(ISEqual(is1[i],is2[i],&flg));
166       if (!flg) {
167         if (rank == 0) {
168           CHKERRQ(ISSort(is1[i]));
169           CHKERRQ(ISSort(is2[i]));
170         }
171         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"i=%" PetscInt_FMT ", is1 != is2",i);
172       }
173     }
174   }
175 
176   /* Test MatCreateSubmatrices */
177   if (TestSubMat) {
178     if (test_sorted) {
179       for (i = 0; i < nd; ++i) {
180         CHKERRQ(ISSort(is1[i]));
181       }
182     }
183     CHKERRQ(MatCreateSubMatrices(A,nd,is1,is1,MAT_INITIAL_MATRIX,&submatA));
184     CHKERRQ(MatCreateSubMatrices(sA,nd,is1,is1,MAT_INITIAL_MATRIX,&submatsA));
185 
186     CHKERRQ(MatMultEqual(A,sA,10,&flg));
187     PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"A != sA");
188 
189     /* Now test MatCreateSubmatrices with MAT_REUSE_MATRIX option */
190     CHKERRQ(MatCreateSubMatrices(A,nd,is1,is1,MAT_REUSE_MATRIX,&submatA));
191     CHKERRQ(MatCreateSubMatrices(sA,nd,is1,is1,MAT_REUSE_MATRIX,&submatsA));
192     CHKERRQ(MatMultEqual(A,sA,10,&flg));
193     PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"MatCreateSubmatrices(): A != sA");
194 
195     CHKERRQ(MatDestroySubMatrices(nd,&submatA));
196     CHKERRQ(MatDestroySubMatrices(nd,&submatsA));
197   }
198 
199   /* Free allocated memory */
200   for (i=0; i<nd; ++i) {
201     CHKERRQ(ISDestroy(&is1[i]));
202     CHKERRQ(ISDestroy(&is2[i]));
203   }
204   CHKERRQ(PetscFree(is1));
205   CHKERRQ(PetscFree(is2));
206   CHKERRQ(PetscFree(idx));
207   CHKERRQ(PetscFree(rows));
208   CHKERRQ(PetscFree(cols));
209   CHKERRQ(PetscFree(vals));
210   CHKERRQ(MatDestroy(&A));
211   CHKERRQ(MatDestroy(&sA));
212   CHKERRQ(PetscRandomDestroy(&rand));
213   ierr = PetscFinalize();
214   return ierr;
215 }
216 
217 /*TEST
218 
219    test:
220       args: -ov {{1 3}} -mat_block_size {{2 8}} -test_overlap -test_submat
221       output_file: output/ex92_1.out
222 
223    test:
224       suffix: 2
225       nsize: {{3 4}}
226       args: -ov {{1 3}} -mat_block_size {{2 8}} -test_overlap -test_submat
227       output_file: output/ex92_1.out
228 
229    test:
230       suffix: 3
231       nsize: {{3 4}}
232       args: -ov {{1 3}} -mat_block_size {{2 8}} -test_overlap -test_allcols
233       output_file: output/ex92_1.out
234 
235    test:
236       suffix: 3_sorted
237       nsize: {{3 4}}
238       args: -ov {{1 3}} -mat_block_size {{2 8}} -test_overlap -test_allcols -test_sorted
239       output_file: output/ex92_1.out
240 
241    test:
242       suffix: 4
243       nsize: {{3 4}}
244       args: -ov {{1 3}} -mat_block_size {{2 8}} -test_submat -test_allcols
245       output_file: output/ex92_1.out
246 
247 TEST*/
248