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; 275aa482453SBarry Smith #if defined(PETSC_USE_LOG) 2763ba16761SJacob Faibussowitsch PetscCall(PetscLogObjectState((PetscObject)mat, "Rows=%" PetscInt_FMT ", Cols=%" PetscInt_FMT ", NZ=%" PetscInt_FMT, mat->rmap->n, mat->cmap->n, a->nz)); 277b97cf49bSBarry Smith #endif 2789566063dSJacob Faibussowitsch PetscCall(PetscFree(a->diag)); 2790400557aSBarry Smith if (a->freeaij) { 28014196391SBarry Smith if (a->freeaijwithfree) { 28114196391SBarry Smith if (a->i) free(a->i); 28214196391SBarry Smith if (a->j) free(a->j); 28314196391SBarry Smith } else { 2849566063dSJacob Faibussowitsch PetscCall(PetscFree(a->i)); 2859566063dSJacob Faibussowitsch PetscCall(PetscFree(a->j)); 2869566063dSJacob Faibussowitsch PetscCall(PetscFree(a->values)); 28714196391SBarry Smith } 2880400557aSBarry Smith } 2899566063dSJacob Faibussowitsch PetscCall(PetscFree(a->rowvalues)); 2909566063dSJacob Faibussowitsch PetscCall(PetscFree(mat->data)); 2919566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)mat, NULL)); 2929566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMPIAdjSetPreallocation_C", NULL)); 2939566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMPIAdjCreateNonemptySubcommMat_C", NULL)); 2942e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMPIAdjToSeq_C", NULL)); 295252a1336SBarry Smith PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMPIAdjToSeqRankZero_C", NULL)); 2963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 297b97cf49bSBarry Smith } 298b97cf49bSBarry Smith 299d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSetOption_MPIAdj(Mat A, MatOption op, PetscBool flg) 300d71ae5a4SJacob Faibussowitsch { 3013eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj *)A->data; 302b97cf49bSBarry Smith 303433994e6SBarry Smith PetscFunctionBegin; 30412c028f9SKris Buschelman switch (op) { 30577e54ba9SKris Buschelman case MAT_SYMMETRIC: 30612c028f9SKris Buschelman case MAT_STRUCTURALLY_SYMMETRIC: 3079a4540c5SBarry Smith case MAT_HERMITIAN: 308d71ae5a4SJacob Faibussowitsch case MAT_SPD: 309d71ae5a4SJacob Faibussowitsch a->symmetric = flg; 310d71ae5a4SJacob Faibussowitsch break; 3119a4540c5SBarry Smith case MAT_SYMMETRY_ETERNAL: 312b94d7dedSBarry Smith case MAT_STRUCTURAL_SYMMETRY_ETERNAL: 313d71ae5a4SJacob Faibussowitsch case MAT_SPD_ETERNAL: 314d71ae5a4SJacob Faibussowitsch break; 315d71ae5a4SJacob Faibussowitsch default: 316d71ae5a4SJacob Faibussowitsch PetscCall(PetscInfo(A, "Option %s ignored\n", MatOptions[op])); 317d71ae5a4SJacob Faibussowitsch break; 318b97cf49bSBarry Smith } 3193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 320b97cf49bSBarry Smith } 321b97cf49bSBarry Smith 322d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatGetRow_MPIAdj(Mat A, PetscInt row, PetscInt *nz, PetscInt **idx, PetscScalar **v) 323d71ae5a4SJacob Faibussowitsch { 3243eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj *)A->data; 325b97cf49bSBarry Smith 326433994e6SBarry Smith PetscFunctionBegin; 327d0f46423SBarry Smith row -= A->rmap->rstart; 328aed4548fSBarry Smith PetscCheck(row >= 0 && row < A->rmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Row out of range"); 329b97cf49bSBarry Smith *nz = a->i[row + 1] - a->i[row]; 3302b1d8763SJed Brown if (v) { 3312b1d8763SJed Brown PetscInt j; 3322b1d8763SJed Brown if (a->rowvalues_alloc < *nz) { 3339566063dSJacob Faibussowitsch PetscCall(PetscFree(a->rowvalues)); 3342b1d8763SJed Brown a->rowvalues_alloc = PetscMax(a->rowvalues_alloc * 2, *nz); 3359566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(a->rowvalues_alloc, &a->rowvalues)); 336b97cf49bSBarry Smith } 337ad540459SPierre Jolivet for (j = 0; j < *nz; j++) a->rowvalues[j] = a->values ? a->values[a->i[row] + j] : 1.0; 3382b1d8763SJed Brown *v = (*nz) ? a->rowvalues : NULL; 339b97cf49bSBarry Smith } 34092b6f2a0SJed Brown if (idx) *idx = (*nz) ? a->j + a->i[row] : NULL; 3413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 342b97cf49bSBarry Smith } 343b97cf49bSBarry Smith 344d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatRestoreRow_MPIAdj(Mat A, PetscInt row, PetscInt *nz, PetscInt **idx, PetscScalar **v) 345d71ae5a4SJacob Faibussowitsch { 346433994e6SBarry Smith PetscFunctionBegin; 3473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 348b97cf49bSBarry Smith } 349b97cf49bSBarry Smith 350d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatEqual_MPIAdj(Mat A, Mat B, PetscBool *flg) 351d71ae5a4SJacob Faibussowitsch { 3523eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj *)A->data, *b = (Mat_MPIAdj *)B->data; 353ace3abfcSBarry Smith PetscBool flag; 354b97cf49bSBarry Smith 355433994e6SBarry Smith PetscFunctionBegin; 356b97cf49bSBarry Smith /* If the matrix dimensions are not equal,or no of nonzeros */ 357ad540459SPierre Jolivet if ((A->rmap->n != B->rmap->n) || (a->nz != b->nz)) flag = PETSC_FALSE; 358b97cf49bSBarry Smith 359b97cf49bSBarry Smith /* if the a->i are the same */ 3609566063dSJacob Faibussowitsch PetscCall(PetscArraycmp(a->i, b->i, A->rmap->n + 1, &flag)); 361b97cf49bSBarry Smith 362b97cf49bSBarry Smith /* if a->j are the same */ 3639566063dSJacob Faibussowitsch PetscCall(PetscMemcmp(a->j, b->j, (a->nz) * sizeof(PetscInt), &flag)); 364b97cf49bSBarry Smith 3651c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(&flag, flg, 1, MPIU_BOOL, MPI_LAND, PetscObjectComm((PetscObject)A))); 3663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 367b97cf49bSBarry Smith } 368b97cf49bSBarry Smith 369d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatGetRowIJ_MPIAdj(Mat A, PetscInt oshift, PetscBool symmetric, PetscBool blockcompressed, PetscInt *m, const PetscInt *inia[], const PetscInt *inja[], PetscBool *done) 370d71ae5a4SJacob Faibussowitsch { 371b24ad042SBarry Smith PetscInt i; 3726cd91f64SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj *)A->data; 3731a83f524SJed Brown PetscInt **ia = (PetscInt **)inia, **ja = (PetscInt **)inja; 3746cd91f64SBarry Smith 3756cd91f64SBarry Smith PetscFunctionBegin; 376d0f46423SBarry Smith *m = A->rmap->n; 3776cd91f64SBarry Smith *ia = a->i; 3786cd91f64SBarry Smith *ja = a->j; 3796cd91f64SBarry Smith *done = PETSC_TRUE; 380d42735a1SBarry Smith if (oshift) { 381ad540459SPierre Jolivet for (i = 0; i < (*ia)[*m]; i++) (*ja)[i]++; 382d42735a1SBarry Smith for (i = 0; i <= (*m); i++) (*ia)[i]++; 383d42735a1SBarry Smith } 3843ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 385d42735a1SBarry Smith } 386d42735a1SBarry Smith 387d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatRestoreRowIJ_MPIAdj(Mat A, PetscInt oshift, PetscBool symmetric, PetscBool blockcompressed, PetscInt *m, const PetscInt *inia[], const PetscInt *inja[], PetscBool *done) 388d71ae5a4SJacob Faibussowitsch { 389b24ad042SBarry Smith PetscInt i; 390d42735a1SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj *)A->data; 3911a83f524SJed Brown PetscInt **ia = (PetscInt **)inia, **ja = (PetscInt **)inja; 392d42735a1SBarry Smith 393d42735a1SBarry Smith PetscFunctionBegin; 394aed4548fSBarry Smith PetscCheck(!ia || a->i == *ia, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "ia passed back is not one obtained with MatGetRowIJ()"); 395aed4548fSBarry Smith PetscCheck(!ja || a->j == *ja, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "ja passed back is not one obtained with MatGetRowIJ()"); 396d42735a1SBarry Smith if (oshift) { 39728b400f6SJacob Faibussowitsch PetscCheck(ia, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "If oshift then you must passed in inia[] argument"); 39828b400f6SJacob Faibussowitsch PetscCheck(ja, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "If oshift then you must passed in inja[] argument"); 399d42735a1SBarry Smith for (i = 0; i <= (*m); i++) (*ia)[i]--; 400ad540459SPierre Jolivet for (i = 0; i < (*ia)[*m]; i++) (*ja)[i]--; 401d42735a1SBarry Smith } 4023ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4036cd91f64SBarry Smith } 404b97cf49bSBarry Smith 405d71ae5a4SJacob Faibussowitsch PetscErrorCode MatConvertFrom_MPIAdj(Mat A, MatType type, MatReuse reuse, Mat *newmat) 406d71ae5a4SJacob Faibussowitsch { 40717667f90SBarry Smith Mat B; 40817667f90SBarry Smith PetscInt i, m, N, nzeros = 0, *ia, *ja, len, rstart, cnt, j, *a; 40917667f90SBarry Smith const PetscInt *rj; 41017667f90SBarry Smith const PetscScalar *ra; 41117667f90SBarry Smith MPI_Comm comm; 41217667f90SBarry Smith 41317667f90SBarry Smith PetscFunctionBegin; 4149566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, NULL, &N)); 4159566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A, &m, NULL)); 4169566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A, &rstart, NULL)); 41717667f90SBarry Smith 41817667f90SBarry Smith /* count the number of nonzeros per row */ 41917667f90SBarry Smith for (i = 0; i < m; i++) { 4209566063dSJacob Faibussowitsch PetscCall(MatGetRow(A, i + rstart, &len, &rj, NULL)); 4215ee9ba1cSJed Brown for (j = 0; j < len; j++) { 4229371c9d4SSatish Balay if (rj[j] == i + rstart) { 4239371c9d4SSatish Balay len--; 4249371c9d4SSatish Balay break; 4259371c9d4SSatish Balay } /* don't count diagonal */ 4265ee9ba1cSJed Brown } 42717667f90SBarry Smith nzeros += len; 4289566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(A, i + rstart, &len, &rj, NULL)); 42917667f90SBarry Smith } 43017667f90SBarry Smith 43117667f90SBarry Smith /* malloc space for nonzeros */ 4329566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nzeros + 1, &a)); 4339566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(N + 1, &ia)); 4349566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nzeros + 1, &ja)); 43517667f90SBarry Smith 43617667f90SBarry Smith nzeros = 0; 43717667f90SBarry Smith ia[0] = 0; 43817667f90SBarry Smith for (i = 0; i < m; i++) { 4399566063dSJacob Faibussowitsch PetscCall(MatGetRow(A, i + rstart, &len, &rj, &ra)); 44017667f90SBarry Smith cnt = 0; 44117667f90SBarry Smith for (j = 0; j < len; j++) { 4425ee9ba1cSJed Brown if (rj[j] != i + rstart) { /* if not diagonal */ 44317667f90SBarry Smith a[nzeros + cnt] = (PetscInt)PetscAbsScalar(ra[j]); 44417667f90SBarry Smith ja[nzeros + cnt++] = rj[j]; 44517667f90SBarry Smith } 4465ee9ba1cSJed Brown } 4479566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(A, i + rstart, &len, &rj, &ra)); 44817667f90SBarry Smith nzeros += cnt; 44917667f90SBarry Smith ia[i + 1] = nzeros; 45017667f90SBarry Smith } 45117667f90SBarry Smith 4529566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)A, &comm)); 4539566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, &B)); 4549566063dSJacob Faibussowitsch PetscCall(MatSetSizes(B, m, PETSC_DETERMINE, PETSC_DETERMINE, N)); 4559566063dSJacob Faibussowitsch PetscCall(MatSetType(B, type)); 4569566063dSJacob Faibussowitsch PetscCall(MatMPIAdjSetPreallocation(B, ia, ja, a)); 45717667f90SBarry Smith 458511c6705SHong Zhang if (reuse == MAT_INPLACE_MATRIX) { 4599566063dSJacob Faibussowitsch PetscCall(MatHeaderReplace(A, &B)); 46017667f90SBarry Smith } else { 46117667f90SBarry Smith *newmat = B; 46217667f90SBarry Smith } 4633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 46417667f90SBarry Smith } 46517667f90SBarry Smith 466d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSetValues_MPIAdj(Mat A, PetscInt m, const PetscInt *rows, PetscInt n, const PetscInt *cols, const PetscScalar *values, InsertMode im) 467d71ae5a4SJacob Faibussowitsch { 4686a09307cSBarry Smith Mat_MPIAdj *adj = (Mat_MPIAdj *)A->data; 4696a09307cSBarry Smith PetscInt rStart, rEnd, cStart, cEnd; 4706a09307cSBarry Smith 4716a09307cSBarry Smith PetscFunctionBegin; 4726a09307cSBarry Smith PetscCheck(!adj->i, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Matrix is already assembled, cannot change its values"); 4736a09307cSBarry Smith PetscCall(MatGetOwnershipRange(A, &rStart, &rEnd)); 4746a09307cSBarry Smith PetscCall(MatGetOwnershipRangeColumn(A, &cStart, &cEnd)); 4756a09307cSBarry Smith if (!adj->ht) { 4766a09307cSBarry Smith PetscCall(PetscHSetIJCreate(&adj->ht)); 4776a09307cSBarry Smith PetscCall(MatStashCreate_Private(PetscObjectComm((PetscObject)A), 1, &A->stash)); 4786a09307cSBarry Smith PetscCall(PetscLayoutSetUp(A->rmap)); 4796a09307cSBarry Smith PetscCall(PetscLayoutSetUp(A->cmap)); 4806a09307cSBarry Smith } 4816a09307cSBarry Smith for (PetscInt r = 0; r < m; ++r) { 4826a09307cSBarry Smith PetscHashIJKey key; 4836a09307cSBarry Smith 4846a09307cSBarry Smith key.i = rows[r]; 4856a09307cSBarry Smith if (key.i < 0) continue; 4866a09307cSBarry Smith if ((key.i < rStart) || (key.i >= rEnd)) { 4876a09307cSBarry Smith PetscCall(MatStashValuesRow_Private(&A->stash, key.i, n, cols, values, PETSC_FALSE)); 4886a09307cSBarry Smith } else { 4896a09307cSBarry Smith for (PetscInt c = 0; c < n; ++c) { 4906a09307cSBarry Smith key.j = cols[c]; 4916a09307cSBarry Smith if (key.j < 0 || key.i == key.j) continue; 4926a09307cSBarry Smith PetscCall(PetscHSetIJAdd(adj->ht, key)); 4936a09307cSBarry Smith } 4946a09307cSBarry Smith } 4956a09307cSBarry Smith } 4963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4976a09307cSBarry Smith } 4986a09307cSBarry Smith 499d71ae5a4SJacob Faibussowitsch PetscErrorCode MatAssemblyBegin_MPIAdj(Mat A, MatAssemblyType type) 500d71ae5a4SJacob Faibussowitsch { 5016a09307cSBarry Smith PetscInt nstash, reallocs; 5026a09307cSBarry Smith Mat_MPIAdj *adj = (Mat_MPIAdj *)A->data; 5036a09307cSBarry Smith 5046a09307cSBarry Smith PetscFunctionBegin; 5056a09307cSBarry Smith if (!adj->ht) { 5066a09307cSBarry Smith PetscCall(PetscHSetIJCreate(&adj->ht)); 5076a09307cSBarry Smith PetscCall(PetscLayoutSetUp(A->rmap)); 5086a09307cSBarry Smith PetscCall(PetscLayoutSetUp(A->cmap)); 5096a09307cSBarry Smith } 5106a09307cSBarry Smith PetscCall(MatStashScatterBegin_Private(A, &A->stash, A->rmap->range)); 5116a09307cSBarry Smith PetscCall(MatStashGetInfo_Private(&A->stash, &nstash, &reallocs)); 5126a09307cSBarry Smith PetscCall(PetscInfo(A, "Stash has %" PetscInt_FMT " entries, uses %" PetscInt_FMT " mallocs.\n", nstash, reallocs)); 5133ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5146a09307cSBarry Smith } 5156a09307cSBarry Smith 516d71ae5a4SJacob Faibussowitsch PetscErrorCode MatAssemblyEnd_MPIAdj(Mat A, MatAssemblyType type) 517d71ae5a4SJacob Faibussowitsch { 5186a09307cSBarry Smith PetscScalar *val; 5196a09307cSBarry Smith PetscInt *row, *col, m, rstart, *rowstarts; 5206a09307cSBarry Smith PetscInt i, j, ncols, flg, nz; 5216a09307cSBarry Smith PetscMPIInt n; 5226a09307cSBarry Smith Mat_MPIAdj *adj = (Mat_MPIAdj *)A->data; 5236a09307cSBarry Smith PetscHashIter hi; 5246a09307cSBarry Smith PetscHashIJKey key; 5256a09307cSBarry Smith PetscHSetIJ ht = adj->ht; 5266a09307cSBarry Smith 5276a09307cSBarry Smith PetscFunctionBegin; 5286a09307cSBarry Smith while (1) { 5296a09307cSBarry Smith PetscCall(MatStashScatterGetMesg_Private(&A->stash, &n, &row, &col, &val, &flg)); 5306a09307cSBarry Smith if (!flg) break; 5316a09307cSBarry Smith 5326a09307cSBarry Smith for (i = 0; i < n;) { 5336a09307cSBarry Smith /* Identify the consecutive vals belonging to the same row */ 5346a09307cSBarry Smith for (j = i, rstart = row[j]; j < n; j++) { 5356a09307cSBarry Smith if (row[j] != rstart) break; 5366a09307cSBarry Smith } 5376a09307cSBarry Smith if (j < n) ncols = j - i; 5386a09307cSBarry Smith else ncols = n - i; 5396a09307cSBarry Smith /* Set all these values with a single function call */ 5406a09307cSBarry Smith PetscCall(MatSetValues_MPIAdj(A, 1, row + i, ncols, col + i, val + i, INSERT_VALUES)); 5416a09307cSBarry Smith i = j; 5426a09307cSBarry Smith } 5436a09307cSBarry Smith } 5446a09307cSBarry Smith PetscCall(MatStashScatterEnd_Private(&A->stash)); 5456a09307cSBarry Smith PetscCall(MatStashDestroy_Private(&A->stash)); 5466a09307cSBarry Smith 5476a09307cSBarry Smith PetscCall(MatGetLocalSize(A, &m, NULL)); 5486a09307cSBarry Smith PetscCall(MatGetOwnershipRange(A, &rstart, NULL)); 5496a09307cSBarry Smith PetscCall(PetscCalloc1(m + 1, &rowstarts)); 5506a09307cSBarry Smith PetscHashIterBegin(ht, hi); 5516a09307cSBarry Smith for (; !PetscHashIterAtEnd(ht, hi);) { 5526a09307cSBarry Smith PetscHashIterGetKey(ht, hi, key); 5536a09307cSBarry Smith rowstarts[key.i - rstart + 1]++; 5546a09307cSBarry Smith PetscHashIterNext(ht, hi); 5556a09307cSBarry Smith } 556ad540459SPierre Jolivet for (i = 1; i < m + 1; i++) rowstarts[i] = rowstarts[i - 1] + rowstarts[i]; 5576a09307cSBarry Smith 5586a09307cSBarry Smith PetscCall(PetscHSetIJGetSize(ht, &nz)); 5596a09307cSBarry Smith PetscCall(PetscMalloc1(nz, &col)); 5606a09307cSBarry Smith PetscHashIterBegin(ht, hi); 5616a09307cSBarry Smith for (; !PetscHashIterAtEnd(ht, hi);) { 5626a09307cSBarry Smith PetscHashIterGetKey(ht, hi, key); 5636a09307cSBarry Smith col[rowstarts[key.i - rstart]++] = key.j; 5646a09307cSBarry Smith PetscHashIterNext(ht, hi); 5656a09307cSBarry Smith } 5666a09307cSBarry Smith PetscCall(PetscHSetIJDestroy(&ht)); 567ad540459SPierre Jolivet for (i = m; i > 0; i--) rowstarts[i] = rowstarts[i - 1]; 5686a09307cSBarry Smith rowstarts[0] = 0; 5696a09307cSBarry Smith 57048a46eb9SPierre Jolivet for (PetscInt i = 0; i < m; i++) PetscCall(PetscSortInt(rowstarts[i + 1] - rowstarts[i], &col[rowstarts[i]])); 5716a09307cSBarry Smith 5726a09307cSBarry Smith adj->i = rowstarts; 5736a09307cSBarry Smith adj->j = col; 574252a1336SBarry Smith adj->nz = rowstarts[m]; 5756a09307cSBarry Smith adj->freeaij = PETSC_TRUE; 5763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5776a09307cSBarry Smith } 5786a09307cSBarry Smith 5796a09307cSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_MPIAdj, 5803eda8832SBarry Smith MatGetRow_MPIAdj, 5813eda8832SBarry Smith MatRestoreRow_MPIAdj, 582f4259b30SLisandro Dalcin NULL, 583f4259b30SLisandro Dalcin /* 4*/ NULL, 584f4259b30SLisandro Dalcin NULL, 585f4259b30SLisandro Dalcin NULL, 586f4259b30SLisandro Dalcin NULL, 587f4259b30SLisandro Dalcin NULL, 588f4259b30SLisandro Dalcin NULL, 589f4259b30SLisandro Dalcin /*10*/ NULL, 590f4259b30SLisandro Dalcin NULL, 591f4259b30SLisandro Dalcin NULL, 592f4259b30SLisandro Dalcin NULL, 593f4259b30SLisandro Dalcin NULL, 594f4259b30SLisandro Dalcin /*15*/ NULL, 5953eda8832SBarry Smith MatEqual_MPIAdj, 596f4259b30SLisandro Dalcin NULL, 597f4259b30SLisandro Dalcin NULL, 598f4259b30SLisandro Dalcin NULL, 5996a09307cSBarry Smith /*20*/ MatAssemblyBegin_MPIAdj, 6006a09307cSBarry Smith MatAssemblyEnd_MPIAdj, 6013eda8832SBarry Smith MatSetOption_MPIAdj, 602f4259b30SLisandro Dalcin NULL, 603f4259b30SLisandro Dalcin /*24*/ NULL, 604f4259b30SLisandro Dalcin NULL, 605f4259b30SLisandro Dalcin NULL, 606f4259b30SLisandro Dalcin NULL, 607f4259b30SLisandro Dalcin NULL, 608f4259b30SLisandro Dalcin /*29*/ NULL, 609f4259b30SLisandro Dalcin NULL, 610f4259b30SLisandro Dalcin NULL, 611f4259b30SLisandro Dalcin NULL, 612f4259b30SLisandro Dalcin NULL, 613f4259b30SLisandro Dalcin /*34*/ NULL, 614f4259b30SLisandro Dalcin NULL, 615f4259b30SLisandro Dalcin NULL, 616f4259b30SLisandro Dalcin NULL, 617f4259b30SLisandro Dalcin NULL, 618f4259b30SLisandro Dalcin /*39*/ NULL, 6197dae84e0SHong Zhang MatCreateSubMatrices_MPIAdj, 620f4259b30SLisandro Dalcin NULL, 621f4259b30SLisandro Dalcin NULL, 622f4259b30SLisandro Dalcin NULL, 623f4259b30SLisandro Dalcin /*44*/ NULL, 624f4259b30SLisandro Dalcin NULL, 6257d68702bSBarry Smith MatShift_Basic, 626f4259b30SLisandro Dalcin NULL, 627f4259b30SLisandro Dalcin NULL, 628f4259b30SLisandro Dalcin /*49*/ NULL, 6296cd91f64SBarry Smith MatGetRowIJ_MPIAdj, 630d42735a1SBarry Smith MatRestoreRowIJ_MPIAdj, 631f4259b30SLisandro Dalcin NULL, 632f4259b30SLisandro Dalcin NULL, 633f4259b30SLisandro Dalcin /*54*/ NULL, 634f4259b30SLisandro Dalcin NULL, 635f4259b30SLisandro Dalcin NULL, 636f4259b30SLisandro Dalcin NULL, 637f4259b30SLisandro Dalcin NULL, 638f4259b30SLisandro Dalcin /*59*/ NULL, 639b9b97703SBarry Smith MatDestroy_MPIAdj, 640b9b97703SBarry Smith MatView_MPIAdj, 64117667f90SBarry Smith MatConvertFrom_MPIAdj, 642f4259b30SLisandro Dalcin NULL, 643f4259b30SLisandro Dalcin /*64*/ NULL, 644f4259b30SLisandro Dalcin NULL, 645f4259b30SLisandro Dalcin NULL, 646f4259b30SLisandro Dalcin NULL, 647f4259b30SLisandro Dalcin NULL, 648f4259b30SLisandro Dalcin /*69*/ NULL, 649f4259b30SLisandro Dalcin NULL, 650f4259b30SLisandro Dalcin NULL, 651f4259b30SLisandro Dalcin NULL, 652f4259b30SLisandro Dalcin NULL, 653f4259b30SLisandro Dalcin /*74*/ NULL, 654f4259b30SLisandro Dalcin NULL, 655f4259b30SLisandro Dalcin NULL, 656f4259b30SLisandro Dalcin NULL, 657f4259b30SLisandro Dalcin NULL, 658f4259b30SLisandro Dalcin /*79*/ NULL, 659f4259b30SLisandro Dalcin NULL, 660f4259b30SLisandro Dalcin NULL, 661f4259b30SLisandro Dalcin NULL, 662f4259b30SLisandro Dalcin NULL, 663f4259b30SLisandro Dalcin /*84*/ NULL, 664f4259b30SLisandro Dalcin NULL, 665f4259b30SLisandro Dalcin NULL, 666f4259b30SLisandro Dalcin NULL, 667f4259b30SLisandro Dalcin NULL, 668f4259b30SLisandro Dalcin /*89*/ NULL, 669f4259b30SLisandro Dalcin NULL, 670f4259b30SLisandro Dalcin NULL, 671f4259b30SLisandro Dalcin NULL, 672f4259b30SLisandro Dalcin NULL, 673f4259b30SLisandro Dalcin /*94*/ NULL, 674f4259b30SLisandro Dalcin NULL, 675f4259b30SLisandro Dalcin NULL, 676f4259b30SLisandro Dalcin NULL, 677f4259b30SLisandro Dalcin NULL, 678f4259b30SLisandro Dalcin /*99*/ NULL, 679f4259b30SLisandro Dalcin NULL, 680f4259b30SLisandro Dalcin NULL, 681f4259b30SLisandro Dalcin NULL, 682f4259b30SLisandro Dalcin NULL, 683f4259b30SLisandro Dalcin /*104*/ NULL, 684f4259b30SLisandro Dalcin NULL, 685f4259b30SLisandro Dalcin NULL, 686f4259b30SLisandro Dalcin NULL, 687f4259b30SLisandro Dalcin NULL, 688f4259b30SLisandro Dalcin /*109*/ NULL, 689f4259b30SLisandro Dalcin NULL, 690f4259b30SLisandro Dalcin NULL, 691f4259b30SLisandro Dalcin NULL, 692f4259b30SLisandro Dalcin NULL, 693f4259b30SLisandro Dalcin /*114*/ NULL, 694f4259b30SLisandro Dalcin NULL, 695f4259b30SLisandro Dalcin NULL, 696f4259b30SLisandro Dalcin NULL, 697f4259b30SLisandro Dalcin NULL, 698f4259b30SLisandro Dalcin /*119*/ NULL, 699f4259b30SLisandro Dalcin NULL, 700f4259b30SLisandro Dalcin NULL, 701f4259b30SLisandro Dalcin NULL, 702f4259b30SLisandro Dalcin NULL, 703f4259b30SLisandro Dalcin /*124*/ NULL, 704f4259b30SLisandro Dalcin NULL, 705f4259b30SLisandro Dalcin NULL, 706f4259b30SLisandro Dalcin NULL, 7077dae84e0SHong Zhang MatCreateSubMatricesMPI_MPIAdj, 708f4259b30SLisandro Dalcin /*129*/ NULL, 709f4259b30SLisandro Dalcin NULL, 710f4259b30SLisandro Dalcin NULL, 711f4259b30SLisandro Dalcin NULL, 712f4259b30SLisandro Dalcin NULL, 713f4259b30SLisandro Dalcin /*134*/ NULL, 714f4259b30SLisandro Dalcin NULL, 715f4259b30SLisandro Dalcin NULL, 716f4259b30SLisandro Dalcin NULL, 717f4259b30SLisandro Dalcin NULL, 718f4259b30SLisandro Dalcin /*139*/ NULL, 719f4259b30SLisandro Dalcin NULL, 720d70f29a3SPierre Jolivet NULL, 721d70f29a3SPierre Jolivet NULL, 722d70f29a3SPierre Jolivet NULL, 723d70f29a3SPierre Jolivet /*144*/ NULL, 724d70f29a3SPierre Jolivet NULL, 725d70f29a3SPierre Jolivet NULL, 72699a7f59eSMark Adams NULL, 72799a7f59eSMark Adams NULL, 7287fb60732SBarry Smith NULL, 729dec0b466SHong Zhang /*150*/ NULL, 730dec0b466SHong Zhang NULL}; 731b97cf49bSBarry Smith 732d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMPIAdjSetPreallocation_MPIAdj(Mat B, PetscInt *i, PetscInt *j, PetscInt *values) 733d71ae5a4SJacob Faibussowitsch { 734a23d5eceSKris Buschelman Mat_MPIAdj *b = (Mat_MPIAdj *)B->data; 735e895ccc0SFande Kong PetscBool useedgeweights; 736a23d5eceSKris Buschelman 737a23d5eceSKris Buschelman PetscFunctionBegin; 7389566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(B->rmap)); 7399566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(B->cmap)); 7409371c9d4SSatish Balay if (values) useedgeweights = PETSC_TRUE; 7419371c9d4SSatish Balay else useedgeweights = PETSC_FALSE; 742e895ccc0SFande Kong /* Make everybody knows if they are using edge weights or not */ 7431c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce((int *)&useedgeweights, (int *)&b->useedgeweights, 1, MPI_INT, MPI_MAX, PetscObjectComm((PetscObject)B))); 74458ec18a5SLisandro Dalcin 74576bd3646SJed Brown if (PetscDefined(USE_DEBUG)) { 74676bd3646SJed Brown PetscInt ii; 74776bd3646SJed Brown 74808401ef6SPierre 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]); 749d0f46423SBarry Smith for (ii = 1; ii < B->rmap->n; ii++) { 750aed4548fSBarry 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]); 751a23d5eceSKris Buschelman } 752ad540459SPierre 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]); 75376bd3646SJed Brown } 754a23d5eceSKris Buschelman b->j = j; 755a23d5eceSKris Buschelman b->i = i; 756a23d5eceSKris Buschelman b->values = values; 757a23d5eceSKris Buschelman 758d0f46423SBarry Smith b->nz = i[B->rmap->n]; 759f4259b30SLisandro Dalcin b->diag = NULL; 760a23d5eceSKris Buschelman b->symmetric = PETSC_FALSE; 761a23d5eceSKris Buschelman b->freeaij = PETSC_TRUE; 762a23d5eceSKris Buschelman 7636a09307cSBarry Smith B->ops->assemblybegin = NULL; 7646a09307cSBarry Smith B->ops->assemblyend = NULL; 7659566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 7669566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 7676a09307cSBarry Smith PetscCall(MatStashDestroy_Private(&B->stash)); 7683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 769a23d5eceSKris Buschelman } 770b97cf49bSBarry Smith 771d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMPIAdjCreateNonemptySubcommMat_MPIAdj(Mat A, Mat *B) 772d71ae5a4SJacob Faibussowitsch { 7739a3665c6SJed Brown Mat_MPIAdj *a = (Mat_MPIAdj *)A->data; 7749a3665c6SJed Brown const PetscInt *ranges; 7759a3665c6SJed Brown MPI_Comm acomm, bcomm; 7769a3665c6SJed Brown MPI_Group agroup, bgroup; 7779a3665c6SJed Brown PetscMPIInt i, rank, size, nranks, *ranks; 7789a3665c6SJed Brown 7799a3665c6SJed Brown PetscFunctionBegin; 7800298fd71SBarry Smith *B = NULL; 7819566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)A, &acomm)); 7829566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(acomm, &size)); 7839566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(acomm, &rank)); 7849566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRanges(A, &ranges)); 7859a3665c6SJed Brown for (i = 0, nranks = 0; i < size; i++) { 7869a3665c6SJed Brown if (ranges[i + 1] - ranges[i] > 0) nranks++; 7879a3665c6SJed Brown } 7889a3665c6SJed Brown if (nranks == size) { /* All ranks have a positive number of rows, so we do not need to create a subcomm; */ 7899566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)A)); 7909a3665c6SJed Brown *B = A; 7913ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7929a3665c6SJed Brown } 7939a3665c6SJed Brown 7949566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nranks, &ranks)); 7959a3665c6SJed Brown for (i = 0, nranks = 0; i < size; i++) { 7969a3665c6SJed Brown if (ranges[i + 1] - ranges[i] > 0) ranks[nranks++] = i; 7979a3665c6SJed Brown } 7989566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_group(acomm, &agroup)); 7999566063dSJacob Faibussowitsch PetscCallMPI(MPI_Group_incl(agroup, nranks, ranks, &bgroup)); 8009566063dSJacob Faibussowitsch PetscCall(PetscFree(ranks)); 8019566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_create(acomm, bgroup, &bcomm)); 8029566063dSJacob Faibussowitsch PetscCallMPI(MPI_Group_free(&agroup)); 8039566063dSJacob Faibussowitsch PetscCallMPI(MPI_Group_free(&bgroup)); 8049a3665c6SJed Brown if (bcomm != MPI_COMM_NULL) { 8059a3665c6SJed Brown PetscInt m, N; 8069a3665c6SJed Brown Mat_MPIAdj *b; 8079566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A, &m, NULL)); 8089566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, NULL, &N)); 8099566063dSJacob Faibussowitsch PetscCall(MatCreateMPIAdj(bcomm, m, N, a->i, a->j, a->values, B)); 8109a3665c6SJed Brown b = (Mat_MPIAdj *)(*B)->data; 8119a3665c6SJed Brown b->freeaij = PETSC_FALSE; 8129566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free(&bcomm)); 8139a3665c6SJed Brown } 8143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8159a3665c6SJed Brown } 8169a3665c6SJed Brown 817d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMPIAdjToSeq_MPIAdj(Mat A, Mat *B) 818d71ae5a4SJacob Faibussowitsch { 819b44e7856SBarry Smith PetscInt M, N, *II, *J, NZ, nz, m, nzstart, i; 820b44e7856SBarry Smith PetscInt *Values = NULL; 821b44e7856SBarry Smith Mat_MPIAdj *adj = (Mat_MPIAdj *)A->data; 822b44e7856SBarry Smith PetscMPIInt mnz, mm, *allnz, *allm, size, *dispnz, *dispm; 823b44e7856SBarry Smith 824b44e7856SBarry Smith PetscFunctionBegin; 8259566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size)); 8269566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &M, &N)); 8279566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A, &m, NULL)); 828b44e7856SBarry Smith nz = adj->nz; 82908401ef6SPierre 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]); 830712fec58SPierre Jolivet PetscCall(MPIU_Allreduce(&nz, &NZ, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)A))); 831d4e69552SBarry Smith 8329566063dSJacob Faibussowitsch PetscCall(PetscMPIIntCast(nz, &mnz)); 8339566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(size, &allnz, size, &dispnz)); 8349566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&mnz, 1, MPI_INT, allnz, 1, MPI_INT, PetscObjectComm((PetscObject)A))); 8359371c9d4SSatish Balay dispnz[0] = 0; 8369371c9d4SSatish Balay for (i = 1; i < size; i++) dispnz[i] = dispnz[i - 1] + allnz[i - 1]; 837b44e7856SBarry Smith if (adj->values) { 8389566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(NZ, &Values)); 8399566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgatherv(adj->values, mnz, MPIU_INT, Values, allnz, dispnz, MPIU_INT, PetscObjectComm((PetscObject)A))); 840b44e7856SBarry Smith } 8419566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(NZ, &J)); 8429566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgatherv(adj->j, mnz, MPIU_INT, J, allnz, dispnz, MPIU_INT, PetscObjectComm((PetscObject)A))); 8439566063dSJacob Faibussowitsch PetscCall(PetscFree2(allnz, dispnz)); 8449566063dSJacob Faibussowitsch PetscCallMPI(MPI_Scan(&nz, &nzstart, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)A))); 845b44e7856SBarry Smith nzstart -= nz; 846b44e7856SBarry Smith /* shift the i[] values so they will be correct after being received */ 847b44e7856SBarry Smith for (i = 0; i < m; i++) adj->i[i] += nzstart; 8489566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(M + 1, &II)); 8499566063dSJacob Faibussowitsch PetscCall(PetscMPIIntCast(m, &mm)); 8509566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(size, &allm, size, &dispm)); 8519566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&mm, 1, MPI_INT, allm, 1, MPI_INT, PetscObjectComm((PetscObject)A))); 8529371c9d4SSatish Balay dispm[0] = 0; 8539371c9d4SSatish Balay for (i = 1; i < size; i++) dispm[i] = dispm[i - 1] + allm[i - 1]; 8549566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgatherv(adj->i, mm, MPIU_INT, II, allm, dispm, MPIU_INT, PetscObjectComm((PetscObject)A))); 8559566063dSJacob Faibussowitsch PetscCall(PetscFree2(allm, dispm)); 856b44e7856SBarry Smith II[M] = NZ; 857b44e7856SBarry Smith /* shift the i[] values back */ 858b44e7856SBarry Smith for (i = 0; i < m; i++) adj->i[i] -= nzstart; 8599566063dSJacob Faibussowitsch PetscCall(MatCreateMPIAdj(PETSC_COMM_SELF, M, N, II, J, Values, B)); 8603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 861b44e7856SBarry Smith } 862b44e7856SBarry Smith 863d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMPIAdjToSeqRankZero_MPIAdj(Mat A, Mat *B) 864d71ae5a4SJacob Faibussowitsch { 865252a1336SBarry Smith PetscInt M, N, *II, *J, NZ, nz, m, nzstart, i; 866252a1336SBarry Smith PetscInt *Values = NULL; 867252a1336SBarry Smith Mat_MPIAdj *adj = (Mat_MPIAdj *)A->data; 868252a1336SBarry Smith PetscMPIInt mnz, mm, *allnz = NULL, *allm, size, *dispnz, *dispm, rank; 869252a1336SBarry Smith 870252a1336SBarry Smith PetscFunctionBegin; 871252a1336SBarry Smith PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size)); 872252a1336SBarry Smith PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)A), &rank)); 873252a1336SBarry Smith PetscCall(MatGetSize(A, &M, &N)); 874252a1336SBarry Smith PetscCall(MatGetLocalSize(A, &m, NULL)); 875252a1336SBarry Smith nz = adj->nz; 876252a1336SBarry 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]); 877712fec58SPierre Jolivet PetscCall(MPIU_Allreduce(&nz, &NZ, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)A))); 878252a1336SBarry Smith 879252a1336SBarry Smith PetscCall(PetscMPIIntCast(nz, &mnz)); 880252a1336SBarry Smith if (!rank) PetscCall(PetscMalloc2(size, &allnz, size, &dispnz)); 881252a1336SBarry Smith PetscCallMPI(MPI_Gather(&mnz, 1, MPI_INT, allnz, 1, MPI_INT, 0, PetscObjectComm((PetscObject)A))); 882252a1336SBarry Smith if (!rank) { 8839371c9d4SSatish Balay dispnz[0] = 0; 8849371c9d4SSatish Balay for (i = 1; i < size; i++) dispnz[i] = dispnz[i - 1] + allnz[i - 1]; 885252a1336SBarry Smith if (adj->values) { 886252a1336SBarry Smith PetscCall(PetscMalloc1(NZ, &Values)); 887252a1336SBarry Smith PetscCallMPI(MPI_Gatherv(adj->values, mnz, MPIU_INT, Values, allnz, dispnz, MPIU_INT, 0, PetscObjectComm((PetscObject)A))); 888252a1336SBarry Smith } 889252a1336SBarry Smith PetscCall(PetscMalloc1(NZ, &J)); 890252a1336SBarry Smith PetscCallMPI(MPI_Gatherv(adj->j, mnz, MPIU_INT, J, allnz, dispnz, MPIU_INT, 0, PetscObjectComm((PetscObject)A))); 891252a1336SBarry Smith PetscCall(PetscFree2(allnz, dispnz)); 892252a1336SBarry Smith } else { 893252a1336SBarry Smith if (adj->values) PetscCallMPI(MPI_Gatherv(adj->values, mnz, MPIU_INT, NULL, NULL, NULL, MPIU_INT, 0, PetscObjectComm((PetscObject)A))); 894252a1336SBarry Smith PetscCallMPI(MPI_Gatherv(adj->j, mnz, MPIU_INT, NULL, NULL, NULL, MPIU_INT, 0, PetscObjectComm((PetscObject)A))); 895252a1336SBarry Smith } 896252a1336SBarry Smith PetscCallMPI(MPI_Scan(&nz, &nzstart, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)A))); 897252a1336SBarry Smith nzstart -= nz; 898252a1336SBarry Smith /* shift the i[] values so they will be correct after being received */ 899252a1336SBarry Smith for (i = 0; i < m; i++) adj->i[i] += nzstart; 900252a1336SBarry Smith PetscCall(PetscMPIIntCast(m, &mm)); 901252a1336SBarry Smith if (!rank) { 902252a1336SBarry Smith PetscCall(PetscMalloc1(M + 1, &II)); 903252a1336SBarry Smith PetscCall(PetscMalloc2(size, &allm, size, &dispm)); 904252a1336SBarry Smith PetscCallMPI(MPI_Gather(&mm, 1, MPI_INT, allm, 1, MPI_INT, 0, PetscObjectComm((PetscObject)A))); 9059371c9d4SSatish Balay dispm[0] = 0; 9069371c9d4SSatish Balay for (i = 1; i < size; i++) dispm[i] = dispm[i - 1] + allm[i - 1]; 907252a1336SBarry Smith PetscCallMPI(MPI_Gatherv(adj->i, mm, MPIU_INT, II, allm, dispm, MPIU_INT, 0, PetscObjectComm((PetscObject)A))); 908252a1336SBarry Smith PetscCall(PetscFree2(allm, dispm)); 909252a1336SBarry Smith II[M] = NZ; 910252a1336SBarry Smith } else { 911252a1336SBarry Smith PetscCallMPI(MPI_Gather(&mm, 1, MPI_INT, NULL, 1, MPI_INT, 0, PetscObjectComm((PetscObject)A))); 912252a1336SBarry Smith PetscCallMPI(MPI_Gatherv(adj->i, mm, MPIU_INT, NULL, NULL, NULL, MPIU_INT, 0, PetscObjectComm((PetscObject)A))); 913252a1336SBarry Smith } 914252a1336SBarry Smith /* shift the i[] values back */ 915252a1336SBarry Smith for (i = 0; i < m; i++) adj->i[i] -= nzstart; 916252a1336SBarry Smith if (!rank) PetscCall(MatCreateMPIAdj(PETSC_COMM_SELF, M, N, II, J, Values, B)); 9173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 918252a1336SBarry Smith } 919252a1336SBarry Smith 9209a3665c6SJed Brown /*@ 92111a5261eSBarry Smith MatMPIAdjCreateNonemptySubcommMat - create the same `MATMPIADJ` matrix on a subcommunicator containing only processes owning a positive number of rows 9229a3665c6SJed Brown 9239a3665c6SJed Brown Collective 9249a3665c6SJed Brown 9254165533cSJose E. Roman Input Parameter: 9262ef1f0ffSBarry Smith . A - original `MATMPIADJ` matrix 9279a3665c6SJed Brown 9284165533cSJose E. Roman Output Parameter: 9292ef1f0ffSBarry Smith . B - matrix on subcommunicator, `NULL` on ranks that owned zero rows of `A` 9309a3665c6SJed Brown 9319a3665c6SJed Brown Level: developer 9329a3665c6SJed Brown 9339a3665c6SJed Brown Note: 9349a3665c6SJed Brown This function is mostly useful for internal use by mesh partitioning packages that require that every process owns at least one row. 9359a3665c6SJed Brown 9362ef1f0ffSBarry Smith The matrix `B` should be destroyed with `MatDestroy()`. The arrays are not copied, so `B` should be destroyed before `A` is destroyed. 9379a3665c6SJed Brown 938*1cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATMPIADJ`, `MatCreateMPIAdj()` 9399a3665c6SJed Brown @*/ 940d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMPIAdjCreateNonemptySubcommMat(Mat A, Mat *B) 941d71ae5a4SJacob Faibussowitsch { 9429a3665c6SJed Brown PetscFunctionBegin; 9439a3665c6SJed Brown PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 944cac4c232SBarry Smith PetscUseMethod(A, "MatMPIAdjCreateNonemptySubcommMat_C", (Mat, Mat *), (A, B)); 9453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9469a3665c6SJed Brown } 9479a3665c6SJed Brown 9480bad9183SKris Buschelman /*MC 949fafad747SKris Buschelman MATMPIADJ - MATMPIADJ = "mpiadj" - A matrix type to be used for distributed adjacency matrices, 9500bad9183SKris Buschelman intended for use constructing orderings and partitionings. 9510bad9183SKris Buschelman 9520bad9183SKris Buschelman Level: beginner 9530bad9183SKris Buschelman 95411a5261eSBarry Smith Note: 95511a5261eSBarry Smith You can provide values to the matrix using `MatMPIAdjSetPreallocation()`, `MatCreateMPIAdj()`, or 95611a5261eSBarry Smith by calling `MatSetValues()` and `MatAssemblyBegin()` followed by `MatAssemblyEnd()` 9576a09307cSBarry Smith 958*1cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatCreateMPIAdj()`, `MatMPIAdjSetPreallocation()`, `MatSetValues()` 9590bad9183SKris Buschelman M*/ 960d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatCreate_MPIAdj(Mat B) 961d71ae5a4SJacob Faibussowitsch { 962273d9f13SBarry Smith Mat_MPIAdj *b; 963273d9f13SBarry Smith 964273d9f13SBarry Smith PetscFunctionBegin; 9654dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&b)); 966b0a32e0cSBarry Smith B->data = (void *)b; 9679566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(B->ops, &MatOps_Values, sizeof(struct _MatOps))); 968273d9f13SBarry Smith B->assembled = PETSC_FALSE; 9696a09307cSBarry Smith B->preallocated = PETSC_TRUE; /* so that MatSetValues() may be used */ 970273d9f13SBarry Smith 9719566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMPIAdjSetPreallocation_C", MatMPIAdjSetPreallocation_MPIAdj)); 9729566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMPIAdjCreateNonemptySubcommMat_C", MatMPIAdjCreateNonemptySubcommMat_MPIAdj)); 9739566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMPIAdjToSeq_C", MatMPIAdjToSeq_MPIAdj)); 974252a1336SBarry Smith PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMPIAdjToSeqRankZero_C", MatMPIAdjToSeqRankZero_MPIAdj)); 9759566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)B, MATMPIADJ)); 9763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 977273d9f13SBarry Smith } 978273d9f13SBarry Smith 979a23d5eceSKris Buschelman /*@C 98011a5261eSBarry Smith MatMPIAdjToSeq - Converts an parallel `MATMPIADJ` matrix to complete `MATMPIADJ` on each process (needed by sequential partitioner) 981b44e7856SBarry Smith 982d083f849SBarry Smith Logically Collective 983b44e7856SBarry Smith 984b44e7856SBarry Smith Input Parameter: 985b44e7856SBarry Smith . A - the matrix 986b44e7856SBarry Smith 987b44e7856SBarry Smith Output Parameter: 988b44e7856SBarry Smith . B - the same matrix on all processes 989b44e7856SBarry Smith 990b44e7856SBarry Smith Level: intermediate 991b44e7856SBarry Smith 992*1cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATMPIADJ`, `MatCreate()`, `MatCreateMPIAdj()`, `MatSetValues()`, `MatMPIAdjToSeqRankZero()` 993b44e7856SBarry Smith @*/ 994d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMPIAdjToSeq(Mat A, Mat *B) 995d71ae5a4SJacob Faibussowitsch { 996b44e7856SBarry Smith PetscFunctionBegin; 997cac4c232SBarry Smith PetscUseMethod(A, "MatMPIAdjToSeq_C", (Mat, Mat *), (A, B)); 9983ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 999b44e7856SBarry Smith } 1000b44e7856SBarry Smith 1001b44e7856SBarry Smith /*@C 100211a5261eSBarry Smith MatMPIAdjToSeqRankZero - Converts an parallel `MATMPIADJ` matrix to complete `MATMPIADJ` on rank zero (needed by sequential partitioner) 1003252a1336SBarry Smith 1004252a1336SBarry Smith Logically Collective 1005252a1336SBarry Smith 1006252a1336SBarry Smith Input Parameter: 1007252a1336SBarry Smith . A - the matrix 1008252a1336SBarry Smith 1009252a1336SBarry Smith Output Parameter: 1010252a1336SBarry Smith . B - the same matrix on rank zero, not set on other ranks 1011252a1336SBarry Smith 1012252a1336SBarry Smith Level: intermediate 1013252a1336SBarry Smith 101411a5261eSBarry Smith Note: 1015252a1336SBarry Smith This routine has the advantage on systems with multiple ranks per node since only one copy of the matrix 1016252a1336SBarry Smith is stored on the first node, instead of the number of ranks copies. This can allow partitioning much larger 1017d5b43468SJose E. Roman parallel graph sequentially. 1018252a1336SBarry Smith 1019*1cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATMPIADJ`, `MatCreate()`, `MatCreateMPIAdj()`, `MatSetValues()`, `MatMPIAdjToSeq()` 1020252a1336SBarry Smith @*/ 1021d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMPIAdjToSeqRankZero(Mat A, Mat *B) 1022d71ae5a4SJacob Faibussowitsch { 1023252a1336SBarry Smith PetscFunctionBegin; 1024252a1336SBarry Smith PetscUseMethod(A, "MatMPIAdjToSeqRankZero_C", (Mat, Mat *), (A, B)); 10253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1026252a1336SBarry Smith } 1027252a1336SBarry Smith 1028252a1336SBarry Smith /*@C 1029a23d5eceSKris Buschelman MatMPIAdjSetPreallocation - Sets the array used for storing the matrix elements 1030a23d5eceSKris Buschelman 1031d083f849SBarry Smith Logically Collective 1032a23d5eceSKris Buschelman 1033a23d5eceSKris Buschelman Input Parameters: 1034a23d5eceSKris Buschelman + A - the matrix 10352ef1f0ffSBarry Smith . i - the indices into `j` for the start of each row 1036a23d5eceSKris Buschelman . j - the column indices for each row (sorted for each row). 10372ef1f0ffSBarry Smith The indices in `i` and `j` start with zero (NOT with one). 10382ef1f0ffSBarry Smith - values - [use `NULL` if not provided] edge weights 1039a23d5eceSKris Buschelman 1040a23d5eceSKris Buschelman Level: intermediate 1041a23d5eceSKris Buschelman 1042*1cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatCreate()`, `MatCreateMPIAdj()`, `MatSetValues()`, `MATMPIADJ`, `MatCreateMPIAdj()` 1043a23d5eceSKris Buschelman @*/ 1044d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMPIAdjSetPreallocation(Mat B, PetscInt *i, PetscInt *j, PetscInt *values) 1045d71ae5a4SJacob Faibussowitsch { 1046273d9f13SBarry Smith PetscFunctionBegin; 1047cac4c232SBarry Smith PetscTryMethod(B, "MatMPIAdjSetPreallocation_C", (Mat, PetscInt *, PetscInt *, PetscInt *), (B, i, j, values)); 10483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1049273d9f13SBarry Smith } 1050273d9f13SBarry Smith 1051b97cf49bSBarry Smith /*@C 10523eda8832SBarry Smith MatCreateMPIAdj - Creates a sparse matrix representing an adjacency list. 1053b97cf49bSBarry Smith The matrix does not have numerical values associated with it, but is 1054b97cf49bSBarry Smith intended for ordering (to reduce bandwidth etc) and partitioning. 1055b97cf49bSBarry Smith 1056d083f849SBarry Smith Collective 1057ef5ee4d1SLois Curfman McInnes 1058b97cf49bSBarry Smith Input Parameters: 1059c2e958feSBarry Smith + comm - MPI communicator 10600752156aSBarry Smith . m - number of local rows 106158ec18a5SLisandro Dalcin . N - number of global columns 10622ef1f0ffSBarry Smith . i - the indices into `j` for the start of each row 1063329f5518SBarry Smith . j - the column indices for each row (sorted for each row). 10642ef1f0ffSBarry Smith The indices in `i` and `j` start with zero (NOT with one). 10652ef1f0ffSBarry Smith - values -[use `NULL` if not provided] edge weights 1066b97cf49bSBarry Smith 1067b97cf49bSBarry Smith Output Parameter: 1068b97cf49bSBarry Smith . A - the matrix 1069b97cf49bSBarry Smith 107015091d37SBarry Smith Level: intermediate 107115091d37SBarry Smith 107295452b02SPatrick Sanan Notes: 10732ef1f0ffSBarry Smith You must NOT free the `i`, `values` and `j` arrays yourself. PETSc will free them 107411a5261eSBarry Smith when the matrix is destroyed; you must allocate them with `PetscMalloc()`. If you 107511a5261eSBarry Smith call from Fortran you need not create the arrays with `PetscMalloc()`. 10766a09307cSBarry Smith 10776a09307cSBarry Smith You should not include the matrix diagonals. 1078b97cf49bSBarry Smith 1079e3f7e1d6Svictorle If you already have a matrix, you can create its adjacency matrix by a call 108011a5261eSBarry Smith to `MatConvert()`, specifying a type of `MATMPIADJ`. 1081ededeb1aSvictorle 108211a5261eSBarry Smith Possible values for `MatSetOption()` - `MAT_STRUCTURALLY_SYMMETRIC` 1083b97cf49bSBarry Smith 1084*1cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatCreate()`, `MatConvert()`, `MatGetOrdering()`, `MATMPIADJ`, `MatMPIAdjSetPreallocation()` 1085b97cf49bSBarry Smith @*/ 1086d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCreateMPIAdj(MPI_Comm comm, PetscInt m, PetscInt N, PetscInt *i, PetscInt *j, PetscInt *values, Mat *A) 1087d71ae5a4SJacob Faibussowitsch { 1088433994e6SBarry Smith PetscFunctionBegin; 10899566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, A)); 10909566063dSJacob Faibussowitsch PetscCall(MatSetSizes(*A, m, PETSC_DETERMINE, PETSC_DETERMINE, N)); 10919566063dSJacob Faibussowitsch PetscCall(MatSetType(*A, MATMPIADJ)); 10929566063dSJacob Faibussowitsch PetscCall(MatMPIAdjSetPreallocation(*A, i, j, values)); 10933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1094b97cf49bSBarry Smith } 1095