1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests MatCreateSubmatrix() in parallel."; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscmat.h> 5c4762a1bSJed Brown #include <../src/mat/impls/aij/mpi/mpiaij.h> 6c4762a1bSJed Brown 7d71ae5a4SJacob Faibussowitsch PetscErrorCode ISGetSeqIS_SameColDist_Private(Mat mat, IS isrow, IS iscol, IS *isrow_d, IS *iscol_d, IS *iscol_o, const PetscInt *garray[]) 8d71ae5a4SJacob Faibussowitsch { 9c4762a1bSJed Brown Vec x, cmap; 10c4762a1bSJed Brown const PetscInt *is_idx; 11c4762a1bSJed Brown PetscScalar *xarray, *cmaparray; 12c4762a1bSJed Brown PetscInt ncols, isstart, *idx, m, rstart, count; 13c4762a1bSJed Brown Mat_MPIAIJ *a = (Mat_MPIAIJ *)mat->data; 14c4762a1bSJed Brown Mat B = a->B; 15c4762a1bSJed Brown Vec lvec = a->lvec, lcmap; 16c4762a1bSJed Brown PetscInt i, cstart, cend, Bn = B->cmap->N; 17c4762a1bSJed Brown MPI_Comm comm; 18c4762a1bSJed Brown PetscMPIInt rank; 19c4762a1bSJed Brown VecScatter Mvctx; 20c4762a1bSJed Brown 21c4762a1bSJed Brown PetscFunctionBegin; 229566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)mat, &comm)); 239566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 249566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(iscol, &ncols)); 25c4762a1bSJed Brown 26c4762a1bSJed Brown /* (1) iscol is a sub-column vector of mat, pad it with '-1.' to form a full vector x */ 279566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(mat, &x, NULL)); 289566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &cmap)); 299566063dSJacob Faibussowitsch PetscCall(VecSet(x, -1.0)); 309566063dSJacob Faibussowitsch PetscCall(VecSet(cmap, -1.0)); 31c4762a1bSJed Brown 329566063dSJacob Faibussowitsch PetscCall(VecDuplicate(lvec, &lcmap)); 33c4762a1bSJed Brown 34c4762a1bSJed Brown /* Get start indices */ 359566063dSJacob Faibussowitsch PetscCallMPI(MPI_Scan(&ncols, &isstart, 1, MPIU_INT, MPI_SUM, comm)); 36c4762a1bSJed Brown isstart -= ncols; 379566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRangeColumn(mat, &cstart, &cend)); 38c4762a1bSJed Brown 399566063dSJacob Faibussowitsch PetscCall(ISGetIndices(iscol, &is_idx)); 409566063dSJacob Faibussowitsch PetscCall(VecGetArray(x, &xarray)); 419566063dSJacob Faibussowitsch PetscCall(VecGetArray(cmap, &cmaparray)); 429566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ncols, &idx)); 43c4762a1bSJed Brown for (i = 0; i < ncols; i++) { 44c4762a1bSJed Brown xarray[is_idx[i] - cstart] = (PetscScalar)is_idx[i]; 45c4762a1bSJed Brown cmaparray[is_idx[i] - cstart] = (PetscScalar)(i + isstart); /* global index of iscol[i] */ 46c4762a1bSJed Brown idx[i] = is_idx[i] - cstart; /* local index of iscol[i] */ 47c4762a1bSJed Brown } 489566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(x, &xarray)); 499566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(cmap, &cmaparray)); 509566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(iscol, &is_idx)); 51c4762a1bSJed Brown 52c4762a1bSJed Brown /* Get iscol_d */ 539566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, ncols, idx, PETSC_OWN_POINTER, iscol_d)); 549566063dSJacob Faibussowitsch PetscCall(ISGetBlockSize(iscol, &i)); 559566063dSJacob Faibussowitsch PetscCall(ISSetBlockSize(*iscol_d, i)); 56c4762a1bSJed Brown 57c4762a1bSJed Brown /* Get isrow_d */ 589566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(isrow, &m)); 59c4762a1bSJed Brown rstart = mat->rmap->rstart; 609566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m, &idx)); 619566063dSJacob Faibussowitsch PetscCall(ISGetIndices(isrow, &is_idx)); 62c4762a1bSJed Brown for (i = 0; i < m; i++) idx[i] = is_idx[i] - rstart; 639566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(isrow, &is_idx)); 64c4762a1bSJed Brown 659566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, m, idx, PETSC_OWN_POINTER, isrow_d)); 669566063dSJacob Faibussowitsch PetscCall(ISGetBlockSize(isrow, &i)); 679566063dSJacob Faibussowitsch PetscCall(ISSetBlockSize(*isrow_d, i)); 68c4762a1bSJed Brown 69c4762a1bSJed Brown /* (2) Scatter x and cmap using aij->Mvctx to get their off-process portions (see MatMult_MPIAIJ) */ 70c4762a1bSJed Brown #if 0 71c4762a1bSJed Brown if (!a->Mvctx_mpi1) { 72c4762a1bSJed Brown /* a->Mvctx causes random 'count' in o-build? See src/mat/tests/runex59_2 */ 73c4762a1bSJed Brown a->Mvctx_mpi1_flg = PETSC_TRUE; 749566063dSJacob Faibussowitsch PetscCall(MatSetUpMultiply_MPIAIJ(mat)); 75c4762a1bSJed Brown } 76c4762a1bSJed Brown Mvctx = a->Mvctx_mpi1; 77c4762a1bSJed Brown #endif 78c4762a1bSJed Brown Mvctx = a->Mvctx; 799566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(Mvctx, x, lvec, INSERT_VALUES, SCATTER_FORWARD)); 809566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(Mvctx, x, lvec, INSERT_VALUES, SCATTER_FORWARD)); 81c4762a1bSJed Brown 829566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(Mvctx, cmap, lcmap, INSERT_VALUES, SCATTER_FORWARD)); 839566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(Mvctx, cmap, lcmap, INSERT_VALUES, SCATTER_FORWARD)); 84c4762a1bSJed Brown 85c4762a1bSJed Brown /* (3) create sequential iscol_o (a subset of iscol) and isgarray */ 86c4762a1bSJed Brown /* off-process column indices */ 87c4762a1bSJed Brown count = 0; 88c4762a1bSJed Brown PetscInt *cmap1; 899566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(Bn, &idx)); 909566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(Bn, &cmap1)); 91c4762a1bSJed Brown 929566063dSJacob Faibussowitsch PetscCall(VecGetArray(lvec, &xarray)); 939566063dSJacob Faibussowitsch PetscCall(VecGetArray(lcmap, &cmaparray)); 94c4762a1bSJed Brown for (i = 0; i < Bn; i++) { 95c4762a1bSJed Brown if (PetscRealPart(xarray[i]) > -1.0) { 96c4762a1bSJed Brown idx[count] = i; /* local column index in off-diagonal part B */ 97c4762a1bSJed Brown cmap1[count] = (PetscInt)(PetscRealPart(cmaparray[i])); /* column index in submat */ 98c4762a1bSJed Brown count++; 99c4762a1bSJed Brown } 100c4762a1bSJed Brown } 101c4762a1bSJed Brown printf("[%d] Bn %d, count %d\n", rank, Bn, count); 1029566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(lvec, &xarray)); 1039566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(lcmap, &cmaparray)); 104c4762a1bSJed Brown if (count != 6) { 105c4762a1bSJed Brown printf("[%d] count %d != 6 lvec:\n", rank, count); 1069566063dSJacob Faibussowitsch PetscCall(VecView(lvec, 0)); 107c4762a1bSJed Brown 108c4762a1bSJed Brown printf("[%d] count %d != 6 lcmap:\n", rank, count); 1099566063dSJacob Faibussowitsch PetscCall(VecView(lcmap, 0)); 11098921bdaSJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "count %d != 6", count); 111c4762a1bSJed Brown } 112c4762a1bSJed Brown 1139566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, count, idx, PETSC_COPY_VALUES, iscol_o)); 114c4762a1bSJed Brown /* cannot ensure iscol_o has same blocksize as iscol! */ 115c4762a1bSJed Brown 1169566063dSJacob Faibussowitsch PetscCall(PetscFree(idx)); 117c4762a1bSJed Brown 118c4762a1bSJed Brown *garray = cmap1; 119c4762a1bSJed Brown 1209566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 1219566063dSJacob Faibussowitsch PetscCall(VecDestroy(&cmap)); 1229566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lcmap)); 123*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 124c4762a1bSJed Brown } 125c4762a1bSJed Brown 126d71ae5a4SJacob Faibussowitsch int main(int argc, char **args) 127d71ae5a4SJacob Faibussowitsch { 128c4762a1bSJed Brown Mat C, A; 129c4762a1bSJed Brown PetscInt i, j, m = 3, n = 2, rstart, rend; 130c4762a1bSJed Brown PetscMPIInt size, rank; 131c4762a1bSJed Brown PetscScalar v; 132c4762a1bSJed Brown IS isrow, iscol; 133c4762a1bSJed Brown 134327415f7SBarry Smith PetscFunctionBeginUser; 1359566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &args, (char *)0, help)); 1369566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 1379566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 138c4762a1bSJed Brown n = 2 * size; 139c4762a1bSJed Brown 1409566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &C)); 1419566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C, PETSC_DECIDE, PETSC_DECIDE, m * n, m * n)); 1429566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(C)); 1439566063dSJacob Faibussowitsch PetscCall(MatSetUp(C)); 144c4762a1bSJed Brown 145c4762a1bSJed Brown /* 146c4762a1bSJed Brown This is JUST to generate a nice test matrix, all processors fill up 147c4762a1bSJed Brown the entire matrix. This is not something one would ever do in practice. 148c4762a1bSJed Brown */ 1499566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(C, &rstart, &rend)); 150c4762a1bSJed Brown for (i = rstart; i < rend; i++) { 151c4762a1bSJed Brown for (j = 0; j < m * n; j++) { 152c4762a1bSJed Brown v = i + j + 1; 1539566063dSJacob Faibussowitsch PetscCall(MatSetValues(C, 1, &i, 1, &j, &v, INSERT_VALUES)); 154c4762a1bSJed Brown } 155c4762a1bSJed Brown } 156c4762a1bSJed Brown 1579566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY)); 1589566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY)); 159c4762a1bSJed Brown 160c4762a1bSJed Brown /* 161c4762a1bSJed Brown Generate a new matrix consisting of every second row and column of 162c4762a1bSJed Brown the original matrix 163c4762a1bSJed Brown */ 1649566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(C, &rstart, &rend)); 165c4762a1bSJed Brown /* Create parallel IS with the rows we want on THIS processor */ 1669566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD, (rend - rstart) / 2, rstart, 2, &isrow)); 167c4762a1bSJed Brown /* Create parallel IS with the rows we want on THIS processor (same as rows for now) */ 1689566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD, (rend - rstart) / 2, rstart, 2, &iscol)); 169c4762a1bSJed Brown 170c4762a1bSJed Brown IS iscol_d, isrow_d, iscol_o; 171c4762a1bSJed Brown const PetscInt *garray; 1729566063dSJacob Faibussowitsch PetscCall(ISGetSeqIS_SameColDist_Private(C, isrow, iscol, &isrow_d, &iscol_d, &iscol_o, &garray)); 173c4762a1bSJed Brown 1749566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isrow_d)); 1759566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iscol_d)); 1769566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iscol_o)); 1779566063dSJacob Faibussowitsch PetscCall(PetscFree(garray)); 178c4762a1bSJed Brown 1799566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(C, isrow, iscol, MAT_INITIAL_MATRIX, &A)); 1809566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(C, isrow, iscol, MAT_REUSE_MATRIX, &A)); 181c4762a1bSJed Brown 1829566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isrow)); 1839566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iscol)); 1849566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 1859566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 1869566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 187b122ec5aSJacob Faibussowitsch return 0; 188c4762a1bSJed Brown } 189