xref: /petsc/src/mat/tests/ex92.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests MatIncreaseOverlap(), MatCreateSubMatrices() for parallel MatSBAIJ format.\n";
3c4762a1bSJed Brown /* Example of usage:
4c4762a1bSJed Brown       mpiexec -n 2 ./ex92 -nd 2 -ov 3 -mat_block_size 2 -view_id 0 -test_overlap -test_submat
5c4762a1bSJed Brown */
6c4762a1bSJed Brown #include <petscmat.h>
7c4762a1bSJed Brown 
8c4762a1bSJed Brown int main(int argc,char **args)
9c4762a1bSJed Brown {
10c4762a1bSJed Brown   Mat            A,Atrans,sA,*submatA,*submatsA;
11c4762a1bSJed Brown   PetscMPIInt    size,rank;
12c4762a1bSJed Brown   PetscInt       bs=1,mbs=10,ov=1,i,j,k,*rows,*cols,nd=2,*idx,rstart,rend,sz,M,N,Mbs;
13c4762a1bSJed Brown   PetscScalar    *vals,rval,one=1.0;
14c4762a1bSJed Brown   IS             *is1,*is2;
15c4762a1bSJed Brown   PetscRandom    rand;
16c4762a1bSJed Brown   PetscBool      flg,TestOverlap,TestSubMat,TestAllcols,test_sorted=PETSC_FALSE;
17c4762a1bSJed Brown   PetscInt       vid = -1;
18c4762a1bSJed Brown #if defined(PETSC_USE_LOG)
19c4762a1bSJed Brown   PetscLogStage  stages[2];
20c4762a1bSJed Brown #endif
21c4762a1bSJed Brown 
22*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help));
235f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
245f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
25c4762a1bSJed Brown 
265f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-mat_block_size",&bs,NULL));
275f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-mat_mbs",&mbs,NULL));
285f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-ov",&ov,NULL));
295f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-nd",&nd,NULL));
305f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-view_id",&vid,NULL));
315f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHasName(NULL,NULL, "-test_overlap", &TestOverlap));
325f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHasName(NULL,NULL, "-test_submat", &TestSubMat));
335f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHasName(NULL,NULL, "-test_allcols", &TestAllcols));
345f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-test_sorted",&test_sorted,NULL));
35c4762a1bSJed Brown 
365f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A));
375f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(A,mbs*bs,mbs*bs,PETSC_DECIDE,PETSC_DECIDE));
385f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(A,MATBAIJ));
395f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSeqBAIJSetPreallocation(A,bs,PETSC_DEFAULT,NULL));
405f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMPIBAIJSetPreallocation(A,bs,PETSC_DEFAULT,NULL,PETSC_DEFAULT,NULL));
41c4762a1bSJed Brown 
425f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomCreate(PETSC_COMM_WORLD,&rand));
435f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomSetFromOptions(rand));
44c4762a1bSJed Brown 
455f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetOwnershipRange(A,&rstart,&rend));
465f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetSize(A,&M,&N));
47c4762a1bSJed Brown   Mbs  = M/bs;
48c4762a1bSJed Brown 
495f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(bs,&rows));
505f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(bs,&cols));
515f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(bs*bs,&vals));
525f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(M,&idx));
53c4762a1bSJed Brown 
54c4762a1bSJed Brown   /* Now set blocks of values */
55c4762a1bSJed Brown   for (j=0; j<bs*bs; j++) vals[j] = 0.0;
56c4762a1bSJed Brown   for (i=0; i<Mbs; i++) {
57c4762a1bSJed Brown     cols[0] = i*bs; rows[0] = i*bs;
58c4762a1bSJed Brown     for (j=1; j<bs; j++) {
59c4762a1bSJed Brown       rows[j] = rows[j-1]+1;
60c4762a1bSJed Brown       cols[j] = cols[j-1]+1;
61c4762a1bSJed Brown     }
625f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(A,bs,rows,bs,cols,vals,ADD_VALUES));
63c4762a1bSJed Brown   }
64c4762a1bSJed Brown   /* second, add random blocks */
65c4762a1bSJed Brown   for (i=0; i<20*bs; i++) {
665f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscRandomGetValue(rand,&rval));
67c4762a1bSJed Brown     cols[0] = bs*(PetscInt)(PetscRealPart(rval)*Mbs);
685f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscRandomGetValue(rand,&rval));
69c4762a1bSJed Brown     rows[0] = rstart + bs*(PetscInt)(PetscRealPart(rval)*mbs);
70c4762a1bSJed Brown     for (j=1; j<bs; j++) {
71c4762a1bSJed Brown       rows[j] = rows[j-1]+1;
72c4762a1bSJed Brown       cols[j] = cols[j-1]+1;
73c4762a1bSJed Brown     }
74c4762a1bSJed Brown 
75c4762a1bSJed Brown     for (j=0; j<bs*bs; j++) {
765f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscRandomGetValue(rand,&rval));
77c4762a1bSJed Brown       vals[j] = rval;
78c4762a1bSJed Brown     }
795f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(A,bs,rows,bs,cols,vals,ADD_VALUES));
80c4762a1bSJed Brown   }
81c4762a1bSJed Brown 
825f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
835f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
84c4762a1bSJed Brown 
85c4762a1bSJed Brown   /* make A a symmetric matrix: A <- A^T + A */
865f80ce2aSJacob Faibussowitsch   CHKERRQ(MatTranspose(A,MAT_INITIAL_MATRIX, &Atrans));
875f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAXPY(A,one,Atrans,DIFFERENT_NONZERO_PATTERN));
885f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&Atrans));
895f80ce2aSJacob Faibussowitsch   CHKERRQ(MatTranspose(A,MAT_INITIAL_MATRIX, &Atrans));
905f80ce2aSJacob Faibussowitsch   CHKERRQ(MatEqual(A, Atrans, &flg));
91c4762a1bSJed Brown   if (flg) {
925f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetOption(A,MAT_SYMMETRIC,PETSC_TRUE));
93c4762a1bSJed Brown   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"A+A^T is non-symmetric");
945f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&Atrans));
95c4762a1bSJed Brown 
96c4762a1bSJed Brown   /* create a SeqSBAIJ matrix sA (= A) */
975f80ce2aSJacob Faibussowitsch   CHKERRQ(MatConvert(A,MATSBAIJ,MAT_INITIAL_MATRIX,&sA));
98c4762a1bSJed Brown   if (vid >= 0 && vid < size) {
995f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"A:\n"));
1005f80ce2aSJacob Faibussowitsch     CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
1015f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"sA:\n"));
1025f80ce2aSJacob Faibussowitsch     CHKERRQ(MatView(sA,PETSC_VIEWER_STDOUT_WORLD));
103c4762a1bSJed Brown   }
104c4762a1bSJed Brown 
105c4762a1bSJed Brown   /* Test sA==A through MatMult() */
1065f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMultEqual(A,sA,10,&flg));
10728b400f6SJacob Faibussowitsch   PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Error in MatConvert(): A != sA");
108c4762a1bSJed Brown 
109c4762a1bSJed Brown   /* Test MatIncreaseOverlap() */
1105f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nd,&is1));
1115f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nd,&is2));
112c4762a1bSJed Brown 
113c4762a1bSJed Brown   for (i=0; i<nd; i++) {
114c4762a1bSJed Brown     if (!TestAllcols) {
1155f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscRandomGetValue(rand,&rval));
116c4762a1bSJed Brown       sz   = (PetscInt)((0.5+0.2*PetscRealPart(rval))*mbs); /* 0.5*mbs < sz < 0.7*mbs */
117c4762a1bSJed Brown 
118c4762a1bSJed Brown       for (j=0; j<sz; j++) {
1195f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscRandomGetValue(rand,&rval));
120c4762a1bSJed Brown         idx[j*bs] = bs*(PetscInt)(PetscRealPart(rval)*Mbs);
121c4762a1bSJed Brown         for (k=1; k<bs; k++) idx[j*bs+k] = idx[j*bs]+k;
122c4762a1bSJed Brown       }
1235f80ce2aSJacob Faibussowitsch       CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF,sz*bs,idx,PETSC_COPY_VALUES,is1+i));
1245f80ce2aSJacob Faibussowitsch       CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF,sz*bs,idx,PETSC_COPY_VALUES,is2+i));
125c4762a1bSJed Brown       if (rank == vid) {
1265f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscPrintf(PETSC_COMM_SELF," [%d] IS sz[%" PetscInt_FMT "]: %" PetscInt_FMT "\n",rank,i,sz));
1275f80ce2aSJacob Faibussowitsch         CHKERRQ(ISView(is2[i],PETSC_VIEWER_STDOUT_SELF));
128c4762a1bSJed Brown       }
129a5b23f4aSJose E. Roman     } else { /* Test all rows and columns */
130c4762a1bSJed Brown       sz   = M;
1315f80ce2aSJacob Faibussowitsch       CHKERRQ(ISCreateStride(PETSC_COMM_SELF,sz,0,1,is1+i));
1325f80ce2aSJacob Faibussowitsch       CHKERRQ(ISCreateStride(PETSC_COMM_SELF,sz,0,1,is2+i));
133c4762a1bSJed Brown 
134c4762a1bSJed Brown       if (rank == vid) {
135c4762a1bSJed Brown         PetscBool colflag;
1365f80ce2aSJacob Faibussowitsch         CHKERRQ(ISIdentity(is2[i],&colflag));
1375f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"[%d] is2[%" PetscInt_FMT "], colflag %d\n",rank,i,colflag));
1385f80ce2aSJacob Faibussowitsch         CHKERRQ(ISView(is2[i],PETSC_VIEWER_STDOUT_SELF));
139c4762a1bSJed Brown       }
140c4762a1bSJed Brown     }
141c4762a1bSJed Brown   }
142c4762a1bSJed Brown 
1435f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogStageRegister("MatOv_SBAIJ",&stages[0]));
1445f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogStageRegister("MatOv_BAIJ",&stages[1]));
145c4762a1bSJed Brown 
146c4762a1bSJed Brown   /* Test MatIncreaseOverlap */
147c4762a1bSJed Brown   if (TestOverlap) {
1485f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscLogStagePush(stages[0]));
1495f80ce2aSJacob Faibussowitsch     CHKERRQ(MatIncreaseOverlap(sA,nd,is2,ov));
1505f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscLogStagePop());
151c4762a1bSJed Brown 
1525f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscLogStagePush(stages[1]));
1535f80ce2aSJacob Faibussowitsch     CHKERRQ(MatIncreaseOverlap(A,nd,is1,ov));
1545f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscLogStagePop());
155c4762a1bSJed Brown 
156c4762a1bSJed Brown     if (rank == vid) {
1575f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"\n[%d] IS from BAIJ:\n",rank));
1585f80ce2aSJacob Faibussowitsch       CHKERRQ(ISView(is1[0],PETSC_VIEWER_STDOUT_SELF));
1595f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"\n[%d] IS from SBAIJ:\n",rank));
1605f80ce2aSJacob Faibussowitsch       CHKERRQ(ISView(is2[0],PETSC_VIEWER_STDOUT_SELF));
161c4762a1bSJed Brown     }
162c4762a1bSJed Brown 
163c4762a1bSJed Brown     for (i=0; i<nd; ++i) {
1645f80ce2aSJacob Faibussowitsch       CHKERRQ(ISEqual(is1[i],is2[i],&flg));
165c4762a1bSJed Brown       if (!flg) {
166dd400576SPatrick Sanan         if (rank == 0) {
1675f80ce2aSJacob Faibussowitsch           CHKERRQ(ISSort(is1[i]));
1685f80ce2aSJacob Faibussowitsch           CHKERRQ(ISSort(is2[i]));
169c4762a1bSJed Brown         }
17098921bdaSJacob Faibussowitsch         SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"i=%" PetscInt_FMT ", is1 != is2",i);
171c4762a1bSJed Brown       }
172c4762a1bSJed Brown     }
173c4762a1bSJed Brown   }
174c4762a1bSJed Brown 
175c4762a1bSJed Brown   /* Test MatCreateSubmatrices */
176c4762a1bSJed Brown   if (TestSubMat) {
177c4762a1bSJed Brown     if (test_sorted) {
178c4762a1bSJed Brown       for (i = 0; i < nd; ++i) {
1795f80ce2aSJacob Faibussowitsch         CHKERRQ(ISSort(is1[i]));
180c4762a1bSJed Brown       }
181c4762a1bSJed Brown     }
1825f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreateSubMatrices(A,nd,is1,is1,MAT_INITIAL_MATRIX,&submatA));
1835f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreateSubMatrices(sA,nd,is1,is1,MAT_INITIAL_MATRIX,&submatsA));
184c4762a1bSJed Brown 
1855f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMultEqual(A,sA,10,&flg));
18628b400f6SJacob Faibussowitsch     PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"A != sA");
187c4762a1bSJed Brown 
188c4762a1bSJed Brown     /* Now test MatCreateSubmatrices with MAT_REUSE_MATRIX option */
1895f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreateSubMatrices(A,nd,is1,is1,MAT_REUSE_MATRIX,&submatA));
1905f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreateSubMatrices(sA,nd,is1,is1,MAT_REUSE_MATRIX,&submatsA));
1915f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMultEqual(A,sA,10,&flg));
19228b400f6SJacob Faibussowitsch     PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"MatCreateSubmatrices(): A != sA");
193c4762a1bSJed Brown 
1945f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroySubMatrices(nd,&submatA));
1955f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroySubMatrices(nd,&submatsA));
196c4762a1bSJed Brown   }
197c4762a1bSJed Brown 
198c4762a1bSJed Brown   /* Free allocated memory */
199c4762a1bSJed Brown   for (i=0; i<nd; ++i) {
2005f80ce2aSJacob Faibussowitsch     CHKERRQ(ISDestroy(&is1[i]));
2015f80ce2aSJacob Faibussowitsch     CHKERRQ(ISDestroy(&is2[i]));
202c4762a1bSJed Brown   }
2035f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(is1));
2045f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(is2));
2055f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(idx));
2065f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(rows));
2075f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(cols));
2085f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(vals));
2095f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
2105f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&sA));
2115f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomDestroy(&rand));
212*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
213*b122ec5aSJacob Faibussowitsch   return 0;
214c4762a1bSJed Brown }
215c4762a1bSJed Brown 
216c4762a1bSJed Brown /*TEST
217c4762a1bSJed Brown 
218c4762a1bSJed Brown    test:
219c4762a1bSJed Brown       args: -ov {{1 3}} -mat_block_size {{2 8}} -test_overlap -test_submat
220c4762a1bSJed Brown       output_file: output/ex92_1.out
221c4762a1bSJed Brown 
222c4762a1bSJed Brown    test:
223c4762a1bSJed Brown       suffix: 2
224c4762a1bSJed Brown       nsize: {{3 4}}
225c4762a1bSJed Brown       args: -ov {{1 3}} -mat_block_size {{2 8}} -test_overlap -test_submat
226c4762a1bSJed Brown       output_file: output/ex92_1.out
227c4762a1bSJed Brown 
228c4762a1bSJed Brown    test:
229c4762a1bSJed Brown       suffix: 3
230c4762a1bSJed Brown       nsize: {{3 4}}
231c4762a1bSJed Brown       args: -ov {{1 3}} -mat_block_size {{2 8}} -test_overlap -test_allcols
232c4762a1bSJed Brown       output_file: output/ex92_1.out
233c4762a1bSJed Brown 
234c4762a1bSJed Brown    test:
235c4762a1bSJed Brown       suffix: 3_sorted
236c4762a1bSJed Brown       nsize: {{3 4}}
237c4762a1bSJed Brown       args: -ov {{1 3}} -mat_block_size {{2 8}} -test_overlap -test_allcols -test_sorted
238c4762a1bSJed Brown       output_file: output/ex92_1.out
239c4762a1bSJed Brown 
240c4762a1bSJed Brown    test:
241c4762a1bSJed Brown       suffix: 4
242c4762a1bSJed Brown       nsize: {{3 4}}
243c4762a1bSJed Brown       args: -ov {{1 3}} -mat_block_size {{2 8}} -test_submat -test_allcols
244c4762a1bSJed Brown       output_file: output/ex92_1.out
245c4762a1bSJed Brown 
246c4762a1bSJed Brown TEST*/
247