xref: /petsc/src/mat/impls/adj/mpi/mpiadj.c (revision 4dfa11a44d5adf2389f1d3acbc8f3c1116dc6c3a)
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 /*
97dae84e0SHong Zhang  * The interface should be easy to use for both MatCreateSubMatrix (parallel sub-matrix) and MatCreateSubMatrices (sequential sub-matrices)
1040244768SBarry Smith  * */
119371c9d4SSatish Balay static PetscErrorCode MatCreateSubMatrix_MPIAdj_data(Mat adj, IS irows, IS icols, PetscInt **sadj_xadj, PetscInt **sadj_adjncy, PetscInt **sadj_values) {
12131c27b5Sprj-   PetscInt        nlrows_is, icols_n, i, j, nroots, nleaves, rlocalindex, *ncols_send, *ncols_recv;
1340244768SBarry Smith   PetscInt        nlrows_mat, *adjncy_recv, Ncols_recv, Ncols_send, *xadj_recv, *values_recv;
14e895ccc0SFande Kong   PetscInt       *ncols_recv_offsets, loc, rnclos, *sadjncy, *sxadj, *svalues;
1540244768SBarry Smith   const PetscInt *irows_indices, *icols_indices, *xadj, *adjncy;
16131c27b5Sprj-   PetscMPIInt     owner;
1740244768SBarry Smith   Mat_MPIAdj     *a = (Mat_MPIAdj *)adj->data;
1840244768SBarry Smith   PetscLayout     rmap;
1940244768SBarry Smith   MPI_Comm        comm;
2040244768SBarry Smith   PetscSF         sf;
2140244768SBarry Smith   PetscSFNode    *iremote;
2240244768SBarry Smith   PetscBool       done;
2340244768SBarry Smith 
2440244768SBarry Smith   PetscFunctionBegin;
259566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)adj, &comm));
269566063dSJacob Faibussowitsch   PetscCall(MatGetLayouts(adj, &rmap, NULL));
279566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(irows, &nlrows_is));
289566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(irows, &irows_indices));
299566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nlrows_is, &iremote));
3040244768SBarry Smith   /* construct sf graph*/
3140244768SBarry Smith   nleaves = nlrows_is;
3240244768SBarry Smith   for (i = 0; i < nlrows_is; i++) {
3340244768SBarry Smith     owner       = -1;
3440244768SBarry Smith     rlocalindex = -1;
359566063dSJacob Faibussowitsch     PetscCall(PetscLayoutFindOwnerIndex(rmap, irows_indices[i], &owner, &rlocalindex));
3640244768SBarry Smith     iremote[i].rank  = owner;
3740244768SBarry Smith     iremote[i].index = rlocalindex;
3840244768SBarry Smith   }
399566063dSJacob Faibussowitsch   PetscCall(MatGetRowIJ(adj, 0, PETSC_FALSE, PETSC_FALSE, &nlrows_mat, &xadj, &adjncy, &done));
409566063dSJacob Faibussowitsch   PetscCall(PetscCalloc4(nlrows_mat, &ncols_send, nlrows_is, &xadj_recv, nlrows_is + 1, &ncols_recv_offsets, nlrows_is, &ncols_recv));
4140244768SBarry Smith   nroots = nlrows_mat;
42ad540459SPierre Jolivet   for (i = 0; i < nlrows_mat; i++) ncols_send[i] = xadj[i + 1] - xadj[i];
439566063dSJacob Faibussowitsch   PetscCall(PetscSFCreate(comm, &sf));
449566063dSJacob Faibussowitsch   PetscCall(PetscSFSetGraph(sf, nroots, nleaves, NULL, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER));
459566063dSJacob Faibussowitsch   PetscCall(PetscSFSetType(sf, PETSCSFBASIC));
469566063dSJacob Faibussowitsch   PetscCall(PetscSFSetFromOptions(sf));
479566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastBegin(sf, MPIU_INT, ncols_send, ncols_recv, MPI_REPLACE));
489566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastEnd(sf, MPIU_INT, ncols_send, ncols_recv, MPI_REPLACE));
499566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastBegin(sf, MPIU_INT, xadj, xadj_recv, MPI_REPLACE));
509566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastEnd(sf, MPIU_INT, xadj, xadj_recv, MPI_REPLACE));
519566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&sf));
5240244768SBarry Smith   Ncols_recv = 0;
5340244768SBarry Smith   for (i = 0; i < nlrows_is; i++) {
5440244768SBarry Smith     Ncols_recv += ncols_recv[i];
5540244768SBarry Smith     ncols_recv_offsets[i + 1] = ncols_recv[i] + ncols_recv_offsets[i];
5640244768SBarry Smith   }
5740244768SBarry Smith   Ncols_send = 0;
58ad540459SPierre Jolivet   for (i = 0; i < nlrows_mat; i++) Ncols_send += ncols_send[i];
599566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(Ncols_recv, &iremote));
609566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(Ncols_recv, &adjncy_recv));
6140244768SBarry Smith   nleaves    = Ncols_recv;
6240244768SBarry Smith   Ncols_recv = 0;
6340244768SBarry Smith   for (i = 0; i < nlrows_is; i++) {
649566063dSJacob Faibussowitsch     PetscCall(PetscLayoutFindOwner(rmap, irows_indices[i], &owner));
6540244768SBarry Smith     for (j = 0; j < ncols_recv[i]; j++) {
6640244768SBarry Smith       iremote[Ncols_recv].rank    = owner;
6740244768SBarry Smith       iremote[Ncols_recv++].index = xadj_recv[i] + j;
6840244768SBarry Smith     }
6940244768SBarry Smith   }
709566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(irows, &irows_indices));
7140244768SBarry Smith   /*if we need to deal with edge weights ???*/
729566063dSJacob Faibussowitsch   if (a->useedgeweights) PetscCall(PetscCalloc1(Ncols_recv, &values_recv));
7340244768SBarry Smith   nroots = Ncols_send;
749566063dSJacob Faibussowitsch   PetscCall(PetscSFCreate(comm, &sf));
759566063dSJacob Faibussowitsch   PetscCall(PetscSFSetGraph(sf, nroots, nleaves, NULL, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER));
769566063dSJacob Faibussowitsch   PetscCall(PetscSFSetType(sf, PETSCSFBASIC));
779566063dSJacob Faibussowitsch   PetscCall(PetscSFSetFromOptions(sf));
789566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastBegin(sf, MPIU_INT, adjncy, adjncy_recv, MPI_REPLACE));
799566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastEnd(sf, MPIU_INT, adjncy, adjncy_recv, MPI_REPLACE));
80e895ccc0SFande Kong   if (a->useedgeweights) {
819566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastBegin(sf, MPIU_INT, a->values, values_recv, MPI_REPLACE));
829566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastEnd(sf, MPIU_INT, a->values, values_recv, MPI_REPLACE));
8340244768SBarry Smith   }
849566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&sf));
859566063dSJacob Faibussowitsch   PetscCall(MatRestoreRowIJ(adj, 0, PETSC_FALSE, PETSC_FALSE, &nlrows_mat, &xadj, &adjncy, &done));
869566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(icols, &icols_n));
879566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(icols, &icols_indices));
8840244768SBarry Smith   rnclos = 0;
8940244768SBarry Smith   for (i = 0; i < nlrows_is; i++) {
9040244768SBarry Smith     for (j = ncols_recv_offsets[i]; j < ncols_recv_offsets[i + 1]; j++) {
919566063dSJacob Faibussowitsch       PetscCall(PetscFindInt(adjncy_recv[j], icols_n, icols_indices, &loc));
9240244768SBarry Smith       if (loc < 0) {
9340244768SBarry Smith         adjncy_recv[j] = -1;
94e895ccc0SFande Kong         if (a->useedgeweights) values_recv[j] = -1;
9540244768SBarry Smith         ncols_recv[i]--;
9640244768SBarry Smith       } else {
9740244768SBarry Smith         rnclos++;
9840244768SBarry Smith       }
9940244768SBarry Smith     }
10040244768SBarry Smith   }
1019566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(icols, &icols_indices));
1029566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(rnclos, &sadjncy));
1039566063dSJacob Faibussowitsch   if (a->useedgeweights) PetscCall(PetscCalloc1(rnclos, &svalues));
1049566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(nlrows_is + 1, &sxadj));
10540244768SBarry Smith   rnclos = 0;
10640244768SBarry Smith   for (i = 0; i < nlrows_is; i++) {
10740244768SBarry Smith     for (j = ncols_recv_offsets[i]; j < ncols_recv_offsets[i + 1]; j++) {
10840244768SBarry Smith       if (adjncy_recv[j] < 0) continue;
10940244768SBarry Smith       sadjncy[rnclos] = adjncy_recv[j];
110e895ccc0SFande Kong       if (a->useedgeweights) svalues[rnclos] = values_recv[j];
11140244768SBarry Smith       rnclos++;
11240244768SBarry Smith     }
11340244768SBarry Smith   }
114ad540459SPierre Jolivet   for (i = 0; i < nlrows_is; i++) sxadj[i + 1] = sxadj[i] + ncols_recv[i];
1159371c9d4SSatish Balay   if (sadj_xadj) {
1169371c9d4SSatish Balay     *sadj_xadj = sxadj;
1179371c9d4SSatish Balay   } else PetscCall(PetscFree(sxadj));
1189371c9d4SSatish Balay   if (sadj_adjncy) {
1199371c9d4SSatish Balay     *sadj_adjncy = sadjncy;
1209371c9d4SSatish Balay   } else PetscCall(PetscFree(sadjncy));
12140244768SBarry Smith   if (sadj_values) {
1229371c9d4SSatish Balay     if (a->useedgeweights) *sadj_values = svalues;
1239371c9d4SSatish Balay     else *sadj_values = NULL;
12440244768SBarry Smith   } else {
1259566063dSJacob Faibussowitsch     if (a->useedgeweights) PetscCall(PetscFree(svalues));
12640244768SBarry Smith   }
1279566063dSJacob Faibussowitsch   PetscCall(PetscFree4(ncols_send, xadj_recv, ncols_recv_offsets, ncols_recv));
1289566063dSJacob Faibussowitsch   PetscCall(PetscFree(adjncy_recv));
1299566063dSJacob Faibussowitsch   if (a->useedgeweights) PetscCall(PetscFree(values_recv));
13040244768SBarry Smith   PetscFunctionReturn(0);
13140244768SBarry Smith }
13240244768SBarry Smith 
1339371c9d4SSatish Balay static PetscErrorCode MatCreateSubMatrices_MPIAdj_Private(Mat mat, PetscInt n, const IS irow[], const IS icol[], PetscBool subcomm, MatReuse scall, Mat *submat[]) {
13440244768SBarry Smith   PetscInt        i, irow_n, icol_n, *sxadj, *sadjncy, *svalues;
13540244768SBarry Smith   PetscInt       *indices, nindx, j, k, loc;
13640244768SBarry Smith   PetscMPIInt     issame;
13740244768SBarry Smith   const PetscInt *irow_indices, *icol_indices;
13840244768SBarry Smith   MPI_Comm        scomm_row, scomm_col, scomm_mat;
13940244768SBarry Smith 
14040244768SBarry Smith   PetscFunctionBegin;
14140244768SBarry Smith   nindx = 0;
14240244768SBarry Smith   /*
14340244768SBarry Smith    * Estimate a maximum number for allocating memory
14440244768SBarry Smith    */
14540244768SBarry Smith   for (i = 0; i < n; i++) {
1469566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(irow[i], &irow_n));
1479566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(icol[i], &icol_n));
14840244768SBarry Smith     nindx = nindx > (irow_n + icol_n) ? nindx : (irow_n + icol_n);
14940244768SBarry Smith   }
1509566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nindx, &indices));
15140244768SBarry Smith   /* construct a submat */
152252a1336SBarry Smith   //  if (scall == MAT_INITIAL_MATRIX) PetscCall(PetscMalloc1(n,submat));
1536a7b62d2SBarry Smith 
15440244768SBarry Smith   for (i = 0; i < n; i++) {
15540244768SBarry Smith     if (subcomm) {
1569566063dSJacob Faibussowitsch       PetscCall(PetscObjectGetComm((PetscObject)irow[i], &scomm_row));
1579566063dSJacob Faibussowitsch       PetscCall(PetscObjectGetComm((PetscObject)icol[i], &scomm_col));
1589566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_compare(scomm_row, scomm_col, &issame));
15908401ef6SPierre 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");
1609566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_compare(scomm_row, PETSC_COMM_SELF, &issame));
16108401ef6SPierre 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");
16240244768SBarry Smith     } else {
16340244768SBarry Smith       scomm_row = PETSC_COMM_SELF;
16440244768SBarry Smith     }
16540244768SBarry Smith     /*get sub-matrix data*/
1669371c9d4SSatish Balay     sxadj   = NULL;
1679371c9d4SSatish Balay     sadjncy = NULL;
1689371c9d4SSatish Balay     svalues = NULL;
1699566063dSJacob Faibussowitsch     PetscCall(MatCreateSubMatrix_MPIAdj_data(mat, irow[i], icol[i], &sxadj, &sadjncy, &svalues));
1709566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(irow[i], &irow_n));
1719566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(icol[i], &icol_n));
1729566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(irow[i], &irow_indices));
1739566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(indices, irow_indices, irow_n));
1749566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(irow[i], &irow_indices));
1759566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(icol[i], &icol_indices));
1769566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(indices + irow_n, icol_indices, icol_n));
1779566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(icol[i], &icol_indices));
17840244768SBarry Smith     nindx = irow_n + icol_n;
1799566063dSJacob Faibussowitsch     PetscCall(PetscSortRemoveDupsInt(&nindx, indices));
18040244768SBarry Smith     /* renumber columns */
18140244768SBarry Smith     for (j = 0; j < irow_n; j++) {
18240244768SBarry Smith       for (k = sxadj[j]; k < sxadj[j + 1]; k++) {
1839566063dSJacob Faibussowitsch         PetscCall(PetscFindInt(sadjncy[k], nindx, indices, &loc));
18408401ef6SPierre Jolivet         PetscCheck(loc >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "can not find col %" PetscInt_FMT, sadjncy[k]);
18540244768SBarry Smith         sadjncy[k] = loc;
18640244768SBarry Smith       }
18740244768SBarry Smith     }
18840244768SBarry Smith     if (scall == MAT_INITIAL_MATRIX) {
1899566063dSJacob Faibussowitsch       PetscCall(MatCreateMPIAdj(scomm_row, irow_n, icol_n, sxadj, sadjncy, svalues, submat[i]));
19040244768SBarry Smith     } else {
19140244768SBarry Smith       Mat         sadj = *(submat[i]);
19240244768SBarry Smith       Mat_MPIAdj *sa   = (Mat_MPIAdj *)((sadj)->data);
1939566063dSJacob Faibussowitsch       PetscCall(PetscObjectGetComm((PetscObject)sadj, &scomm_mat));
1949566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_compare(scomm_row, scomm_mat, &issame));
19508401ef6SPierre Jolivet       PetscCheck(issame == MPI_IDENT, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "submatrix  must have the same comm as the col index set");
1969566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(sa->i, sxadj, irow_n + 1));
1979566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(sa->j, sadjncy, sxadj[irow_n]));
1989566063dSJacob Faibussowitsch       if (svalues) PetscCall(PetscArraycpy(sa->values, svalues, sxadj[irow_n]));
1999566063dSJacob Faibussowitsch       PetscCall(PetscFree(sxadj));
2009566063dSJacob Faibussowitsch       PetscCall(PetscFree(sadjncy));
201252a1336SBarry Smith       PetscCall(PetscFree(svalues));
20240244768SBarry Smith     }
20340244768SBarry Smith   }
2049566063dSJacob Faibussowitsch   PetscCall(PetscFree(indices));
20540244768SBarry Smith   PetscFunctionReturn(0);
20640244768SBarry Smith }
20740244768SBarry Smith 
2089371c9d4SSatish Balay static PetscErrorCode MatCreateSubMatricesMPI_MPIAdj(Mat mat, PetscInt n, const IS irow[], const IS icol[], MatReuse scall, Mat *submat[]) {
20940244768SBarry Smith   /*get sub-matrices across a sub communicator */
21040244768SBarry Smith   PetscFunctionBegin;
2119566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices_MPIAdj_Private(mat, n, irow, icol, PETSC_TRUE, scall, submat));
21240244768SBarry Smith   PetscFunctionReturn(0);
21340244768SBarry Smith }
21440244768SBarry Smith 
2159371c9d4SSatish Balay static PetscErrorCode MatCreateSubMatrices_MPIAdj(Mat mat, PetscInt n, const IS irow[], const IS icol[], MatReuse scall, Mat *submat[]) {
21640244768SBarry Smith   PetscFunctionBegin;
21740244768SBarry Smith   /*get sub-matrices based on PETSC_COMM_SELF */
2189566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices_MPIAdj_Private(mat, n, irow, icol, PETSC_FALSE, scall, submat));
21940244768SBarry Smith   PetscFunctionReturn(0);
22040244768SBarry Smith }
22140244768SBarry Smith 
2229371c9d4SSatish Balay static PetscErrorCode MatView_MPIAdj_ASCII(Mat A, PetscViewer viewer) {
2233eda8832SBarry Smith   Mat_MPIAdj       *a = (Mat_MPIAdj *)A->data;
224d0f46423SBarry Smith   PetscInt          i, j, m = A->rmap->n;
2252dcb1b2aSMatthew Knepley   const char       *name;
226ce1411ecSBarry Smith   PetscViewerFormat format;
227b97cf49bSBarry Smith 
228433994e6SBarry Smith   PetscFunctionBegin;
2299566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetName((PetscObject)A, &name));
2309566063dSJacob Faibussowitsch   PetscCall(PetscViewerGetFormat(viewer, &format));
231fb9695e5SSatish Balay   if (format == PETSC_VIEWER_ASCII_INFO) {
2323a40ed3dSBarry Smith     PetscFunctionReturn(0);
233f7d195e4SLawrence Mitchell   } else {
234f7d195e4SLawrence Mitchell     PetscCheck(format != PETSC_VIEWER_ASCII_MATLAB, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "MATLAB format not supported");
2359566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE));
2369566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushSynchronized(viewer));
237b97cf49bSBarry Smith     for (i = 0; i < m; i++) {
2389566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "row %" PetscInt_FMT ":", i + A->rmap->rstart));
239b97cf49bSBarry Smith       for (j = a->i[i]; j < a->i[i + 1]; j++) {
2400170919eSFande Kong         if (a->values) {
2419566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " (%" PetscInt_FMT ", %" PetscInt_FMT ") ", a->j[j], a->values[j]));
2420170919eSFande Kong         } else {
2439566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " %" PetscInt_FMT " ", a->j[j]));
2440752156aSBarry Smith         }
2450170919eSFande Kong       }
2469566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "\n"));
247b97cf49bSBarry Smith     }
2489566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE));
2499566063dSJacob Faibussowitsch     PetscCall(PetscViewerFlush(viewer));
2509566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopSynchronized(viewer));
2517b23a99aSBarry Smith   }
2523a40ed3dSBarry Smith   PetscFunctionReturn(0);
253b97cf49bSBarry Smith }
254b97cf49bSBarry Smith 
2559371c9d4SSatish Balay static PetscErrorCode MatView_MPIAdj(Mat A, PetscViewer viewer) {
256ace3abfcSBarry Smith   PetscBool iascii;
257b97cf49bSBarry Smith 
258433994e6SBarry Smith   PetscFunctionBegin;
2599566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
2601baa6e33SBarry Smith   if (iascii) PetscCall(MatView_MPIAdj_ASCII(A, viewer));
2613a40ed3dSBarry Smith   PetscFunctionReturn(0);
262b97cf49bSBarry Smith }
263b97cf49bSBarry Smith 
2649371c9d4SSatish Balay static PetscErrorCode MatDestroy_MPIAdj(Mat mat) {
2653eda8832SBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj *)mat->data;
26694d884c6SBarry Smith 
26794d884c6SBarry Smith   PetscFunctionBegin;
268aa482453SBarry Smith #if defined(PETSC_USE_LOG)
269c0aa6a63SJacob Faibussowitsch   PetscLogObjectState((PetscObject)mat, "Rows=%" PetscInt_FMT ", Cols=%" PetscInt_FMT ", NZ=%" PetscInt_FMT, mat->rmap->n, mat->cmap->n, a->nz);
270b97cf49bSBarry Smith #endif
2719566063dSJacob Faibussowitsch   PetscCall(PetscFree(a->diag));
2720400557aSBarry Smith   if (a->freeaij) {
27314196391SBarry Smith     if (a->freeaijwithfree) {
27414196391SBarry Smith       if (a->i) free(a->i);
27514196391SBarry Smith       if (a->j) free(a->j);
27614196391SBarry Smith     } else {
2779566063dSJacob Faibussowitsch       PetscCall(PetscFree(a->i));
2789566063dSJacob Faibussowitsch       PetscCall(PetscFree(a->j));
2799566063dSJacob Faibussowitsch       PetscCall(PetscFree(a->values));
28014196391SBarry Smith     }
2810400557aSBarry Smith   }
2829566063dSJacob Faibussowitsch   PetscCall(PetscFree(a->rowvalues));
2839566063dSJacob Faibussowitsch   PetscCall(PetscFree(mat->data));
2849566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)mat, NULL));
2859566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMPIAdjSetPreallocation_C", NULL));
2869566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMPIAdjCreateNonemptySubcommMat_C", NULL));
2872e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMPIAdjToSeq_C", NULL));
288252a1336SBarry Smith   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMPIAdjToSeqRankZero_C", NULL));
2893a40ed3dSBarry Smith   PetscFunctionReturn(0);
290b97cf49bSBarry Smith }
291b97cf49bSBarry Smith 
2929371c9d4SSatish Balay static PetscErrorCode MatSetOption_MPIAdj(Mat A, MatOption op, PetscBool flg) {
2933eda8832SBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj *)A->data;
294b97cf49bSBarry Smith 
295433994e6SBarry Smith   PetscFunctionBegin;
29612c028f9SKris Buschelman   switch (op) {
29777e54ba9SKris Buschelman   case MAT_SYMMETRIC:
29812c028f9SKris Buschelman   case MAT_STRUCTURALLY_SYMMETRIC:
2999a4540c5SBarry Smith   case MAT_HERMITIAN:
3009371c9d4SSatish Balay   case MAT_SPD: a->symmetric = flg; break;
3019a4540c5SBarry Smith   case MAT_SYMMETRY_ETERNAL:
302b94d7dedSBarry Smith   case MAT_STRUCTURAL_SYMMETRY_ETERNAL:
3039371c9d4SSatish Balay   case MAT_SPD_ETERNAL: break;
3049371c9d4SSatish Balay   default: PetscCall(PetscInfo(A, "Option %s ignored\n", MatOptions[op])); break;
305b97cf49bSBarry Smith   }
3063a40ed3dSBarry Smith   PetscFunctionReturn(0);
307b97cf49bSBarry Smith }
308b97cf49bSBarry Smith 
3099371c9d4SSatish Balay static PetscErrorCode MatGetRow_MPIAdj(Mat A, PetscInt row, PetscInt *nz, PetscInt **idx, PetscScalar **v) {
3103eda8832SBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj *)A->data;
311b97cf49bSBarry Smith 
312433994e6SBarry Smith   PetscFunctionBegin;
313d0f46423SBarry Smith   row -= A->rmap->rstart;
314aed4548fSBarry Smith   PetscCheck(row >= 0 && row < A->rmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Row out of range");
315b97cf49bSBarry Smith   *nz = a->i[row + 1] - a->i[row];
3162b1d8763SJed Brown   if (v) {
3172b1d8763SJed Brown     PetscInt j;
3182b1d8763SJed Brown     if (a->rowvalues_alloc < *nz) {
3199566063dSJacob Faibussowitsch       PetscCall(PetscFree(a->rowvalues));
3202b1d8763SJed Brown       a->rowvalues_alloc = PetscMax(a->rowvalues_alloc * 2, *nz);
3219566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(a->rowvalues_alloc, &a->rowvalues));
322b97cf49bSBarry Smith     }
323ad540459SPierre Jolivet     for (j = 0; j < *nz; j++) a->rowvalues[j] = a->values ? a->values[a->i[row] + j] : 1.0;
3242b1d8763SJed Brown     *v = (*nz) ? a->rowvalues : NULL;
325b97cf49bSBarry Smith   }
32692b6f2a0SJed Brown   if (idx) *idx = (*nz) ? a->j + a->i[row] : NULL;
3273a40ed3dSBarry Smith   PetscFunctionReturn(0);
328b97cf49bSBarry Smith }
329b97cf49bSBarry Smith 
3309371c9d4SSatish Balay static PetscErrorCode MatRestoreRow_MPIAdj(Mat A, PetscInt row, PetscInt *nz, PetscInt **idx, PetscScalar **v) {
331433994e6SBarry Smith   PetscFunctionBegin;
3323a40ed3dSBarry Smith   PetscFunctionReturn(0);
333b97cf49bSBarry Smith }
334b97cf49bSBarry Smith 
3359371c9d4SSatish Balay static PetscErrorCode MatEqual_MPIAdj(Mat A, Mat B, PetscBool *flg) {
3363eda8832SBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj *)A->data, *b = (Mat_MPIAdj *)B->data;
337ace3abfcSBarry Smith   PetscBool   flag;
338b97cf49bSBarry Smith 
339433994e6SBarry Smith   PetscFunctionBegin;
340b97cf49bSBarry Smith   /* If the  matrix dimensions are not equal,or no of nonzeros */
341ad540459SPierre Jolivet   if ((A->rmap->n != B->rmap->n) || (a->nz != b->nz)) flag = PETSC_FALSE;
342b97cf49bSBarry Smith 
343b97cf49bSBarry Smith   /* if the a->i are the same */
3449566063dSJacob Faibussowitsch   PetscCall(PetscArraycmp(a->i, b->i, A->rmap->n + 1, &flag));
345b97cf49bSBarry Smith 
346b97cf49bSBarry Smith   /* if a->j are the same */
3479566063dSJacob Faibussowitsch   PetscCall(PetscMemcmp(a->j, b->j, (a->nz) * sizeof(PetscInt), &flag));
348b97cf49bSBarry Smith 
3491c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(&flag, flg, 1, MPIU_BOOL, MPI_LAND, PetscObjectComm((PetscObject)A)));
3503a40ed3dSBarry Smith   PetscFunctionReturn(0);
351b97cf49bSBarry Smith }
352b97cf49bSBarry Smith 
3539371c9d4SSatish Balay static PetscErrorCode MatGetRowIJ_MPIAdj(Mat A, PetscInt oshift, PetscBool symmetric, PetscBool blockcompressed, PetscInt *m, const PetscInt *inia[], const PetscInt *inja[], PetscBool *done) {
354b24ad042SBarry Smith   PetscInt    i;
3556cd91f64SBarry Smith   Mat_MPIAdj *a  = (Mat_MPIAdj *)A->data;
3561a83f524SJed Brown   PetscInt  **ia = (PetscInt **)inia, **ja = (PetscInt **)inja;
3576cd91f64SBarry Smith 
3586cd91f64SBarry Smith   PetscFunctionBegin;
359d0f46423SBarry Smith   *m    = A->rmap->n;
3606cd91f64SBarry Smith   *ia   = a->i;
3616cd91f64SBarry Smith   *ja   = a->j;
3626cd91f64SBarry Smith   *done = PETSC_TRUE;
363d42735a1SBarry Smith   if (oshift) {
364ad540459SPierre Jolivet     for (i = 0; i < (*ia)[*m]; i++) (*ja)[i]++;
365d42735a1SBarry Smith     for (i = 0; i <= (*m); i++) (*ia)[i]++;
366d42735a1SBarry Smith   }
367d42735a1SBarry Smith   PetscFunctionReturn(0);
368d42735a1SBarry Smith }
369d42735a1SBarry Smith 
3709371c9d4SSatish Balay static PetscErrorCode MatRestoreRowIJ_MPIAdj(Mat A, PetscInt oshift, PetscBool symmetric, PetscBool blockcompressed, PetscInt *m, const PetscInt *inia[], const PetscInt *inja[], PetscBool *done) {
371b24ad042SBarry Smith   PetscInt    i;
372d42735a1SBarry Smith   Mat_MPIAdj *a  = (Mat_MPIAdj *)A->data;
3731a83f524SJed Brown   PetscInt  **ia = (PetscInt **)inia, **ja = (PetscInt **)inja;
374d42735a1SBarry Smith 
375d42735a1SBarry Smith   PetscFunctionBegin;
376aed4548fSBarry Smith   PetscCheck(!ia || a->i == *ia, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "ia passed back is not one obtained with MatGetRowIJ()");
377aed4548fSBarry Smith   PetscCheck(!ja || a->j == *ja, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "ja passed back is not one obtained with MatGetRowIJ()");
378d42735a1SBarry Smith   if (oshift) {
37928b400f6SJacob Faibussowitsch     PetscCheck(ia, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "If oshift then you must passed in inia[] argument");
38028b400f6SJacob Faibussowitsch     PetscCheck(ja, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "If oshift then you must passed in inja[] argument");
381d42735a1SBarry Smith     for (i = 0; i <= (*m); i++) (*ia)[i]--;
382ad540459SPierre Jolivet     for (i = 0; i < (*ia)[*m]; i++) (*ja)[i]--;
383d42735a1SBarry Smith   }
3846cd91f64SBarry Smith   PetscFunctionReturn(0);
3856cd91f64SBarry Smith }
386b97cf49bSBarry Smith 
3879371c9d4SSatish Balay PetscErrorCode MatConvertFrom_MPIAdj(Mat A, MatType type, MatReuse reuse, Mat *newmat) {
38817667f90SBarry Smith   Mat                B;
38917667f90SBarry Smith   PetscInt           i, m, N, nzeros = 0, *ia, *ja, len, rstart, cnt, j, *a;
39017667f90SBarry Smith   const PetscInt    *rj;
39117667f90SBarry Smith   const PetscScalar *ra;
39217667f90SBarry Smith   MPI_Comm           comm;
39317667f90SBarry Smith 
39417667f90SBarry Smith   PetscFunctionBegin;
3959566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A, NULL, &N));
3969566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A, &m, NULL));
3979566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(A, &rstart, NULL));
39817667f90SBarry Smith 
39917667f90SBarry Smith   /* count the number of nonzeros per row */
40017667f90SBarry Smith   for (i = 0; i < m; i++) {
4019566063dSJacob Faibussowitsch     PetscCall(MatGetRow(A, i + rstart, &len, &rj, NULL));
4025ee9ba1cSJed Brown     for (j = 0; j < len; j++) {
4039371c9d4SSatish Balay       if (rj[j] == i + rstart) {
4049371c9d4SSatish Balay         len--;
4059371c9d4SSatish Balay         break;
4069371c9d4SSatish Balay       } /* don't count diagonal */
4075ee9ba1cSJed Brown     }
40817667f90SBarry Smith     nzeros += len;
4099566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(A, i + rstart, &len, &rj, NULL));
41017667f90SBarry Smith   }
41117667f90SBarry Smith 
41217667f90SBarry Smith   /* malloc space for nonzeros */
4139566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nzeros + 1, &a));
4149566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(N + 1, &ia));
4159566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nzeros + 1, &ja));
41617667f90SBarry Smith 
41717667f90SBarry Smith   nzeros = 0;
41817667f90SBarry Smith   ia[0]  = 0;
41917667f90SBarry Smith   for (i = 0; i < m; i++) {
4209566063dSJacob Faibussowitsch     PetscCall(MatGetRow(A, i + rstart, &len, &rj, &ra));
42117667f90SBarry Smith     cnt = 0;
42217667f90SBarry Smith     for (j = 0; j < len; j++) {
4235ee9ba1cSJed Brown       if (rj[j] != i + rstart) { /* if not diagonal */
42417667f90SBarry Smith         a[nzeros + cnt]    = (PetscInt)PetscAbsScalar(ra[j]);
42517667f90SBarry Smith         ja[nzeros + cnt++] = rj[j];
42617667f90SBarry Smith       }
4275ee9ba1cSJed Brown     }
4289566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(A, i + rstart, &len, &rj, &ra));
42917667f90SBarry Smith     nzeros += cnt;
43017667f90SBarry Smith     ia[i + 1] = nzeros;
43117667f90SBarry Smith   }
43217667f90SBarry Smith 
4339566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)A, &comm));
4349566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, &B));
4359566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(B, m, PETSC_DETERMINE, PETSC_DETERMINE, N));
4369566063dSJacob Faibussowitsch   PetscCall(MatSetType(B, type));
4379566063dSJacob Faibussowitsch   PetscCall(MatMPIAdjSetPreallocation(B, ia, ja, a));
43817667f90SBarry Smith 
439511c6705SHong Zhang   if (reuse == MAT_INPLACE_MATRIX) {
4409566063dSJacob Faibussowitsch     PetscCall(MatHeaderReplace(A, &B));
44117667f90SBarry Smith   } else {
44217667f90SBarry Smith     *newmat = B;
44317667f90SBarry Smith   }
44417667f90SBarry Smith   PetscFunctionReturn(0);
44517667f90SBarry Smith }
44617667f90SBarry Smith 
4479371c9d4SSatish Balay PetscErrorCode MatSetValues_MPIAdj(Mat A, PetscInt m, const PetscInt *rows, PetscInt n, const PetscInt *cols, const PetscScalar *values, InsertMode im) {
4486a09307cSBarry Smith   Mat_MPIAdj *adj = (Mat_MPIAdj *)A->data;
4496a09307cSBarry Smith   PetscInt    rStart, rEnd, cStart, cEnd;
4506a09307cSBarry Smith 
4516a09307cSBarry Smith   PetscFunctionBegin;
4526a09307cSBarry Smith   PetscCheck(!adj->i, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Matrix is already assembled, cannot change its values");
4536a09307cSBarry Smith   PetscCall(MatGetOwnershipRange(A, &rStart, &rEnd));
4546a09307cSBarry Smith   PetscCall(MatGetOwnershipRangeColumn(A, &cStart, &cEnd));
4556a09307cSBarry Smith   if (!adj->ht) {
4566a09307cSBarry Smith     PetscCall(PetscHSetIJCreate(&adj->ht));
4576a09307cSBarry Smith     PetscCall(MatStashCreate_Private(PetscObjectComm((PetscObject)A), 1, &A->stash));
4586a09307cSBarry Smith     PetscCall(PetscLayoutSetUp(A->rmap));
4596a09307cSBarry Smith     PetscCall(PetscLayoutSetUp(A->cmap));
4606a09307cSBarry Smith   }
4616a09307cSBarry Smith   for (PetscInt r = 0; r < m; ++r) {
4626a09307cSBarry Smith     PetscHashIJKey key;
4636a09307cSBarry Smith 
4646a09307cSBarry Smith     key.i = rows[r];
4656a09307cSBarry Smith     if (key.i < 0) continue;
4666a09307cSBarry Smith     if ((key.i < rStart) || (key.i >= rEnd)) {
4676a09307cSBarry Smith       PetscCall(MatStashValuesRow_Private(&A->stash, key.i, n, cols, values, PETSC_FALSE));
4686a09307cSBarry Smith     } else {
4696a09307cSBarry Smith       for (PetscInt c = 0; c < n; ++c) {
4706a09307cSBarry Smith         key.j = cols[c];
4716a09307cSBarry Smith         if (key.j < 0 || key.i == key.j) continue;
4726a09307cSBarry Smith         PetscCall(PetscHSetIJAdd(adj->ht, key));
4736a09307cSBarry Smith       }
4746a09307cSBarry Smith     }
4756a09307cSBarry Smith   }
4766a09307cSBarry Smith   PetscFunctionReturn(0);
4776a09307cSBarry Smith }
4786a09307cSBarry Smith 
4799371c9d4SSatish Balay PetscErrorCode MatAssemblyBegin_MPIAdj(Mat A, MatAssemblyType type) {
4806a09307cSBarry Smith   PetscInt    nstash, reallocs;
4816a09307cSBarry Smith   Mat_MPIAdj *adj = (Mat_MPIAdj *)A->data;
4826a09307cSBarry Smith 
4836a09307cSBarry Smith   PetscFunctionBegin;
4846a09307cSBarry Smith   if (!adj->ht) {
4856a09307cSBarry Smith     PetscCall(PetscHSetIJCreate(&adj->ht));
4866a09307cSBarry Smith     PetscCall(PetscLayoutSetUp(A->rmap));
4876a09307cSBarry Smith     PetscCall(PetscLayoutSetUp(A->cmap));
4886a09307cSBarry Smith   }
4896a09307cSBarry Smith   PetscCall(MatStashScatterBegin_Private(A, &A->stash, A->rmap->range));
4906a09307cSBarry Smith   PetscCall(MatStashGetInfo_Private(&A->stash, &nstash, &reallocs));
4916a09307cSBarry Smith   PetscCall(PetscInfo(A, "Stash has %" PetscInt_FMT " entries, uses %" PetscInt_FMT " mallocs.\n", nstash, reallocs));
4926a09307cSBarry Smith   PetscFunctionReturn(0);
4936a09307cSBarry Smith }
4946a09307cSBarry Smith 
4959371c9d4SSatish Balay PetscErrorCode MatAssemblyEnd_MPIAdj(Mat A, MatAssemblyType type) {
4966a09307cSBarry Smith   PetscScalar   *val;
4976a09307cSBarry Smith   PetscInt      *row, *col, m, rstart, *rowstarts;
4986a09307cSBarry Smith   PetscInt       i, j, ncols, flg, nz;
4996a09307cSBarry Smith   PetscMPIInt    n;
5006a09307cSBarry Smith   Mat_MPIAdj    *adj = (Mat_MPIAdj *)A->data;
5016a09307cSBarry Smith   PetscHashIter  hi;
5026a09307cSBarry Smith   PetscHashIJKey key;
5036a09307cSBarry Smith   PetscHSetIJ    ht = adj->ht;
5046a09307cSBarry Smith 
5056a09307cSBarry Smith   PetscFunctionBegin;
5066a09307cSBarry Smith   while (1) {
5076a09307cSBarry Smith     PetscCall(MatStashScatterGetMesg_Private(&A->stash, &n, &row, &col, &val, &flg));
5086a09307cSBarry Smith     if (!flg) break;
5096a09307cSBarry Smith 
5106a09307cSBarry Smith     for (i = 0; i < n;) {
5116a09307cSBarry Smith       /* Identify the consecutive vals belonging to the same row */
5126a09307cSBarry Smith       for (j = i, rstart = row[j]; j < n; j++) {
5136a09307cSBarry Smith         if (row[j] != rstart) break;
5146a09307cSBarry Smith       }
5156a09307cSBarry Smith       if (j < n) ncols = j - i;
5166a09307cSBarry Smith       else ncols = n - i;
5176a09307cSBarry Smith       /* Set all these values with a single function call */
5186a09307cSBarry Smith       PetscCall(MatSetValues_MPIAdj(A, 1, row + i, ncols, col + i, val + i, INSERT_VALUES));
5196a09307cSBarry Smith       i = j;
5206a09307cSBarry Smith     }
5216a09307cSBarry Smith   }
5226a09307cSBarry Smith   PetscCall(MatStashScatterEnd_Private(&A->stash));
5236a09307cSBarry Smith   PetscCall(MatStashDestroy_Private(&A->stash));
5246a09307cSBarry Smith 
5256a09307cSBarry Smith   PetscCall(MatGetLocalSize(A, &m, NULL));
5266a09307cSBarry Smith   PetscCall(MatGetOwnershipRange(A, &rstart, NULL));
5276a09307cSBarry Smith   PetscCall(PetscCalloc1(m + 1, &rowstarts));
5286a09307cSBarry Smith   PetscHashIterBegin(ht, hi);
5296a09307cSBarry Smith   for (; !PetscHashIterAtEnd(ht, hi);) {
5306a09307cSBarry Smith     PetscHashIterGetKey(ht, hi, key);
5316a09307cSBarry Smith     rowstarts[key.i - rstart + 1]++;
5326a09307cSBarry Smith     PetscHashIterNext(ht, hi);
5336a09307cSBarry Smith   }
534ad540459SPierre Jolivet   for (i = 1; i < m + 1; i++) rowstarts[i] = rowstarts[i - 1] + rowstarts[i];
5356a09307cSBarry Smith 
5366a09307cSBarry Smith   PetscCall(PetscHSetIJGetSize(ht, &nz));
5376a09307cSBarry Smith   PetscCall(PetscMalloc1(nz, &col));
5386a09307cSBarry Smith   PetscHashIterBegin(ht, hi);
5396a09307cSBarry Smith   for (; !PetscHashIterAtEnd(ht, hi);) {
5406a09307cSBarry Smith     PetscHashIterGetKey(ht, hi, key);
5416a09307cSBarry Smith     col[rowstarts[key.i - rstart]++] = key.j;
5426a09307cSBarry Smith     PetscHashIterNext(ht, hi);
5436a09307cSBarry Smith   }
5446a09307cSBarry Smith   PetscCall(PetscHSetIJDestroy(&ht));
545ad540459SPierre Jolivet   for (i = m; i > 0; i--) rowstarts[i] = rowstarts[i - 1];
5466a09307cSBarry Smith   rowstarts[0] = 0;
5476a09307cSBarry Smith 
54848a46eb9SPierre Jolivet   for (PetscInt i = 0; i < m; i++) PetscCall(PetscSortInt(rowstarts[i + 1] - rowstarts[i], &col[rowstarts[i]]));
5496a09307cSBarry Smith 
5506a09307cSBarry Smith   adj->i       = rowstarts;
5516a09307cSBarry Smith   adj->j       = col;
552252a1336SBarry Smith   adj->nz      = rowstarts[m];
5536a09307cSBarry Smith   adj->freeaij = PETSC_TRUE;
5546a09307cSBarry Smith   PetscFunctionReturn(0);
5556a09307cSBarry Smith }
5566a09307cSBarry Smith 
557b97cf49bSBarry Smith /* -------------------------------------------------------------------*/
5586a09307cSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_MPIAdj,
5593eda8832SBarry Smith                                        MatGetRow_MPIAdj,
5603eda8832SBarry Smith                                        MatRestoreRow_MPIAdj,
561f4259b30SLisandro Dalcin                                        NULL,
562f4259b30SLisandro Dalcin                                        /* 4*/ NULL,
563f4259b30SLisandro Dalcin                                        NULL,
564f4259b30SLisandro Dalcin                                        NULL,
565f4259b30SLisandro Dalcin                                        NULL,
566f4259b30SLisandro Dalcin                                        NULL,
567f4259b30SLisandro Dalcin                                        NULL,
568f4259b30SLisandro Dalcin                                        /*10*/ NULL,
569f4259b30SLisandro Dalcin                                        NULL,
570f4259b30SLisandro Dalcin                                        NULL,
571f4259b30SLisandro Dalcin                                        NULL,
572f4259b30SLisandro Dalcin                                        NULL,
573f4259b30SLisandro Dalcin                                        /*15*/ NULL,
5743eda8832SBarry Smith                                        MatEqual_MPIAdj,
575f4259b30SLisandro Dalcin                                        NULL,
576f4259b30SLisandro Dalcin                                        NULL,
577f4259b30SLisandro Dalcin                                        NULL,
5786a09307cSBarry Smith                                        /*20*/ MatAssemblyBegin_MPIAdj,
5796a09307cSBarry Smith                                        MatAssemblyEnd_MPIAdj,
5803eda8832SBarry Smith                                        MatSetOption_MPIAdj,
581f4259b30SLisandro Dalcin                                        NULL,
582f4259b30SLisandro Dalcin                                        /*24*/ NULL,
583f4259b30SLisandro Dalcin                                        NULL,
584f4259b30SLisandro Dalcin                                        NULL,
585f4259b30SLisandro Dalcin                                        NULL,
586f4259b30SLisandro Dalcin                                        NULL,
587f4259b30SLisandro Dalcin                                        /*29*/ NULL,
588f4259b30SLisandro Dalcin                                        NULL,
589f4259b30SLisandro Dalcin                                        NULL,
590f4259b30SLisandro Dalcin                                        NULL,
591f4259b30SLisandro Dalcin                                        NULL,
592f4259b30SLisandro Dalcin                                        /*34*/ NULL,
593f4259b30SLisandro Dalcin                                        NULL,
594f4259b30SLisandro Dalcin                                        NULL,
595f4259b30SLisandro Dalcin                                        NULL,
596f4259b30SLisandro Dalcin                                        NULL,
597f4259b30SLisandro Dalcin                                        /*39*/ NULL,
5987dae84e0SHong Zhang                                        MatCreateSubMatrices_MPIAdj,
599f4259b30SLisandro Dalcin                                        NULL,
600f4259b30SLisandro Dalcin                                        NULL,
601f4259b30SLisandro Dalcin                                        NULL,
602f4259b30SLisandro Dalcin                                        /*44*/ NULL,
603f4259b30SLisandro Dalcin                                        NULL,
6047d68702bSBarry Smith                                        MatShift_Basic,
605f4259b30SLisandro Dalcin                                        NULL,
606f4259b30SLisandro Dalcin                                        NULL,
607f4259b30SLisandro Dalcin                                        /*49*/ NULL,
6086cd91f64SBarry Smith                                        MatGetRowIJ_MPIAdj,
609d42735a1SBarry Smith                                        MatRestoreRowIJ_MPIAdj,
610f4259b30SLisandro Dalcin                                        NULL,
611f4259b30SLisandro Dalcin                                        NULL,
612f4259b30SLisandro Dalcin                                        /*54*/ NULL,
613f4259b30SLisandro Dalcin                                        NULL,
614f4259b30SLisandro Dalcin                                        NULL,
615f4259b30SLisandro Dalcin                                        NULL,
616f4259b30SLisandro Dalcin                                        NULL,
617f4259b30SLisandro Dalcin                                        /*59*/ NULL,
618b9b97703SBarry Smith                                        MatDestroy_MPIAdj,
619b9b97703SBarry Smith                                        MatView_MPIAdj,
62017667f90SBarry Smith                                        MatConvertFrom_MPIAdj,
621f4259b30SLisandro Dalcin                                        NULL,
622f4259b30SLisandro Dalcin                                        /*64*/ NULL,
623f4259b30SLisandro Dalcin                                        NULL,
624f4259b30SLisandro Dalcin                                        NULL,
625f4259b30SLisandro Dalcin                                        NULL,
626f4259b30SLisandro Dalcin                                        NULL,
627f4259b30SLisandro Dalcin                                        /*69*/ NULL,
628f4259b30SLisandro Dalcin                                        NULL,
629f4259b30SLisandro Dalcin                                        NULL,
630f4259b30SLisandro Dalcin                                        NULL,
631f4259b30SLisandro Dalcin                                        NULL,
632f4259b30SLisandro Dalcin                                        /*74*/ NULL,
633f4259b30SLisandro Dalcin                                        NULL,
634f4259b30SLisandro Dalcin                                        NULL,
635f4259b30SLisandro Dalcin                                        NULL,
636f4259b30SLisandro Dalcin                                        NULL,
637f4259b30SLisandro Dalcin                                        /*79*/ NULL,
638f4259b30SLisandro Dalcin                                        NULL,
639f4259b30SLisandro Dalcin                                        NULL,
640f4259b30SLisandro Dalcin                                        NULL,
641f4259b30SLisandro Dalcin                                        NULL,
642f4259b30SLisandro Dalcin                                        /*84*/ NULL,
643f4259b30SLisandro Dalcin                                        NULL,
644f4259b30SLisandro Dalcin                                        NULL,
645f4259b30SLisandro Dalcin                                        NULL,
646f4259b30SLisandro Dalcin                                        NULL,
647f4259b30SLisandro Dalcin                                        /*89*/ NULL,
648f4259b30SLisandro Dalcin                                        NULL,
649f4259b30SLisandro Dalcin                                        NULL,
650f4259b30SLisandro Dalcin                                        NULL,
651f4259b30SLisandro Dalcin                                        NULL,
652f4259b30SLisandro Dalcin                                        /*94*/ NULL,
653f4259b30SLisandro Dalcin                                        NULL,
654f4259b30SLisandro Dalcin                                        NULL,
655f4259b30SLisandro Dalcin                                        NULL,
656f4259b30SLisandro Dalcin                                        NULL,
657f4259b30SLisandro Dalcin                                        /*99*/ NULL,
658f4259b30SLisandro Dalcin                                        NULL,
659f4259b30SLisandro Dalcin                                        NULL,
660f4259b30SLisandro Dalcin                                        NULL,
661f4259b30SLisandro Dalcin                                        NULL,
662f4259b30SLisandro Dalcin                                        /*104*/ NULL,
663f4259b30SLisandro Dalcin                                        NULL,
664f4259b30SLisandro Dalcin                                        NULL,
665f4259b30SLisandro Dalcin                                        NULL,
666f4259b30SLisandro Dalcin                                        NULL,
667f4259b30SLisandro Dalcin                                        /*109*/ NULL,
668f4259b30SLisandro Dalcin                                        NULL,
669f4259b30SLisandro Dalcin                                        NULL,
670f4259b30SLisandro Dalcin                                        NULL,
671f4259b30SLisandro Dalcin                                        NULL,
672f4259b30SLisandro Dalcin                                        /*114*/ NULL,
673f4259b30SLisandro Dalcin                                        NULL,
674f4259b30SLisandro Dalcin                                        NULL,
675f4259b30SLisandro Dalcin                                        NULL,
676f4259b30SLisandro Dalcin                                        NULL,
677f4259b30SLisandro Dalcin                                        /*119*/ NULL,
678f4259b30SLisandro Dalcin                                        NULL,
679f4259b30SLisandro Dalcin                                        NULL,
680f4259b30SLisandro Dalcin                                        NULL,
681f4259b30SLisandro Dalcin                                        NULL,
682f4259b30SLisandro Dalcin                                        /*124*/ NULL,
683f4259b30SLisandro Dalcin                                        NULL,
684f4259b30SLisandro Dalcin                                        NULL,
685f4259b30SLisandro Dalcin                                        NULL,
6867dae84e0SHong Zhang                                        MatCreateSubMatricesMPI_MPIAdj,
687f4259b30SLisandro Dalcin                                        /*129*/ NULL,
688f4259b30SLisandro Dalcin                                        NULL,
689f4259b30SLisandro Dalcin                                        NULL,
690f4259b30SLisandro Dalcin                                        NULL,
691f4259b30SLisandro Dalcin                                        NULL,
692f4259b30SLisandro Dalcin                                        /*134*/ NULL,
693f4259b30SLisandro Dalcin                                        NULL,
694f4259b30SLisandro Dalcin                                        NULL,
695f4259b30SLisandro Dalcin                                        NULL,
696f4259b30SLisandro Dalcin                                        NULL,
697f4259b30SLisandro Dalcin                                        /*139*/ NULL,
698f4259b30SLisandro Dalcin                                        NULL,
699d70f29a3SPierre Jolivet                                        NULL,
700d70f29a3SPierre Jolivet                                        NULL,
701d70f29a3SPierre Jolivet                                        NULL,
702d70f29a3SPierre Jolivet                                        /*144*/ NULL,
703d70f29a3SPierre Jolivet                                        NULL,
704d70f29a3SPierre Jolivet                                        NULL,
70599a7f59eSMark Adams                                        NULL,
70699a7f59eSMark Adams                                        NULL,
7077fb60732SBarry Smith                                        NULL,
7089371c9d4SSatish Balay                                        /*150*/ NULL};
709b97cf49bSBarry Smith 
7109371c9d4SSatish Balay static PetscErrorCode MatMPIAdjSetPreallocation_MPIAdj(Mat B, PetscInt *i, PetscInt *j, PetscInt *values) {
711a23d5eceSKris Buschelman   Mat_MPIAdj *b = (Mat_MPIAdj *)B->data;
712e895ccc0SFande Kong   PetscBool   useedgeweights;
713a23d5eceSKris Buschelman 
714a23d5eceSKris Buschelman   PetscFunctionBegin;
7159566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(B->rmap));
7169566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(B->cmap));
7179371c9d4SSatish Balay   if (values) useedgeweights = PETSC_TRUE;
7189371c9d4SSatish Balay   else useedgeweights = PETSC_FALSE;
719e895ccc0SFande Kong   /* Make everybody knows if they are using edge weights or not */
7201c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce((int *)&useedgeweights, (int *)&b->useedgeweights, 1, MPI_INT, MPI_MAX, PetscObjectComm((PetscObject)B)));
72158ec18a5SLisandro Dalcin 
72276bd3646SJed Brown   if (PetscDefined(USE_DEBUG)) {
72376bd3646SJed Brown     PetscInt ii;
72476bd3646SJed Brown 
72508401ef6SPierre 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]);
726d0f46423SBarry Smith     for (ii = 1; ii < B->rmap->n; ii++) {
727aed4548fSBarry 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]);
728a23d5eceSKris Buschelman     }
729ad540459SPierre 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]);
73076bd3646SJed Brown   }
731a23d5eceSKris Buschelman   b->j      = j;
732a23d5eceSKris Buschelman   b->i      = i;
733a23d5eceSKris Buschelman   b->values = values;
734a23d5eceSKris Buschelman 
735d0f46423SBarry Smith   b->nz        = i[B->rmap->n];
736f4259b30SLisandro Dalcin   b->diag      = NULL;
737a23d5eceSKris Buschelman   b->symmetric = PETSC_FALSE;
738a23d5eceSKris Buschelman   b->freeaij   = PETSC_TRUE;
739a23d5eceSKris Buschelman 
7406a09307cSBarry Smith   B->ops->assemblybegin = NULL;
7416a09307cSBarry Smith   B->ops->assemblyend   = NULL;
7429566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
7439566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
7446a09307cSBarry Smith   PetscCall(MatStashDestroy_Private(&B->stash));
745a23d5eceSKris Buschelman   PetscFunctionReturn(0);
746a23d5eceSKris Buschelman }
747b97cf49bSBarry Smith 
7489371c9d4SSatish Balay static PetscErrorCode MatMPIAdjCreateNonemptySubcommMat_MPIAdj(Mat A, Mat *B) {
7499a3665c6SJed Brown   Mat_MPIAdj     *a = (Mat_MPIAdj *)A->data;
7509a3665c6SJed Brown   const PetscInt *ranges;
7519a3665c6SJed Brown   MPI_Comm        acomm, bcomm;
7529a3665c6SJed Brown   MPI_Group       agroup, bgroup;
7539a3665c6SJed Brown   PetscMPIInt     i, rank, size, nranks, *ranks;
7549a3665c6SJed Brown 
7559a3665c6SJed Brown   PetscFunctionBegin;
7560298fd71SBarry Smith   *B = NULL;
7579566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)A, &acomm));
7589566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(acomm, &size));
7599566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(acomm, &rank));
7609566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRanges(A, &ranges));
7619a3665c6SJed Brown   for (i = 0, nranks = 0; i < size; i++) {
7629a3665c6SJed Brown     if (ranges[i + 1] - ranges[i] > 0) nranks++;
7639a3665c6SJed Brown   }
7649a3665c6SJed Brown   if (nranks == size) { /* All ranks have a positive number of rows, so we do not need to create a subcomm; */
7659566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)A));
7669a3665c6SJed Brown     *B = A;
7679a3665c6SJed Brown     PetscFunctionReturn(0);
7689a3665c6SJed Brown   }
7699a3665c6SJed Brown 
7709566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nranks, &ranks));
7719a3665c6SJed Brown   for (i = 0, nranks = 0; i < size; i++) {
7729a3665c6SJed Brown     if (ranges[i + 1] - ranges[i] > 0) ranks[nranks++] = i;
7739a3665c6SJed Brown   }
7749566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_group(acomm, &agroup));
7759566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Group_incl(agroup, nranks, ranks, &bgroup));
7769566063dSJacob Faibussowitsch   PetscCall(PetscFree(ranks));
7779566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_create(acomm, bgroup, &bcomm));
7789566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Group_free(&agroup));
7799566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Group_free(&bgroup));
7809a3665c6SJed Brown   if (bcomm != MPI_COMM_NULL) {
7819a3665c6SJed Brown     PetscInt    m, N;
7829a3665c6SJed Brown     Mat_MPIAdj *b;
7839566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(A, &m, NULL));
7849566063dSJacob Faibussowitsch     PetscCall(MatGetSize(A, NULL, &N));
7859566063dSJacob Faibussowitsch     PetscCall(MatCreateMPIAdj(bcomm, m, N, a->i, a->j, a->values, B));
7869a3665c6SJed Brown     b          = (Mat_MPIAdj *)(*B)->data;
7879a3665c6SJed Brown     b->freeaij = PETSC_FALSE;
7889566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_free(&bcomm));
7899a3665c6SJed Brown   }
7909a3665c6SJed Brown   PetscFunctionReturn(0);
7919a3665c6SJed Brown }
7929a3665c6SJed Brown 
7939371c9d4SSatish Balay PetscErrorCode MatMPIAdjToSeq_MPIAdj(Mat A, Mat *B) {
794b44e7856SBarry Smith   PetscInt    M, N, *II, *J, NZ, nz, m, nzstart, i;
795b44e7856SBarry Smith   PetscInt   *Values = NULL;
796b44e7856SBarry Smith   Mat_MPIAdj *adj    = (Mat_MPIAdj *)A->data;
797b44e7856SBarry Smith   PetscMPIInt mnz, mm, *allnz, *allm, size, *dispnz, *dispm;
798b44e7856SBarry Smith 
799b44e7856SBarry Smith   PetscFunctionBegin;
8009566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
8019566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A, &M, &N));
8029566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A, &m, NULL));
803b44e7856SBarry Smith   nz = adj->nz;
80408401ef6SPierre 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]);
8059566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allreduce(&nz, &NZ, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)A)));
806d4e69552SBarry Smith 
8079566063dSJacob Faibussowitsch   PetscCall(PetscMPIIntCast(nz, &mnz));
8089566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(size, &allnz, size, &dispnz));
8099566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allgather(&mnz, 1, MPI_INT, allnz, 1, MPI_INT, PetscObjectComm((PetscObject)A)));
8109371c9d4SSatish Balay   dispnz[0] = 0;
8119371c9d4SSatish Balay   for (i = 1; i < size; i++) dispnz[i] = dispnz[i - 1] + allnz[i - 1];
812b44e7856SBarry Smith   if (adj->values) {
8139566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(NZ, &Values));
8149566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Allgatherv(adj->values, mnz, MPIU_INT, Values, allnz, dispnz, MPIU_INT, PetscObjectComm((PetscObject)A)));
815b44e7856SBarry Smith   }
8169566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(NZ, &J));
8179566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allgatherv(adj->j, mnz, MPIU_INT, J, allnz, dispnz, MPIU_INT, PetscObjectComm((PetscObject)A)));
8189566063dSJacob Faibussowitsch   PetscCall(PetscFree2(allnz, dispnz));
8199566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Scan(&nz, &nzstart, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)A)));
820b44e7856SBarry Smith   nzstart -= nz;
821b44e7856SBarry Smith   /* shift the i[] values so they will be correct after being received */
822b44e7856SBarry Smith   for (i = 0; i < m; i++) adj->i[i] += nzstart;
8239566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(M + 1, &II));
8249566063dSJacob Faibussowitsch   PetscCall(PetscMPIIntCast(m, &mm));
8259566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(size, &allm, size, &dispm));
8269566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allgather(&mm, 1, MPI_INT, allm, 1, MPI_INT, PetscObjectComm((PetscObject)A)));
8279371c9d4SSatish Balay   dispm[0] = 0;
8289371c9d4SSatish Balay   for (i = 1; i < size; i++) dispm[i] = dispm[i - 1] + allm[i - 1];
8299566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allgatherv(adj->i, mm, MPIU_INT, II, allm, dispm, MPIU_INT, PetscObjectComm((PetscObject)A)));
8309566063dSJacob Faibussowitsch   PetscCall(PetscFree2(allm, dispm));
831b44e7856SBarry Smith   II[M] = NZ;
832b44e7856SBarry Smith   /* shift the i[] values back */
833b44e7856SBarry Smith   for (i = 0; i < m; i++) adj->i[i] -= nzstart;
8349566063dSJacob Faibussowitsch   PetscCall(MatCreateMPIAdj(PETSC_COMM_SELF, M, N, II, J, Values, B));
835b44e7856SBarry Smith   PetscFunctionReturn(0);
836b44e7856SBarry Smith }
837b44e7856SBarry Smith 
8389371c9d4SSatish Balay PetscErrorCode MatMPIAdjToSeqRankZero_MPIAdj(Mat A, Mat *B) {
839252a1336SBarry Smith   PetscInt    M, N, *II, *J, NZ, nz, m, nzstart, i;
840252a1336SBarry Smith   PetscInt   *Values = NULL;
841252a1336SBarry Smith   Mat_MPIAdj *adj    = (Mat_MPIAdj *)A->data;
842252a1336SBarry Smith   PetscMPIInt mnz, mm, *allnz = NULL, *allm, size, *dispnz, *dispm, rank;
843252a1336SBarry Smith 
844252a1336SBarry Smith   PetscFunctionBegin;
845252a1336SBarry Smith   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
846252a1336SBarry Smith   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)A), &rank));
847252a1336SBarry Smith   PetscCall(MatGetSize(A, &M, &N));
848252a1336SBarry Smith   PetscCall(MatGetLocalSize(A, &m, NULL));
849252a1336SBarry Smith   nz = adj->nz;
850252a1336SBarry 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]);
851252a1336SBarry Smith   PetscCallMPI(MPI_Allreduce(&nz, &NZ, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)A)));
852252a1336SBarry Smith 
853252a1336SBarry Smith   PetscCall(PetscMPIIntCast(nz, &mnz));
854252a1336SBarry Smith   if (!rank) PetscCall(PetscMalloc2(size, &allnz, size, &dispnz));
855252a1336SBarry Smith   PetscCallMPI(MPI_Gather(&mnz, 1, MPI_INT, allnz, 1, MPI_INT, 0, PetscObjectComm((PetscObject)A)));
856252a1336SBarry Smith   if (!rank) {
8579371c9d4SSatish Balay     dispnz[0] = 0;
8589371c9d4SSatish Balay     for (i = 1; i < size; i++) dispnz[i] = dispnz[i - 1] + allnz[i - 1];
859252a1336SBarry Smith     if (adj->values) {
860252a1336SBarry Smith       PetscCall(PetscMalloc1(NZ, &Values));
861252a1336SBarry Smith       PetscCallMPI(MPI_Gatherv(adj->values, mnz, MPIU_INT, Values, allnz, dispnz, MPIU_INT, 0, PetscObjectComm((PetscObject)A)));
862252a1336SBarry Smith     }
863252a1336SBarry Smith     PetscCall(PetscMalloc1(NZ, &J));
864252a1336SBarry Smith     PetscCallMPI(MPI_Gatherv(adj->j, mnz, MPIU_INT, J, allnz, dispnz, MPIU_INT, 0, PetscObjectComm((PetscObject)A)));
865252a1336SBarry Smith     PetscCall(PetscFree2(allnz, dispnz));
866252a1336SBarry Smith   } else {
867252a1336SBarry Smith     if (adj->values) PetscCallMPI(MPI_Gatherv(adj->values, mnz, MPIU_INT, NULL, NULL, NULL, MPIU_INT, 0, PetscObjectComm((PetscObject)A)));
868252a1336SBarry Smith     PetscCallMPI(MPI_Gatherv(adj->j, mnz, MPIU_INT, NULL, NULL, NULL, MPIU_INT, 0, PetscObjectComm((PetscObject)A)));
869252a1336SBarry Smith   }
870252a1336SBarry Smith   PetscCallMPI(MPI_Scan(&nz, &nzstart, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)A)));
871252a1336SBarry Smith   nzstart -= nz;
872252a1336SBarry Smith   /* shift the i[] values so they will be correct after being received */
873252a1336SBarry Smith   for (i = 0; i < m; i++) adj->i[i] += nzstart;
874252a1336SBarry Smith   PetscCall(PetscMPIIntCast(m, &mm));
875252a1336SBarry Smith   if (!rank) {
876252a1336SBarry Smith     PetscCall(PetscMalloc1(M + 1, &II));
877252a1336SBarry Smith     PetscCall(PetscMalloc2(size, &allm, size, &dispm));
878252a1336SBarry Smith     PetscCallMPI(MPI_Gather(&mm, 1, MPI_INT, allm, 1, MPI_INT, 0, PetscObjectComm((PetscObject)A)));
8799371c9d4SSatish Balay     dispm[0] = 0;
8809371c9d4SSatish Balay     for (i = 1; i < size; i++) dispm[i] = dispm[i - 1] + allm[i - 1];
881252a1336SBarry Smith     PetscCallMPI(MPI_Gatherv(adj->i, mm, MPIU_INT, II, allm, dispm, MPIU_INT, 0, PetscObjectComm((PetscObject)A)));
882252a1336SBarry Smith     PetscCall(PetscFree2(allm, dispm));
883252a1336SBarry Smith     II[M] = NZ;
884252a1336SBarry Smith   } else {
885252a1336SBarry Smith     PetscCallMPI(MPI_Gather(&mm, 1, MPI_INT, NULL, 1, MPI_INT, 0, PetscObjectComm((PetscObject)A)));
886252a1336SBarry Smith     PetscCallMPI(MPI_Gatherv(adj->i, mm, MPIU_INT, NULL, NULL, NULL, MPIU_INT, 0, PetscObjectComm((PetscObject)A)));
887252a1336SBarry Smith   }
888252a1336SBarry Smith   /* shift the i[] values back */
889252a1336SBarry Smith   for (i = 0; i < m; i++) adj->i[i] -= nzstart;
890252a1336SBarry Smith   if (!rank) PetscCall(MatCreateMPIAdj(PETSC_COMM_SELF, M, N, II, J, Values, B));
891252a1336SBarry Smith   PetscFunctionReturn(0);
892252a1336SBarry Smith }
893252a1336SBarry Smith 
8949a3665c6SJed Brown /*@
89511a5261eSBarry Smith    MatMPIAdjCreateNonemptySubcommMat - create the same `MATMPIADJ` matrix on a subcommunicator containing only processes owning a positive number of rows
8969a3665c6SJed Brown 
8979a3665c6SJed Brown    Collective
8989a3665c6SJed Brown 
8994165533cSJose E. Roman    Input Parameter:
9009a3665c6SJed Brown .  A - original MPIAdj matrix
9019a3665c6SJed Brown 
9024165533cSJose E. Roman    Output Parameter:
9030298fd71SBarry Smith .  B - matrix on subcommunicator, NULL on ranks that owned zero rows of A
9049a3665c6SJed Brown 
9059a3665c6SJed Brown    Level: developer
9069a3665c6SJed Brown 
9079a3665c6SJed Brown    Note:
9089a3665c6SJed Brown    This function is mostly useful for internal use by mesh partitioning packages that require that every process owns at least one row.
9099a3665c6SJed Brown 
91011a5261eSBarry Smith    The matrix B should be destroyed with `MatDestroy()`. The arrays are not copied, so B should be destroyed before A is destroyed.
9119a3665c6SJed Brown 
91211a5261eSBarry Smith .seealso: `MATMPIADJ`, `MatCreateMPIAdj()`
9139a3665c6SJed Brown @*/
9149371c9d4SSatish Balay PetscErrorCode MatMPIAdjCreateNonemptySubcommMat(Mat A, Mat *B) {
9159a3665c6SJed Brown   PetscFunctionBegin;
9169a3665c6SJed Brown   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
917cac4c232SBarry Smith   PetscUseMethod(A, "MatMPIAdjCreateNonemptySubcommMat_C", (Mat, Mat *), (A, B));
9189a3665c6SJed Brown   PetscFunctionReturn(0);
9199a3665c6SJed Brown }
9209a3665c6SJed Brown 
9210bad9183SKris Buschelman /*MC
922fafad747SKris Buschelman    MATMPIADJ - MATMPIADJ = "mpiadj" - A matrix type to be used for distributed adjacency matrices,
9230bad9183SKris Buschelman    intended for use constructing orderings and partitionings.
9240bad9183SKris Buschelman 
9250bad9183SKris Buschelman   Level: beginner
9260bad9183SKris Buschelman 
92711a5261eSBarry Smith   Note:
92811a5261eSBarry Smith     You can provide values to the matrix using `MatMPIAdjSetPreallocation()`, `MatCreateMPIAdj()`, or
92911a5261eSBarry Smith     by calling `MatSetValues()` and `MatAssemblyBegin()` followed by `MatAssemblyEnd()`
9306a09307cSBarry Smith 
9316a09307cSBarry Smith .seealso: `MatCreateMPIAdj()`, `MatMPIAdjSetPreallocation()`, `MatSetValues()`
9320bad9183SKris Buschelman M*/
9330bad9183SKris Buschelman 
9349371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode MatCreate_MPIAdj(Mat B) {
935273d9f13SBarry Smith   Mat_MPIAdj *b;
936273d9f13SBarry Smith 
937273d9f13SBarry Smith   PetscFunctionBegin;
938*4dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&b));
939b0a32e0cSBarry Smith   B->data = (void *)b;
9409566063dSJacob Faibussowitsch   PetscCall(PetscMemcpy(B->ops, &MatOps_Values, sizeof(struct _MatOps)));
941273d9f13SBarry Smith   B->assembled    = PETSC_FALSE;
9426a09307cSBarry Smith   B->preallocated = PETSC_TRUE; /* so that MatSetValues() may be used */
943273d9f13SBarry Smith 
9449566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMPIAdjSetPreallocation_C", MatMPIAdjSetPreallocation_MPIAdj));
9459566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMPIAdjCreateNonemptySubcommMat_C", MatMPIAdjCreateNonemptySubcommMat_MPIAdj));
9469566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMPIAdjToSeq_C", MatMPIAdjToSeq_MPIAdj));
947252a1336SBarry Smith   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMPIAdjToSeqRankZero_C", MatMPIAdjToSeqRankZero_MPIAdj));
9489566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)B, MATMPIADJ));
949273d9f13SBarry Smith   PetscFunctionReturn(0);
950273d9f13SBarry Smith }
951273d9f13SBarry Smith 
952a23d5eceSKris Buschelman /*@C
95311a5261eSBarry Smith    MatMPIAdjToSeq - Converts an parallel `MATMPIADJ` matrix to complete `MATMPIADJ` on each process (needed by sequential partitioner)
954b44e7856SBarry Smith 
955d083f849SBarry Smith    Logically Collective
956b44e7856SBarry Smith 
957b44e7856SBarry Smith    Input Parameter:
958b44e7856SBarry Smith .  A - the matrix
959b44e7856SBarry Smith 
960b44e7856SBarry Smith    Output Parameter:
961b44e7856SBarry Smith .  B - the same matrix on all processes
962b44e7856SBarry Smith 
963b44e7856SBarry Smith    Level: intermediate
964b44e7856SBarry Smith 
96511a5261eSBarry Smith .seealso: `MATMPIADJ`, `MatCreate()`, `MatCreateMPIAdj()`, `MatSetValues()`, `MatMPIAdjToSeqRankZero()`
966b44e7856SBarry Smith @*/
9679371c9d4SSatish Balay PetscErrorCode MatMPIAdjToSeq(Mat A, Mat *B) {
968b44e7856SBarry Smith   PetscFunctionBegin;
969cac4c232SBarry Smith   PetscUseMethod(A, "MatMPIAdjToSeq_C", (Mat, Mat *), (A, B));
970b44e7856SBarry Smith   PetscFunctionReturn(0);
971b44e7856SBarry Smith }
972b44e7856SBarry Smith 
973b44e7856SBarry Smith /*@C
97411a5261eSBarry Smith    MatMPIAdjToSeqRankZero - Converts an parallel `MATMPIADJ` matrix to complete `MATMPIADJ` on rank zero (needed by sequential partitioner)
975252a1336SBarry Smith 
976252a1336SBarry Smith    Logically Collective
977252a1336SBarry Smith 
978252a1336SBarry Smith    Input Parameter:
979252a1336SBarry Smith .  A - the matrix
980252a1336SBarry Smith 
981252a1336SBarry Smith    Output Parameter:
982252a1336SBarry Smith .  B - the same matrix on rank zero, not set on other ranks
983252a1336SBarry Smith 
984252a1336SBarry Smith    Level: intermediate
985252a1336SBarry Smith 
98611a5261eSBarry Smith    Note:
987252a1336SBarry Smith      This routine has the advantage on systems with multiple ranks per node since only one copy of the matrix
988252a1336SBarry Smith      is stored on the first node, instead of the number of ranks copies. This can allow partitioning much larger
989252a1336SBarry Smith      paralllel graph sequentially.
990252a1336SBarry Smith 
99111a5261eSBarry Smith .seealso: `MATMPIADJ`, `MatCreate()`, `MatCreateMPIAdj()`, `MatSetValues()`, `MatMPIAdjToSeq()`
992252a1336SBarry Smith @*/
9939371c9d4SSatish Balay PetscErrorCode MatMPIAdjToSeqRankZero(Mat A, Mat *B) {
994252a1336SBarry Smith   PetscFunctionBegin;
995252a1336SBarry Smith   PetscUseMethod(A, "MatMPIAdjToSeqRankZero_C", (Mat, Mat *), (A, B));
996252a1336SBarry Smith   PetscFunctionReturn(0);
997252a1336SBarry Smith }
998252a1336SBarry Smith 
999252a1336SBarry Smith /*@C
1000a23d5eceSKris Buschelman    MatMPIAdjSetPreallocation - Sets the array used for storing the matrix elements
1001a23d5eceSKris Buschelman 
1002d083f849SBarry Smith    Logically Collective
1003a23d5eceSKris Buschelman 
1004a23d5eceSKris Buschelman    Input Parameters:
1005a23d5eceSKris Buschelman +  A - the matrix
1006a23d5eceSKris Buschelman .  i - the indices into j for the start of each row
1007a23d5eceSKris Buschelman .  j - the column indices for each row (sorted for each row).
1008a23d5eceSKris Buschelman        The indices in i and j start with zero (NOT with one).
1009a23d5eceSKris Buschelman -  values - [optional] edge weights
1010a23d5eceSKris Buschelman 
1011a23d5eceSKris Buschelman    Level: intermediate
1012a23d5eceSKris Buschelman 
10136a09307cSBarry Smith .seealso: `MatCreate()`, `MatCreateMPIAdj()`, `MatSetValues()`, `MATMPIADJ`, `MatCreateMPIAdj()`
1014a23d5eceSKris Buschelman @*/
10159371c9d4SSatish Balay PetscErrorCode MatMPIAdjSetPreallocation(Mat B, PetscInt *i, PetscInt *j, PetscInt *values) {
1016273d9f13SBarry Smith   PetscFunctionBegin;
1017cac4c232SBarry Smith   PetscTryMethod(B, "MatMPIAdjSetPreallocation_C", (Mat, PetscInt *, PetscInt *, PetscInt *), (B, i, j, values));
1018273d9f13SBarry Smith   PetscFunctionReturn(0);
1019273d9f13SBarry Smith }
1020273d9f13SBarry Smith 
1021b97cf49bSBarry Smith /*@C
10223eda8832SBarry Smith    MatCreateMPIAdj - Creates a sparse matrix representing an adjacency list.
1023b97cf49bSBarry Smith    The matrix does not have numerical values associated with it, but is
1024b97cf49bSBarry Smith    intended for ordering (to reduce bandwidth etc) and partitioning.
1025b97cf49bSBarry Smith 
1026d083f849SBarry Smith    Collective
1027ef5ee4d1SLois Curfman McInnes 
1028b97cf49bSBarry Smith    Input Parameters:
1029c2e958feSBarry Smith +  comm - MPI communicator
10300752156aSBarry Smith .  m - number of local rows
103158ec18a5SLisandro Dalcin .  N - number of global columns
1032b97cf49bSBarry Smith .  i - the indices into j for the start of each row
1033329f5518SBarry Smith .  j - the column indices for each row (sorted for each row).
1034ef5ee4d1SLois Curfman McInnes        The indices in i and j start with zero (NOT with one).
1035329f5518SBarry Smith -  values -[optional] edge weights
1036b97cf49bSBarry Smith 
1037b97cf49bSBarry Smith    Output Parameter:
1038b97cf49bSBarry Smith .  A - the matrix
1039b97cf49bSBarry Smith 
104015091d37SBarry Smith    Level: intermediate
104115091d37SBarry Smith 
104295452b02SPatrick Sanan    Notes:
1043c2e958feSBarry Smith    You must NOT free the ii, values and jj arrays yourself. PETSc will free them
104411a5261eSBarry Smith    when the matrix is destroyed; you must allocate them with `PetscMalloc()`. If you
104511a5261eSBarry Smith    call from Fortran you need not create the arrays with `PetscMalloc()`.
10466a09307cSBarry Smith 
10476a09307cSBarry Smith    You should not include the matrix diagonals.
1048b97cf49bSBarry Smith 
1049e3f7e1d6Svictorle    If you already have a matrix, you can create its adjacency matrix by a call
105011a5261eSBarry Smith    to `MatConvert()`, specifying a type of `MATMPIADJ`.
1051ededeb1aSvictorle 
105211a5261eSBarry Smith    Possible values for `MatSetOption()` - `MAT_STRUCTURALLY_SYMMETRIC`
1053b97cf49bSBarry Smith 
10546a09307cSBarry Smith .seealso: `MatCreate()`, `MatConvert()`, `MatGetOrdering()`, `MATMPIADJ`, `MatMPIAdjSetPreallocation()`
1055b97cf49bSBarry Smith @*/
10569371c9d4SSatish Balay PetscErrorCode MatCreateMPIAdj(MPI_Comm comm, PetscInt m, PetscInt N, PetscInt *i, PetscInt *j, PetscInt *values, Mat *A) {
1057433994e6SBarry Smith   PetscFunctionBegin;
10589566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, A));
10599566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*A, m, PETSC_DETERMINE, PETSC_DETERMINE, N));
10609566063dSJacob Faibussowitsch   PetscCall(MatSetType(*A, MATMPIADJ));
10619566063dSJacob Faibussowitsch   PetscCall(MatMPIAdjSetPreallocation(*A, i, j, values));
10623a40ed3dSBarry Smith   PetscFunctionReturn(0);
1063b97cf49bSBarry Smith }
1064