1c4762a1bSJed Brown static char help[] = "Example of extracting an array of MPI submatrices from a given MPI matrix.\n"
2c4762a1bSJed Brown "This test can only be run in parallel.\n"
3c4762a1bSJed Brown "\n";
4c4762a1bSJed Brown
5c4762a1bSJed Brown #include <petscmat.h>
6c4762a1bSJed Brown
MyISView(IS * rowis,IS * colis,PetscInt gs,PetscInt ss,PetscViewer viewer)7fe8fb074SBarry Smith PetscErrorCode MyISView(IS *rowis, IS *colis, PetscInt gs, PetscInt ss, PetscViewer viewer)
8fe8fb074SBarry Smith {
9fe8fb074SBarry Smith PetscViewer subviewer = NULL;
10fe8fb074SBarry Smith
11fe8fb074SBarry Smith PetscFunctionBeginUser;
12b4025f61SBarry Smith PetscCheck(ss <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "ss must be less than or equal to 1");
13fe8fb074SBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, "Row IS %" PetscInt_FMT "\n", gs));
14b4025f61SBarry Smith PetscCall(PetscViewerGetSubViewer(viewer, ss > -1 ? PetscObjectComm((PetscObject)rowis[ss]) : PETSC_COMM_SELF, &subviewer));
15b4025f61SBarry Smith if (ss > -1) PetscCall(ISView(rowis[ss], subviewer));
16b4025f61SBarry Smith PetscCall(PetscViewerRestoreSubViewer(viewer, ss > -1 ? PetscObjectComm((PetscObject)rowis[ss]) : PETSC_COMM_SELF, &subviewer));
17fe8fb074SBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, "Col IS %" PetscInt_FMT "\n", gs));
18b4025f61SBarry Smith PetscCall(PetscViewerGetSubViewer(viewer, ss > -1 ? PetscObjectComm((PetscObject)rowis[ss]) : PETSC_COMM_SELF, &subviewer));
19b4025f61SBarry Smith if (ss > -1) PetscCall(ISView(colis[ss], subviewer));
20b4025f61SBarry Smith PetscCall(PetscViewerRestoreSubViewer(viewer, ss > -1 ? PetscObjectComm((PetscObject)rowis[ss]) : PETSC_COMM_SELF, &subviewer));
21fe8fb074SBarry Smith PetscFunctionReturn(PETSC_SUCCESS);
22fe8fb074SBarry Smith }
23fe8fb074SBarry Smith
main(int argc,char ** args)24d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
25d71ae5a4SJacob Faibussowitsch {
26c4762a1bSJed Brown Mat A, *submats;
27c4762a1bSJed Brown MPI_Comm subcomm;
28c4762a1bSJed Brown PetscMPIInt rank, size, subrank, subsize, color;
29c4762a1bSJed Brown PetscInt m, n, N, bs, rstart, rend, i, j, k, total_subdomains, hash, nsubdomains = 1;
30c4762a1bSJed Brown PetscInt nis, *cols, gnsubdomains, gsubdomainnums[1], gsubdomainperm[1], s, gs;
31c4762a1bSJed Brown PetscInt *rowindices, *colindices, idx, rep;
32c4762a1bSJed Brown PetscScalar *vals;
33c4762a1bSJed Brown IS rowis[1], colis[1];
34c4762a1bSJed Brown PetscViewer viewer;
35c4762a1bSJed Brown PetscBool permute_indices, flg;
36c4762a1bSJed Brown
37327415f7SBarry Smith PetscFunctionBeginUser;
38*c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &args, NULL, help));
399566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
409566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
41c4762a1bSJed Brown
42d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "ex183", "Mat");
43c4762a1bSJed Brown m = 5;
449566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-m", "Local matrix size", "MatSetSizes", m, &m, &flg));
45c4762a1bSJed Brown total_subdomains = size - 1;
4652ce0ab5SPierre Jolivet PetscCall(PetscOptionsRangeInt("-total_subdomains", "Number of submatrices where 0 < n < comm size", "MatCreateSubMatricesMPI", total_subdomains, &total_subdomains, &flg, 1, size));
47c4762a1bSJed Brown permute_indices = PETSC_FALSE;
489566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-permute_indices", "Whether to permute indices before breaking them into subdomains", "ISCreateGeneral", permute_indices, &permute_indices, &flg));
49c4762a1bSJed Brown hash = 7;
509566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-hash", "Permutation factor, which has to be relatively prime to M = size*m (total matrix size)", "ISCreateGeneral", hash, &hash, &flg));
51c4762a1bSJed Brown rep = 2;
5252ce0ab5SPierre Jolivet PetscCall(PetscOptionsRangeInt("-rep", "Number of times to carry out submatrix extractions; currently only 1 & 2 are supported", NULL, rep, &rep, &flg, 1, 2));
53d0609cedSBarry Smith PetscOptionsEnd();
54c4762a1bSJed Brown
55c4762a1bSJed Brown viewer = PETSC_VIEWER_STDOUT_WORLD;
56c4762a1bSJed Brown /* Create logically sparse, but effectively dense matrix for easy verification of submatrix extraction correctness. */
579566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
589566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, m, m, PETSC_DECIDE, PETSC_DECIDE));
599566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A));
609566063dSJacob Faibussowitsch PetscCall(MatSetUp(A));
619566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, NULL, &N));
629566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A, NULL, &n));
639566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(A, &bs));
649566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(A, n, NULL));
659566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(A, n, NULL, N - n, NULL));
669566063dSJacob Faibussowitsch PetscCall(MatSeqBAIJSetPreallocation(A, bs, n / bs, NULL));
679566063dSJacob Faibussowitsch PetscCall(MatMPIBAIJSetPreallocation(A, bs, n / bs, NULL, (N - n) / bs, NULL));
689566063dSJacob Faibussowitsch PetscCall(MatSeqSBAIJSetPreallocation(A, bs, n / bs, NULL));
699566063dSJacob Faibussowitsch PetscCall(MatMPISBAIJSetPreallocation(A, bs, n / bs, NULL, (N - n) / bs, NULL));
70c4762a1bSJed Brown
719566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(N, &cols, N, &vals));
729566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A, &rstart, &rend));
73c4762a1bSJed Brown for (j = 0; j < N; ++j) cols[j] = j;
74c4762a1bSJed Brown for (i = rstart; i < rend; i++) {
75ad540459SPierre Jolivet for (j = 0; j < N; ++j) vals[j] = i * 10000 + j;
769566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &i, N, cols, vals, INSERT_VALUES));
77c4762a1bSJed Brown }
789566063dSJacob Faibussowitsch PetscCall(PetscFree2(cols, vals));
799566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
809566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
81c4762a1bSJed Brown
829566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Initial matrix:\n"));
839566063dSJacob Faibussowitsch PetscCall(MatView(A, viewer));
84c4762a1bSJed Brown
85c4762a1bSJed Brown /*
86c4762a1bSJed Brown Create subcomms and ISs so that each rank participates in one IS.
87c4762a1bSJed Brown The IS either coalesces adjacent rank indices (contiguous),
88c4762a1bSJed Brown or selects indices by scrambling them using a hash.
89c4762a1bSJed Brown */
90c4762a1bSJed Brown k = size / total_subdomains + (size % total_subdomains > 0); /* There are up to k ranks to a color */
91c4762a1bSJed Brown color = rank / k;
929566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_split(PETSC_COMM_WORLD, color, rank, &subcomm));
939566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(subcomm, &subsize));
949566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(subcomm, &subrank));
959566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A, &rstart, &rend));
96c4762a1bSJed Brown nis = 1;
979566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(rend - rstart, &rowindices, rend - rstart, &colindices));
98c4762a1bSJed Brown
99c4762a1bSJed Brown for (j = rstart; j < rend; ++j) {
100c4762a1bSJed Brown if (permute_indices) {
101c4762a1bSJed Brown idx = (j * hash);
102c4762a1bSJed Brown } else {
103c4762a1bSJed Brown idx = j;
104c4762a1bSJed Brown }
105c4762a1bSJed Brown rowindices[j - rstart] = idx % N;
106c4762a1bSJed Brown colindices[j - rstart] = (idx + m) % N;
107c4762a1bSJed Brown }
1089566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(subcomm, rend - rstart, rowindices, PETSC_COPY_VALUES, &rowis[0]));
1099566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(subcomm, rend - rstart, colindices, PETSC_COPY_VALUES, &colis[0]));
1109566063dSJacob Faibussowitsch PetscCall(ISSort(rowis[0]));
1119566063dSJacob Faibussowitsch PetscCall(ISSort(colis[0]));
1129566063dSJacob Faibussowitsch PetscCall(PetscFree2(rowindices, colindices));
113c4762a1bSJed Brown /*
114c4762a1bSJed Brown Now view the ISs. To avoid deadlock when viewing a list of objects on different subcomms,
115c4762a1bSJed Brown we need to obtain the global numbers of our local objects and wait for the corresponding global
116c4762a1bSJed Brown number to be viewed.
117c4762a1bSJed Brown */
1189566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Subdomains"));
11948a46eb9SPierre Jolivet if (permute_indices) PetscCall(PetscViewerASCIIPrintf(viewer, " (hash=%" PetscInt_FMT ")", hash));
1209566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, ":\n"));
1219566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer));
122c4762a1bSJed Brown
123c4762a1bSJed Brown nsubdomains = 1;
124c4762a1bSJed Brown for (s = 0; s < nsubdomains; ++s) gsubdomainperm[s] = s;
1259566063dSJacob Faibussowitsch PetscCall(PetscObjectsListGetGlobalNumbering(PETSC_COMM_WORLD, 1, (PetscObject *)rowis, &gnsubdomains, gsubdomainnums));
1269566063dSJacob Faibussowitsch PetscCall(PetscSortIntWithPermutation(nsubdomains, gsubdomainnums, gsubdomainperm));
127c4762a1bSJed Brown for (gs = 0, s = 0; gs < gnsubdomains; ++gs) {
128c4762a1bSJed Brown PetscInt ss;
129fe8fb074SBarry Smith if (s < nsubdomains) {
130c4762a1bSJed Brown ss = gsubdomainperm[s];
131c4762a1bSJed Brown if (gs == gsubdomainnums[ss]) { /* Global subdomain gs being viewed is my subdomain with local number ss. */
132c4762a1bSJed Brown ++s;
133fe8fb074SBarry Smith } else ss = -1;
134fe8fb074SBarry Smith } else ss = -1;
135fe8fb074SBarry Smith PetscCall(MyISView(rowis, colis, gs, ss, viewer));
136c4762a1bSJed Brown }
1379566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer));
1389566063dSJacob Faibussowitsch PetscCall(ISSort(rowis[0]));
1399566063dSJacob Faibussowitsch PetscCall(ISSort(colis[0]));
140c4762a1bSJed Brown nsubdomains = 1;
1419566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatricesMPI(A, nsubdomains, rowis, colis, MAT_INITIAL_MATRIX, &submats));
142c4762a1bSJed Brown /*
143c4762a1bSJed Brown Now view the matrices. To avoid deadlock when viewing a list of objects on different subcomms,
144c4762a1bSJed Brown we need to obtain the global numbers of our local objects and wait for the corresponding global
14514c49948SBarry Smith number to be viewed. Also all MPI processes need to call PetscViewerGetSubViewer() the same number of times
146c4762a1bSJed Brown */
1479566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Submatrices (repetition 1):\n"));
148c4762a1bSJed Brown for (s = 0; s < nsubdomains; ++s) gsubdomainperm[s] = s;
1499566063dSJacob Faibussowitsch PetscCall(PetscObjectsListGetGlobalNumbering(PETSC_COMM_WORLD, 1, (PetscObject *)submats, &gnsubdomains, gsubdomainnums));
1509566063dSJacob Faibussowitsch PetscCall(PetscSortIntWithPermutation(nsubdomains, gsubdomainnums, gsubdomainperm));
151c4762a1bSJed Brown for (gs = 0, s = 0; gs < gnsubdomains; ++gs) {
15214c49948SBarry Smith PetscViewer subviewer = NULL;
153c4762a1bSJed Brown if (s < nsubdomains) {
154c4762a1bSJed Brown PetscInt ss;
155c4762a1bSJed Brown ss = gsubdomainperm[s];
156c4762a1bSJed Brown if (gs == gsubdomainnums[ss]) { /* Global subdomain gs being viewed is my subdomain with local number ss. */
1579566063dSJacob Faibussowitsch PetscCall(PetscViewerGetSubViewer(viewer, PetscObjectComm((PetscObject)submats[ss]), &subviewer));
1589566063dSJacob Faibussowitsch PetscCall(MatView(submats[ss], subviewer));
1599566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(viewer, PetscObjectComm((PetscObject)submats[ss]), &subviewer));
160c4762a1bSJed Brown ++s;
16114c49948SBarry Smith } else {
16214c49948SBarry Smith PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &subviewer));
16314c49948SBarry Smith PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &subviewer));
16414c49948SBarry Smith }
16514c49948SBarry Smith } else {
16614c49948SBarry Smith PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &subviewer));
16714c49948SBarry Smith PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &subviewer));
168c4762a1bSJed Brown }
169c4762a1bSJed Brown }
170c4762a1bSJed Brown if (rep == 1) goto cleanup;
171c4762a1bSJed Brown nsubdomains = 1;
1729566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatricesMPI(A, nsubdomains, rowis, colis, MAT_REUSE_MATRIX, &submats));
173c4762a1bSJed Brown /*
174c4762a1bSJed Brown Now view the matrices. To avoid deadlock when viewing a list of objects on different subcomms,
175c4762a1bSJed Brown we need to obtain the global numbers of our local objects and wait for the corresponding global
176c4762a1bSJed Brown number to be viewed.
177c4762a1bSJed Brown */
1789566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Submatrices (repetition 2):\n"));
179c4762a1bSJed Brown for (s = 0; s < nsubdomains; ++s) gsubdomainperm[s] = s;
1809566063dSJacob Faibussowitsch PetscCall(PetscObjectsListGetGlobalNumbering(PETSC_COMM_WORLD, 1, (PetscObject *)submats, &gnsubdomains, gsubdomainnums));
1819566063dSJacob Faibussowitsch PetscCall(PetscSortIntWithPermutation(nsubdomains, gsubdomainnums, gsubdomainperm));
182c4762a1bSJed Brown for (gs = 0, s = 0; gs < gnsubdomains; ++gs) {
18314c49948SBarry Smith PetscViewer subviewer = NULL;
184c4762a1bSJed Brown if (s < nsubdomains) {
185c4762a1bSJed Brown PetscInt ss;
186c4762a1bSJed Brown ss = gsubdomainperm[s];
187c4762a1bSJed Brown if (gs == gsubdomainnums[ss]) { /* Global subdomain gs being viewed is my subdomain with local number ss. */
1889566063dSJacob Faibussowitsch PetscCall(PetscViewerGetSubViewer(viewer, PetscObjectComm((PetscObject)submats[ss]), &subviewer));
1899566063dSJacob Faibussowitsch PetscCall(MatView(submats[ss], subviewer));
1909566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(viewer, PetscObjectComm((PetscObject)submats[ss]), &subviewer));
191c4762a1bSJed Brown ++s;
19214c49948SBarry Smith } else {
19314c49948SBarry Smith PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &subviewer));
19414c49948SBarry Smith PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &subviewer));
19514c49948SBarry Smith }
19614c49948SBarry Smith } else {
19714c49948SBarry Smith PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &subviewer));
19814c49948SBarry Smith PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &subviewer));
199c4762a1bSJed Brown }
200c4762a1bSJed Brown }
201c4762a1bSJed Brown cleanup:
20248a46eb9SPierre Jolivet for (k = 0; k < nsubdomains; ++k) PetscCall(MatDestroy(submats + k));
2039566063dSJacob Faibussowitsch PetscCall(PetscFree(submats));
204c4762a1bSJed Brown for (k = 0; k < nis; ++k) {
2059566063dSJacob Faibussowitsch PetscCall(ISDestroy(rowis + k));
2069566063dSJacob Faibussowitsch PetscCall(ISDestroy(colis + k));
207c4762a1bSJed Brown }
2089566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A));
2099566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free(&subcomm));
2109566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
211b122ec5aSJacob Faibussowitsch return 0;
212c4762a1bSJed Brown }
213c4762a1bSJed Brown
214c4762a1bSJed Brown /*TEST
215c4762a1bSJed Brown
216c4762a1bSJed Brown test:
217c4762a1bSJed Brown nsize: 2
218c4762a1bSJed Brown args: -total_subdomains 1
219c4762a1bSJed Brown output_file: output/ex183_2_1.out
220c4762a1bSJed Brown
221c4762a1bSJed Brown test:
222c4762a1bSJed Brown suffix: 2
223c4762a1bSJed Brown nsize: 3
224c4762a1bSJed Brown args: -total_subdomains 2
225c4762a1bSJed Brown output_file: output/ex183_3_2.out
226c4762a1bSJed Brown
227c4762a1bSJed Brown test:
228c4762a1bSJed Brown suffix: 3
229c4762a1bSJed Brown nsize: 4
230c4762a1bSJed Brown args: -total_subdomains 2
231c4762a1bSJed Brown output_file: output/ex183_4_2.out
232c4762a1bSJed Brown
233c4762a1bSJed Brown test:
234c4762a1bSJed Brown suffix: 4
235c4762a1bSJed Brown nsize: 6
236c4762a1bSJed Brown args: -total_subdomains 2
237c4762a1bSJed Brown output_file: output/ex183_6_2.out
238c4762a1bSJed Brown
239c4762a1bSJed Brown TEST*/
240