1be1d678aSKris Buschelman 2b97cf49bSBarry Smith /* 3b97cf49bSBarry Smith Defines the basic matrix operations for the ADJ adjacency list matrix data-structure. 4b97cf49bSBarry Smith */ 5e37d9405SJed Brown #include <../src/mat/impls/adj/mpi/mpiadj.h> /*I "petscmat.h" I*/ 6841d17a1SFande Kong #include <petscsf.h> 7b97cf49bSBarry Smith 840244768SBarry Smith /* 92ef1f0ffSBarry Smith The interface should be easy to use for both MatCreateSubMatrix (parallel sub-matrix) and MatCreateSubMatrices (sequential sub-matrices) 102ef1f0ffSBarry Smith */ 11d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatCreateSubMatrix_MPIAdj_data(Mat adj, IS irows, IS icols, PetscInt **sadj_xadj, PetscInt **sadj_adjncy, PetscInt **sadj_values) 12d71ae5a4SJacob Faibussowitsch { 13131c27b5Sprj- PetscInt nlrows_is, icols_n, i, j, nroots, nleaves, rlocalindex, *ncols_send, *ncols_recv; 1440244768SBarry Smith PetscInt nlrows_mat, *adjncy_recv, Ncols_recv, Ncols_send, *xadj_recv, *values_recv; 15e895ccc0SFande Kong PetscInt *ncols_recv_offsets, loc, rnclos, *sadjncy, *sxadj, *svalues; 1640244768SBarry Smith const PetscInt *irows_indices, *icols_indices, *xadj, *adjncy; 17131c27b5Sprj- PetscMPIInt owner; 1840244768SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj *)adj->data; 1940244768SBarry Smith PetscLayout rmap; 2040244768SBarry Smith MPI_Comm comm; 2140244768SBarry Smith PetscSF sf; 2240244768SBarry Smith PetscSFNode *iremote; 2340244768SBarry Smith PetscBool done; 2440244768SBarry Smith 2540244768SBarry Smith PetscFunctionBegin; 269566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)adj, &comm)); 279566063dSJacob Faibussowitsch PetscCall(MatGetLayouts(adj, &rmap, NULL)); 289566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(irows, &nlrows_is)); 299566063dSJacob Faibussowitsch PetscCall(ISGetIndices(irows, &irows_indices)); 309566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nlrows_is, &iremote)); 3140244768SBarry Smith /* construct sf graph*/ 3240244768SBarry Smith nleaves = nlrows_is; 3340244768SBarry Smith for (i = 0; i < nlrows_is; i++) { 3440244768SBarry Smith owner = -1; 3540244768SBarry Smith rlocalindex = -1; 369566063dSJacob Faibussowitsch PetscCall(PetscLayoutFindOwnerIndex(rmap, irows_indices[i], &owner, &rlocalindex)); 3740244768SBarry Smith iremote[i].rank = owner; 3840244768SBarry Smith iremote[i].index = rlocalindex; 3940244768SBarry Smith } 409566063dSJacob Faibussowitsch PetscCall(MatGetRowIJ(adj, 0, PETSC_FALSE, PETSC_FALSE, &nlrows_mat, &xadj, &adjncy, &done)); 419566063dSJacob Faibussowitsch PetscCall(PetscCalloc4(nlrows_mat, &ncols_send, nlrows_is, &xadj_recv, nlrows_is + 1, &ncols_recv_offsets, nlrows_is, &ncols_recv)); 4240244768SBarry Smith nroots = nlrows_mat; 43ad540459SPierre Jolivet for (i = 0; i < nlrows_mat; i++) ncols_send[i] = xadj[i + 1] - xadj[i]; 449566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(comm, &sf)); 459566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(sf, nroots, nleaves, NULL, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER)); 469566063dSJacob Faibussowitsch PetscCall(PetscSFSetType(sf, PETSCSFBASIC)); 479566063dSJacob Faibussowitsch PetscCall(PetscSFSetFromOptions(sf)); 489566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(sf, MPIU_INT, ncols_send, ncols_recv, MPI_REPLACE)); 499566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(sf, MPIU_INT, ncols_send, ncols_recv, MPI_REPLACE)); 509566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(sf, MPIU_INT, xadj, xadj_recv, MPI_REPLACE)); 519566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(sf, MPIU_INT, xadj, xadj_recv, MPI_REPLACE)); 529566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sf)); 5340244768SBarry Smith Ncols_recv = 0; 5440244768SBarry Smith for (i = 0; i < nlrows_is; i++) { 5540244768SBarry Smith Ncols_recv += ncols_recv[i]; 5640244768SBarry Smith ncols_recv_offsets[i + 1] = ncols_recv[i] + ncols_recv_offsets[i]; 5740244768SBarry Smith } 5840244768SBarry Smith Ncols_send = 0; 59ad540459SPierre Jolivet for (i = 0; i < nlrows_mat; i++) Ncols_send += ncols_send[i]; 609566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(Ncols_recv, &iremote)); 619566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(Ncols_recv, &adjncy_recv)); 6240244768SBarry Smith nleaves = Ncols_recv; 6340244768SBarry Smith Ncols_recv = 0; 6440244768SBarry Smith for (i = 0; i < nlrows_is; i++) { 659566063dSJacob Faibussowitsch PetscCall(PetscLayoutFindOwner(rmap, irows_indices[i], &owner)); 6640244768SBarry Smith for (j = 0; j < ncols_recv[i]; j++) { 6740244768SBarry Smith iremote[Ncols_recv].rank = owner; 6840244768SBarry Smith iremote[Ncols_recv++].index = xadj_recv[i] + j; 6940244768SBarry Smith } 7040244768SBarry Smith } 719566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(irows, &irows_indices)); 7240244768SBarry Smith /*if we need to deal with edge weights ???*/ 739566063dSJacob Faibussowitsch if (a->useedgeweights) PetscCall(PetscCalloc1(Ncols_recv, &values_recv)); 7440244768SBarry Smith nroots = Ncols_send; 759566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(comm, &sf)); 769566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(sf, nroots, nleaves, NULL, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER)); 779566063dSJacob Faibussowitsch PetscCall(PetscSFSetType(sf, PETSCSFBASIC)); 789566063dSJacob Faibussowitsch PetscCall(PetscSFSetFromOptions(sf)); 799566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(sf, MPIU_INT, adjncy, adjncy_recv, MPI_REPLACE)); 809566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(sf, MPIU_INT, adjncy, adjncy_recv, MPI_REPLACE)); 81e895ccc0SFande Kong if (a->useedgeweights) { 829566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(sf, MPIU_INT, a->values, values_recv, MPI_REPLACE)); 839566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(sf, MPIU_INT, a->values, values_recv, MPI_REPLACE)); 8440244768SBarry Smith } 859566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sf)); 869566063dSJacob Faibussowitsch PetscCall(MatRestoreRowIJ(adj, 0, PETSC_FALSE, PETSC_FALSE, &nlrows_mat, &xadj, &adjncy, &done)); 879566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(icols, &icols_n)); 889566063dSJacob Faibussowitsch PetscCall(ISGetIndices(icols, &icols_indices)); 8940244768SBarry Smith rnclos = 0; 9040244768SBarry Smith for (i = 0; i < nlrows_is; i++) { 9140244768SBarry Smith for (j = ncols_recv_offsets[i]; j < ncols_recv_offsets[i + 1]; j++) { 929566063dSJacob Faibussowitsch PetscCall(PetscFindInt(adjncy_recv[j], icols_n, icols_indices, &loc)); 9340244768SBarry Smith if (loc < 0) { 9440244768SBarry Smith adjncy_recv[j] = -1; 95e895ccc0SFande Kong if (a->useedgeweights) values_recv[j] = -1; 9640244768SBarry Smith ncols_recv[i]--; 9740244768SBarry Smith } else { 9840244768SBarry Smith rnclos++; 9940244768SBarry Smith } 10040244768SBarry Smith } 10140244768SBarry Smith } 1029566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(icols, &icols_indices)); 1039566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(rnclos, &sadjncy)); 1049566063dSJacob Faibussowitsch if (a->useedgeweights) PetscCall(PetscCalloc1(rnclos, &svalues)); 1059566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(nlrows_is + 1, &sxadj)); 10640244768SBarry Smith rnclos = 0; 10740244768SBarry Smith for (i = 0; i < nlrows_is; i++) { 10840244768SBarry Smith for (j = ncols_recv_offsets[i]; j < ncols_recv_offsets[i + 1]; j++) { 10940244768SBarry Smith if (adjncy_recv[j] < 0) continue; 11040244768SBarry Smith sadjncy[rnclos] = adjncy_recv[j]; 111e895ccc0SFande Kong if (a->useedgeweights) svalues[rnclos] = values_recv[j]; 11240244768SBarry Smith rnclos++; 11340244768SBarry Smith } 11440244768SBarry Smith } 115ad540459SPierre Jolivet for (i = 0; i < nlrows_is; i++) sxadj[i + 1] = sxadj[i] + ncols_recv[i]; 1169371c9d4SSatish Balay if (sadj_xadj) { 1179371c9d4SSatish Balay *sadj_xadj = sxadj; 1189371c9d4SSatish Balay } else PetscCall(PetscFree(sxadj)); 1199371c9d4SSatish Balay if (sadj_adjncy) { 1209371c9d4SSatish Balay *sadj_adjncy = sadjncy; 1219371c9d4SSatish Balay } else PetscCall(PetscFree(sadjncy)); 12240244768SBarry Smith if (sadj_values) { 1239371c9d4SSatish Balay if (a->useedgeweights) *sadj_values = svalues; 1249371c9d4SSatish Balay else *sadj_values = NULL; 12540244768SBarry Smith } else { 1269566063dSJacob Faibussowitsch if (a->useedgeweights) PetscCall(PetscFree(svalues)); 12740244768SBarry Smith } 1289566063dSJacob Faibussowitsch PetscCall(PetscFree4(ncols_send, xadj_recv, ncols_recv_offsets, ncols_recv)); 1299566063dSJacob Faibussowitsch PetscCall(PetscFree(adjncy_recv)); 1309566063dSJacob Faibussowitsch if (a->useedgeweights) PetscCall(PetscFree(values_recv)); 1313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13240244768SBarry Smith } 13340244768SBarry Smith 134d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatCreateSubMatrices_MPIAdj_Private(Mat mat, PetscInt n, const IS irow[], const IS icol[], PetscBool subcomm, MatReuse scall, Mat *submat[]) 135d71ae5a4SJacob Faibussowitsch { 13640244768SBarry Smith PetscInt i, irow_n, icol_n, *sxadj, *sadjncy, *svalues; 13740244768SBarry Smith PetscInt *indices, nindx, j, k, loc; 13840244768SBarry Smith PetscMPIInt issame; 13940244768SBarry Smith const PetscInt *irow_indices, *icol_indices; 14040244768SBarry Smith MPI_Comm scomm_row, scomm_col, scomm_mat; 14140244768SBarry Smith 14240244768SBarry Smith PetscFunctionBegin; 14340244768SBarry Smith nindx = 0; 14440244768SBarry Smith /* 14540244768SBarry Smith * Estimate a maximum number for allocating memory 14640244768SBarry Smith */ 14740244768SBarry Smith for (i = 0; i < n; i++) { 1489566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(irow[i], &irow_n)); 1499566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(icol[i], &icol_n)); 15040244768SBarry Smith nindx = nindx > (irow_n + icol_n) ? nindx : (irow_n + icol_n); 15140244768SBarry Smith } 1529566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nindx, &indices)); 15340244768SBarry Smith /* construct a submat */ 154252a1336SBarry Smith // if (scall == MAT_INITIAL_MATRIX) PetscCall(PetscMalloc1(n,submat)); 1556a7b62d2SBarry Smith 15640244768SBarry Smith for (i = 0; i < n; i++) { 15740244768SBarry Smith if (subcomm) { 1589566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)irow[i], &scomm_row)); 1599566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)icol[i], &scomm_col)); 1609566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_compare(scomm_row, scomm_col, &issame)); 16108401ef6SPierre Jolivet PetscCheck(issame == MPI_IDENT, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "row index set must have the same comm as the col index set"); 1629566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_compare(scomm_row, PETSC_COMM_SELF, &issame)); 16308401ef6SPierre Jolivet PetscCheck(issame != MPI_IDENT, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, " can not use PETSC_COMM_SELF as comm when extracting a parallel submatrix"); 16440244768SBarry Smith } else { 16540244768SBarry Smith scomm_row = PETSC_COMM_SELF; 16640244768SBarry Smith } 16740244768SBarry Smith /*get sub-matrix data*/ 1689371c9d4SSatish Balay sxadj = NULL; 1699371c9d4SSatish Balay sadjncy = NULL; 1709371c9d4SSatish Balay svalues = NULL; 1719566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix_MPIAdj_data(mat, irow[i], icol[i], &sxadj, &sadjncy, &svalues)); 1729566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(irow[i], &irow_n)); 1739566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(icol[i], &icol_n)); 1749566063dSJacob Faibussowitsch PetscCall(ISGetIndices(irow[i], &irow_indices)); 1759566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(indices, irow_indices, irow_n)); 1769566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(irow[i], &irow_indices)); 1779566063dSJacob Faibussowitsch PetscCall(ISGetIndices(icol[i], &icol_indices)); 1789566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(indices + irow_n, icol_indices, icol_n)); 1799566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(icol[i], &icol_indices)); 18040244768SBarry Smith nindx = irow_n + icol_n; 1819566063dSJacob Faibussowitsch PetscCall(PetscSortRemoveDupsInt(&nindx, indices)); 18240244768SBarry Smith /* renumber columns */ 18340244768SBarry Smith for (j = 0; j < irow_n; j++) { 18440244768SBarry Smith for (k = sxadj[j]; k < sxadj[j + 1]; k++) { 1859566063dSJacob Faibussowitsch PetscCall(PetscFindInt(sadjncy[k], nindx, indices, &loc)); 18608401ef6SPierre Jolivet PetscCheck(loc >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "can not find col %" PetscInt_FMT, sadjncy[k]); 18740244768SBarry Smith sadjncy[k] = loc; 18840244768SBarry Smith } 18940244768SBarry Smith } 19040244768SBarry Smith if (scall == MAT_INITIAL_MATRIX) { 1919566063dSJacob Faibussowitsch PetscCall(MatCreateMPIAdj(scomm_row, irow_n, icol_n, sxadj, sadjncy, svalues, submat[i])); 19240244768SBarry Smith } else { 19340244768SBarry Smith Mat sadj = *(submat[i]); 19440244768SBarry Smith Mat_MPIAdj *sa = (Mat_MPIAdj *)((sadj)->data); 1959566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)sadj, &scomm_mat)); 1969566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_compare(scomm_row, scomm_mat, &issame)); 19708401ef6SPierre Jolivet PetscCheck(issame == MPI_IDENT, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "submatrix must have the same comm as the col index set"); 1989566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(sa->i, sxadj, irow_n + 1)); 1999566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(sa->j, sadjncy, sxadj[irow_n])); 2009566063dSJacob Faibussowitsch if (svalues) PetscCall(PetscArraycpy(sa->values, svalues, sxadj[irow_n])); 2019566063dSJacob Faibussowitsch PetscCall(PetscFree(sxadj)); 2029566063dSJacob Faibussowitsch PetscCall(PetscFree(sadjncy)); 203252a1336SBarry Smith PetscCall(PetscFree(svalues)); 20440244768SBarry Smith } 20540244768SBarry Smith } 2069566063dSJacob Faibussowitsch PetscCall(PetscFree(indices)); 2073ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 20840244768SBarry Smith } 20940244768SBarry Smith 210d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatCreateSubMatricesMPI_MPIAdj(Mat mat, PetscInt n, const IS irow[], const IS icol[], MatReuse scall, Mat *submat[]) 211d71ae5a4SJacob Faibussowitsch { 21240244768SBarry Smith /*get sub-matrices across a sub communicator */ 21340244768SBarry Smith PetscFunctionBegin; 2149566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices_MPIAdj_Private(mat, n, irow, icol, PETSC_TRUE, scall, submat)); 2153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 21640244768SBarry Smith } 21740244768SBarry Smith 218d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatCreateSubMatrices_MPIAdj(Mat mat, PetscInt n, const IS irow[], const IS icol[], MatReuse scall, Mat *submat[]) 219d71ae5a4SJacob Faibussowitsch { 22040244768SBarry Smith PetscFunctionBegin; 22140244768SBarry Smith /*get sub-matrices based on PETSC_COMM_SELF */ 2229566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices_MPIAdj_Private(mat, n, irow, icol, PETSC_FALSE, scall, submat)); 2233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 22440244768SBarry Smith } 22540244768SBarry Smith 226d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatView_MPIAdj_ASCII(Mat A, PetscViewer viewer) 227d71ae5a4SJacob Faibussowitsch { 2283eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj *)A->data; 229d0f46423SBarry Smith PetscInt i, j, m = A->rmap->n; 2302dcb1b2aSMatthew Knepley const char *name; 231ce1411ecSBarry Smith PetscViewerFormat format; 232b97cf49bSBarry Smith 233433994e6SBarry Smith PetscFunctionBegin; 2349566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName((PetscObject)A, &name)); 2359566063dSJacob Faibussowitsch PetscCall(PetscViewerGetFormat(viewer, &format)); 236fb9695e5SSatish Balay if (format == PETSC_VIEWER_ASCII_INFO) { 2373ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 238f7d195e4SLawrence Mitchell } else { 239f7d195e4SLawrence Mitchell PetscCheck(format != PETSC_VIEWER_ASCII_MATLAB, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "MATLAB format not supported"); 2409566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE)); 2419566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 242b97cf49bSBarry Smith for (i = 0; i < m; i++) { 2439566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "row %" PetscInt_FMT ":", i + A->rmap->rstart)); 244b97cf49bSBarry Smith for (j = a->i[i]; j < a->i[i + 1]; j++) { 2450170919eSFande Kong if (a->values) { 2469566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " (%" PetscInt_FMT ", %" PetscInt_FMT ") ", a->j[j], a->values[j])); 2470170919eSFande Kong } else { 2489566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " %" PetscInt_FMT " ", a->j[j])); 2490752156aSBarry Smith } 2500170919eSFande Kong } 2519566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "\n")); 252b97cf49bSBarry Smith } 2539566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE)); 2549566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 2559566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 2567b23a99aSBarry Smith } 2573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 258b97cf49bSBarry Smith } 259b97cf49bSBarry Smith 260d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatView_MPIAdj(Mat A, PetscViewer viewer) 261d71ae5a4SJacob Faibussowitsch { 262ace3abfcSBarry Smith PetscBool iascii; 263b97cf49bSBarry Smith 264433994e6SBarry Smith PetscFunctionBegin; 2659566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 2661baa6e33SBarry Smith if (iascii) PetscCall(MatView_MPIAdj_ASCII(A, viewer)); 2673ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 268b97cf49bSBarry Smith } 269b97cf49bSBarry Smith 270d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatDestroy_MPIAdj(Mat mat) 271d71ae5a4SJacob Faibussowitsch { 2723eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj *)mat->data; 27394d884c6SBarry Smith 27494d884c6SBarry Smith PetscFunctionBegin; 2753ba16761SJacob Faibussowitsch PetscCall(PetscLogObjectState((PetscObject)mat, "Rows=%" PetscInt_FMT ", Cols=%" PetscInt_FMT ", NZ=%" PetscInt_FMT, mat->rmap->n, mat->cmap->n, a->nz)); 2769566063dSJacob Faibussowitsch PetscCall(PetscFree(a->diag)); 2770400557aSBarry Smith if (a->freeaij) { 27814196391SBarry Smith if (a->freeaijwithfree) { 27914196391SBarry Smith if (a->i) free(a->i); 28014196391SBarry Smith if (a->j) free(a->j); 28114196391SBarry Smith } else { 2829566063dSJacob Faibussowitsch PetscCall(PetscFree(a->i)); 2839566063dSJacob Faibussowitsch PetscCall(PetscFree(a->j)); 2849566063dSJacob Faibussowitsch PetscCall(PetscFree(a->values)); 28514196391SBarry Smith } 2860400557aSBarry Smith } 2879566063dSJacob Faibussowitsch PetscCall(PetscFree(a->rowvalues)); 2889566063dSJacob Faibussowitsch PetscCall(PetscFree(mat->data)); 2899566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)mat, NULL)); 2909566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMPIAdjSetPreallocation_C", NULL)); 2919566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMPIAdjCreateNonemptySubcommMat_C", NULL)); 2922e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMPIAdjToSeq_C", NULL)); 293252a1336SBarry Smith PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMPIAdjToSeqRankZero_C", NULL)); 2943ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 295b97cf49bSBarry Smith } 296b97cf49bSBarry Smith 297d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSetOption_MPIAdj(Mat A, MatOption op, PetscBool flg) 298d71ae5a4SJacob Faibussowitsch { 2993eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj *)A->data; 300b97cf49bSBarry Smith 301433994e6SBarry Smith PetscFunctionBegin; 30212c028f9SKris Buschelman switch (op) { 30377e54ba9SKris Buschelman case MAT_SYMMETRIC: 30412c028f9SKris Buschelman case MAT_STRUCTURALLY_SYMMETRIC: 3059a4540c5SBarry Smith case MAT_HERMITIAN: 306d71ae5a4SJacob Faibussowitsch case MAT_SPD: 307d71ae5a4SJacob Faibussowitsch a->symmetric = flg; 308d71ae5a4SJacob Faibussowitsch break; 3099a4540c5SBarry Smith case MAT_SYMMETRY_ETERNAL: 310b94d7dedSBarry Smith case MAT_STRUCTURAL_SYMMETRY_ETERNAL: 311d71ae5a4SJacob Faibussowitsch case MAT_SPD_ETERNAL: 312d71ae5a4SJacob Faibussowitsch break; 313d71ae5a4SJacob Faibussowitsch default: 314d71ae5a4SJacob Faibussowitsch PetscCall(PetscInfo(A, "Option %s ignored\n", MatOptions[op])); 315d71ae5a4SJacob Faibussowitsch break; 316b97cf49bSBarry Smith } 3173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 318b97cf49bSBarry Smith } 319b97cf49bSBarry Smith 320d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatGetRow_MPIAdj(Mat A, PetscInt row, PetscInt *nz, PetscInt **idx, PetscScalar **v) 321d71ae5a4SJacob Faibussowitsch { 3223eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj *)A->data; 323b97cf49bSBarry Smith 324433994e6SBarry Smith PetscFunctionBegin; 325d0f46423SBarry Smith row -= A->rmap->rstart; 326aed4548fSBarry Smith PetscCheck(row >= 0 && row < A->rmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Row out of range"); 327b97cf49bSBarry Smith *nz = a->i[row + 1] - a->i[row]; 3282b1d8763SJed Brown if (v) { 3292b1d8763SJed Brown PetscInt j; 3302b1d8763SJed Brown if (a->rowvalues_alloc < *nz) { 3319566063dSJacob Faibussowitsch PetscCall(PetscFree(a->rowvalues)); 3322b1d8763SJed Brown a->rowvalues_alloc = PetscMax(a->rowvalues_alloc * 2, *nz); 3339566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(a->rowvalues_alloc, &a->rowvalues)); 334b97cf49bSBarry Smith } 335ad540459SPierre Jolivet for (j = 0; j < *nz; j++) a->rowvalues[j] = a->values ? a->values[a->i[row] + j] : 1.0; 3362b1d8763SJed Brown *v = (*nz) ? a->rowvalues : NULL; 337b97cf49bSBarry Smith } 33892b6f2a0SJed Brown if (idx) *idx = (*nz) ? a->j + a->i[row] : NULL; 3393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 340b97cf49bSBarry Smith } 341b97cf49bSBarry Smith 342d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatRestoreRow_MPIAdj(Mat A, PetscInt row, PetscInt *nz, PetscInt **idx, PetscScalar **v) 343d71ae5a4SJacob Faibussowitsch { 344433994e6SBarry Smith PetscFunctionBegin; 3453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 346b97cf49bSBarry Smith } 347b97cf49bSBarry Smith 348d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatEqual_MPIAdj(Mat A, Mat B, PetscBool *flg) 349d71ae5a4SJacob Faibussowitsch { 3503eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj *)A->data, *b = (Mat_MPIAdj *)B->data; 351ace3abfcSBarry Smith PetscBool flag; 352b97cf49bSBarry Smith 353433994e6SBarry Smith PetscFunctionBegin; 354b97cf49bSBarry Smith /* If the matrix dimensions are not equal,or no of nonzeros */ 355ad540459SPierre Jolivet if ((A->rmap->n != B->rmap->n) || (a->nz != b->nz)) flag = PETSC_FALSE; 356b97cf49bSBarry Smith 357b97cf49bSBarry Smith /* if the a->i are the same */ 3589566063dSJacob Faibussowitsch PetscCall(PetscArraycmp(a->i, b->i, A->rmap->n + 1, &flag)); 359b97cf49bSBarry Smith 360b97cf49bSBarry Smith /* if a->j are the same */ 3619566063dSJacob Faibussowitsch PetscCall(PetscMemcmp(a->j, b->j, (a->nz) * sizeof(PetscInt), &flag)); 362b97cf49bSBarry Smith 3631c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(&flag, flg, 1, MPIU_BOOL, MPI_LAND, PetscObjectComm((PetscObject)A))); 3643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 365b97cf49bSBarry Smith } 366b97cf49bSBarry Smith 367d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatGetRowIJ_MPIAdj(Mat A, PetscInt oshift, PetscBool symmetric, PetscBool blockcompressed, PetscInt *m, const PetscInt *inia[], const PetscInt *inja[], PetscBool *done) 368d71ae5a4SJacob Faibussowitsch { 369b24ad042SBarry Smith PetscInt i; 3706cd91f64SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj *)A->data; 3711a83f524SJed Brown PetscInt **ia = (PetscInt **)inia, **ja = (PetscInt **)inja; 3726cd91f64SBarry Smith 3736cd91f64SBarry Smith PetscFunctionBegin; 374d0f46423SBarry Smith *m = A->rmap->n; 3756cd91f64SBarry Smith *ia = a->i; 3766cd91f64SBarry Smith *ja = a->j; 3776cd91f64SBarry Smith *done = PETSC_TRUE; 378d42735a1SBarry Smith if (oshift) { 379ad540459SPierre Jolivet for (i = 0; i < (*ia)[*m]; i++) (*ja)[i]++; 380d42735a1SBarry Smith for (i = 0; i <= (*m); i++) (*ia)[i]++; 381d42735a1SBarry Smith } 3823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 383d42735a1SBarry Smith } 384d42735a1SBarry Smith 385d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatRestoreRowIJ_MPIAdj(Mat A, PetscInt oshift, PetscBool symmetric, PetscBool blockcompressed, PetscInt *m, const PetscInt *inia[], const PetscInt *inja[], PetscBool *done) 386d71ae5a4SJacob Faibussowitsch { 387b24ad042SBarry Smith PetscInt i; 388d42735a1SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj *)A->data; 3891a83f524SJed Brown PetscInt **ia = (PetscInt **)inia, **ja = (PetscInt **)inja; 390d42735a1SBarry Smith 391d42735a1SBarry Smith PetscFunctionBegin; 392aed4548fSBarry Smith PetscCheck(!ia || a->i == *ia, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "ia passed back is not one obtained with MatGetRowIJ()"); 393aed4548fSBarry Smith PetscCheck(!ja || a->j == *ja, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "ja passed back is not one obtained with MatGetRowIJ()"); 394d42735a1SBarry Smith if (oshift) { 39528b400f6SJacob Faibussowitsch PetscCheck(ia, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "If oshift then you must passed in inia[] argument"); 39628b400f6SJacob Faibussowitsch PetscCheck(ja, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "If oshift then you must passed in inja[] argument"); 397d42735a1SBarry Smith for (i = 0; i <= (*m); i++) (*ia)[i]--; 398ad540459SPierre Jolivet for (i = 0; i < (*ia)[*m]; i++) (*ja)[i]--; 399d42735a1SBarry Smith } 4003ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4016cd91f64SBarry Smith } 402b97cf49bSBarry Smith 403*66976f2fSJacob Faibussowitsch static PetscErrorCode MatConvertFrom_MPIAdj(Mat A, MatType type, MatReuse reuse, Mat *newmat) 404d71ae5a4SJacob Faibussowitsch { 40517667f90SBarry Smith Mat B; 40617667f90SBarry Smith PetscInt i, m, N, nzeros = 0, *ia, *ja, len, rstart, cnt, j, *a; 40717667f90SBarry Smith const PetscInt *rj; 40817667f90SBarry Smith const PetscScalar *ra; 40917667f90SBarry Smith MPI_Comm comm; 41017667f90SBarry Smith 41117667f90SBarry Smith PetscFunctionBegin; 4129566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, NULL, &N)); 4139566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A, &m, NULL)); 4149566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A, &rstart, NULL)); 41517667f90SBarry Smith 41617667f90SBarry Smith /* count the number of nonzeros per row */ 41717667f90SBarry Smith for (i = 0; i < m; i++) { 4189566063dSJacob Faibussowitsch PetscCall(MatGetRow(A, i + rstart, &len, &rj, NULL)); 4195ee9ba1cSJed Brown for (j = 0; j < len; j++) { 4209371c9d4SSatish Balay if (rj[j] == i + rstart) { 4219371c9d4SSatish Balay len--; 4229371c9d4SSatish Balay break; 4239371c9d4SSatish Balay } /* don't count diagonal */ 4245ee9ba1cSJed Brown } 42517667f90SBarry Smith nzeros += len; 4269566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(A, i + rstart, &len, &rj, NULL)); 42717667f90SBarry Smith } 42817667f90SBarry Smith 42917667f90SBarry Smith /* malloc space for nonzeros */ 4309566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nzeros + 1, &a)); 4319566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(N + 1, &ia)); 4329566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nzeros + 1, &ja)); 43317667f90SBarry Smith 43417667f90SBarry Smith nzeros = 0; 43517667f90SBarry Smith ia[0] = 0; 43617667f90SBarry Smith for (i = 0; i < m; i++) { 4379566063dSJacob Faibussowitsch PetscCall(MatGetRow(A, i + rstart, &len, &rj, &ra)); 43817667f90SBarry Smith cnt = 0; 43917667f90SBarry Smith for (j = 0; j < len; j++) { 4405ee9ba1cSJed Brown if (rj[j] != i + rstart) { /* if not diagonal */ 44117667f90SBarry Smith a[nzeros + cnt] = (PetscInt)PetscAbsScalar(ra[j]); 44217667f90SBarry Smith ja[nzeros + cnt++] = rj[j]; 44317667f90SBarry Smith } 4445ee9ba1cSJed Brown } 4459566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(A, i + rstart, &len, &rj, &ra)); 44617667f90SBarry Smith nzeros += cnt; 44717667f90SBarry Smith ia[i + 1] = nzeros; 44817667f90SBarry Smith } 44917667f90SBarry Smith 4509566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)A, &comm)); 4519566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, &B)); 4529566063dSJacob Faibussowitsch PetscCall(MatSetSizes(B, m, PETSC_DETERMINE, PETSC_DETERMINE, N)); 4539566063dSJacob Faibussowitsch PetscCall(MatSetType(B, type)); 4549566063dSJacob Faibussowitsch PetscCall(MatMPIAdjSetPreallocation(B, ia, ja, a)); 45517667f90SBarry Smith 456511c6705SHong Zhang if (reuse == MAT_INPLACE_MATRIX) { 4579566063dSJacob Faibussowitsch PetscCall(MatHeaderReplace(A, &B)); 45817667f90SBarry Smith } else { 45917667f90SBarry Smith *newmat = B; 46017667f90SBarry Smith } 4613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 46217667f90SBarry Smith } 46317667f90SBarry Smith 464*66976f2fSJacob Faibussowitsch static PetscErrorCode MatSetValues_MPIAdj(Mat A, PetscInt m, const PetscInt *rows, PetscInt n, const PetscInt *cols, const PetscScalar *values, InsertMode im) 465d71ae5a4SJacob Faibussowitsch { 4666a09307cSBarry Smith Mat_MPIAdj *adj = (Mat_MPIAdj *)A->data; 4676a09307cSBarry Smith PetscInt rStart, rEnd, cStart, cEnd; 4686a09307cSBarry Smith 4696a09307cSBarry Smith PetscFunctionBegin; 4706a09307cSBarry Smith PetscCheck(!adj->i, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Matrix is already assembled, cannot change its values"); 4716a09307cSBarry Smith PetscCall(MatGetOwnershipRange(A, &rStart, &rEnd)); 4726a09307cSBarry Smith PetscCall(MatGetOwnershipRangeColumn(A, &cStart, &cEnd)); 4736a09307cSBarry Smith if (!adj->ht) { 4746a09307cSBarry Smith PetscCall(PetscHSetIJCreate(&adj->ht)); 4756a09307cSBarry Smith PetscCall(MatStashCreate_Private(PetscObjectComm((PetscObject)A), 1, &A->stash)); 4766a09307cSBarry Smith PetscCall(PetscLayoutSetUp(A->rmap)); 4776a09307cSBarry Smith PetscCall(PetscLayoutSetUp(A->cmap)); 4786a09307cSBarry Smith } 4796a09307cSBarry Smith for (PetscInt r = 0; r < m; ++r) { 4806a09307cSBarry Smith PetscHashIJKey key; 4816a09307cSBarry Smith 4826a09307cSBarry Smith key.i = rows[r]; 4836a09307cSBarry Smith if (key.i < 0) continue; 4846a09307cSBarry Smith if ((key.i < rStart) || (key.i >= rEnd)) { 4856a09307cSBarry Smith PetscCall(MatStashValuesRow_Private(&A->stash, key.i, n, cols, values, PETSC_FALSE)); 4866a09307cSBarry Smith } else { 4876a09307cSBarry Smith for (PetscInt c = 0; c < n; ++c) { 4886a09307cSBarry Smith key.j = cols[c]; 4896a09307cSBarry Smith if (key.j < 0 || key.i == key.j) continue; 4906a09307cSBarry Smith PetscCall(PetscHSetIJAdd(adj->ht, key)); 4916a09307cSBarry Smith } 4926a09307cSBarry Smith } 4936a09307cSBarry Smith } 4943ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4956a09307cSBarry Smith } 4966a09307cSBarry Smith 497*66976f2fSJacob Faibussowitsch static PetscErrorCode MatAssemblyBegin_MPIAdj(Mat A, MatAssemblyType type) 498d71ae5a4SJacob Faibussowitsch { 4996a09307cSBarry Smith PetscInt nstash, reallocs; 5006a09307cSBarry Smith Mat_MPIAdj *adj = (Mat_MPIAdj *)A->data; 5016a09307cSBarry Smith 5026a09307cSBarry Smith PetscFunctionBegin; 5036a09307cSBarry Smith if (!adj->ht) { 5046a09307cSBarry Smith PetscCall(PetscHSetIJCreate(&adj->ht)); 5056a09307cSBarry Smith PetscCall(PetscLayoutSetUp(A->rmap)); 5066a09307cSBarry Smith PetscCall(PetscLayoutSetUp(A->cmap)); 5076a09307cSBarry Smith } 5086a09307cSBarry Smith PetscCall(MatStashScatterBegin_Private(A, &A->stash, A->rmap->range)); 5096a09307cSBarry Smith PetscCall(MatStashGetInfo_Private(&A->stash, &nstash, &reallocs)); 5106a09307cSBarry Smith PetscCall(PetscInfo(A, "Stash has %" PetscInt_FMT " entries, uses %" PetscInt_FMT " mallocs.\n", nstash, reallocs)); 5113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5126a09307cSBarry Smith } 5136a09307cSBarry Smith 514*66976f2fSJacob Faibussowitsch static PetscErrorCode MatAssemblyEnd_MPIAdj(Mat A, MatAssemblyType type) 515d71ae5a4SJacob Faibussowitsch { 5166a09307cSBarry Smith PetscScalar *val; 5176a09307cSBarry Smith PetscInt *row, *col, m, rstart, *rowstarts; 5186a09307cSBarry Smith PetscInt i, j, ncols, flg, nz; 5196a09307cSBarry Smith PetscMPIInt n; 5206a09307cSBarry Smith Mat_MPIAdj *adj = (Mat_MPIAdj *)A->data; 5216a09307cSBarry Smith PetscHashIter hi; 5226a09307cSBarry Smith PetscHashIJKey key; 5236a09307cSBarry Smith PetscHSetIJ ht = adj->ht; 5246a09307cSBarry Smith 5256a09307cSBarry Smith PetscFunctionBegin; 5266a09307cSBarry Smith while (1) { 5276a09307cSBarry Smith PetscCall(MatStashScatterGetMesg_Private(&A->stash, &n, &row, &col, &val, &flg)); 5286a09307cSBarry Smith if (!flg) break; 5296a09307cSBarry Smith 5306a09307cSBarry Smith for (i = 0; i < n;) { 5316a09307cSBarry Smith /* Identify the consecutive vals belonging to the same row */ 5326a09307cSBarry Smith for (j = i, rstart = row[j]; j < n; j++) { 5336a09307cSBarry Smith if (row[j] != rstart) break; 5346a09307cSBarry Smith } 5356a09307cSBarry Smith if (j < n) ncols = j - i; 5366a09307cSBarry Smith else ncols = n - i; 5376a09307cSBarry Smith /* Set all these values with a single function call */ 5386a09307cSBarry Smith PetscCall(MatSetValues_MPIAdj(A, 1, row + i, ncols, col + i, val + i, INSERT_VALUES)); 5396a09307cSBarry Smith i = j; 5406a09307cSBarry Smith } 5416a09307cSBarry Smith } 5426a09307cSBarry Smith PetscCall(MatStashScatterEnd_Private(&A->stash)); 5436a09307cSBarry Smith PetscCall(MatStashDestroy_Private(&A->stash)); 5446a09307cSBarry Smith 5456a09307cSBarry Smith PetscCall(MatGetLocalSize(A, &m, NULL)); 5466a09307cSBarry Smith PetscCall(MatGetOwnershipRange(A, &rstart, NULL)); 5476a09307cSBarry Smith PetscCall(PetscCalloc1(m + 1, &rowstarts)); 5486a09307cSBarry Smith PetscHashIterBegin(ht, hi); 5496a09307cSBarry Smith for (; !PetscHashIterAtEnd(ht, hi);) { 5506a09307cSBarry Smith PetscHashIterGetKey(ht, hi, key); 5516a09307cSBarry Smith rowstarts[key.i - rstart + 1]++; 5526a09307cSBarry Smith PetscHashIterNext(ht, hi); 5536a09307cSBarry Smith } 554ad540459SPierre Jolivet for (i = 1; i < m + 1; i++) rowstarts[i] = rowstarts[i - 1] + rowstarts[i]; 5556a09307cSBarry Smith 5566a09307cSBarry Smith PetscCall(PetscHSetIJGetSize(ht, &nz)); 5576a09307cSBarry Smith PetscCall(PetscMalloc1(nz, &col)); 5586a09307cSBarry Smith PetscHashIterBegin(ht, hi); 5596a09307cSBarry Smith for (; !PetscHashIterAtEnd(ht, hi);) { 5606a09307cSBarry Smith PetscHashIterGetKey(ht, hi, key); 5616a09307cSBarry Smith col[rowstarts[key.i - rstart]++] = key.j; 5626a09307cSBarry Smith PetscHashIterNext(ht, hi); 5636a09307cSBarry Smith } 5646a09307cSBarry Smith PetscCall(PetscHSetIJDestroy(&ht)); 565ad540459SPierre Jolivet for (i = m; i > 0; i--) rowstarts[i] = rowstarts[i - 1]; 5666a09307cSBarry Smith rowstarts[0] = 0; 5676a09307cSBarry Smith 56848a46eb9SPierre Jolivet for (PetscInt i = 0; i < m; i++) PetscCall(PetscSortInt(rowstarts[i + 1] - rowstarts[i], &col[rowstarts[i]])); 5696a09307cSBarry Smith 5706a09307cSBarry Smith adj->i = rowstarts; 5716a09307cSBarry Smith adj->j = col; 572252a1336SBarry Smith adj->nz = rowstarts[m]; 5736a09307cSBarry Smith adj->freeaij = PETSC_TRUE; 5743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5756a09307cSBarry Smith } 5766a09307cSBarry Smith 5776a09307cSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_MPIAdj, 5783eda8832SBarry Smith MatGetRow_MPIAdj, 5793eda8832SBarry Smith MatRestoreRow_MPIAdj, 580f4259b30SLisandro Dalcin NULL, 581f4259b30SLisandro Dalcin /* 4*/ NULL, 582f4259b30SLisandro Dalcin NULL, 583f4259b30SLisandro Dalcin NULL, 584f4259b30SLisandro Dalcin NULL, 585f4259b30SLisandro Dalcin NULL, 586f4259b30SLisandro Dalcin NULL, 587f4259b30SLisandro Dalcin /*10*/ NULL, 588f4259b30SLisandro Dalcin NULL, 589f4259b30SLisandro Dalcin NULL, 590f4259b30SLisandro Dalcin NULL, 591f4259b30SLisandro Dalcin NULL, 592f4259b30SLisandro Dalcin /*15*/ NULL, 5933eda8832SBarry Smith MatEqual_MPIAdj, 594f4259b30SLisandro Dalcin NULL, 595f4259b30SLisandro Dalcin NULL, 596f4259b30SLisandro Dalcin NULL, 5976a09307cSBarry Smith /*20*/ MatAssemblyBegin_MPIAdj, 5986a09307cSBarry Smith MatAssemblyEnd_MPIAdj, 5993eda8832SBarry Smith MatSetOption_MPIAdj, 600f4259b30SLisandro Dalcin NULL, 601f4259b30SLisandro Dalcin /*24*/ NULL, 602f4259b30SLisandro Dalcin NULL, 603f4259b30SLisandro Dalcin NULL, 604f4259b30SLisandro Dalcin NULL, 605f4259b30SLisandro Dalcin NULL, 606f4259b30SLisandro Dalcin /*29*/ NULL, 607f4259b30SLisandro Dalcin NULL, 608f4259b30SLisandro Dalcin NULL, 609f4259b30SLisandro Dalcin NULL, 610f4259b30SLisandro Dalcin NULL, 611f4259b30SLisandro Dalcin /*34*/ NULL, 612f4259b30SLisandro Dalcin NULL, 613f4259b30SLisandro Dalcin NULL, 614f4259b30SLisandro Dalcin NULL, 615f4259b30SLisandro Dalcin NULL, 616f4259b30SLisandro Dalcin /*39*/ NULL, 6177dae84e0SHong Zhang MatCreateSubMatrices_MPIAdj, 618f4259b30SLisandro Dalcin NULL, 619f4259b30SLisandro Dalcin NULL, 620f4259b30SLisandro Dalcin NULL, 621f4259b30SLisandro Dalcin /*44*/ NULL, 622f4259b30SLisandro Dalcin NULL, 6237d68702bSBarry Smith MatShift_Basic, 624f4259b30SLisandro Dalcin NULL, 625f4259b30SLisandro Dalcin NULL, 626f4259b30SLisandro Dalcin /*49*/ NULL, 6276cd91f64SBarry Smith MatGetRowIJ_MPIAdj, 628d42735a1SBarry Smith MatRestoreRowIJ_MPIAdj, 629f4259b30SLisandro Dalcin NULL, 630f4259b30SLisandro Dalcin NULL, 631f4259b30SLisandro Dalcin /*54*/ NULL, 632f4259b30SLisandro Dalcin NULL, 633f4259b30SLisandro Dalcin NULL, 634f4259b30SLisandro Dalcin NULL, 635f4259b30SLisandro Dalcin NULL, 636f4259b30SLisandro Dalcin /*59*/ NULL, 637b9b97703SBarry Smith MatDestroy_MPIAdj, 638b9b97703SBarry Smith MatView_MPIAdj, 63917667f90SBarry Smith MatConvertFrom_MPIAdj, 640f4259b30SLisandro Dalcin NULL, 641f4259b30SLisandro Dalcin /*64*/ NULL, 642f4259b30SLisandro Dalcin NULL, 643f4259b30SLisandro Dalcin NULL, 644f4259b30SLisandro Dalcin NULL, 645f4259b30SLisandro Dalcin NULL, 646f4259b30SLisandro Dalcin /*69*/ NULL, 647f4259b30SLisandro Dalcin NULL, 648f4259b30SLisandro Dalcin NULL, 649f4259b30SLisandro Dalcin NULL, 650f4259b30SLisandro Dalcin NULL, 651f4259b30SLisandro Dalcin /*74*/ NULL, 652f4259b30SLisandro Dalcin NULL, 653f4259b30SLisandro Dalcin NULL, 654f4259b30SLisandro Dalcin NULL, 655f4259b30SLisandro Dalcin NULL, 656f4259b30SLisandro Dalcin /*79*/ NULL, 657f4259b30SLisandro Dalcin NULL, 658f4259b30SLisandro Dalcin NULL, 659f4259b30SLisandro Dalcin NULL, 660f4259b30SLisandro Dalcin NULL, 661f4259b30SLisandro Dalcin /*84*/ NULL, 662f4259b30SLisandro Dalcin NULL, 663f4259b30SLisandro Dalcin NULL, 664f4259b30SLisandro Dalcin NULL, 665f4259b30SLisandro Dalcin NULL, 666f4259b30SLisandro Dalcin /*89*/ NULL, 667f4259b30SLisandro Dalcin NULL, 668f4259b30SLisandro Dalcin NULL, 669f4259b30SLisandro Dalcin NULL, 670f4259b30SLisandro Dalcin NULL, 671f4259b30SLisandro Dalcin /*94*/ NULL, 672f4259b30SLisandro Dalcin NULL, 673f4259b30SLisandro Dalcin NULL, 674f4259b30SLisandro Dalcin NULL, 675f4259b30SLisandro Dalcin NULL, 676f4259b30SLisandro Dalcin /*99*/ NULL, 677f4259b30SLisandro Dalcin NULL, 678f4259b30SLisandro Dalcin NULL, 679f4259b30SLisandro Dalcin NULL, 680f4259b30SLisandro Dalcin NULL, 681f4259b30SLisandro Dalcin /*104*/ NULL, 682f4259b30SLisandro Dalcin NULL, 683f4259b30SLisandro Dalcin NULL, 684f4259b30SLisandro Dalcin NULL, 685f4259b30SLisandro Dalcin NULL, 686f4259b30SLisandro Dalcin /*109*/ NULL, 687f4259b30SLisandro Dalcin NULL, 688f4259b30SLisandro Dalcin NULL, 689f4259b30SLisandro Dalcin NULL, 690f4259b30SLisandro Dalcin NULL, 691f4259b30SLisandro Dalcin /*114*/ NULL, 692f4259b30SLisandro Dalcin NULL, 693f4259b30SLisandro Dalcin NULL, 694f4259b30SLisandro Dalcin NULL, 695f4259b30SLisandro Dalcin NULL, 696f4259b30SLisandro Dalcin /*119*/ NULL, 697f4259b30SLisandro Dalcin NULL, 698f4259b30SLisandro Dalcin NULL, 699f4259b30SLisandro Dalcin NULL, 700f4259b30SLisandro Dalcin NULL, 701f4259b30SLisandro Dalcin /*124*/ NULL, 702f4259b30SLisandro Dalcin NULL, 703f4259b30SLisandro Dalcin NULL, 704f4259b30SLisandro Dalcin NULL, 7057dae84e0SHong Zhang MatCreateSubMatricesMPI_MPIAdj, 706f4259b30SLisandro Dalcin /*129*/ NULL, 707f4259b30SLisandro Dalcin NULL, 708f4259b30SLisandro Dalcin NULL, 709f4259b30SLisandro Dalcin NULL, 710f4259b30SLisandro Dalcin NULL, 711f4259b30SLisandro Dalcin /*134*/ NULL, 712f4259b30SLisandro Dalcin NULL, 713f4259b30SLisandro Dalcin NULL, 714f4259b30SLisandro Dalcin NULL, 715f4259b30SLisandro Dalcin NULL, 716f4259b30SLisandro Dalcin /*139*/ NULL, 717f4259b30SLisandro Dalcin NULL, 718d70f29a3SPierre Jolivet NULL, 719d70f29a3SPierre Jolivet NULL, 720d70f29a3SPierre Jolivet NULL, 721d70f29a3SPierre Jolivet /*144*/ NULL, 722d70f29a3SPierre Jolivet NULL, 723d70f29a3SPierre Jolivet NULL, 72499a7f59eSMark Adams NULL, 72599a7f59eSMark Adams NULL, 7267fb60732SBarry Smith NULL, 727dec0b466SHong Zhang /*150*/ NULL, 728dec0b466SHong Zhang NULL}; 729b97cf49bSBarry Smith 730d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMPIAdjSetPreallocation_MPIAdj(Mat B, PetscInt *i, PetscInt *j, PetscInt *values) 731d71ae5a4SJacob Faibussowitsch { 732a23d5eceSKris Buschelman Mat_MPIAdj *b = (Mat_MPIAdj *)B->data; 733e895ccc0SFande Kong PetscBool useedgeweights; 734a23d5eceSKris Buschelman 735a23d5eceSKris Buschelman PetscFunctionBegin; 7369566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(B->rmap)); 7379566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(B->cmap)); 7389371c9d4SSatish Balay if (values) useedgeweights = PETSC_TRUE; 7399371c9d4SSatish Balay else useedgeweights = PETSC_FALSE; 740e895ccc0SFande Kong /* Make everybody knows if they are using edge weights or not */ 7411c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce((int *)&useedgeweights, (int *)&b->useedgeweights, 1, MPI_INT, MPI_MAX, PetscObjectComm((PetscObject)B))); 74258ec18a5SLisandro Dalcin 74376bd3646SJed Brown if (PetscDefined(USE_DEBUG)) { 74476bd3646SJed Brown PetscInt ii; 74576bd3646SJed Brown 74608401ef6SPierre Jolivet PetscCheck(i[0] == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "First i[] index must be zero, instead it is %" PetscInt_FMT, i[0]); 747d0f46423SBarry Smith for (ii = 1; ii < B->rmap->n; ii++) { 748aed4548fSBarry Smith PetscCheck(i[ii] >= 0 && i[ii] >= i[ii - 1], PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "i[%" PetscInt_FMT "]=%" PetscInt_FMT " index is out of range: i[%" PetscInt_FMT "]=%" PetscInt_FMT, ii, i[ii], ii - 1, i[ii - 1]); 749a23d5eceSKris Buschelman } 750ad540459SPierre Jolivet for (ii = 0; ii < i[B->rmap->n]; ii++) PetscCheck(j[ii] >= 0 && j[ii] < B->cmap->N, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Column index %" PetscInt_FMT " out of range %" PetscInt_FMT, ii, j[ii]); 75176bd3646SJed Brown } 752a23d5eceSKris Buschelman b->j = j; 753a23d5eceSKris Buschelman b->i = i; 754a23d5eceSKris Buschelman b->values = values; 755a23d5eceSKris Buschelman 756d0f46423SBarry Smith b->nz = i[B->rmap->n]; 757f4259b30SLisandro Dalcin b->diag = NULL; 758a23d5eceSKris Buschelman b->symmetric = PETSC_FALSE; 759a23d5eceSKris Buschelman b->freeaij = PETSC_TRUE; 760a23d5eceSKris Buschelman 7616a09307cSBarry Smith B->ops->assemblybegin = NULL; 7626a09307cSBarry Smith B->ops->assemblyend = NULL; 7639566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 7649566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 7656a09307cSBarry Smith PetscCall(MatStashDestroy_Private(&B->stash)); 7663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 767a23d5eceSKris Buschelman } 768b97cf49bSBarry Smith 769d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMPIAdjCreateNonemptySubcommMat_MPIAdj(Mat A, Mat *B) 770d71ae5a4SJacob Faibussowitsch { 7719a3665c6SJed Brown Mat_MPIAdj *a = (Mat_MPIAdj *)A->data; 7729a3665c6SJed Brown const PetscInt *ranges; 7739a3665c6SJed Brown MPI_Comm acomm, bcomm; 7749a3665c6SJed Brown MPI_Group agroup, bgroup; 7759a3665c6SJed Brown PetscMPIInt i, rank, size, nranks, *ranks; 7769a3665c6SJed Brown 7779a3665c6SJed Brown PetscFunctionBegin; 7780298fd71SBarry Smith *B = NULL; 7799566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)A, &acomm)); 7809566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(acomm, &size)); 7819566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(acomm, &rank)); 7829566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRanges(A, &ranges)); 7839a3665c6SJed Brown for (i = 0, nranks = 0; i < size; i++) { 7849a3665c6SJed Brown if (ranges[i + 1] - ranges[i] > 0) nranks++; 7859a3665c6SJed Brown } 7869a3665c6SJed Brown if (nranks == size) { /* All ranks have a positive number of rows, so we do not need to create a subcomm; */ 7879566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)A)); 7889a3665c6SJed Brown *B = A; 7893ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7909a3665c6SJed Brown } 7919a3665c6SJed Brown 7929566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nranks, &ranks)); 7939a3665c6SJed Brown for (i = 0, nranks = 0; i < size; i++) { 7949a3665c6SJed Brown if (ranges[i + 1] - ranges[i] > 0) ranks[nranks++] = i; 7959a3665c6SJed Brown } 7969566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_group(acomm, &agroup)); 7979566063dSJacob Faibussowitsch PetscCallMPI(MPI_Group_incl(agroup, nranks, ranks, &bgroup)); 7989566063dSJacob Faibussowitsch PetscCall(PetscFree(ranks)); 7999566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_create(acomm, bgroup, &bcomm)); 8009566063dSJacob Faibussowitsch PetscCallMPI(MPI_Group_free(&agroup)); 8019566063dSJacob Faibussowitsch PetscCallMPI(MPI_Group_free(&bgroup)); 8029a3665c6SJed Brown if (bcomm != MPI_COMM_NULL) { 8039a3665c6SJed Brown PetscInt m, N; 8049a3665c6SJed Brown Mat_MPIAdj *b; 8059566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A, &m, NULL)); 8069566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, NULL, &N)); 8079566063dSJacob Faibussowitsch PetscCall(MatCreateMPIAdj(bcomm, m, N, a->i, a->j, a->values, B)); 8089a3665c6SJed Brown b = (Mat_MPIAdj *)(*B)->data; 8099a3665c6SJed Brown b->freeaij = PETSC_FALSE; 8109566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free(&bcomm)); 8119a3665c6SJed Brown } 8123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8139a3665c6SJed Brown } 8149a3665c6SJed Brown 815*66976f2fSJacob Faibussowitsch static PetscErrorCode MatMPIAdjToSeq_MPIAdj(Mat A, Mat *B) 816d71ae5a4SJacob Faibussowitsch { 817b44e7856SBarry Smith PetscInt M, N, *II, *J, NZ, nz, m, nzstart, i; 818b44e7856SBarry Smith PetscInt *Values = NULL; 819b44e7856SBarry Smith Mat_MPIAdj *adj = (Mat_MPIAdj *)A->data; 820b44e7856SBarry Smith PetscMPIInt mnz, mm, *allnz, *allm, size, *dispnz, *dispm; 821b44e7856SBarry Smith 822b44e7856SBarry Smith PetscFunctionBegin; 8239566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size)); 8249566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &M, &N)); 8259566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A, &m, NULL)); 826b44e7856SBarry Smith nz = adj->nz; 82708401ef6SPierre Jolivet PetscCheck(adj->i[m] == nz, PETSC_COMM_SELF, PETSC_ERR_PLIB, "nz %" PetscInt_FMT " not correct i[m] %" PetscInt_FMT, nz, adj->i[m]); 828712fec58SPierre Jolivet PetscCall(MPIU_Allreduce(&nz, &NZ, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)A))); 829d4e69552SBarry Smith 8309566063dSJacob Faibussowitsch PetscCall(PetscMPIIntCast(nz, &mnz)); 8319566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(size, &allnz, size, &dispnz)); 8329566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&mnz, 1, MPI_INT, allnz, 1, MPI_INT, PetscObjectComm((PetscObject)A))); 8339371c9d4SSatish Balay dispnz[0] = 0; 8349371c9d4SSatish Balay for (i = 1; i < size; i++) dispnz[i] = dispnz[i - 1] + allnz[i - 1]; 835b44e7856SBarry Smith if (adj->values) { 8369566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(NZ, &Values)); 8379566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgatherv(adj->values, mnz, MPIU_INT, Values, allnz, dispnz, MPIU_INT, PetscObjectComm((PetscObject)A))); 838b44e7856SBarry Smith } 8399566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(NZ, &J)); 8409566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgatherv(adj->j, mnz, MPIU_INT, J, allnz, dispnz, MPIU_INT, PetscObjectComm((PetscObject)A))); 8419566063dSJacob Faibussowitsch PetscCall(PetscFree2(allnz, dispnz)); 8429566063dSJacob Faibussowitsch PetscCallMPI(MPI_Scan(&nz, &nzstart, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)A))); 843b44e7856SBarry Smith nzstart -= nz; 844b44e7856SBarry Smith /* shift the i[] values so they will be correct after being received */ 845b44e7856SBarry Smith for (i = 0; i < m; i++) adj->i[i] += nzstart; 8469566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(M + 1, &II)); 8479566063dSJacob Faibussowitsch PetscCall(PetscMPIIntCast(m, &mm)); 8489566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(size, &allm, size, &dispm)); 8499566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&mm, 1, MPI_INT, allm, 1, MPI_INT, PetscObjectComm((PetscObject)A))); 8509371c9d4SSatish Balay dispm[0] = 0; 8519371c9d4SSatish Balay for (i = 1; i < size; i++) dispm[i] = dispm[i - 1] + allm[i - 1]; 8529566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgatherv(adj->i, mm, MPIU_INT, II, allm, dispm, MPIU_INT, PetscObjectComm((PetscObject)A))); 8539566063dSJacob Faibussowitsch PetscCall(PetscFree2(allm, dispm)); 854b44e7856SBarry Smith II[M] = NZ; 855b44e7856SBarry Smith /* shift the i[] values back */ 856b44e7856SBarry Smith for (i = 0; i < m; i++) adj->i[i] -= nzstart; 8579566063dSJacob Faibussowitsch PetscCall(MatCreateMPIAdj(PETSC_COMM_SELF, M, N, II, J, Values, B)); 8583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 859b44e7856SBarry Smith } 860b44e7856SBarry Smith 861*66976f2fSJacob Faibussowitsch static PetscErrorCode MatMPIAdjToSeqRankZero_MPIAdj(Mat A, Mat *B) 862d71ae5a4SJacob Faibussowitsch { 863252a1336SBarry Smith PetscInt M, N, *II, *J, NZ, nz, m, nzstart, i; 864252a1336SBarry Smith PetscInt *Values = NULL; 865252a1336SBarry Smith Mat_MPIAdj *adj = (Mat_MPIAdj *)A->data; 866252a1336SBarry Smith PetscMPIInt mnz, mm, *allnz = NULL, *allm, size, *dispnz, *dispm, rank; 867252a1336SBarry Smith 868252a1336SBarry Smith PetscFunctionBegin; 869252a1336SBarry Smith PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size)); 870252a1336SBarry Smith PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)A), &rank)); 871252a1336SBarry Smith PetscCall(MatGetSize(A, &M, &N)); 872252a1336SBarry Smith PetscCall(MatGetLocalSize(A, &m, NULL)); 873252a1336SBarry Smith nz = adj->nz; 874252a1336SBarry Smith PetscCheck(adj->i[m] == nz, PETSC_COMM_SELF, PETSC_ERR_PLIB, "nz %" PetscInt_FMT " not correct i[m] %" PetscInt_FMT, nz, adj->i[m]); 875712fec58SPierre Jolivet PetscCall(MPIU_Allreduce(&nz, &NZ, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)A))); 876252a1336SBarry Smith 877252a1336SBarry Smith PetscCall(PetscMPIIntCast(nz, &mnz)); 878252a1336SBarry Smith if (!rank) PetscCall(PetscMalloc2(size, &allnz, size, &dispnz)); 879252a1336SBarry Smith PetscCallMPI(MPI_Gather(&mnz, 1, MPI_INT, allnz, 1, MPI_INT, 0, PetscObjectComm((PetscObject)A))); 880252a1336SBarry Smith if (!rank) { 8819371c9d4SSatish Balay dispnz[0] = 0; 8829371c9d4SSatish Balay for (i = 1; i < size; i++) dispnz[i] = dispnz[i - 1] + allnz[i - 1]; 883252a1336SBarry Smith if (adj->values) { 884252a1336SBarry Smith PetscCall(PetscMalloc1(NZ, &Values)); 885252a1336SBarry Smith PetscCallMPI(MPI_Gatherv(adj->values, mnz, MPIU_INT, Values, allnz, dispnz, MPIU_INT, 0, PetscObjectComm((PetscObject)A))); 886252a1336SBarry Smith } 887252a1336SBarry Smith PetscCall(PetscMalloc1(NZ, &J)); 888252a1336SBarry Smith PetscCallMPI(MPI_Gatherv(adj->j, mnz, MPIU_INT, J, allnz, dispnz, MPIU_INT, 0, PetscObjectComm((PetscObject)A))); 889252a1336SBarry Smith PetscCall(PetscFree2(allnz, dispnz)); 890252a1336SBarry Smith } else { 891252a1336SBarry Smith if (adj->values) PetscCallMPI(MPI_Gatherv(adj->values, mnz, MPIU_INT, NULL, NULL, NULL, MPIU_INT, 0, PetscObjectComm((PetscObject)A))); 892252a1336SBarry Smith PetscCallMPI(MPI_Gatherv(adj->j, mnz, MPIU_INT, NULL, NULL, NULL, MPIU_INT, 0, PetscObjectComm((PetscObject)A))); 893252a1336SBarry Smith } 894252a1336SBarry Smith PetscCallMPI(MPI_Scan(&nz, &nzstart, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)A))); 895252a1336SBarry Smith nzstart -= nz; 896252a1336SBarry Smith /* shift the i[] values so they will be correct after being received */ 897252a1336SBarry Smith for (i = 0; i < m; i++) adj->i[i] += nzstart; 898252a1336SBarry Smith PetscCall(PetscMPIIntCast(m, &mm)); 899252a1336SBarry Smith if (!rank) { 900252a1336SBarry Smith PetscCall(PetscMalloc1(M + 1, &II)); 901252a1336SBarry Smith PetscCall(PetscMalloc2(size, &allm, size, &dispm)); 902252a1336SBarry Smith PetscCallMPI(MPI_Gather(&mm, 1, MPI_INT, allm, 1, MPI_INT, 0, PetscObjectComm((PetscObject)A))); 9039371c9d4SSatish Balay dispm[0] = 0; 9049371c9d4SSatish Balay for (i = 1; i < size; i++) dispm[i] = dispm[i - 1] + allm[i - 1]; 905252a1336SBarry Smith PetscCallMPI(MPI_Gatherv(adj->i, mm, MPIU_INT, II, allm, dispm, MPIU_INT, 0, PetscObjectComm((PetscObject)A))); 906252a1336SBarry Smith PetscCall(PetscFree2(allm, dispm)); 907252a1336SBarry Smith II[M] = NZ; 908252a1336SBarry Smith } else { 909252a1336SBarry Smith PetscCallMPI(MPI_Gather(&mm, 1, MPI_INT, NULL, 1, MPI_INT, 0, PetscObjectComm((PetscObject)A))); 910252a1336SBarry Smith PetscCallMPI(MPI_Gatherv(adj->i, mm, MPIU_INT, NULL, NULL, NULL, MPIU_INT, 0, PetscObjectComm((PetscObject)A))); 911252a1336SBarry Smith } 912252a1336SBarry Smith /* shift the i[] values back */ 913252a1336SBarry Smith for (i = 0; i < m; i++) adj->i[i] -= nzstart; 914252a1336SBarry Smith if (!rank) PetscCall(MatCreateMPIAdj(PETSC_COMM_SELF, M, N, II, J, Values, B)); 9153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 916252a1336SBarry Smith } 917252a1336SBarry Smith 9189a3665c6SJed Brown /*@ 91911a5261eSBarry Smith MatMPIAdjCreateNonemptySubcommMat - create the same `MATMPIADJ` matrix on a subcommunicator containing only processes owning a positive number of rows 9209a3665c6SJed Brown 9219a3665c6SJed Brown Collective 9229a3665c6SJed Brown 9234165533cSJose E. Roman Input Parameter: 9242ef1f0ffSBarry Smith . A - original `MATMPIADJ` matrix 9259a3665c6SJed Brown 9264165533cSJose E. Roman Output Parameter: 9272ef1f0ffSBarry Smith . B - matrix on subcommunicator, `NULL` on ranks that owned zero rows of `A` 9289a3665c6SJed Brown 9299a3665c6SJed Brown Level: developer 9309a3665c6SJed Brown 9319a3665c6SJed Brown Note: 9329a3665c6SJed Brown This function is mostly useful for internal use by mesh partitioning packages that require that every process owns at least one row. 9339a3665c6SJed Brown 9342ef1f0ffSBarry Smith The matrix `B` should be destroyed with `MatDestroy()`. The arrays are not copied, so `B` should be destroyed before `A` is destroyed. 9359a3665c6SJed Brown 9361cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATMPIADJ`, `MatCreateMPIAdj()` 9379a3665c6SJed Brown @*/ 938d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMPIAdjCreateNonemptySubcommMat(Mat A, Mat *B) 939d71ae5a4SJacob Faibussowitsch { 9409a3665c6SJed Brown PetscFunctionBegin; 9419a3665c6SJed Brown PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 942cac4c232SBarry Smith PetscUseMethod(A, "MatMPIAdjCreateNonemptySubcommMat_C", (Mat, Mat *), (A, B)); 9433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9449a3665c6SJed Brown } 9459a3665c6SJed Brown 9460bad9183SKris Buschelman /*MC 947fafad747SKris Buschelman MATMPIADJ - MATMPIADJ = "mpiadj" - A matrix type to be used for distributed adjacency matrices, 9480bad9183SKris Buschelman intended for use constructing orderings and partitionings. 9490bad9183SKris Buschelman 9500bad9183SKris Buschelman Level: beginner 9510bad9183SKris Buschelman 95211a5261eSBarry Smith Note: 95311a5261eSBarry Smith You can provide values to the matrix using `MatMPIAdjSetPreallocation()`, `MatCreateMPIAdj()`, or 95411a5261eSBarry Smith by calling `MatSetValues()` and `MatAssemblyBegin()` followed by `MatAssemblyEnd()` 9556a09307cSBarry Smith 9561cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatCreateMPIAdj()`, `MatMPIAdjSetPreallocation()`, `MatSetValues()` 9570bad9183SKris Buschelman M*/ 958d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatCreate_MPIAdj(Mat B) 959d71ae5a4SJacob Faibussowitsch { 960273d9f13SBarry Smith Mat_MPIAdj *b; 961273d9f13SBarry Smith 962273d9f13SBarry Smith PetscFunctionBegin; 9634dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&b)); 964b0a32e0cSBarry Smith B->data = (void *)b; 965aea10558SJacob Faibussowitsch B->ops[0] = MatOps_Values; 966273d9f13SBarry Smith B->assembled = PETSC_FALSE; 9676a09307cSBarry Smith B->preallocated = PETSC_TRUE; /* so that MatSetValues() may be used */ 968273d9f13SBarry Smith 9699566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMPIAdjSetPreallocation_C", MatMPIAdjSetPreallocation_MPIAdj)); 9709566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMPIAdjCreateNonemptySubcommMat_C", MatMPIAdjCreateNonemptySubcommMat_MPIAdj)); 9719566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMPIAdjToSeq_C", MatMPIAdjToSeq_MPIAdj)); 972252a1336SBarry Smith PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMPIAdjToSeqRankZero_C", MatMPIAdjToSeqRankZero_MPIAdj)); 9739566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)B, MATMPIADJ)); 9743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 975273d9f13SBarry Smith } 976273d9f13SBarry Smith 977a23d5eceSKris Buschelman /*@C 97811a5261eSBarry Smith MatMPIAdjToSeq - Converts an parallel `MATMPIADJ` matrix to complete `MATMPIADJ` on each process (needed by sequential partitioner) 979b44e7856SBarry Smith 980d083f849SBarry Smith Logically Collective 981b44e7856SBarry Smith 982b44e7856SBarry Smith Input Parameter: 983b44e7856SBarry Smith . A - the matrix 984b44e7856SBarry Smith 985b44e7856SBarry Smith Output Parameter: 986b44e7856SBarry Smith . B - the same matrix on all processes 987b44e7856SBarry Smith 988b44e7856SBarry Smith Level: intermediate 989b44e7856SBarry Smith 9901cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATMPIADJ`, `MatCreate()`, `MatCreateMPIAdj()`, `MatSetValues()`, `MatMPIAdjToSeqRankZero()` 991b44e7856SBarry Smith @*/ 992d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMPIAdjToSeq(Mat A, Mat *B) 993d71ae5a4SJacob Faibussowitsch { 994b44e7856SBarry Smith PetscFunctionBegin; 995cac4c232SBarry Smith PetscUseMethod(A, "MatMPIAdjToSeq_C", (Mat, Mat *), (A, B)); 9963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 997b44e7856SBarry Smith } 998b44e7856SBarry Smith 999b44e7856SBarry Smith /*@C 100011a5261eSBarry Smith MatMPIAdjToSeqRankZero - Converts an parallel `MATMPIADJ` matrix to complete `MATMPIADJ` on rank zero (needed by sequential partitioner) 1001252a1336SBarry Smith 1002252a1336SBarry Smith Logically Collective 1003252a1336SBarry Smith 1004252a1336SBarry Smith Input Parameter: 1005252a1336SBarry Smith . A - the matrix 1006252a1336SBarry Smith 1007252a1336SBarry Smith Output Parameter: 1008252a1336SBarry Smith . B - the same matrix on rank zero, not set on other ranks 1009252a1336SBarry Smith 1010252a1336SBarry Smith Level: intermediate 1011252a1336SBarry Smith 101211a5261eSBarry Smith Note: 1013252a1336SBarry Smith This routine has the advantage on systems with multiple ranks per node since only one copy of the matrix 1014252a1336SBarry Smith is stored on the first node, instead of the number of ranks copies. This can allow partitioning much larger 1015d5b43468SJose E. Roman parallel graph sequentially. 1016252a1336SBarry Smith 10171cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATMPIADJ`, `MatCreate()`, `MatCreateMPIAdj()`, `MatSetValues()`, `MatMPIAdjToSeq()` 1018252a1336SBarry Smith @*/ 1019d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMPIAdjToSeqRankZero(Mat A, Mat *B) 1020d71ae5a4SJacob Faibussowitsch { 1021252a1336SBarry Smith PetscFunctionBegin; 1022252a1336SBarry Smith PetscUseMethod(A, "MatMPIAdjToSeqRankZero_C", (Mat, Mat *), (A, B)); 10233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1024252a1336SBarry Smith } 1025252a1336SBarry Smith 1026252a1336SBarry Smith /*@C 1027a23d5eceSKris Buschelman MatMPIAdjSetPreallocation - Sets the array used for storing the matrix elements 1028a23d5eceSKris Buschelman 1029d083f849SBarry Smith Logically Collective 1030a23d5eceSKris Buschelman 1031a23d5eceSKris Buschelman Input Parameters: 1032fe59aa6dSJacob Faibussowitsch + B - the matrix 10332ef1f0ffSBarry Smith . i - the indices into `j` for the start of each row 1034a23d5eceSKris Buschelman . j - the column indices for each row (sorted for each row). 10352ef1f0ffSBarry Smith The indices in `i` and `j` start with zero (NOT with one). 10362ef1f0ffSBarry Smith - values - [use `NULL` if not provided] edge weights 1037a23d5eceSKris Buschelman 1038a23d5eceSKris Buschelman Level: intermediate 1039a23d5eceSKris Buschelman 1040fe59aa6dSJacob Faibussowitsch .seealso: [](ch_matrices), `Mat`, `MatCreate()`, `MatCreateMPIAdj()`, `MatSetValues()`, `MATMPIADJ` 1041a23d5eceSKris Buschelman @*/ 1042d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMPIAdjSetPreallocation(Mat B, PetscInt *i, PetscInt *j, PetscInt *values) 1043d71ae5a4SJacob Faibussowitsch { 1044273d9f13SBarry Smith PetscFunctionBegin; 1045cac4c232SBarry Smith PetscTryMethod(B, "MatMPIAdjSetPreallocation_C", (Mat, PetscInt *, PetscInt *, PetscInt *), (B, i, j, values)); 10463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1047273d9f13SBarry Smith } 1048273d9f13SBarry Smith 1049b97cf49bSBarry Smith /*@C 10503eda8832SBarry Smith MatCreateMPIAdj - Creates a sparse matrix representing an adjacency list. 1051b97cf49bSBarry Smith The matrix does not have numerical values associated with it, but is 1052b97cf49bSBarry Smith intended for ordering (to reduce bandwidth etc) and partitioning. 1053b97cf49bSBarry Smith 1054d083f849SBarry Smith Collective 1055ef5ee4d1SLois Curfman McInnes 1056b97cf49bSBarry Smith Input Parameters: 1057c2e958feSBarry Smith + comm - MPI communicator 10580752156aSBarry Smith . m - number of local rows 105958ec18a5SLisandro Dalcin . N - number of global columns 10602ef1f0ffSBarry Smith . i - the indices into `j` for the start of each row 1061329f5518SBarry Smith . j - the column indices for each row (sorted for each row). 1062fe59aa6dSJacob Faibussowitsch - values - the values 1063b97cf49bSBarry Smith 1064b97cf49bSBarry Smith Output Parameter: 1065b97cf49bSBarry Smith . A - the matrix 1066b97cf49bSBarry Smith 106715091d37SBarry Smith Level: intermediate 106815091d37SBarry Smith 106995452b02SPatrick Sanan Notes: 1070fe59aa6dSJacob Faibussowitsch The indices in `i` and `j` start with zero (NOT with one). 1071fe59aa6dSJacob Faibussowitsch 10722ef1f0ffSBarry Smith You must NOT free the `i`, `values` and `j` arrays yourself. PETSc will free them 107311a5261eSBarry Smith when the matrix is destroyed; you must allocate them with `PetscMalloc()`. If you 107411a5261eSBarry Smith call from Fortran you need not create the arrays with `PetscMalloc()`. 10756a09307cSBarry Smith 10766a09307cSBarry Smith You should not include the matrix diagonals. 1077b97cf49bSBarry Smith 1078e3f7e1d6Svictorle If you already have a matrix, you can create its adjacency matrix by a call 107911a5261eSBarry Smith to `MatConvert()`, specifying a type of `MATMPIADJ`. 1080ededeb1aSvictorle 108111a5261eSBarry Smith Possible values for `MatSetOption()` - `MAT_STRUCTURALLY_SYMMETRIC` 1082b97cf49bSBarry Smith 10831cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatCreate()`, `MatConvert()`, `MatGetOrdering()`, `MATMPIADJ`, `MatMPIAdjSetPreallocation()` 1084b97cf49bSBarry Smith @*/ 1085d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCreateMPIAdj(MPI_Comm comm, PetscInt m, PetscInt N, PetscInt *i, PetscInt *j, PetscInt *values, Mat *A) 1086d71ae5a4SJacob Faibussowitsch { 1087433994e6SBarry Smith PetscFunctionBegin; 10889566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, A)); 10899566063dSJacob Faibussowitsch PetscCall(MatSetSizes(*A, m, PETSC_DETERMINE, PETSC_DETERMINE, N)); 10909566063dSJacob Faibussowitsch PetscCall(MatSetType(*A, MATMPIADJ)); 10919566063dSJacob Faibussowitsch PetscCall(MatMPIAdjSetPreallocation(*A, i, j, values)); 10923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1093b97cf49bSBarry Smith } 1094