xref: /petsc/src/mat/impls/adj/mpi/mpiadj.c (revision d8a51d2aac522cef49bfbd84fd4cf473dc3fa6ee)
1b97cf49bSBarry Smith /*
2b97cf49bSBarry Smith     Defines the basic matrix operations for the ADJ adjacency list matrix data-structure.
3b97cf49bSBarry Smith */
4e37d9405SJed Brown #include <../src/mat/impls/adj/mpi/mpiadj.h> /*I "petscmat.h" I*/
5841d17a1SFande Kong #include <petscsf.h>
6b97cf49bSBarry Smith 
740244768SBarry Smith /*
82ef1f0ffSBarry Smith   The interface should be easy to use for both MatCreateSubMatrix (parallel sub-matrix) and MatCreateSubMatrices (sequential sub-matrices)
92ef1f0ffSBarry Smith */
10d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatCreateSubMatrix_MPIAdj_data(Mat adj, IS irows, IS icols, PetscInt **sadj_xadj, PetscInt **sadj_adjncy, PetscInt **sadj_values)
11d71ae5a4SJacob Faibussowitsch {
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));
1303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13140244768SBarry Smith }
13240244768SBarry Smith 
133d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatCreateSubMatrices_MPIAdj_Private(Mat mat, PetscInt n, const IS irow[], const IS icol[], PetscBool subcomm, MatReuse scall, Mat *submat[])
134d71ae5a4SJacob Faibussowitsch {
13540244768SBarry Smith   PetscInt        i, irow_n, icol_n, *sxadj, *sadjncy, *svalues;
13640244768SBarry Smith   PetscInt       *indices, nindx, j, k, loc;
13740244768SBarry Smith   PetscMPIInt     issame;
13840244768SBarry Smith   const PetscInt *irow_indices, *icol_indices;
13940244768SBarry Smith   MPI_Comm        scomm_row, scomm_col, scomm_mat;
14040244768SBarry Smith 
14140244768SBarry Smith   PetscFunctionBegin;
14240244768SBarry Smith   nindx = 0;
14340244768SBarry Smith   /*
14440244768SBarry Smith    * Estimate a maximum number for allocating memory
14540244768SBarry Smith    */
14640244768SBarry Smith   for (i = 0; i < n; i++) {
1479566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(irow[i], &irow_n));
1489566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(icol[i], &icol_n));
14940244768SBarry Smith     nindx = nindx > (irow_n + icol_n) ? nindx : (irow_n + icol_n);
15040244768SBarry Smith   }
1519566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nindx, &indices));
15240244768SBarry Smith   /* construct a submat */
153252a1336SBarry Smith   //  if (scall == MAT_INITIAL_MATRIX) PetscCall(PetscMalloc1(n,submat));
1546a7b62d2SBarry Smith 
15540244768SBarry Smith   for (i = 0; i < n; i++) {
15640244768SBarry Smith     if (subcomm) {
1579566063dSJacob Faibussowitsch       PetscCall(PetscObjectGetComm((PetscObject)irow[i], &scomm_row));
1589566063dSJacob Faibussowitsch       PetscCall(PetscObjectGetComm((PetscObject)icol[i], &scomm_col));
1599566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_compare(scomm_row, scomm_col, &issame));
16008401ef6SPierre 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");
1619566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_compare(scomm_row, PETSC_COMM_SELF, &issame));
16208401ef6SPierre 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");
16340244768SBarry Smith     } else {
16440244768SBarry Smith       scomm_row = PETSC_COMM_SELF;
16540244768SBarry Smith     }
16640244768SBarry Smith     /*get sub-matrix data*/
1679371c9d4SSatish Balay     sxadj   = NULL;
1689371c9d4SSatish Balay     sadjncy = NULL;
1699371c9d4SSatish Balay     svalues = NULL;
1709566063dSJacob Faibussowitsch     PetscCall(MatCreateSubMatrix_MPIAdj_data(mat, irow[i], icol[i], &sxadj, &sadjncy, &svalues));
1719566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(irow[i], &irow_n));
1729566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(icol[i], &icol_n));
1739566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(irow[i], &irow_indices));
1749566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(indices, irow_indices, irow_n));
1759566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(irow[i], &irow_indices));
1769566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(icol[i], &icol_indices));
1778e3a54c0SPierre Jolivet     PetscCall(PetscArraycpy(PetscSafePointerPlusOffset(indices, irow_n), icol_indices, icol_n));
1789566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(icol[i], &icol_indices));
17940244768SBarry Smith     nindx = irow_n + icol_n;
1809566063dSJacob Faibussowitsch     PetscCall(PetscSortRemoveDupsInt(&nindx, indices));
18140244768SBarry Smith     /* renumber columns */
18240244768SBarry Smith     for (j = 0; j < irow_n; j++) {
18340244768SBarry Smith       for (k = sxadj[j]; k < sxadj[j + 1]; k++) {
1849566063dSJacob Faibussowitsch         PetscCall(PetscFindInt(sadjncy[k], nindx, indices, &loc));
18508401ef6SPierre Jolivet         PetscCheck(loc >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "can not find col %" PetscInt_FMT, sadjncy[k]);
18640244768SBarry Smith         sadjncy[k] = loc;
18740244768SBarry Smith       }
18840244768SBarry Smith     }
18940244768SBarry Smith     if (scall == MAT_INITIAL_MATRIX) {
1909566063dSJacob Faibussowitsch       PetscCall(MatCreateMPIAdj(scomm_row, irow_n, icol_n, sxadj, sadjncy, svalues, submat[i]));
19140244768SBarry Smith     } else {
19240244768SBarry Smith       Mat         sadj = *(submat[i]);
19340244768SBarry Smith       Mat_MPIAdj *sa   = (Mat_MPIAdj *)((sadj)->data);
1949566063dSJacob Faibussowitsch       PetscCall(PetscObjectGetComm((PetscObject)sadj, &scomm_mat));
1959566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_compare(scomm_row, scomm_mat, &issame));
19608401ef6SPierre Jolivet       PetscCheck(issame == MPI_IDENT, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "submatrix  must have the same comm as the col index set");
1979566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(sa->i, sxadj, irow_n + 1));
1989566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(sa->j, sadjncy, sxadj[irow_n]));
1999566063dSJacob Faibussowitsch       if (svalues) PetscCall(PetscArraycpy(sa->values, svalues, sxadj[irow_n]));
2009566063dSJacob Faibussowitsch       PetscCall(PetscFree(sxadj));
2019566063dSJacob Faibussowitsch       PetscCall(PetscFree(sadjncy));
202252a1336SBarry Smith       PetscCall(PetscFree(svalues));
20340244768SBarry Smith     }
20440244768SBarry Smith   }
2059566063dSJacob Faibussowitsch   PetscCall(PetscFree(indices));
2063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
20740244768SBarry Smith }
20840244768SBarry Smith 
209d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatCreateSubMatricesMPI_MPIAdj(Mat mat, PetscInt n, const IS irow[], const IS icol[], MatReuse scall, Mat *submat[])
210d71ae5a4SJacob Faibussowitsch {
21140244768SBarry Smith   /*get sub-matrices across a sub communicator */
21240244768SBarry Smith   PetscFunctionBegin;
2139566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices_MPIAdj_Private(mat, n, irow, icol, PETSC_TRUE, scall, submat));
2143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
21540244768SBarry Smith }
21640244768SBarry Smith 
217d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatCreateSubMatrices_MPIAdj(Mat mat, PetscInt n, const IS irow[], const IS icol[], MatReuse scall, Mat *submat[])
218d71ae5a4SJacob Faibussowitsch {
21940244768SBarry Smith   PetscFunctionBegin;
22040244768SBarry Smith   /*get sub-matrices based on PETSC_COMM_SELF */
2219566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices_MPIAdj_Private(mat, n, irow, icol, PETSC_FALSE, scall, submat));
2223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
22340244768SBarry Smith }
22440244768SBarry Smith 
225d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatView_MPIAdj_ASCII(Mat A, PetscViewer viewer)
226d71ae5a4SJacob Faibussowitsch {
2273eda8832SBarry Smith   Mat_MPIAdj       *a = (Mat_MPIAdj *)A->data;
228d0f46423SBarry Smith   PetscInt          i, j, m = A->rmap->n;
2292dcb1b2aSMatthew Knepley   const char       *name;
230ce1411ecSBarry Smith   PetscViewerFormat format;
231b97cf49bSBarry Smith 
232433994e6SBarry Smith   PetscFunctionBegin;
2339566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetName((PetscObject)A, &name));
2349566063dSJacob Faibussowitsch   PetscCall(PetscViewerGetFormat(viewer, &format));
235fb9695e5SSatish Balay   if (format == PETSC_VIEWER_ASCII_INFO) {
2363ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
237f7d195e4SLawrence Mitchell   } else {
238f7d195e4SLawrence Mitchell     PetscCheck(format != PETSC_VIEWER_ASCII_MATLAB, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "MATLAB format not supported");
2399566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE));
2409566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushSynchronized(viewer));
241b97cf49bSBarry Smith     for (i = 0; i < m; i++) {
2429566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "row %" PetscInt_FMT ":", i + A->rmap->rstart));
243b97cf49bSBarry Smith       for (j = a->i[i]; j < a->i[i + 1]; j++) {
2440170919eSFande Kong         if (a->values) {
2459566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " (%" PetscInt_FMT ", %" PetscInt_FMT ") ", a->j[j], a->values[j]));
2460170919eSFande Kong         } else {
2479566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " %" PetscInt_FMT " ", a->j[j]));
2480752156aSBarry Smith         }
2490170919eSFande Kong       }
2509566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "\n"));
251b97cf49bSBarry Smith     }
2529566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE));
2539566063dSJacob Faibussowitsch     PetscCall(PetscViewerFlush(viewer));
2549566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopSynchronized(viewer));
2557b23a99aSBarry Smith   }
2563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
257b97cf49bSBarry Smith }
258b97cf49bSBarry Smith 
259d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatView_MPIAdj(Mat A, PetscViewer viewer)
260d71ae5a4SJacob Faibussowitsch {
261ace3abfcSBarry Smith   PetscBool iascii;
262b97cf49bSBarry Smith 
263433994e6SBarry Smith   PetscFunctionBegin;
2649566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
2651baa6e33SBarry Smith   if (iascii) PetscCall(MatView_MPIAdj_ASCII(A, viewer));
2663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
267b97cf49bSBarry Smith }
268b97cf49bSBarry Smith 
269d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatDestroy_MPIAdj(Mat mat)
270d71ae5a4SJacob Faibussowitsch {
2713eda8832SBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj *)mat->data;
27294d884c6SBarry Smith 
27394d884c6SBarry Smith   PetscFunctionBegin;
2743ba16761SJacob Faibussowitsch   PetscCall(PetscLogObjectState((PetscObject)mat, "Rows=%" PetscInt_FMT ", Cols=%" PetscInt_FMT ", NZ=%" PetscInt_FMT, mat->rmap->n, mat->cmap->n, a->nz));
2759566063dSJacob Faibussowitsch   PetscCall(PetscFree(a->diag));
2760400557aSBarry Smith   if (a->freeaij) {
27714196391SBarry Smith     if (a->freeaijwithfree) {
27814196391SBarry Smith       if (a->i) free(a->i);
27914196391SBarry Smith       if (a->j) free(a->j);
28014196391SBarry Smith     } else {
2819566063dSJacob Faibussowitsch       PetscCall(PetscFree(a->i));
2829566063dSJacob Faibussowitsch       PetscCall(PetscFree(a->j));
2839566063dSJacob Faibussowitsch       PetscCall(PetscFree(a->values));
28414196391SBarry Smith     }
2850400557aSBarry Smith   }
2869566063dSJacob Faibussowitsch   PetscCall(PetscFree(a->rowvalues));
2879566063dSJacob Faibussowitsch   PetscCall(PetscFree(mat->data));
2889566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)mat, NULL));
2899566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMPIAdjSetPreallocation_C", NULL));
2909566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMPIAdjCreateNonemptySubcommMat_C", NULL));
2912e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMPIAdjToSeq_C", NULL));
292252a1336SBarry Smith   PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatMPIAdjToSeqRankZero_C", NULL));
2933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
294b97cf49bSBarry Smith }
295b97cf49bSBarry Smith 
296d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSetOption_MPIAdj(Mat A, MatOption op, PetscBool flg)
297d71ae5a4SJacob Faibussowitsch {
2983eda8832SBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj *)A->data;
299b97cf49bSBarry Smith 
300433994e6SBarry Smith   PetscFunctionBegin;
30112c028f9SKris Buschelman   switch (op) {
30277e54ba9SKris Buschelman   case MAT_SYMMETRIC:
30312c028f9SKris Buschelman   case MAT_STRUCTURALLY_SYMMETRIC:
3049a4540c5SBarry Smith   case MAT_HERMITIAN:
305d71ae5a4SJacob Faibussowitsch   case MAT_SPD:
306d71ae5a4SJacob Faibussowitsch     a->symmetric = flg;
307d71ae5a4SJacob Faibussowitsch     break;
3089a4540c5SBarry Smith   case MAT_SYMMETRY_ETERNAL:
309b94d7dedSBarry Smith   case MAT_STRUCTURAL_SYMMETRY_ETERNAL:
310d71ae5a4SJacob Faibussowitsch   case MAT_SPD_ETERNAL:
311d71ae5a4SJacob Faibussowitsch     break;
312d71ae5a4SJacob Faibussowitsch   default:
313d71ae5a4SJacob Faibussowitsch     PetscCall(PetscInfo(A, "Option %s ignored\n", MatOptions[op]));
314d71ae5a4SJacob Faibussowitsch     break;
315b97cf49bSBarry Smith   }
3163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
317b97cf49bSBarry Smith }
318b97cf49bSBarry Smith 
319d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatGetRow_MPIAdj(Mat A, PetscInt row, PetscInt *nz, PetscInt **idx, PetscScalar **v)
320d71ae5a4SJacob Faibussowitsch {
3213eda8832SBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj *)A->data;
322b97cf49bSBarry Smith 
323433994e6SBarry Smith   PetscFunctionBegin;
324d0f46423SBarry Smith   row -= A->rmap->rstart;
325aed4548fSBarry Smith   PetscCheck(row >= 0 && row < A->rmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Row out of range");
326b97cf49bSBarry Smith   *nz = a->i[row + 1] - a->i[row];
3272b1d8763SJed Brown   if (v) {
3282b1d8763SJed Brown     PetscInt j;
3292b1d8763SJed Brown     if (a->rowvalues_alloc < *nz) {
3309566063dSJacob Faibussowitsch       PetscCall(PetscFree(a->rowvalues));
3312b1d8763SJed Brown       a->rowvalues_alloc = PetscMax(a->rowvalues_alloc * 2, *nz);
3329566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(a->rowvalues_alloc, &a->rowvalues));
333b97cf49bSBarry Smith     }
334ad540459SPierre Jolivet     for (j = 0; j < *nz; j++) a->rowvalues[j] = a->values ? a->values[a->i[row] + j] : 1.0;
3352b1d8763SJed Brown     *v = (*nz) ? a->rowvalues : NULL;
336b97cf49bSBarry Smith   }
33792b6f2a0SJed Brown   if (idx) *idx = (*nz) ? a->j + a->i[row] : NULL;
3383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
339b97cf49bSBarry Smith }
340b97cf49bSBarry Smith 
341d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatRestoreRow_MPIAdj(Mat A, PetscInt row, PetscInt *nz, PetscInt **idx, PetscScalar **v)
342d71ae5a4SJacob Faibussowitsch {
343433994e6SBarry Smith   PetscFunctionBegin;
3443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
345b97cf49bSBarry Smith }
346b97cf49bSBarry Smith 
347d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatEqual_MPIAdj(Mat A, Mat B, PetscBool *flg)
348d71ae5a4SJacob Faibussowitsch {
3493eda8832SBarry Smith   Mat_MPIAdj *a = (Mat_MPIAdj *)A->data, *b = (Mat_MPIAdj *)B->data;
350ace3abfcSBarry Smith   PetscBool   flag;
351b97cf49bSBarry Smith 
352433994e6SBarry Smith   PetscFunctionBegin;
353b97cf49bSBarry Smith   /* If the  matrix dimensions are not equal,or no of nonzeros */
354ad540459SPierre Jolivet   if ((A->rmap->n != B->rmap->n) || (a->nz != b->nz)) flag = PETSC_FALSE;
355b97cf49bSBarry Smith 
356b97cf49bSBarry Smith   /* if the a->i are the same */
3579566063dSJacob Faibussowitsch   PetscCall(PetscArraycmp(a->i, b->i, A->rmap->n + 1, &flag));
358b97cf49bSBarry Smith 
359b97cf49bSBarry Smith   /* if a->j are the same */
3609566063dSJacob Faibussowitsch   PetscCall(PetscMemcmp(a->j, b->j, (a->nz) * sizeof(PetscInt), &flag));
361b97cf49bSBarry Smith 
3621c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(&flag, flg, 1, MPIU_BOOL, MPI_LAND, PetscObjectComm((PetscObject)A)));
3633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
364b97cf49bSBarry Smith }
365b97cf49bSBarry Smith 
366d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatGetRowIJ_MPIAdj(Mat A, PetscInt oshift, PetscBool symmetric, PetscBool blockcompressed, PetscInt *m, const PetscInt *inia[], const PetscInt *inja[], PetscBool *done)
367d71ae5a4SJacob Faibussowitsch {
368b24ad042SBarry Smith   PetscInt    i;
3696cd91f64SBarry Smith   Mat_MPIAdj *a  = (Mat_MPIAdj *)A->data;
3701a83f524SJed Brown   PetscInt  **ia = (PetscInt **)inia, **ja = (PetscInt **)inja;
3716cd91f64SBarry Smith 
3726cd91f64SBarry Smith   PetscFunctionBegin;
373d0f46423SBarry Smith   *m    = A->rmap->n;
3746cd91f64SBarry Smith   *ia   = a->i;
3756cd91f64SBarry Smith   *ja   = a->j;
3766cd91f64SBarry Smith   *done = PETSC_TRUE;
377d42735a1SBarry Smith   if (oshift) {
378ad540459SPierre Jolivet     for (i = 0; i < (*ia)[*m]; i++) (*ja)[i]++;
379d42735a1SBarry Smith     for (i = 0; i <= (*m); i++) (*ia)[i]++;
380d42735a1SBarry Smith   }
3813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
382d42735a1SBarry Smith }
383d42735a1SBarry Smith 
384d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatRestoreRowIJ_MPIAdj(Mat A, PetscInt oshift, PetscBool symmetric, PetscBool blockcompressed, PetscInt *m, const PetscInt *inia[], const PetscInt *inja[], PetscBool *done)
385d71ae5a4SJacob Faibussowitsch {
386b24ad042SBarry Smith   PetscInt    i;
387d42735a1SBarry Smith   Mat_MPIAdj *a  = (Mat_MPIAdj *)A->data;
3881a83f524SJed Brown   PetscInt  **ia = (PetscInt **)inia, **ja = (PetscInt **)inja;
389d42735a1SBarry Smith 
390d42735a1SBarry Smith   PetscFunctionBegin;
391aed4548fSBarry Smith   PetscCheck(!ia || a->i == *ia, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "ia passed back is not one obtained with MatGetRowIJ()");
392aed4548fSBarry Smith   PetscCheck(!ja || a->j == *ja, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "ja passed back is not one obtained with MatGetRowIJ()");
393d42735a1SBarry Smith   if (oshift) {
39428b400f6SJacob Faibussowitsch     PetscCheck(ia, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "If oshift then you must passed in inia[] argument");
39528b400f6SJacob Faibussowitsch     PetscCheck(ja, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "If oshift then you must passed in inja[] argument");
396d42735a1SBarry Smith     for (i = 0; i <= (*m); i++) (*ia)[i]--;
397ad540459SPierre Jolivet     for (i = 0; i < (*ia)[*m]; i++) (*ja)[i]--;
398d42735a1SBarry Smith   }
3993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4006cd91f64SBarry Smith }
401b97cf49bSBarry Smith 
40266976f2fSJacob Faibussowitsch static PetscErrorCode MatConvertFrom_MPIAdj(Mat A, MatType type, MatReuse reuse, Mat *newmat)
403d71ae5a4SJacob Faibussowitsch {
40417667f90SBarry Smith   Mat                B;
40517667f90SBarry Smith   PetscInt           i, m, N, nzeros = 0, *ia, *ja, len, rstart, cnt, j, *a;
40617667f90SBarry Smith   const PetscInt    *rj;
40717667f90SBarry Smith   const PetscScalar *ra;
40817667f90SBarry Smith   MPI_Comm           comm;
40917667f90SBarry Smith 
41017667f90SBarry Smith   PetscFunctionBegin;
4119566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A, NULL, &N));
4129566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A, &m, NULL));
4139566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(A, &rstart, NULL));
41417667f90SBarry Smith 
41517667f90SBarry Smith   /* count the number of nonzeros per row */
41617667f90SBarry Smith   for (i = 0; i < m; i++) {
4179566063dSJacob Faibussowitsch     PetscCall(MatGetRow(A, i + rstart, &len, &rj, NULL));
4185ee9ba1cSJed Brown     for (j = 0; j < len; j++) {
4199371c9d4SSatish Balay       if (rj[j] == i + rstart) {
4209371c9d4SSatish Balay         len--;
4219371c9d4SSatish Balay         break;
4229371c9d4SSatish Balay       } /* don't count diagonal */
4235ee9ba1cSJed Brown     }
42417667f90SBarry Smith     nzeros += len;
4259566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(A, i + rstart, &len, &rj, NULL));
42617667f90SBarry Smith   }
42717667f90SBarry Smith 
42817667f90SBarry Smith   /* malloc space for nonzeros */
4299566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nzeros + 1, &a));
4309566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(N + 1, &ia));
4319566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nzeros + 1, &ja));
43217667f90SBarry Smith 
43317667f90SBarry Smith   nzeros = 0;
43417667f90SBarry Smith   ia[0]  = 0;
43517667f90SBarry Smith   for (i = 0; i < m; i++) {
4369566063dSJacob Faibussowitsch     PetscCall(MatGetRow(A, i + rstart, &len, &rj, &ra));
43717667f90SBarry Smith     cnt = 0;
43817667f90SBarry Smith     for (j = 0; j < len; j++) {
4395ee9ba1cSJed Brown       if (rj[j] != i + rstart) { /* if not diagonal */
44017667f90SBarry Smith         a[nzeros + cnt]    = (PetscInt)PetscAbsScalar(ra[j]);
44117667f90SBarry Smith         ja[nzeros + cnt++] = rj[j];
44217667f90SBarry Smith       }
4435ee9ba1cSJed Brown     }
4449566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(A, i + rstart, &len, &rj, &ra));
44517667f90SBarry Smith     nzeros += cnt;
44617667f90SBarry Smith     ia[i + 1] = nzeros;
44717667f90SBarry Smith   }
44817667f90SBarry Smith 
4499566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)A, &comm));
4509566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, &B));
4519566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(B, m, PETSC_DETERMINE, PETSC_DETERMINE, N));
4529566063dSJacob Faibussowitsch   PetscCall(MatSetType(B, type));
4539566063dSJacob Faibussowitsch   PetscCall(MatMPIAdjSetPreallocation(B, ia, ja, a));
45417667f90SBarry Smith 
455511c6705SHong Zhang   if (reuse == MAT_INPLACE_MATRIX) {
4569566063dSJacob Faibussowitsch     PetscCall(MatHeaderReplace(A, &B));
45717667f90SBarry Smith   } else {
45817667f90SBarry Smith     *newmat = B;
45917667f90SBarry Smith   }
4603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
46117667f90SBarry Smith }
46217667f90SBarry Smith 
46366976f2fSJacob Faibussowitsch static PetscErrorCode MatSetValues_MPIAdj(Mat A, PetscInt m, const PetscInt *rows, PetscInt n, const PetscInt *cols, const PetscScalar *values, InsertMode im)
464d71ae5a4SJacob Faibussowitsch {
4656a09307cSBarry Smith   Mat_MPIAdj *adj = (Mat_MPIAdj *)A->data;
4666a09307cSBarry Smith   PetscInt    rStart, rEnd, cStart, cEnd;
4676a09307cSBarry Smith 
4686a09307cSBarry Smith   PetscFunctionBegin;
4696a09307cSBarry Smith   PetscCheck(!adj->i, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Matrix is already assembled, cannot change its values");
4706a09307cSBarry Smith   PetscCall(MatGetOwnershipRange(A, &rStart, &rEnd));
4716a09307cSBarry Smith   PetscCall(MatGetOwnershipRangeColumn(A, &cStart, &cEnd));
4726a09307cSBarry Smith   if (!adj->ht) {
4736a09307cSBarry Smith     PetscCall(PetscHSetIJCreate(&adj->ht));
4746a09307cSBarry Smith     PetscCall(MatStashCreate_Private(PetscObjectComm((PetscObject)A), 1, &A->stash));
4756a09307cSBarry Smith     PetscCall(PetscLayoutSetUp(A->rmap));
4766a09307cSBarry Smith     PetscCall(PetscLayoutSetUp(A->cmap));
4776a09307cSBarry Smith   }
4786a09307cSBarry Smith   for (PetscInt r = 0; r < m; ++r) {
4796a09307cSBarry Smith     PetscHashIJKey key;
4806a09307cSBarry Smith 
4816a09307cSBarry Smith     key.i = rows[r];
4826a09307cSBarry Smith     if (key.i < 0) continue;
4836a09307cSBarry Smith     if ((key.i < rStart) || (key.i >= rEnd)) {
4846a09307cSBarry Smith       PetscCall(MatStashValuesRow_Private(&A->stash, key.i, n, cols, values, PETSC_FALSE));
4856a09307cSBarry Smith     } else {
4866a09307cSBarry Smith       for (PetscInt c = 0; c < n; ++c) {
4876a09307cSBarry Smith         key.j = cols[c];
4886a09307cSBarry Smith         if (key.j < 0 || key.i == key.j) continue;
4896a09307cSBarry Smith         PetscCall(PetscHSetIJAdd(adj->ht, key));
4906a09307cSBarry Smith       }
4916a09307cSBarry Smith     }
4926a09307cSBarry Smith   }
4933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4946a09307cSBarry Smith }
4956a09307cSBarry Smith 
49666976f2fSJacob Faibussowitsch static PetscErrorCode MatAssemblyBegin_MPIAdj(Mat A, MatAssemblyType type)
497d71ae5a4SJacob Faibussowitsch {
4986a09307cSBarry Smith   PetscInt    nstash, reallocs;
4996a09307cSBarry Smith   Mat_MPIAdj *adj = (Mat_MPIAdj *)A->data;
5006a09307cSBarry Smith 
5016a09307cSBarry Smith   PetscFunctionBegin;
5026a09307cSBarry Smith   if (!adj->ht) {
5036a09307cSBarry Smith     PetscCall(PetscHSetIJCreate(&adj->ht));
5046a09307cSBarry Smith     PetscCall(PetscLayoutSetUp(A->rmap));
5056a09307cSBarry Smith     PetscCall(PetscLayoutSetUp(A->cmap));
5066a09307cSBarry Smith   }
5076a09307cSBarry Smith   PetscCall(MatStashScatterBegin_Private(A, &A->stash, A->rmap->range));
5086a09307cSBarry Smith   PetscCall(MatStashGetInfo_Private(&A->stash, &nstash, &reallocs));
5096a09307cSBarry Smith   PetscCall(PetscInfo(A, "Stash has %" PetscInt_FMT " entries, uses %" PetscInt_FMT " mallocs.\n", nstash, reallocs));
5103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5116a09307cSBarry Smith }
5126a09307cSBarry Smith 
51366976f2fSJacob Faibussowitsch static PetscErrorCode MatAssemblyEnd_MPIAdj(Mat A, MatAssemblyType type)
514d71ae5a4SJacob Faibussowitsch {
5156a09307cSBarry Smith   PetscScalar   *val;
5166a09307cSBarry Smith   PetscInt      *row, *col, m, rstart, *rowstarts;
5176a09307cSBarry Smith   PetscInt       i, j, ncols, flg, nz;
5186a09307cSBarry Smith   PetscMPIInt    n;
5196a09307cSBarry Smith   Mat_MPIAdj    *adj = (Mat_MPIAdj *)A->data;
5206a09307cSBarry Smith   PetscHashIter  hi;
5216a09307cSBarry Smith   PetscHashIJKey key;
5226a09307cSBarry Smith   PetscHSetIJ    ht = adj->ht;
5236a09307cSBarry Smith 
5246a09307cSBarry Smith   PetscFunctionBegin;
5256a09307cSBarry Smith   while (1) {
5266a09307cSBarry Smith     PetscCall(MatStashScatterGetMesg_Private(&A->stash, &n, &row, &col, &val, &flg));
5276a09307cSBarry Smith     if (!flg) break;
5286a09307cSBarry Smith 
5296a09307cSBarry Smith     for (i = 0; i < n;) {
5306a09307cSBarry Smith       /* Identify the consecutive vals belonging to the same row */
5316a09307cSBarry Smith       for (j = i, rstart = row[j]; j < n; j++) {
5326a09307cSBarry Smith         if (row[j] != rstart) break;
5336a09307cSBarry Smith       }
5346a09307cSBarry Smith       if (j < n) ncols = j - i;
5356a09307cSBarry Smith       else ncols = n - i;
5366a09307cSBarry Smith       /* Set all these values with a single function call */
5376a09307cSBarry Smith       PetscCall(MatSetValues_MPIAdj(A, 1, row + i, ncols, col + i, val + i, INSERT_VALUES));
5386a09307cSBarry Smith       i = j;
5396a09307cSBarry Smith     }
5406a09307cSBarry Smith   }
5416a09307cSBarry Smith   PetscCall(MatStashScatterEnd_Private(&A->stash));
5426a09307cSBarry Smith   PetscCall(MatStashDestroy_Private(&A->stash));
5436a09307cSBarry Smith 
5446a09307cSBarry Smith   PetscCall(MatGetLocalSize(A, &m, NULL));
5456a09307cSBarry Smith   PetscCall(MatGetOwnershipRange(A, &rstart, NULL));
5466a09307cSBarry Smith   PetscCall(PetscCalloc1(m + 1, &rowstarts));
5476a09307cSBarry Smith   PetscHashIterBegin(ht, hi);
5486a09307cSBarry Smith   for (; !PetscHashIterAtEnd(ht, hi);) {
5496a09307cSBarry Smith     PetscHashIterGetKey(ht, hi, key);
5506a09307cSBarry Smith     rowstarts[key.i - rstart + 1]++;
5516a09307cSBarry Smith     PetscHashIterNext(ht, hi);
5526a09307cSBarry Smith   }
553ad540459SPierre Jolivet   for (i = 1; i < m + 1; i++) rowstarts[i] = rowstarts[i - 1] + rowstarts[i];
5546a09307cSBarry Smith 
5556a09307cSBarry Smith   PetscCall(PetscHSetIJGetSize(ht, &nz));
5566a09307cSBarry Smith   PetscCall(PetscMalloc1(nz, &col));
5576a09307cSBarry Smith   PetscHashIterBegin(ht, hi);
5586a09307cSBarry Smith   for (; !PetscHashIterAtEnd(ht, hi);) {
5596a09307cSBarry Smith     PetscHashIterGetKey(ht, hi, key);
5606a09307cSBarry Smith     col[rowstarts[key.i - rstart]++] = key.j;
5616a09307cSBarry Smith     PetscHashIterNext(ht, hi);
5626a09307cSBarry Smith   }
5636a09307cSBarry Smith   PetscCall(PetscHSetIJDestroy(&ht));
564ad540459SPierre Jolivet   for (i = m; i > 0; i--) rowstarts[i] = rowstarts[i - 1];
5656a09307cSBarry Smith   rowstarts[0] = 0;
5666a09307cSBarry Smith 
56748a46eb9SPierre Jolivet   for (PetscInt i = 0; i < m; i++) PetscCall(PetscSortInt(rowstarts[i + 1] - rowstarts[i], &col[rowstarts[i]]));
5686a09307cSBarry Smith 
5696a09307cSBarry Smith   adj->i       = rowstarts;
5706a09307cSBarry Smith   adj->j       = col;
571252a1336SBarry Smith   adj->nz      = rowstarts[m];
5726a09307cSBarry Smith   adj->freeaij = PETSC_TRUE;
5733ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5746a09307cSBarry Smith }
5756a09307cSBarry Smith 
5766a09307cSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_MPIAdj,
5773eda8832SBarry Smith                                        MatGetRow_MPIAdj,
5783eda8832SBarry Smith                                        MatRestoreRow_MPIAdj,
579f4259b30SLisandro Dalcin                                        NULL,
580f4259b30SLisandro Dalcin                                        /* 4*/ NULL,
581f4259b30SLisandro Dalcin                                        NULL,
582f4259b30SLisandro Dalcin                                        NULL,
583f4259b30SLisandro Dalcin                                        NULL,
584f4259b30SLisandro Dalcin                                        NULL,
585f4259b30SLisandro Dalcin                                        NULL,
586f4259b30SLisandro Dalcin                                        /*10*/ NULL,
587f4259b30SLisandro Dalcin                                        NULL,
588f4259b30SLisandro Dalcin                                        NULL,
589f4259b30SLisandro Dalcin                                        NULL,
590f4259b30SLisandro Dalcin                                        NULL,
591f4259b30SLisandro Dalcin                                        /*15*/ NULL,
5923eda8832SBarry Smith                                        MatEqual_MPIAdj,
593f4259b30SLisandro Dalcin                                        NULL,
594f4259b30SLisandro Dalcin                                        NULL,
595f4259b30SLisandro Dalcin                                        NULL,
5966a09307cSBarry Smith                                        /*20*/ MatAssemblyBegin_MPIAdj,
5976a09307cSBarry Smith                                        MatAssemblyEnd_MPIAdj,
5983eda8832SBarry Smith                                        MatSetOption_MPIAdj,
599f4259b30SLisandro Dalcin                                        NULL,
600f4259b30SLisandro Dalcin                                        /*24*/ NULL,
601f4259b30SLisandro Dalcin                                        NULL,
602f4259b30SLisandro Dalcin                                        NULL,
603f4259b30SLisandro Dalcin                                        NULL,
604f4259b30SLisandro Dalcin                                        NULL,
605f4259b30SLisandro Dalcin                                        /*29*/ NULL,
606f4259b30SLisandro Dalcin                                        NULL,
607f4259b30SLisandro Dalcin                                        NULL,
608f4259b30SLisandro Dalcin                                        NULL,
609f4259b30SLisandro Dalcin                                        NULL,
610f4259b30SLisandro Dalcin                                        /*34*/ NULL,
611f4259b30SLisandro Dalcin                                        NULL,
612f4259b30SLisandro Dalcin                                        NULL,
613f4259b30SLisandro Dalcin                                        NULL,
614f4259b30SLisandro Dalcin                                        NULL,
615f4259b30SLisandro Dalcin                                        /*39*/ NULL,
6167dae84e0SHong Zhang                                        MatCreateSubMatrices_MPIAdj,
617f4259b30SLisandro Dalcin                                        NULL,
618f4259b30SLisandro Dalcin                                        NULL,
619f4259b30SLisandro Dalcin                                        NULL,
620f4259b30SLisandro Dalcin                                        /*44*/ NULL,
621f4259b30SLisandro Dalcin                                        NULL,
6227d68702bSBarry Smith                                        MatShift_Basic,
623f4259b30SLisandro Dalcin                                        NULL,
624f4259b30SLisandro Dalcin                                        NULL,
625f4259b30SLisandro Dalcin                                        /*49*/ NULL,
6266cd91f64SBarry Smith                                        MatGetRowIJ_MPIAdj,
627d42735a1SBarry Smith                                        MatRestoreRowIJ_MPIAdj,
628f4259b30SLisandro Dalcin                                        NULL,
629f4259b30SLisandro Dalcin                                        NULL,
630f4259b30SLisandro Dalcin                                        /*54*/ NULL,
631f4259b30SLisandro Dalcin                                        NULL,
632f4259b30SLisandro Dalcin                                        NULL,
633f4259b30SLisandro Dalcin                                        NULL,
634f4259b30SLisandro Dalcin                                        NULL,
635f4259b30SLisandro Dalcin                                        /*59*/ NULL,
636b9b97703SBarry Smith                                        MatDestroy_MPIAdj,
637b9b97703SBarry Smith                                        MatView_MPIAdj,
63817667f90SBarry Smith                                        MatConvertFrom_MPIAdj,
639f4259b30SLisandro Dalcin                                        NULL,
640f4259b30SLisandro Dalcin                                        /*64*/ NULL,
641f4259b30SLisandro Dalcin                                        NULL,
642f4259b30SLisandro Dalcin                                        NULL,
643f4259b30SLisandro Dalcin                                        NULL,
644f4259b30SLisandro Dalcin                                        NULL,
645f4259b30SLisandro Dalcin                                        /*69*/ NULL,
646f4259b30SLisandro Dalcin                                        NULL,
647f4259b30SLisandro Dalcin                                        NULL,
648f4259b30SLisandro Dalcin                                        NULL,
649f4259b30SLisandro Dalcin                                        NULL,
650f4259b30SLisandro Dalcin                                        /*74*/ NULL,
651f4259b30SLisandro Dalcin                                        NULL,
652f4259b30SLisandro Dalcin                                        NULL,
653f4259b30SLisandro Dalcin                                        NULL,
654f4259b30SLisandro Dalcin                                        NULL,
655f4259b30SLisandro Dalcin                                        /*79*/ NULL,
656f4259b30SLisandro Dalcin                                        NULL,
657f4259b30SLisandro Dalcin                                        NULL,
658f4259b30SLisandro Dalcin                                        NULL,
659f4259b30SLisandro Dalcin                                        NULL,
660f4259b30SLisandro Dalcin                                        /*84*/ NULL,
661f4259b30SLisandro Dalcin                                        NULL,
662f4259b30SLisandro Dalcin                                        NULL,
663f4259b30SLisandro Dalcin                                        NULL,
664f4259b30SLisandro Dalcin                                        NULL,
665f4259b30SLisandro Dalcin                                        /*89*/ NULL,
666f4259b30SLisandro Dalcin                                        NULL,
667f4259b30SLisandro Dalcin                                        NULL,
668f4259b30SLisandro Dalcin                                        NULL,
669f4259b30SLisandro Dalcin                                        NULL,
670f4259b30SLisandro Dalcin                                        /*94*/ NULL,
671f4259b30SLisandro Dalcin                                        NULL,
672f4259b30SLisandro Dalcin                                        NULL,
673f4259b30SLisandro Dalcin                                        NULL,
674f4259b30SLisandro Dalcin                                        NULL,
675f4259b30SLisandro Dalcin                                        /*99*/ NULL,
676f4259b30SLisandro Dalcin                                        NULL,
677f4259b30SLisandro Dalcin                                        NULL,
678f4259b30SLisandro Dalcin                                        NULL,
679f4259b30SLisandro Dalcin                                        NULL,
680f4259b30SLisandro Dalcin                                        /*104*/ NULL,
681f4259b30SLisandro Dalcin                                        NULL,
682f4259b30SLisandro Dalcin                                        NULL,
683f4259b30SLisandro Dalcin                                        NULL,
684f4259b30SLisandro Dalcin                                        NULL,
685f4259b30SLisandro Dalcin                                        /*109*/ NULL,
686f4259b30SLisandro Dalcin                                        NULL,
687f4259b30SLisandro Dalcin                                        NULL,
688f4259b30SLisandro Dalcin                                        NULL,
689f4259b30SLisandro Dalcin                                        NULL,
690f4259b30SLisandro Dalcin                                        /*114*/ NULL,
691f4259b30SLisandro Dalcin                                        NULL,
692f4259b30SLisandro Dalcin                                        NULL,
693f4259b30SLisandro Dalcin                                        NULL,
694f4259b30SLisandro Dalcin                                        NULL,
695f4259b30SLisandro Dalcin                                        /*119*/ NULL,
696f4259b30SLisandro Dalcin                                        NULL,
697f4259b30SLisandro Dalcin                                        NULL,
698f4259b30SLisandro Dalcin                                        NULL,
699f4259b30SLisandro Dalcin                                        NULL,
700f4259b30SLisandro Dalcin                                        /*124*/ NULL,
701f4259b30SLisandro Dalcin                                        NULL,
702f4259b30SLisandro Dalcin                                        NULL,
703f4259b30SLisandro Dalcin                                        NULL,
7047dae84e0SHong Zhang                                        MatCreateSubMatricesMPI_MPIAdj,
705f4259b30SLisandro Dalcin                                        /*129*/ NULL,
706f4259b30SLisandro Dalcin                                        NULL,
707f4259b30SLisandro Dalcin                                        NULL,
708f4259b30SLisandro Dalcin                                        NULL,
709f4259b30SLisandro Dalcin                                        NULL,
710f4259b30SLisandro Dalcin                                        /*134*/ NULL,
711f4259b30SLisandro Dalcin                                        NULL,
712f4259b30SLisandro Dalcin                                        NULL,
713f4259b30SLisandro Dalcin                                        NULL,
714f4259b30SLisandro Dalcin                                        NULL,
715f4259b30SLisandro Dalcin                                        /*139*/ NULL,
716f4259b30SLisandro Dalcin                                        NULL,
717d70f29a3SPierre Jolivet                                        NULL,
718d70f29a3SPierre Jolivet                                        NULL,
719d70f29a3SPierre Jolivet                                        NULL,
720d70f29a3SPierre Jolivet                                        /*144*/ NULL,
721d70f29a3SPierre Jolivet                                        NULL,
722d70f29a3SPierre Jolivet                                        NULL,
72399a7f59eSMark Adams                                        NULL,
72499a7f59eSMark Adams                                        NULL,
7257fb60732SBarry Smith                                        NULL,
726dec0b466SHong Zhang                                        /*150*/ NULL,
727eede4a3fSMark Adams                                        NULL,
728dec0b466SHong Zhang                                        NULL};
729b97cf49bSBarry Smith 
730d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMPIAdjSetPreallocation_MPIAdj(Mat B, PetscInt *i, PetscInt *j, PetscInt *values)
731d71ae5a4SJacob Faibussowitsch {
732a23d5eceSKris Buschelman   Mat_MPIAdj *b = (Mat_MPIAdj *)B->data;
733e895ccc0SFande Kong   PetscBool   useedgeweights;
734a23d5eceSKris Buschelman 
735a23d5eceSKris Buschelman   PetscFunctionBegin;
7369566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(B->rmap));
7379566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(B->cmap));
7389371c9d4SSatish Balay   if (values) useedgeweights = PETSC_TRUE;
7399371c9d4SSatish Balay   else useedgeweights = PETSC_FALSE;
740e895ccc0SFande Kong   /* Make everybody knows if they are using edge weights or not */
7411c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce((int *)&useedgeweights, (int *)&b->useedgeweights, 1, MPI_INT, MPI_MAX, PetscObjectComm((PetscObject)B)));
74258ec18a5SLisandro Dalcin 
74376bd3646SJed Brown   if (PetscDefined(USE_DEBUG)) {
74476bd3646SJed Brown     PetscInt ii;
74576bd3646SJed Brown 
74608401ef6SPierre Jolivet     PetscCheck(i[0] == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "First i[] index must be zero, instead it is %" PetscInt_FMT, i[0]);
747d0f46423SBarry Smith     for (ii = 1; ii < B->rmap->n; ii++) {
748aed4548fSBarry Smith       PetscCheck(i[ii] >= 0 && i[ii] >= i[ii - 1], PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "i[%" PetscInt_FMT "]=%" PetscInt_FMT " index is out of range: i[%" PetscInt_FMT "]=%" PetscInt_FMT, ii, i[ii], ii - 1, i[ii - 1]);
749a23d5eceSKris Buschelman     }
750ad540459SPierre Jolivet     for (ii = 0; ii < i[B->rmap->n]; ii++) PetscCheck(j[ii] >= 0 && j[ii] < B->cmap->N, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Column index %" PetscInt_FMT " out of range %" PetscInt_FMT, ii, j[ii]);
75176bd3646SJed Brown   }
752a23d5eceSKris Buschelman   b->j      = j;
753a23d5eceSKris Buschelman   b->i      = i;
754a23d5eceSKris Buschelman   b->values = values;
755a23d5eceSKris Buschelman 
756d0f46423SBarry Smith   b->nz        = i[B->rmap->n];
757f4259b30SLisandro Dalcin   b->diag      = NULL;
758a23d5eceSKris Buschelman   b->symmetric = PETSC_FALSE;
759a23d5eceSKris Buschelman   b->freeaij   = PETSC_TRUE;
760a23d5eceSKris Buschelman 
7616a09307cSBarry Smith   B->ops->assemblybegin = NULL;
7626a09307cSBarry Smith   B->ops->assemblyend   = NULL;
7639566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
7649566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
7656a09307cSBarry Smith   PetscCall(MatStashDestroy_Private(&B->stash));
7663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
767a23d5eceSKris Buschelman }
768b97cf49bSBarry Smith 
769d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMPIAdjCreateNonemptySubcommMat_MPIAdj(Mat A, Mat *B)
770d71ae5a4SJacob Faibussowitsch {
7719a3665c6SJed Brown   Mat_MPIAdj     *a = (Mat_MPIAdj *)A->data;
7729a3665c6SJed Brown   const PetscInt *ranges;
7739a3665c6SJed Brown   MPI_Comm        acomm, bcomm;
7749a3665c6SJed Brown   MPI_Group       agroup, bgroup;
7759a3665c6SJed Brown   PetscMPIInt     i, rank, size, nranks, *ranks;
7769a3665c6SJed Brown 
7779a3665c6SJed Brown   PetscFunctionBegin;
7780298fd71SBarry Smith   *B = NULL;
7799566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)A, &acomm));
7809566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(acomm, &size));
7819566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(acomm, &rank));
7829566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRanges(A, &ranges));
7839a3665c6SJed Brown   for (i = 0, nranks = 0; i < size; i++) {
7849a3665c6SJed Brown     if (ranges[i + 1] - ranges[i] > 0) nranks++;
7859a3665c6SJed Brown   }
7869a3665c6SJed Brown   if (nranks == size) { /* All ranks have a positive number of rows, so we do not need to create a subcomm; */
7879566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)A));
7889a3665c6SJed Brown     *B = A;
7893ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
7909a3665c6SJed Brown   }
7919a3665c6SJed Brown 
7929566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nranks, &ranks));
7939a3665c6SJed Brown   for (i = 0, nranks = 0; i < size; i++) {
7949a3665c6SJed Brown     if (ranges[i + 1] - ranges[i] > 0) ranks[nranks++] = i;
7959a3665c6SJed Brown   }
7969566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_group(acomm, &agroup));
7979566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Group_incl(agroup, nranks, ranks, &bgroup));
7989566063dSJacob Faibussowitsch   PetscCall(PetscFree(ranks));
7999566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_create(acomm, bgroup, &bcomm));
8009566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Group_free(&agroup));
8019566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Group_free(&bgroup));
8029a3665c6SJed Brown   if (bcomm != MPI_COMM_NULL) {
8039a3665c6SJed Brown     PetscInt    m, N;
8049a3665c6SJed Brown     Mat_MPIAdj *b;
8059566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(A, &m, NULL));
8069566063dSJacob Faibussowitsch     PetscCall(MatGetSize(A, NULL, &N));
8079566063dSJacob Faibussowitsch     PetscCall(MatCreateMPIAdj(bcomm, m, N, a->i, a->j, a->values, B));
8089a3665c6SJed Brown     b          = (Mat_MPIAdj *)(*B)->data;
8099a3665c6SJed Brown     b->freeaij = PETSC_FALSE;
8109566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_free(&bcomm));
8119a3665c6SJed Brown   }
8123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8139a3665c6SJed Brown }
8149a3665c6SJed Brown 
81566976f2fSJacob Faibussowitsch static PetscErrorCode MatMPIAdjToSeq_MPIAdj(Mat A, Mat *B)
816d71ae5a4SJacob Faibussowitsch {
817b44e7856SBarry Smith   PetscInt    M, N, *II, *J, NZ, nz, m, nzstart, i;
818b44e7856SBarry Smith   PetscInt   *Values = NULL;
819b44e7856SBarry Smith   Mat_MPIAdj *adj    = (Mat_MPIAdj *)A->data;
820b44e7856SBarry Smith   PetscMPIInt mnz, mm, *allnz, *allm, size, *dispnz, *dispm;
821b44e7856SBarry Smith 
822b44e7856SBarry Smith   PetscFunctionBegin;
8239566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
8249566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A, &M, &N));
8259566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A, &m, NULL));
826b44e7856SBarry Smith   nz = adj->nz;
82708401ef6SPierre Jolivet   PetscCheck(adj->i[m] == nz, PETSC_COMM_SELF, PETSC_ERR_PLIB, "nz %" PetscInt_FMT " not correct i[m] %" PetscInt_FMT, nz, adj->i[m]);
828712fec58SPierre Jolivet   PetscCall(MPIU_Allreduce(&nz, &NZ, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)A)));
829d4e69552SBarry Smith 
8309566063dSJacob Faibussowitsch   PetscCall(PetscMPIIntCast(nz, &mnz));
8319566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(size, &allnz, size, &dispnz));
8329566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allgather(&mnz, 1, MPI_INT, allnz, 1, MPI_INT, PetscObjectComm((PetscObject)A)));
8339371c9d4SSatish Balay   dispnz[0] = 0;
8349371c9d4SSatish Balay   for (i = 1; i < size; i++) dispnz[i] = dispnz[i - 1] + allnz[i - 1];
835b44e7856SBarry Smith   if (adj->values) {
8369566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(NZ, &Values));
8379566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Allgatherv(adj->values, mnz, MPIU_INT, Values, allnz, dispnz, MPIU_INT, PetscObjectComm((PetscObject)A)));
838b44e7856SBarry Smith   }
8399566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(NZ, &J));
8409566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allgatherv(adj->j, mnz, MPIU_INT, J, allnz, dispnz, MPIU_INT, PetscObjectComm((PetscObject)A)));
8419566063dSJacob Faibussowitsch   PetscCall(PetscFree2(allnz, dispnz));
8429566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Scan(&nz, &nzstart, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)A)));
843b44e7856SBarry Smith   nzstart -= nz;
844b44e7856SBarry Smith   /* shift the i[] values so they will be correct after being received */
845b44e7856SBarry Smith   for (i = 0; i < m; i++) adj->i[i] += nzstart;
8469566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(M + 1, &II));
8479566063dSJacob Faibussowitsch   PetscCall(PetscMPIIntCast(m, &mm));
8489566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(size, &allm, size, &dispm));
8499566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allgather(&mm, 1, MPI_INT, allm, 1, MPI_INT, PetscObjectComm((PetscObject)A)));
8509371c9d4SSatish Balay   dispm[0] = 0;
8519371c9d4SSatish Balay   for (i = 1; i < size; i++) dispm[i] = dispm[i - 1] + allm[i - 1];
8529566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allgatherv(adj->i, mm, MPIU_INT, II, allm, dispm, MPIU_INT, PetscObjectComm((PetscObject)A)));
8539566063dSJacob Faibussowitsch   PetscCall(PetscFree2(allm, dispm));
854b44e7856SBarry Smith   II[M] = NZ;
855b44e7856SBarry Smith   /* shift the i[] values back */
856b44e7856SBarry Smith   for (i = 0; i < m; i++) adj->i[i] -= nzstart;
8579566063dSJacob Faibussowitsch   PetscCall(MatCreateMPIAdj(PETSC_COMM_SELF, M, N, II, J, Values, B));
8583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
859b44e7856SBarry Smith }
860b44e7856SBarry Smith 
86166976f2fSJacob Faibussowitsch static PetscErrorCode MatMPIAdjToSeqRankZero_MPIAdj(Mat A, Mat *B)
862d71ae5a4SJacob Faibussowitsch {
863252a1336SBarry Smith   PetscInt    M, N, *II, *J, NZ, nz, m, nzstart, i;
864252a1336SBarry Smith   PetscInt   *Values = NULL;
865252a1336SBarry Smith   Mat_MPIAdj *adj    = (Mat_MPIAdj *)A->data;
866252a1336SBarry Smith   PetscMPIInt mnz, mm, *allnz = NULL, *allm, size, *dispnz, *dispm, rank;
867252a1336SBarry Smith 
868252a1336SBarry Smith   PetscFunctionBegin;
869252a1336SBarry Smith   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
870252a1336SBarry Smith   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)A), &rank));
871252a1336SBarry Smith   PetscCall(MatGetSize(A, &M, &N));
872252a1336SBarry Smith   PetscCall(MatGetLocalSize(A, &m, NULL));
873252a1336SBarry Smith   nz = adj->nz;
874252a1336SBarry Smith   PetscCheck(adj->i[m] == nz, PETSC_COMM_SELF, PETSC_ERR_PLIB, "nz %" PetscInt_FMT " not correct i[m] %" PetscInt_FMT, nz, adj->i[m]);
875712fec58SPierre Jolivet   PetscCall(MPIU_Allreduce(&nz, &NZ, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)A)));
876252a1336SBarry Smith 
877252a1336SBarry Smith   PetscCall(PetscMPIIntCast(nz, &mnz));
878252a1336SBarry Smith   if (!rank) PetscCall(PetscMalloc2(size, &allnz, size, &dispnz));
879252a1336SBarry Smith   PetscCallMPI(MPI_Gather(&mnz, 1, MPI_INT, allnz, 1, MPI_INT, 0, PetscObjectComm((PetscObject)A)));
880252a1336SBarry Smith   if (!rank) {
8819371c9d4SSatish Balay     dispnz[0] = 0;
8829371c9d4SSatish Balay     for (i = 1; i < size; i++) dispnz[i] = dispnz[i - 1] + allnz[i - 1];
883252a1336SBarry Smith     if (adj->values) {
884252a1336SBarry Smith       PetscCall(PetscMalloc1(NZ, &Values));
885252a1336SBarry Smith       PetscCallMPI(MPI_Gatherv(adj->values, mnz, MPIU_INT, Values, allnz, dispnz, MPIU_INT, 0, PetscObjectComm((PetscObject)A)));
886252a1336SBarry Smith     }
887252a1336SBarry Smith     PetscCall(PetscMalloc1(NZ, &J));
888252a1336SBarry Smith     PetscCallMPI(MPI_Gatherv(adj->j, mnz, MPIU_INT, J, allnz, dispnz, MPIU_INT, 0, PetscObjectComm((PetscObject)A)));
889252a1336SBarry Smith     PetscCall(PetscFree2(allnz, dispnz));
890252a1336SBarry Smith   } else {
891252a1336SBarry Smith     if (adj->values) PetscCallMPI(MPI_Gatherv(adj->values, mnz, MPIU_INT, NULL, NULL, NULL, MPIU_INT, 0, PetscObjectComm((PetscObject)A)));
892252a1336SBarry Smith     PetscCallMPI(MPI_Gatherv(adj->j, mnz, MPIU_INT, NULL, NULL, NULL, MPIU_INT, 0, PetscObjectComm((PetscObject)A)));
893252a1336SBarry Smith   }
894252a1336SBarry Smith   PetscCallMPI(MPI_Scan(&nz, &nzstart, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)A)));
895252a1336SBarry Smith   nzstart -= nz;
896252a1336SBarry Smith   /* shift the i[] values so they will be correct after being received */
897252a1336SBarry Smith   for (i = 0; i < m; i++) adj->i[i] += nzstart;
898252a1336SBarry Smith   PetscCall(PetscMPIIntCast(m, &mm));
899252a1336SBarry Smith   if (!rank) {
900252a1336SBarry Smith     PetscCall(PetscMalloc1(M + 1, &II));
901252a1336SBarry Smith     PetscCall(PetscMalloc2(size, &allm, size, &dispm));
902252a1336SBarry Smith     PetscCallMPI(MPI_Gather(&mm, 1, MPI_INT, allm, 1, MPI_INT, 0, PetscObjectComm((PetscObject)A)));
9039371c9d4SSatish Balay     dispm[0] = 0;
9049371c9d4SSatish Balay     for (i = 1; i < size; i++) dispm[i] = dispm[i - 1] + allm[i - 1];
905252a1336SBarry Smith     PetscCallMPI(MPI_Gatherv(adj->i, mm, MPIU_INT, II, allm, dispm, MPIU_INT, 0, PetscObjectComm((PetscObject)A)));
906252a1336SBarry Smith     PetscCall(PetscFree2(allm, dispm));
907252a1336SBarry Smith     II[M] = NZ;
908252a1336SBarry Smith   } else {
909252a1336SBarry Smith     PetscCallMPI(MPI_Gather(&mm, 1, MPI_INT, NULL, 1, MPI_INT, 0, PetscObjectComm((PetscObject)A)));
910252a1336SBarry Smith     PetscCallMPI(MPI_Gatherv(adj->i, mm, MPIU_INT, NULL, NULL, NULL, MPIU_INT, 0, PetscObjectComm((PetscObject)A)));
911252a1336SBarry Smith   }
912252a1336SBarry Smith   /* shift the i[] values back */
913252a1336SBarry Smith   for (i = 0; i < m; i++) adj->i[i] -= nzstart;
914252a1336SBarry Smith   if (!rank) PetscCall(MatCreateMPIAdj(PETSC_COMM_SELF, M, N, II, J, Values, B));
9153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
916252a1336SBarry Smith }
917252a1336SBarry Smith 
9189a3665c6SJed Brown /*@
91911a5261eSBarry Smith   MatMPIAdjCreateNonemptySubcommMat - create the same `MATMPIADJ` matrix on a subcommunicator containing only processes owning a positive number of rows
9209a3665c6SJed Brown 
9219a3665c6SJed Brown   Collective
9229a3665c6SJed Brown 
9234165533cSJose E. Roman   Input Parameter:
9242ef1f0ffSBarry Smith . A - original `MATMPIADJ` matrix
9259a3665c6SJed Brown 
9264165533cSJose E. Roman   Output Parameter:
927*d8a51d2aSBarry Smith . B - matrix on subcommunicator, `NULL` on MPI processes that own zero rows of `A`
9289a3665c6SJed Brown 
9299a3665c6SJed Brown   Level: developer
9309a3665c6SJed Brown 
9319a3665c6SJed Brown   Note:
9322ef1f0ffSBarry Smith   The matrix `B` should be destroyed with `MatDestroy()`. The arrays are not copied, so `B` should be destroyed before `A` is destroyed.
9339a3665c6SJed Brown 
934*d8a51d2aSBarry Smith   Developer Note:
935*d8a51d2aSBarry Smith   This function is mostly useful for internal use by mesh partitioning packages, such as ParMETIS that require that every process owns at least one row.
936*d8a51d2aSBarry Smith 
9371cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATMPIADJ`, `MatCreateMPIAdj()`
9389a3665c6SJed Brown @*/
939d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMPIAdjCreateNonemptySubcommMat(Mat A, Mat *B)
940d71ae5a4SJacob Faibussowitsch {
9419a3665c6SJed Brown   PetscFunctionBegin;
9429a3665c6SJed Brown   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
943cac4c232SBarry Smith   PetscUseMethod(A, "MatMPIAdjCreateNonemptySubcommMat_C", (Mat, Mat *), (A, B));
9443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9459a3665c6SJed Brown }
9469a3665c6SJed Brown 
9470bad9183SKris Buschelman /*MC
948fafad747SKris Buschelman   MATMPIADJ - MATMPIADJ = "mpiadj" - A matrix type to be used for distributed adjacency matrices,
9490bad9183SKris Buschelman   intended for use constructing orderings and partitionings.
9500bad9183SKris Buschelman 
9510bad9183SKris Buschelman   Level: beginner
9520bad9183SKris Buschelman 
95311a5261eSBarry Smith   Note:
95411a5261eSBarry Smith   You can provide values to the matrix using `MatMPIAdjSetPreallocation()`, `MatCreateMPIAdj()`, or
95511a5261eSBarry Smith   by calling `MatSetValues()` and `MatAssemblyBegin()` followed by `MatAssemblyEnd()`
9566a09307cSBarry Smith 
9571cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatCreateMPIAdj()`, `MatMPIAdjSetPreallocation()`, `MatSetValues()`
9580bad9183SKris Buschelman M*/
959d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatCreate_MPIAdj(Mat B)
960d71ae5a4SJacob Faibussowitsch {
961273d9f13SBarry Smith   Mat_MPIAdj *b;
962273d9f13SBarry Smith 
963273d9f13SBarry Smith   PetscFunctionBegin;
9644dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&b));
965b0a32e0cSBarry Smith   B->data         = (void *)b;
966aea10558SJacob Faibussowitsch   B->ops[0]       = MatOps_Values;
967273d9f13SBarry Smith   B->assembled    = PETSC_FALSE;
9686a09307cSBarry Smith   B->preallocated = PETSC_TRUE; /* so that MatSetValues() may be used */
969273d9f13SBarry Smith 
9709566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMPIAdjSetPreallocation_C", MatMPIAdjSetPreallocation_MPIAdj));
9719566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMPIAdjCreateNonemptySubcommMat_C", MatMPIAdjCreateNonemptySubcommMat_MPIAdj));
9729566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMPIAdjToSeq_C", MatMPIAdjToSeq_MPIAdj));
973252a1336SBarry Smith   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMPIAdjToSeqRankZero_C", MatMPIAdjToSeqRankZero_MPIAdj));
9749566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)B, MATMPIADJ));
9753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
976273d9f13SBarry Smith }
977273d9f13SBarry Smith 
978a23d5eceSKris Buschelman /*@C
979*d8a51d2aSBarry Smith   MatMPIAdjToSeq - Converts an parallel `MATMPIADJ` matrix to complete `MATMPIADJ` on each process (needed by sequential partitioners)
980b44e7856SBarry Smith 
981d083f849SBarry Smith   Logically Collective
982b44e7856SBarry Smith 
983b44e7856SBarry Smith   Input Parameter:
984b44e7856SBarry Smith . A - the matrix
985b44e7856SBarry Smith 
986b44e7856SBarry Smith   Output Parameter:
987b44e7856SBarry Smith . B - the same matrix on all processes
988b44e7856SBarry Smith 
989b44e7856SBarry Smith   Level: intermediate
990b44e7856SBarry Smith 
9911cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATMPIADJ`, `MatCreate()`, `MatCreateMPIAdj()`, `MatSetValues()`, `MatMPIAdjToSeqRankZero()`
992b44e7856SBarry Smith @*/
993d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMPIAdjToSeq(Mat A, Mat *B)
994d71ae5a4SJacob Faibussowitsch {
995b44e7856SBarry Smith   PetscFunctionBegin;
996cac4c232SBarry Smith   PetscUseMethod(A, "MatMPIAdjToSeq_C", (Mat, Mat *), (A, B));
9973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
998b44e7856SBarry Smith }
999b44e7856SBarry Smith 
1000b44e7856SBarry Smith /*@C
1001*d8a51d2aSBarry Smith   MatMPIAdjToSeqRankZero - Converts an parallel `MATMPIADJ` matrix to complete `MATMPIADJ` on rank zero (needed by sequential partitioners)
1002252a1336SBarry Smith 
1003252a1336SBarry Smith   Logically Collective
1004252a1336SBarry Smith 
1005252a1336SBarry Smith   Input Parameter:
1006252a1336SBarry Smith . A - the matrix
1007252a1336SBarry Smith 
1008252a1336SBarry Smith   Output Parameter:
1009252a1336SBarry Smith . B - the same matrix on rank zero, not set on other ranks
1010252a1336SBarry Smith 
1011252a1336SBarry Smith   Level: intermediate
1012252a1336SBarry Smith 
101311a5261eSBarry Smith   Note:
1014252a1336SBarry Smith   This routine has the advantage on systems with multiple ranks per node since only one copy of the matrix
1015252a1336SBarry Smith   is stored on the first node, instead of the number of ranks copies. This can allow partitioning much larger
1016d5b43468SJose E. Roman   parallel graph sequentially.
1017252a1336SBarry Smith 
10181cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATMPIADJ`, `MatCreate()`, `MatCreateMPIAdj()`, `MatSetValues()`, `MatMPIAdjToSeq()`
1019252a1336SBarry Smith @*/
1020d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMPIAdjToSeqRankZero(Mat A, Mat *B)
1021d71ae5a4SJacob Faibussowitsch {
1022252a1336SBarry Smith   PetscFunctionBegin;
1023252a1336SBarry Smith   PetscUseMethod(A, "MatMPIAdjToSeqRankZero_C", (Mat, Mat *), (A, B));
10243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1025252a1336SBarry Smith }
1026252a1336SBarry Smith 
1027252a1336SBarry Smith /*@C
1028a23d5eceSKris Buschelman   MatMPIAdjSetPreallocation - Sets the array used for storing the matrix elements
1029a23d5eceSKris Buschelman 
1030d083f849SBarry Smith   Logically Collective
1031a23d5eceSKris Buschelman 
1032a23d5eceSKris Buschelman   Input Parameters:
1033fe59aa6dSJacob Faibussowitsch + B      - the matrix
10342ef1f0ffSBarry Smith . i      - the indices into `j` for the start of each row
1035a23d5eceSKris Buschelman . j      - the column indices for each row (sorted for each row).
10362ef1f0ffSBarry Smith            The indices in `i` and `j` start with zero (NOT with one).
10372ef1f0ffSBarry Smith - values - [use `NULL` if not provided] edge weights
1038a23d5eceSKris Buschelman 
1039a23d5eceSKris Buschelman   Level: intermediate
1040a23d5eceSKris Buschelman 
1041*d8a51d2aSBarry Smith   Notes:
1042*d8a51d2aSBarry Smith   The indices in `i` and `j` start with zero (NOT with one).
1043*d8a51d2aSBarry Smith 
1044*d8a51d2aSBarry Smith   You must NOT free the `i`, `values` and `j` arrays yourself. PETSc will free them
1045*d8a51d2aSBarry Smith   when the matrix is destroyed; you must allocate them with `PetscMalloc()`.
1046*d8a51d2aSBarry Smith 
1047*d8a51d2aSBarry Smith   You should not include the matrix diagonal elements.
1048*d8a51d2aSBarry Smith 
1049*d8a51d2aSBarry Smith   If you already have a matrix, you can create its adjacency matrix by a call
1050*d8a51d2aSBarry Smith   to `MatConvert()`, specifying a type of `MATMPIADJ`.
1051*d8a51d2aSBarry Smith 
1052*d8a51d2aSBarry Smith   Possible values for `MatSetOption()` - `MAT_STRUCTURALLY_SYMMETRIC`
1053*d8a51d2aSBarry Smith 
1054*d8a51d2aSBarry Smith   Fortran Note:
1055*d8a51d2aSBarry Smith   From Fortran the indices and values are copied so the array space need not be provided with `PetscMalloc()`.
1056*d8a51d2aSBarry Smith 
1057fe59aa6dSJacob Faibussowitsch .seealso: [](ch_matrices), `Mat`, `MatCreate()`, `MatCreateMPIAdj()`, `MatSetValues()`, `MATMPIADJ`
1058a23d5eceSKris Buschelman @*/
1059d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMPIAdjSetPreallocation(Mat B, PetscInt *i, PetscInt *j, PetscInt *values)
1060d71ae5a4SJacob Faibussowitsch {
1061273d9f13SBarry Smith   PetscFunctionBegin;
1062cac4c232SBarry Smith   PetscTryMethod(B, "MatMPIAdjSetPreallocation_C", (Mat, PetscInt *, PetscInt *, PetscInt *), (B, i, j, values));
10633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1064273d9f13SBarry Smith }
1065273d9f13SBarry Smith 
1066b97cf49bSBarry Smith /*@C
10673eda8832SBarry Smith   MatCreateMPIAdj - Creates a sparse matrix representing an adjacency list.
1068*d8a51d2aSBarry Smith   The matrix need not have numerical values associated with it, it is
1069b97cf49bSBarry Smith   intended for ordering (to reduce bandwidth etc) and partitioning.
1070b97cf49bSBarry Smith 
1071d083f849SBarry Smith   Collective
1072ef5ee4d1SLois Curfman McInnes 
1073b97cf49bSBarry Smith   Input Parameters:
1074c2e958feSBarry Smith + comm   - MPI communicator
10750752156aSBarry Smith . m      - number of local rows
107658ec18a5SLisandro Dalcin . N      - number of global columns
10772ef1f0ffSBarry Smith . i      - the indices into `j` for the start of each row
1078329f5518SBarry Smith . j      - the column indices for each row (sorted for each row).
1079*d8a51d2aSBarry Smith - values - the values, optional, use `NULL` if not provided
1080b97cf49bSBarry Smith 
1081b97cf49bSBarry Smith   Output Parameter:
1082b97cf49bSBarry Smith . A - the matrix
1083b97cf49bSBarry Smith 
108415091d37SBarry Smith   Level: intermediate
108515091d37SBarry Smith 
108695452b02SPatrick Sanan   Notes:
1087fe59aa6dSJacob Faibussowitsch   The indices in `i` and `j` start with zero (NOT with one).
1088fe59aa6dSJacob Faibussowitsch 
10892ef1f0ffSBarry Smith   You must NOT free the `i`, `values` and `j` arrays yourself. PETSc will free them
1090*d8a51d2aSBarry Smith   when the matrix is destroyed; you must allocate them with `PetscMalloc()`.
10916a09307cSBarry Smith 
10926a09307cSBarry Smith   You should not include the matrix diagonals.
1093b97cf49bSBarry Smith 
1094e3f7e1d6Svictorle   If you already have a matrix, you can create its adjacency matrix by a call
109511a5261eSBarry Smith   to `MatConvert()`, specifying a type of `MATMPIADJ`.
1096ededeb1aSvictorle 
109711a5261eSBarry Smith   Possible values for `MatSetOption()` - `MAT_STRUCTURALLY_SYMMETRIC`
1098b97cf49bSBarry Smith 
1099*d8a51d2aSBarry Smith   Fortran Note:
1100*d8a51d2aSBarry Smith   From Fortran the indices and values are copied so the array space need not be provided with `PetscMalloc()`.
1101*d8a51d2aSBarry Smith 
11021cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatCreate()`, `MatConvert()`, `MatGetOrdering()`, `MATMPIADJ`, `MatMPIAdjSetPreallocation()`
1103b97cf49bSBarry Smith @*/
1104d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCreateMPIAdj(MPI_Comm comm, PetscInt m, PetscInt N, PetscInt *i, PetscInt *j, PetscInt *values, Mat *A)
1105d71ae5a4SJacob Faibussowitsch {
1106433994e6SBarry Smith   PetscFunctionBegin;
11079566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, A));
11089566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*A, m, PETSC_DETERMINE, PETSC_DETERMINE, N));
11099566063dSJacob Faibussowitsch   PetscCall(MatSetType(*A, MATMPIADJ));
11109566063dSJacob Faibussowitsch   PetscCall(MatMPIAdjSetPreallocation(*A, i, j, values));
11113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1112b97cf49bSBarry Smith }
1113