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