xref: /petsc/src/mat/tests/ex92.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
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 
8*d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
9*d71ae5a4SJacob Faibussowitsch {
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 
22327415f7SBarry Smith   PetscFunctionBeginUser;
239566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
249566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
259566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
26c4762a1bSJed Brown 
279566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-mat_block_size", &bs, NULL));
289566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-mat_mbs", &mbs, NULL));
299566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-ov", &ov, NULL));
309566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-nd", &nd, NULL));
319566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-view_id", &vid, NULL));
329566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-test_overlap", &TestOverlap));
339566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-test_submat", &TestSubMat));
349566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-test_allcols", &TestAllcols));
359566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_sorted", &test_sorted, NULL));
36c4762a1bSJed Brown 
379566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
389566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A, mbs * bs, mbs * bs, PETSC_DECIDE, PETSC_DECIDE));
399566063dSJacob Faibussowitsch   PetscCall(MatSetType(A, MATBAIJ));
409566063dSJacob Faibussowitsch   PetscCall(MatSeqBAIJSetPreallocation(A, bs, PETSC_DEFAULT, NULL));
419566063dSJacob Faibussowitsch   PetscCall(MatMPIBAIJSetPreallocation(A, bs, PETSC_DEFAULT, NULL, PETSC_DEFAULT, NULL));
42c4762a1bSJed Brown 
439566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rand));
449566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rand));
45c4762a1bSJed Brown 
469566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(A, &rstart, &rend));
479566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A, &M, &N));
48c4762a1bSJed Brown   Mbs = M / bs;
49c4762a1bSJed Brown 
509566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(bs, &rows));
519566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(bs, &cols));
529566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(bs * bs, &vals));
539566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(M, &idx));
54c4762a1bSJed Brown 
55c4762a1bSJed Brown   /* Now set blocks of values */
56c4762a1bSJed Brown   for (j = 0; j < bs * bs; j++) vals[j] = 0.0;
57c4762a1bSJed Brown   for (i = 0; i < Mbs; i++) {
589371c9d4SSatish Balay     cols[0] = i * bs;
599371c9d4SSatish Balay     rows[0] = i * bs;
60c4762a1bSJed Brown     for (j = 1; j < bs; j++) {
61c4762a1bSJed Brown       rows[j] = rows[j - 1] + 1;
62c4762a1bSJed Brown       cols[j] = cols[j - 1] + 1;
63c4762a1bSJed Brown     }
649566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A, bs, rows, bs, cols, vals, ADD_VALUES));
65c4762a1bSJed Brown   }
66c4762a1bSJed Brown   /* second, add random blocks */
67c4762a1bSJed Brown   for (i = 0; i < 20 * bs; i++) {
689566063dSJacob Faibussowitsch     PetscCall(PetscRandomGetValue(rand, &rval));
69c4762a1bSJed Brown     cols[0] = bs * (PetscInt)(PetscRealPart(rval) * Mbs);
709566063dSJacob Faibussowitsch     PetscCall(PetscRandomGetValue(rand, &rval));
71c4762a1bSJed Brown     rows[0] = rstart + bs * (PetscInt)(PetscRealPart(rval) * mbs);
72c4762a1bSJed Brown     for (j = 1; j < bs; j++) {
73c4762a1bSJed Brown       rows[j] = rows[j - 1] + 1;
74c4762a1bSJed Brown       cols[j] = cols[j - 1] + 1;
75c4762a1bSJed Brown     }
76c4762a1bSJed Brown 
77c4762a1bSJed Brown     for (j = 0; j < bs * bs; j++) {
789566063dSJacob Faibussowitsch       PetscCall(PetscRandomGetValue(rand, &rval));
79c4762a1bSJed Brown       vals[j] = rval;
80c4762a1bSJed Brown     }
819566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A, bs, rows, bs, cols, vals, ADD_VALUES));
82c4762a1bSJed Brown   }
83c4762a1bSJed Brown 
849566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
859566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
86c4762a1bSJed Brown 
87c4762a1bSJed Brown   /* make A a symmetric matrix: A <- A^T + A */
889566063dSJacob Faibussowitsch   PetscCall(MatTranspose(A, MAT_INITIAL_MATRIX, &Atrans));
899566063dSJacob Faibussowitsch   PetscCall(MatAXPY(A, one, Atrans, DIFFERENT_NONZERO_PATTERN));
909566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Atrans));
919566063dSJacob Faibussowitsch   PetscCall(MatTranspose(A, MAT_INITIAL_MATRIX, &Atrans));
929566063dSJacob Faibussowitsch   PetscCall(MatEqual(A, Atrans, &flg));
93c4762a1bSJed Brown   if (flg) {
949566063dSJacob Faibussowitsch     PetscCall(MatSetOption(A, MAT_SYMMETRIC, PETSC_TRUE));
95c4762a1bSJed Brown   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "A+A^T is non-symmetric");
969566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Atrans));
97c4762a1bSJed Brown 
98c4762a1bSJed Brown   /* create a SeqSBAIJ matrix sA (= A) */
999566063dSJacob Faibussowitsch   PetscCall(MatConvert(A, MATSBAIJ, MAT_INITIAL_MATRIX, &sA));
100c4762a1bSJed Brown   if (vid >= 0 && vid < size) {
1019566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "A:\n"));
1029566063dSJacob Faibussowitsch     PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD));
1039566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "sA:\n"));
1049566063dSJacob Faibussowitsch     PetscCall(MatView(sA, PETSC_VIEWER_STDOUT_WORLD));
105c4762a1bSJed Brown   }
106c4762a1bSJed Brown 
107c4762a1bSJed Brown   /* Test sA==A through MatMult() */
1089566063dSJacob Faibussowitsch   PetscCall(MatMultEqual(A, sA, 10, &flg));
10928b400f6SJacob Faibussowitsch   PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Error in MatConvert(): A != sA");
110c4762a1bSJed Brown 
111c4762a1bSJed Brown   /* Test MatIncreaseOverlap() */
1129566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nd, &is1));
1139566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nd, &is2));
114c4762a1bSJed Brown 
115c4762a1bSJed Brown   for (i = 0; i < nd; i++) {
116c4762a1bSJed Brown     if (!TestAllcols) {
1179566063dSJacob Faibussowitsch       PetscCall(PetscRandomGetValue(rand, &rval));
118c4762a1bSJed Brown       sz = (PetscInt)((0.5 + 0.2 * PetscRealPart(rval)) * mbs); /* 0.5*mbs < sz < 0.7*mbs */
119c4762a1bSJed Brown 
120c4762a1bSJed Brown       for (j = 0; j < sz; j++) {
1219566063dSJacob Faibussowitsch         PetscCall(PetscRandomGetValue(rand, &rval));
122c4762a1bSJed Brown         idx[j * bs] = bs * (PetscInt)(PetscRealPart(rval) * Mbs);
123c4762a1bSJed Brown         for (k = 1; k < bs; k++) idx[j * bs + k] = idx[j * bs] + k;
124c4762a1bSJed Brown       }
1259566063dSJacob Faibussowitsch       PetscCall(ISCreateGeneral(PETSC_COMM_SELF, sz * bs, idx, PETSC_COPY_VALUES, is1 + i));
1269566063dSJacob Faibussowitsch       PetscCall(ISCreateGeneral(PETSC_COMM_SELF, sz * bs, idx, PETSC_COPY_VALUES, is2 + i));
127c4762a1bSJed Brown       if (rank == vid) {
1289566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, " [%d] IS sz[%" PetscInt_FMT "]: %" PetscInt_FMT "\n", rank, i, sz));
1299566063dSJacob Faibussowitsch         PetscCall(ISView(is2[i], PETSC_VIEWER_STDOUT_SELF));
130c4762a1bSJed Brown       }
131a5b23f4aSJose E. Roman     } else { /* Test all rows and columns */
132c4762a1bSJed Brown       sz = M;
1339566063dSJacob Faibussowitsch       PetscCall(ISCreateStride(PETSC_COMM_SELF, sz, 0, 1, is1 + i));
1349566063dSJacob Faibussowitsch       PetscCall(ISCreateStride(PETSC_COMM_SELF, sz, 0, 1, is2 + i));
135c4762a1bSJed Brown 
136c4762a1bSJed Brown       if (rank == vid) {
137c4762a1bSJed Brown         PetscBool colflag;
1389566063dSJacob Faibussowitsch         PetscCall(ISIdentity(is2[i], &colflag));
1399566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "[%d] is2[%" PetscInt_FMT "], colflag %d\n", rank, i, colflag));
1409566063dSJacob Faibussowitsch         PetscCall(ISView(is2[i], PETSC_VIEWER_STDOUT_SELF));
141c4762a1bSJed Brown       }
142c4762a1bSJed Brown     }
143c4762a1bSJed Brown   }
144c4762a1bSJed Brown 
1459566063dSJacob Faibussowitsch   PetscCall(PetscLogStageRegister("MatOv_SBAIJ", &stages[0]));
1469566063dSJacob Faibussowitsch   PetscCall(PetscLogStageRegister("MatOv_BAIJ", &stages[1]));
147c4762a1bSJed Brown 
148c4762a1bSJed Brown   /* Test MatIncreaseOverlap */
149c4762a1bSJed Brown   if (TestOverlap) {
1509566063dSJacob Faibussowitsch     PetscCall(PetscLogStagePush(stages[0]));
1519566063dSJacob Faibussowitsch     PetscCall(MatIncreaseOverlap(sA, nd, is2, ov));
1529566063dSJacob Faibussowitsch     PetscCall(PetscLogStagePop());
153c4762a1bSJed Brown 
1549566063dSJacob Faibussowitsch     PetscCall(PetscLogStagePush(stages[1]));
1559566063dSJacob Faibussowitsch     PetscCall(MatIncreaseOverlap(A, nd, is1, ov));
1569566063dSJacob Faibussowitsch     PetscCall(PetscLogStagePop());
157c4762a1bSJed Brown 
158c4762a1bSJed Brown     if (rank == vid) {
1599566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n[%d] IS from BAIJ:\n", rank));
1609566063dSJacob Faibussowitsch       PetscCall(ISView(is1[0], PETSC_VIEWER_STDOUT_SELF));
1619566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n[%d] IS from SBAIJ:\n", rank));
1629566063dSJacob Faibussowitsch       PetscCall(ISView(is2[0], PETSC_VIEWER_STDOUT_SELF));
163c4762a1bSJed Brown     }
164c4762a1bSJed Brown 
165c4762a1bSJed Brown     for (i = 0; i < nd; ++i) {
1669566063dSJacob Faibussowitsch       PetscCall(ISEqual(is1[i], is2[i], &flg));
167c4762a1bSJed Brown       if (!flg) {
168dd400576SPatrick Sanan         if (rank == 0) {
1699566063dSJacob Faibussowitsch           PetscCall(ISSort(is1[i]));
1709566063dSJacob Faibussowitsch           PetscCall(ISSort(is2[i]));
171c4762a1bSJed Brown         }
17298921bdaSJacob Faibussowitsch         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "i=%" PetscInt_FMT ", is1 != is2", i);
173c4762a1bSJed Brown       }
174c4762a1bSJed Brown     }
175c4762a1bSJed Brown   }
176c4762a1bSJed Brown 
177c4762a1bSJed Brown   /* Test MatCreateSubmatrices */
178c4762a1bSJed Brown   if (TestSubMat) {
179c4762a1bSJed Brown     if (test_sorted) {
18048a46eb9SPierre Jolivet       for (i = 0; i < nd; ++i) PetscCall(ISSort(is1[i]));
181c4762a1bSJed Brown     }
1829566063dSJacob Faibussowitsch     PetscCall(MatCreateSubMatrices(A, nd, is1, is1, MAT_INITIAL_MATRIX, &submatA));
1839566063dSJacob Faibussowitsch     PetscCall(MatCreateSubMatrices(sA, nd, is1, is1, MAT_INITIAL_MATRIX, &submatsA));
184c4762a1bSJed Brown 
1859566063dSJacob Faibussowitsch     PetscCall(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 */
1899566063dSJacob Faibussowitsch     PetscCall(MatCreateSubMatrices(A, nd, is1, is1, MAT_REUSE_MATRIX, &submatA));
1909566063dSJacob Faibussowitsch     PetscCall(MatCreateSubMatrices(sA, nd, is1, is1, MAT_REUSE_MATRIX, &submatsA));
1919566063dSJacob Faibussowitsch     PetscCall(MatMultEqual(A, sA, 10, &flg));
19228b400f6SJacob Faibussowitsch     PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "MatCreateSubmatrices(): A != sA");
193c4762a1bSJed Brown 
1949566063dSJacob Faibussowitsch     PetscCall(MatDestroySubMatrices(nd, &submatA));
1959566063dSJacob Faibussowitsch     PetscCall(MatDestroySubMatrices(nd, &submatsA));
196c4762a1bSJed Brown   }
197c4762a1bSJed Brown 
198c4762a1bSJed Brown   /* Free allocated memory */
199c4762a1bSJed Brown   for (i = 0; i < nd; ++i) {
2009566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&is1[i]));
2019566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&is2[i]));
202c4762a1bSJed Brown   }
2039566063dSJacob Faibussowitsch   PetscCall(PetscFree(is1));
2049566063dSJacob Faibussowitsch   PetscCall(PetscFree(is2));
2059566063dSJacob Faibussowitsch   PetscCall(PetscFree(idx));
2069566063dSJacob Faibussowitsch   PetscCall(PetscFree(rows));
2079566063dSJacob Faibussowitsch   PetscCall(PetscFree(cols));
2089566063dSJacob Faibussowitsch   PetscCall(PetscFree(vals));
2099566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
2109566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&sA));
2119566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rand));
2129566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
213b122ec5aSJacob 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