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 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; 12fe8fb074SBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, "Row IS %" PetscInt_FMT "\n", gs)); 13fe8fb074SBarry Smith if (ss > -1) { 14fe8fb074SBarry Smith PetscCall(PetscViewerGetSubViewer(viewer, PetscObjectComm((PetscObject)rowis[ss]), &subviewer)); 15fe8fb074SBarry Smith PetscCall(ISView(rowis[ss], subviewer)); 16fe8fb074SBarry Smith PetscCall(PetscViewerRestoreSubViewer(viewer, PetscObjectComm((PetscObject)rowis[ss]), &subviewer)); 17fe8fb074SBarry Smith } 18fe8fb074SBarry Smith PetscCall(PetscViewerFlush(viewer)); 19fe8fb074SBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, "Col IS %" PetscInt_FMT "\n", gs)); 20fe8fb074SBarry Smith if (ss > -1) { 21fe8fb074SBarry Smith PetscCall(PetscViewerGetSubViewer(viewer, PetscObjectComm((PetscObject)rowis[ss]), &subviewer)); 22fe8fb074SBarry Smith PetscCall(ISView(colis[ss], subviewer)); 23fe8fb074SBarry Smith PetscCall(PetscViewerRestoreSubViewer(viewer, PetscObjectComm((PetscObject)rowis[ss]), &subviewer)); 24fe8fb074SBarry Smith } 25fe8fb074SBarry Smith PetscCall(PetscViewerFlush(viewer)); 26fe8fb074SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 27fe8fb074SBarry Smith } 28fe8fb074SBarry Smith 29d71ae5a4SJacob Faibussowitsch int main(int argc, char **args) 30d71ae5a4SJacob Faibussowitsch { 31c4762a1bSJed Brown Mat A, *submats; 32c4762a1bSJed Brown MPI_Comm subcomm; 33c4762a1bSJed Brown PetscMPIInt rank, size, subrank, subsize, color; 34c4762a1bSJed Brown PetscInt m, n, N, bs, rstart, rend, i, j, k, total_subdomains, hash, nsubdomains = 1; 35c4762a1bSJed Brown PetscInt nis, *cols, gnsubdomains, gsubdomainnums[1], gsubdomainperm[1], s, gs; 36c4762a1bSJed Brown PetscInt *rowindices, *colindices, idx, rep; 37c4762a1bSJed Brown PetscScalar *vals; 38c4762a1bSJed Brown IS rowis[1], colis[1]; 39c4762a1bSJed Brown PetscViewer viewer; 40c4762a1bSJed Brown PetscBool permute_indices, flg; 41c4762a1bSJed Brown 42327415f7SBarry Smith PetscFunctionBeginUser; 439566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &args, (char *)0, help)); 449566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 459566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 46c4762a1bSJed Brown 47d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "ex183", "Mat"); 48c4762a1bSJed Brown m = 5; 499566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-m", "Local matrix size", "MatSetSizes", m, &m, &flg)); 50c4762a1bSJed Brown total_subdomains = size - 1; 51*52ce0ab5SPierre Jolivet PetscCall(PetscOptionsRangeInt("-total_subdomains", "Number of submatrices where 0 < n < comm size", "MatCreateSubMatricesMPI", total_subdomains, &total_subdomains, &flg, 1, size)); 52c4762a1bSJed Brown permute_indices = PETSC_FALSE; 539566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-permute_indices", "Whether to permute indices before breaking them into subdomains", "ISCreateGeneral", permute_indices, &permute_indices, &flg)); 54c4762a1bSJed Brown hash = 7; 559566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-hash", "Permutation factor, which has to be relatively prime to M = size*m (total matrix size)", "ISCreateGeneral", hash, &hash, &flg)); 56c4762a1bSJed Brown rep = 2; 57*52ce0ab5SPierre Jolivet PetscCall(PetscOptionsRangeInt("-rep", "Number of times to carry out submatrix extractions; currently only 1 & 2 are supported", NULL, rep, &rep, &flg, 1, 2)); 58d0609cedSBarry Smith PetscOptionsEnd(); 59c4762a1bSJed Brown 60c4762a1bSJed Brown viewer = PETSC_VIEWER_STDOUT_WORLD; 61c4762a1bSJed Brown /* Create logically sparse, but effectively dense matrix for easy verification of submatrix extraction correctness. */ 629566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 639566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, m, m, PETSC_DECIDE, PETSC_DECIDE)); 649566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 659566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 669566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, NULL, &N)); 679566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A, NULL, &n)); 689566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(A, &bs)); 699566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(A, n, NULL)); 709566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(A, n, NULL, N - n, NULL)); 719566063dSJacob Faibussowitsch PetscCall(MatSeqBAIJSetPreallocation(A, bs, n / bs, NULL)); 729566063dSJacob Faibussowitsch PetscCall(MatMPIBAIJSetPreallocation(A, bs, n / bs, NULL, (N - n) / bs, NULL)); 739566063dSJacob Faibussowitsch PetscCall(MatSeqSBAIJSetPreallocation(A, bs, n / bs, NULL)); 749566063dSJacob Faibussowitsch PetscCall(MatMPISBAIJSetPreallocation(A, bs, n / bs, NULL, (N - n) / bs, NULL)); 75c4762a1bSJed Brown 769566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(N, &cols, N, &vals)); 779566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A, &rstart, &rend)); 78c4762a1bSJed Brown for (j = 0; j < N; ++j) cols[j] = j; 79c4762a1bSJed Brown for (i = rstart; i < rend; i++) { 80ad540459SPierre Jolivet for (j = 0; j < N; ++j) vals[j] = i * 10000 + j; 819566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &i, N, cols, vals, INSERT_VALUES)); 82c4762a1bSJed Brown } 839566063dSJacob Faibussowitsch PetscCall(PetscFree2(cols, vals)); 849566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 859566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 86c4762a1bSJed Brown 879566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Initial matrix:\n")); 889566063dSJacob Faibussowitsch PetscCall(MatView(A, viewer)); 89c4762a1bSJed Brown 90c4762a1bSJed Brown /* 91c4762a1bSJed Brown Create subcomms and ISs so that each rank participates in one IS. 92c4762a1bSJed Brown The IS either coalesces adjacent rank indices (contiguous), 93c4762a1bSJed Brown or selects indices by scrambling them using a hash. 94c4762a1bSJed Brown */ 95c4762a1bSJed Brown k = size / total_subdomains + (size % total_subdomains > 0); /* There are up to k ranks to a color */ 96c4762a1bSJed Brown color = rank / k; 979566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_split(PETSC_COMM_WORLD, color, rank, &subcomm)); 989566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(subcomm, &subsize)); 999566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(subcomm, &subrank)); 1009566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A, &rstart, &rend)); 101c4762a1bSJed Brown nis = 1; 1029566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(rend - rstart, &rowindices, rend - rstart, &colindices)); 103c4762a1bSJed Brown 104c4762a1bSJed Brown for (j = rstart; j < rend; ++j) { 105c4762a1bSJed Brown if (permute_indices) { 106c4762a1bSJed Brown idx = (j * hash); 107c4762a1bSJed Brown } else { 108c4762a1bSJed Brown idx = j; 109c4762a1bSJed Brown } 110c4762a1bSJed Brown rowindices[j - rstart] = idx % N; 111c4762a1bSJed Brown colindices[j - rstart] = (idx + m) % N; 112c4762a1bSJed Brown } 1139566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(subcomm, rend - rstart, rowindices, PETSC_COPY_VALUES, &rowis[0])); 1149566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(subcomm, rend - rstart, colindices, PETSC_COPY_VALUES, &colis[0])); 1159566063dSJacob Faibussowitsch PetscCall(ISSort(rowis[0])); 1169566063dSJacob Faibussowitsch PetscCall(ISSort(colis[0])); 1179566063dSJacob Faibussowitsch PetscCall(PetscFree2(rowindices, colindices)); 118c4762a1bSJed Brown /* 119c4762a1bSJed Brown Now view the ISs. To avoid deadlock when viewing a list of objects on different subcomms, 120c4762a1bSJed Brown we need to obtain the global numbers of our local objects and wait for the corresponding global 121c4762a1bSJed Brown number to be viewed. 122c4762a1bSJed Brown */ 1239566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Subdomains")); 12448a46eb9SPierre Jolivet if (permute_indices) PetscCall(PetscViewerASCIIPrintf(viewer, " (hash=%" PetscInt_FMT ")", hash)); 1259566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, ":\n")); 1269566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 127c4762a1bSJed Brown 128c4762a1bSJed Brown nsubdomains = 1; 129c4762a1bSJed Brown for (s = 0; s < nsubdomains; ++s) gsubdomainperm[s] = s; 1309566063dSJacob Faibussowitsch PetscCall(PetscObjectsListGetGlobalNumbering(PETSC_COMM_WORLD, 1, (PetscObject *)rowis, &gnsubdomains, gsubdomainnums)); 1319566063dSJacob Faibussowitsch PetscCall(PetscSortIntWithPermutation(nsubdomains, gsubdomainnums, gsubdomainperm)); 132c4762a1bSJed Brown for (gs = 0, s = 0; gs < gnsubdomains; ++gs) { 133c4762a1bSJed Brown PetscInt ss; 134fe8fb074SBarry Smith if (s < nsubdomains) { 135c4762a1bSJed Brown ss = gsubdomainperm[s]; 136c4762a1bSJed Brown if (gs == gsubdomainnums[ss]) { /* Global subdomain gs being viewed is my subdomain with local number ss. */ 137c4762a1bSJed Brown ++s; 138fe8fb074SBarry Smith } else ss = -1; 139fe8fb074SBarry Smith } else ss = -1; 140fe8fb074SBarry Smith PetscCall(MyISView(rowis, colis, gs, ss, viewer)); 141c4762a1bSJed Brown } 1429566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 1439566063dSJacob Faibussowitsch PetscCall(ISSort(rowis[0])); 1449566063dSJacob Faibussowitsch PetscCall(ISSort(colis[0])); 145c4762a1bSJed Brown nsubdomains = 1; 1469566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatricesMPI(A, nsubdomains, rowis, colis, MAT_INITIAL_MATRIX, &submats)); 147c4762a1bSJed Brown /* 148c4762a1bSJed Brown Now view the matrices. To avoid deadlock when viewing a list of objects on different subcomms, 149c4762a1bSJed Brown we need to obtain the global numbers of our local objects and wait for the corresponding global 15014c49948SBarry Smith number to be viewed. Also all MPI processes need to call PetscViewerGetSubViewer() the same number of times 151c4762a1bSJed Brown */ 1529566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Submatrices (repetition 1):\n")); 153c4762a1bSJed Brown for (s = 0; s < nsubdomains; ++s) gsubdomainperm[s] = s; 1549566063dSJacob Faibussowitsch PetscCall(PetscObjectsListGetGlobalNumbering(PETSC_COMM_WORLD, 1, (PetscObject *)submats, &gnsubdomains, gsubdomainnums)); 1559566063dSJacob Faibussowitsch PetscCall(PetscSortIntWithPermutation(nsubdomains, gsubdomainnums, gsubdomainperm)); 156c4762a1bSJed Brown for (gs = 0, s = 0; gs < gnsubdomains; ++gs) { 15714c49948SBarry Smith PetscViewer subviewer = NULL; 158c4762a1bSJed Brown if (s < nsubdomains) { 159c4762a1bSJed Brown PetscInt ss; 160c4762a1bSJed Brown ss = gsubdomainperm[s]; 161c4762a1bSJed Brown if (gs == gsubdomainnums[ss]) { /* Global subdomain gs being viewed is my subdomain with local number ss. */ 1629566063dSJacob Faibussowitsch PetscCall(PetscViewerGetSubViewer(viewer, PetscObjectComm((PetscObject)submats[ss]), &subviewer)); 1639566063dSJacob Faibussowitsch PetscCall(MatView(submats[ss], subviewer)); 1649566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(viewer, PetscObjectComm((PetscObject)submats[ss]), &subviewer)); 165c4762a1bSJed Brown ++s; 16614c49948SBarry Smith } else { 16714c49948SBarry Smith PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &subviewer)); 16814c49948SBarry Smith PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &subviewer)); 16914c49948SBarry Smith } 17014c49948SBarry Smith } else { 17114c49948SBarry Smith PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &subviewer)); 17214c49948SBarry Smith PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &subviewer)); 173c4762a1bSJed Brown } 174c4762a1bSJed Brown } 175c4762a1bSJed Brown if (rep == 1) goto cleanup; 176c4762a1bSJed Brown nsubdomains = 1; 1779566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatricesMPI(A, nsubdomains, rowis, colis, MAT_REUSE_MATRIX, &submats)); 178c4762a1bSJed Brown /* 179c4762a1bSJed Brown Now view the matrices. To avoid deadlock when viewing a list of objects on different subcomms, 180c4762a1bSJed Brown we need to obtain the global numbers of our local objects and wait for the corresponding global 181c4762a1bSJed Brown number to be viewed. 182c4762a1bSJed Brown */ 1839566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Submatrices (repetition 2):\n")); 184c4762a1bSJed Brown for (s = 0; s < nsubdomains; ++s) gsubdomainperm[s] = s; 1859566063dSJacob Faibussowitsch PetscCall(PetscObjectsListGetGlobalNumbering(PETSC_COMM_WORLD, 1, (PetscObject *)submats, &gnsubdomains, gsubdomainnums)); 1869566063dSJacob Faibussowitsch PetscCall(PetscSortIntWithPermutation(nsubdomains, gsubdomainnums, gsubdomainperm)); 187c4762a1bSJed Brown for (gs = 0, s = 0; gs < gnsubdomains; ++gs) { 18814c49948SBarry Smith PetscViewer subviewer = NULL; 189c4762a1bSJed Brown if (s < nsubdomains) { 190c4762a1bSJed Brown PetscInt ss; 191c4762a1bSJed Brown ss = gsubdomainperm[s]; 192c4762a1bSJed Brown if (gs == gsubdomainnums[ss]) { /* Global subdomain gs being viewed is my subdomain with local number ss. */ 1939566063dSJacob Faibussowitsch PetscCall(PetscViewerGetSubViewer(viewer, PetscObjectComm((PetscObject)submats[ss]), &subviewer)); 1949566063dSJacob Faibussowitsch PetscCall(MatView(submats[ss], subviewer)); 1959566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(viewer, PetscObjectComm((PetscObject)submats[ss]), &subviewer)); 196c4762a1bSJed Brown ++s; 19714c49948SBarry Smith } else { 19814c49948SBarry Smith PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &subviewer)); 19914c49948SBarry Smith PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &subviewer)); 20014c49948SBarry Smith } 20114c49948SBarry Smith } else { 20214c49948SBarry Smith PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &subviewer)); 20314c49948SBarry Smith PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &subviewer)); 204c4762a1bSJed Brown } 205c4762a1bSJed Brown } 206c4762a1bSJed Brown cleanup: 20748a46eb9SPierre Jolivet for (k = 0; k < nsubdomains; ++k) PetscCall(MatDestroy(submats + k)); 2089566063dSJacob Faibussowitsch PetscCall(PetscFree(submats)); 209c4762a1bSJed Brown for (k = 0; k < nis; ++k) { 2109566063dSJacob Faibussowitsch PetscCall(ISDestroy(rowis + k)); 2119566063dSJacob Faibussowitsch PetscCall(ISDestroy(colis + k)); 212c4762a1bSJed Brown } 2139566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 2149566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free(&subcomm)); 2159566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 216b122ec5aSJacob Faibussowitsch return 0; 217c4762a1bSJed Brown } 218c4762a1bSJed Brown 219c4762a1bSJed Brown /*TEST 220c4762a1bSJed Brown 221c4762a1bSJed Brown test: 222c4762a1bSJed Brown nsize: 2 223c4762a1bSJed Brown args: -total_subdomains 1 224c4762a1bSJed Brown output_file: output/ex183_2_1.out 225c4762a1bSJed Brown 226c4762a1bSJed Brown test: 227c4762a1bSJed Brown suffix: 2 228c4762a1bSJed Brown nsize: 3 229c4762a1bSJed Brown args: -total_subdomains 2 230c4762a1bSJed Brown output_file: output/ex183_3_2.out 231c4762a1bSJed Brown 232c4762a1bSJed Brown test: 233c4762a1bSJed Brown suffix: 3 234c4762a1bSJed Brown nsize: 4 235c4762a1bSJed Brown args: -total_subdomains 2 236c4762a1bSJed Brown output_file: output/ex183_4_2.out 237c4762a1bSJed Brown 238c4762a1bSJed Brown test: 239c4762a1bSJed Brown suffix: 4 240c4762a1bSJed Brown nsize: 6 241c4762a1bSJed Brown args: -total_subdomains 2 242c4762a1bSJed Brown output_file: output/ex183_6_2.out 243c4762a1bSJed Brown 244c4762a1bSJed Brown TEST*/ 245