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 */
MatCreateSubMatrix_MPIAdj_data(Mat adj,IS irows,IS icols,PetscInt ** sadj_xadj,PetscInt ** sadj_adjncy,PetscInt ** sadj_values)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
MatCreateSubMatrices_MPIAdj_Private(Mat mat,PetscInt n,const IS irow[],const IS icol[],PetscBool subcomm,MatReuse scall,Mat * submat[])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 {
192f4f49eeaSPierre Jolivet 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
MatCreateSubMatricesMPI_MPIAdj(Mat mat,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat * submat[])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
MatCreateSubMatrices_MPIAdj(Mat mat,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat * submat[])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
MatView_MPIAdj_ASCII(Mat A,PetscViewer viewer)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
MatView_MPIAdj(Mat A,PetscViewer viewer)259d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatView_MPIAdj(Mat A, PetscViewer viewer)
260d71ae5a4SJacob Faibussowitsch {
2619f196a02SMartin Diehl PetscBool isascii;
262b97cf49bSBarry Smith
263433994e6SBarry Smith PetscFunctionBegin;
2649f196a02SMartin Diehl PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
2659f196a02SMartin Diehl if (isascii) PetscCall(MatView_MPIAdj_ASCII(A, viewer));
2663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
267b97cf49bSBarry Smith }
268b97cf49bSBarry Smith
MatDestroy_MPIAdj(Mat mat)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
MatSetOption_MPIAdj(Mat A,MatOption op,PetscBool flg)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;
308d71ae5a4SJacob Faibussowitsch default:
309d71ae5a4SJacob Faibussowitsch break;
310b97cf49bSBarry Smith }
3113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
312b97cf49bSBarry Smith }
313b97cf49bSBarry Smith
MatGetRow_MPIAdj(Mat A,PetscInt row,PetscInt * nz,PetscInt ** idx,PetscScalar ** v)314d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatGetRow_MPIAdj(Mat A, PetscInt row, PetscInt *nz, PetscInt **idx, PetscScalar **v)
315d71ae5a4SJacob Faibussowitsch {
3163eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj *)A->data;
317b97cf49bSBarry Smith
318433994e6SBarry Smith PetscFunctionBegin;
319d0f46423SBarry Smith row -= A->rmap->rstart;
320aed4548fSBarry Smith PetscCheck(row >= 0 && row < A->rmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Row out of range");
321b97cf49bSBarry Smith *nz = a->i[row + 1] - a->i[row];
3222b1d8763SJed Brown if (v) {
3232b1d8763SJed Brown PetscInt j;
3242b1d8763SJed Brown if (a->rowvalues_alloc < *nz) {
3259566063dSJacob Faibussowitsch PetscCall(PetscFree(a->rowvalues));
3262b1d8763SJed Brown a->rowvalues_alloc = PetscMax(a->rowvalues_alloc * 2, *nz);
3279566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(a->rowvalues_alloc, &a->rowvalues));
328b97cf49bSBarry Smith }
329ad540459SPierre Jolivet for (j = 0; j < *nz; j++) a->rowvalues[j] = a->values ? a->values[a->i[row] + j] : 1.0;
3302b1d8763SJed Brown *v = (*nz) ? a->rowvalues : NULL;
331b97cf49bSBarry Smith }
33292b6f2a0SJed Brown if (idx) *idx = (*nz) ? a->j + a->i[row] : NULL;
3333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
334b97cf49bSBarry Smith }
335b97cf49bSBarry Smith
MatEqual_MPIAdj(Mat A,Mat B,PetscBool * flg)336d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatEqual_MPIAdj(Mat A, Mat B, PetscBool *flg)
337d71ae5a4SJacob Faibussowitsch {
3383eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj *)A->data, *b = (Mat_MPIAdj *)B->data;
339ace3abfcSBarry Smith PetscBool flag;
340b97cf49bSBarry Smith
341433994e6SBarry Smith PetscFunctionBegin;
342b97cf49bSBarry Smith /* If the matrix dimensions are not equal,or no of nonzeros */
343ad540459SPierre Jolivet if ((A->rmap->n != B->rmap->n) || (a->nz != b->nz)) flag = PETSC_FALSE;
344b97cf49bSBarry Smith
345b97cf49bSBarry Smith /* if the a->i are the same */
3469566063dSJacob Faibussowitsch PetscCall(PetscArraycmp(a->i, b->i, A->rmap->n + 1, &flag));
347b97cf49bSBarry Smith
348b97cf49bSBarry Smith /* if a->j are the same */
349418fb43bSPierre Jolivet PetscCall(PetscArraycmp(a->j, b->j, a->nz, &flag));
350b97cf49bSBarry Smith
3515440e5dcSBarry Smith PetscCallMPI(MPIU_Allreduce(&flag, flg, 1, MPI_C_BOOL, MPI_LAND, PetscObjectComm((PetscObject)A)));
3523ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
353b97cf49bSBarry Smith }
354b97cf49bSBarry Smith
MatGetRowIJ_MPIAdj(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool blockcompressed,PetscInt * m,const PetscInt * inia[],const PetscInt * inja[],PetscBool * done)355d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatGetRowIJ_MPIAdj(Mat A, PetscInt oshift, PetscBool symmetric, PetscBool blockcompressed, PetscInt *m, const PetscInt *inia[], const PetscInt *inja[], PetscBool *done)
356d71ae5a4SJacob Faibussowitsch {
357b24ad042SBarry Smith PetscInt i;
3586cd91f64SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj *)A->data;
3591a83f524SJed Brown PetscInt **ia = (PetscInt **)inia, **ja = (PetscInt **)inja;
3606cd91f64SBarry Smith
3616cd91f64SBarry Smith PetscFunctionBegin;
362d0f46423SBarry Smith *m = A->rmap->n;
3636cd91f64SBarry Smith *ia = a->i;
3646cd91f64SBarry Smith *ja = a->j;
3656cd91f64SBarry Smith *done = PETSC_TRUE;
366d42735a1SBarry Smith if (oshift) {
367ad540459SPierre Jolivet for (i = 0; i < (*ia)[*m]; i++) (*ja)[i]++;
368d42735a1SBarry Smith for (i = 0; i <= (*m); i++) (*ia)[i]++;
369d42735a1SBarry Smith }
3703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
371d42735a1SBarry Smith }
372d42735a1SBarry Smith
MatRestoreRowIJ_MPIAdj(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool blockcompressed,PetscInt * m,const PetscInt * inia[],const PetscInt * inja[],PetscBool * done)373d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatRestoreRowIJ_MPIAdj(Mat A, PetscInt oshift, PetscBool symmetric, PetscBool blockcompressed, PetscInt *m, const PetscInt *inia[], const PetscInt *inja[], PetscBool *done)
374d71ae5a4SJacob Faibussowitsch {
375b24ad042SBarry Smith PetscInt i;
376d42735a1SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj *)A->data;
3771a83f524SJed Brown PetscInt **ia = (PetscInt **)inia, **ja = (PetscInt **)inja;
378d42735a1SBarry Smith
379d42735a1SBarry Smith PetscFunctionBegin;
380aed4548fSBarry Smith PetscCheck(!ia || a->i == *ia, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "ia passed back is not one obtained with MatGetRowIJ()");
381aed4548fSBarry Smith PetscCheck(!ja || a->j == *ja, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "ja passed back is not one obtained with MatGetRowIJ()");
382d42735a1SBarry Smith if (oshift) {
38328b400f6SJacob Faibussowitsch PetscCheck(ia, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "If oshift then you must passed in inia[] argument");
38428b400f6SJacob Faibussowitsch PetscCheck(ja, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "If oshift then you must passed in inja[] argument");
385d42735a1SBarry Smith for (i = 0; i <= (*m); i++) (*ia)[i]--;
386ad540459SPierre Jolivet for (i = 0; i < (*ia)[*m]; i++) (*ja)[i]--;
387d42735a1SBarry Smith }
3883ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
3896cd91f64SBarry Smith }
390b97cf49bSBarry Smith
MatConvertFrom_MPIAdj(Mat A,MatType type,MatReuse reuse,Mat * newmat)39166976f2fSJacob Faibussowitsch static PetscErrorCode MatConvertFrom_MPIAdj(Mat A, MatType type, MatReuse reuse, Mat *newmat)
392d71ae5a4SJacob Faibussowitsch {
39317667f90SBarry Smith Mat B;
39417667f90SBarry Smith PetscInt i, m, N, nzeros = 0, *ia, *ja, len, rstart, cnt, j, *a;
39517667f90SBarry Smith const PetscInt *rj;
39617667f90SBarry Smith const PetscScalar *ra;
39717667f90SBarry Smith MPI_Comm comm;
39817667f90SBarry Smith
39917667f90SBarry Smith PetscFunctionBegin;
4009566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, NULL, &N));
4019566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A, &m, NULL));
4029566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A, &rstart, NULL));
40317667f90SBarry Smith
40417667f90SBarry Smith /* count the number of nonzeros per row */
40517667f90SBarry Smith for (i = 0; i < m; i++) {
4069566063dSJacob Faibussowitsch PetscCall(MatGetRow(A, i + rstart, &len, &rj, NULL));
4075ee9ba1cSJed Brown for (j = 0; j < len; j++) {
4089371c9d4SSatish Balay if (rj[j] == i + rstart) {
4099371c9d4SSatish Balay len--;
4109371c9d4SSatish Balay break;
4119371c9d4SSatish Balay } /* don't count diagonal */
4125ee9ba1cSJed Brown }
41317667f90SBarry Smith nzeros += len;
4149566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(A, i + rstart, &len, &rj, NULL));
41517667f90SBarry Smith }
41617667f90SBarry Smith
41717667f90SBarry Smith /* malloc space for nonzeros */
4189566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nzeros + 1, &a));
4199566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(N + 1, &ia));
4209566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nzeros + 1, &ja));
42117667f90SBarry Smith
42217667f90SBarry Smith nzeros = 0;
42317667f90SBarry Smith ia[0] = 0;
42417667f90SBarry Smith for (i = 0; i < m; i++) {
4259566063dSJacob Faibussowitsch PetscCall(MatGetRow(A, i + rstart, &len, &rj, &ra));
42617667f90SBarry Smith cnt = 0;
42717667f90SBarry Smith for (j = 0; j < len; j++) {
4285ee9ba1cSJed Brown if (rj[j] != i + rstart) { /* if not diagonal */
42917667f90SBarry Smith a[nzeros + cnt] = (PetscInt)PetscAbsScalar(ra[j]);
43017667f90SBarry Smith ja[nzeros + cnt++] = rj[j];
43117667f90SBarry Smith }
4325ee9ba1cSJed Brown }
4339566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(A, i + rstart, &len, &rj, &ra));
43417667f90SBarry Smith nzeros += cnt;
43517667f90SBarry Smith ia[i + 1] = nzeros;
43617667f90SBarry Smith }
43717667f90SBarry Smith
4389566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)A, &comm));
4399566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, &B));
4409566063dSJacob Faibussowitsch PetscCall(MatSetSizes(B, m, PETSC_DETERMINE, PETSC_DETERMINE, N));
4419566063dSJacob Faibussowitsch PetscCall(MatSetType(B, type));
4429566063dSJacob Faibussowitsch PetscCall(MatMPIAdjSetPreallocation(B, ia, ja, a));
44317667f90SBarry Smith
444511c6705SHong Zhang if (reuse == MAT_INPLACE_MATRIX) {
4459566063dSJacob Faibussowitsch PetscCall(MatHeaderReplace(A, &B));
44617667f90SBarry Smith } else {
44717667f90SBarry Smith *newmat = B;
44817667f90SBarry Smith }
4493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
45017667f90SBarry Smith }
45117667f90SBarry Smith
MatSetValues_MPIAdj(Mat A,PetscInt m,const PetscInt * rows,PetscInt n,const PetscInt * cols,const PetscScalar * values,InsertMode im)45266976f2fSJacob Faibussowitsch static PetscErrorCode MatSetValues_MPIAdj(Mat A, PetscInt m, const PetscInt *rows, PetscInt n, const PetscInt *cols, const PetscScalar *values, InsertMode im)
453d71ae5a4SJacob Faibussowitsch {
4546a09307cSBarry Smith Mat_MPIAdj *adj = (Mat_MPIAdj *)A->data;
4556a09307cSBarry Smith PetscInt rStart, rEnd, cStart, cEnd;
4566a09307cSBarry Smith
4576a09307cSBarry Smith PetscFunctionBegin;
4586a09307cSBarry Smith PetscCheck(!adj->i, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Matrix is already assembled, cannot change its values");
4596a09307cSBarry Smith PetscCall(MatGetOwnershipRange(A, &rStart, &rEnd));
4606a09307cSBarry Smith PetscCall(MatGetOwnershipRangeColumn(A, &cStart, &cEnd));
4616a09307cSBarry Smith if (!adj->ht) {
4626a09307cSBarry Smith PetscCall(PetscHSetIJCreate(&adj->ht));
4636a09307cSBarry Smith PetscCall(MatStashCreate_Private(PetscObjectComm((PetscObject)A), 1, &A->stash));
4646a09307cSBarry Smith PetscCall(PetscLayoutSetUp(A->rmap));
4656a09307cSBarry Smith PetscCall(PetscLayoutSetUp(A->cmap));
4666a09307cSBarry Smith }
4676a09307cSBarry Smith for (PetscInt r = 0; r < m; ++r) {
4686a09307cSBarry Smith PetscHashIJKey key;
4696a09307cSBarry Smith
4706a09307cSBarry Smith key.i = rows[r];
4716a09307cSBarry Smith if (key.i < 0) continue;
4726a09307cSBarry Smith if ((key.i < rStart) || (key.i >= rEnd)) {
4736a09307cSBarry Smith PetscCall(MatStashValuesRow_Private(&A->stash, key.i, n, cols, values, PETSC_FALSE));
4746a09307cSBarry Smith } else {
4756a09307cSBarry Smith for (PetscInt c = 0; c < n; ++c) {
4766a09307cSBarry Smith key.j = cols[c];
4776a09307cSBarry Smith if (key.j < 0 || key.i == key.j) continue;
4786a09307cSBarry Smith PetscCall(PetscHSetIJAdd(adj->ht, key));
4796a09307cSBarry Smith }
4806a09307cSBarry Smith }
4816a09307cSBarry Smith }
4823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
4836a09307cSBarry Smith }
4846a09307cSBarry Smith
MatAssemblyBegin_MPIAdj(Mat A,MatAssemblyType type)48566976f2fSJacob Faibussowitsch static PetscErrorCode MatAssemblyBegin_MPIAdj(Mat A, MatAssemblyType type)
486d71ae5a4SJacob Faibussowitsch {
4876a09307cSBarry Smith PetscInt nstash, reallocs;
4886a09307cSBarry Smith Mat_MPIAdj *adj = (Mat_MPIAdj *)A->data;
4896a09307cSBarry Smith
4906a09307cSBarry Smith PetscFunctionBegin;
4916a09307cSBarry Smith if (!adj->ht) {
4926a09307cSBarry Smith PetscCall(PetscHSetIJCreate(&adj->ht));
4936a09307cSBarry Smith PetscCall(PetscLayoutSetUp(A->rmap));
4946a09307cSBarry Smith PetscCall(PetscLayoutSetUp(A->cmap));
4956a09307cSBarry Smith }
4966a09307cSBarry Smith PetscCall(MatStashScatterBegin_Private(A, &A->stash, A->rmap->range));
4976a09307cSBarry Smith PetscCall(MatStashGetInfo_Private(&A->stash, &nstash, &reallocs));
4986a09307cSBarry Smith PetscCall(PetscInfo(A, "Stash has %" PetscInt_FMT " entries, uses %" PetscInt_FMT " mallocs.\n", nstash, reallocs));
4993ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
5006a09307cSBarry Smith }
5016a09307cSBarry Smith
MatAssemblyEnd_MPIAdj(Mat A,MatAssemblyType type)50266976f2fSJacob Faibussowitsch static PetscErrorCode MatAssemblyEnd_MPIAdj(Mat A, MatAssemblyType type)
503d71ae5a4SJacob Faibussowitsch {
5046a09307cSBarry Smith PetscScalar *val;
5056a09307cSBarry Smith PetscInt *row, *col, m, rstart, *rowstarts;
5066a09307cSBarry Smith PetscInt i, j, ncols, flg, nz;
5076a09307cSBarry Smith PetscMPIInt n;
5086a09307cSBarry Smith Mat_MPIAdj *adj = (Mat_MPIAdj *)A->data;
5096a09307cSBarry Smith PetscHashIter hi;
5106a09307cSBarry Smith PetscHashIJKey key;
5116a09307cSBarry Smith PetscHSetIJ ht = adj->ht;
5126a09307cSBarry Smith
5136a09307cSBarry Smith PetscFunctionBegin;
5146a09307cSBarry Smith while (1) {
5156a09307cSBarry Smith PetscCall(MatStashScatterGetMesg_Private(&A->stash, &n, &row, &col, &val, &flg));
5166a09307cSBarry Smith if (!flg) break;
5176a09307cSBarry Smith
5186a09307cSBarry Smith for (i = 0; i < n;) {
5196a09307cSBarry Smith /* Identify the consecutive vals belonging to the same row */
5206a09307cSBarry Smith for (j = i, rstart = row[j]; j < n; j++) {
5216a09307cSBarry Smith if (row[j] != rstart) break;
5226a09307cSBarry Smith }
5236a09307cSBarry Smith if (j < n) ncols = j - i;
5246a09307cSBarry Smith else ncols = n - i;
5256a09307cSBarry Smith /* Set all these values with a single function call */
5266a09307cSBarry Smith PetscCall(MatSetValues_MPIAdj(A, 1, row + i, ncols, col + i, val + i, INSERT_VALUES));
5276a09307cSBarry Smith i = j;
5286a09307cSBarry Smith }
5296a09307cSBarry Smith }
5306a09307cSBarry Smith PetscCall(MatStashScatterEnd_Private(&A->stash));
5316a09307cSBarry Smith PetscCall(MatStashDestroy_Private(&A->stash));
5326a09307cSBarry Smith
5336a09307cSBarry Smith PetscCall(MatGetLocalSize(A, &m, NULL));
5346a09307cSBarry Smith PetscCall(MatGetOwnershipRange(A, &rstart, NULL));
5356a09307cSBarry Smith PetscCall(PetscCalloc1(m + 1, &rowstarts));
5366a09307cSBarry Smith PetscHashIterBegin(ht, hi);
5376a09307cSBarry Smith for (; !PetscHashIterAtEnd(ht, hi);) {
5386a09307cSBarry Smith PetscHashIterGetKey(ht, hi, key);
5396a09307cSBarry Smith rowstarts[key.i - rstart + 1]++;
5406a09307cSBarry Smith PetscHashIterNext(ht, hi);
5416a09307cSBarry Smith }
542ad540459SPierre Jolivet for (i = 1; i < m + 1; i++) rowstarts[i] = rowstarts[i - 1] + rowstarts[i];
5436a09307cSBarry Smith
5446a09307cSBarry Smith PetscCall(PetscHSetIJGetSize(ht, &nz));
5456a09307cSBarry Smith PetscCall(PetscMalloc1(nz, &col));
5466a09307cSBarry Smith PetscHashIterBegin(ht, hi);
5476a09307cSBarry Smith for (; !PetscHashIterAtEnd(ht, hi);) {
5486a09307cSBarry Smith PetscHashIterGetKey(ht, hi, key);
5496a09307cSBarry Smith col[rowstarts[key.i - rstart]++] = key.j;
5506a09307cSBarry Smith PetscHashIterNext(ht, hi);
5516a09307cSBarry Smith }
5526a09307cSBarry Smith PetscCall(PetscHSetIJDestroy(&ht));
553ad540459SPierre Jolivet for (i = m; i > 0; i--) rowstarts[i] = rowstarts[i - 1];
5546a09307cSBarry Smith rowstarts[0] = 0;
5556a09307cSBarry Smith
55648a46eb9SPierre Jolivet for (PetscInt i = 0; i < m; i++) PetscCall(PetscSortInt(rowstarts[i + 1] - rowstarts[i], &col[rowstarts[i]]));
5576a09307cSBarry Smith
5586a09307cSBarry Smith adj->i = rowstarts;
5596a09307cSBarry Smith adj->j = col;
560252a1336SBarry Smith adj->nz = rowstarts[m];
5616a09307cSBarry Smith adj->freeaij = PETSC_TRUE;
5623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
5636a09307cSBarry Smith }
5646a09307cSBarry Smith
5656a09307cSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_MPIAdj,
5663eda8832SBarry Smith MatGetRow_MPIAdj,
567c67b119cSPierre Jolivet NULL,
568f4259b30SLisandro Dalcin NULL,
569f4259b30SLisandro Dalcin /* 4*/ NULL,
570f4259b30SLisandro Dalcin NULL,
571f4259b30SLisandro Dalcin NULL,
572f4259b30SLisandro Dalcin NULL,
573f4259b30SLisandro Dalcin NULL,
574f4259b30SLisandro Dalcin NULL,
575f4259b30SLisandro Dalcin /*10*/ NULL,
576f4259b30SLisandro Dalcin NULL,
577f4259b30SLisandro Dalcin NULL,
578f4259b30SLisandro Dalcin NULL,
579f4259b30SLisandro Dalcin NULL,
580f4259b30SLisandro Dalcin /*15*/ NULL,
5813eda8832SBarry Smith MatEqual_MPIAdj,
582f4259b30SLisandro Dalcin NULL,
583f4259b30SLisandro Dalcin NULL,
584f4259b30SLisandro Dalcin NULL,
5856a09307cSBarry Smith /*20*/ MatAssemblyBegin_MPIAdj,
5866a09307cSBarry Smith MatAssemblyEnd_MPIAdj,
5873eda8832SBarry Smith MatSetOption_MPIAdj,
588f4259b30SLisandro Dalcin NULL,
589f4259b30SLisandro Dalcin /*24*/ NULL,
590f4259b30SLisandro Dalcin NULL,
591f4259b30SLisandro Dalcin NULL,
592f4259b30SLisandro Dalcin NULL,
593f4259b30SLisandro Dalcin NULL,
594f4259b30SLisandro Dalcin /*29*/ NULL,
595f4259b30SLisandro Dalcin NULL,
596f4259b30SLisandro Dalcin NULL,
597f4259b30SLisandro Dalcin NULL,
598f4259b30SLisandro Dalcin NULL,
599f4259b30SLisandro Dalcin /*34*/ NULL,
600f4259b30SLisandro Dalcin NULL,
601f4259b30SLisandro Dalcin NULL,
602f4259b30SLisandro Dalcin NULL,
603f4259b30SLisandro Dalcin NULL,
604f4259b30SLisandro Dalcin /*39*/ NULL,
6057dae84e0SHong Zhang MatCreateSubMatrices_MPIAdj,
606f4259b30SLisandro Dalcin NULL,
607f4259b30SLisandro Dalcin NULL,
608f4259b30SLisandro Dalcin NULL,
609f4259b30SLisandro Dalcin /*44*/ NULL,
610f4259b30SLisandro Dalcin NULL,
6117d68702bSBarry Smith MatShift_Basic,
612f4259b30SLisandro Dalcin NULL,
613f4259b30SLisandro Dalcin NULL,
614f4259b30SLisandro Dalcin /*49*/ NULL,
6156cd91f64SBarry Smith MatGetRowIJ_MPIAdj,
616d42735a1SBarry Smith MatRestoreRowIJ_MPIAdj,
617f4259b30SLisandro Dalcin NULL,
618f4259b30SLisandro Dalcin NULL,
619f4259b30SLisandro Dalcin /*54*/ NULL,
620f4259b30SLisandro Dalcin NULL,
621f4259b30SLisandro Dalcin NULL,
622f4259b30SLisandro Dalcin NULL,
623f4259b30SLisandro Dalcin NULL,
624f4259b30SLisandro Dalcin /*59*/ NULL,
625b9b97703SBarry Smith MatDestroy_MPIAdj,
626b9b97703SBarry Smith MatView_MPIAdj,
62717667f90SBarry Smith MatConvertFrom_MPIAdj,
628f4259b30SLisandro Dalcin NULL,
629f4259b30SLisandro Dalcin /*64*/ NULL,
630f4259b30SLisandro Dalcin NULL,
631f4259b30SLisandro Dalcin NULL,
632f4259b30SLisandro Dalcin NULL,
633f4259b30SLisandro Dalcin NULL,
634f4259b30SLisandro Dalcin /*69*/ NULL,
635f4259b30SLisandro Dalcin NULL,
636f4259b30SLisandro Dalcin NULL,
637f4259b30SLisandro Dalcin NULL,
638f4259b30SLisandro Dalcin NULL,
639f4259b30SLisandro Dalcin /*74*/ NULL,
640f4259b30SLisandro Dalcin NULL,
641f4259b30SLisandro Dalcin NULL,
642f4259b30SLisandro Dalcin NULL,
643f4259b30SLisandro Dalcin NULL,
644f4259b30SLisandro Dalcin /*79*/ NULL,
645f4259b30SLisandro Dalcin NULL,
646f4259b30SLisandro Dalcin NULL,
647f4259b30SLisandro Dalcin NULL,
648f4259b30SLisandro Dalcin NULL,
649f4259b30SLisandro Dalcin /*84*/ NULL,
650f4259b30SLisandro Dalcin NULL,
651f4259b30SLisandro Dalcin NULL,
652f4259b30SLisandro Dalcin NULL,
653f4259b30SLisandro Dalcin NULL,
654f4259b30SLisandro Dalcin /*89*/ NULL,
655f4259b30SLisandro Dalcin NULL,
656f4259b30SLisandro Dalcin NULL,
657f4259b30SLisandro Dalcin NULL,
658f4259b30SLisandro Dalcin NULL,
659f4259b30SLisandro Dalcin /*94*/ NULL,
660f4259b30SLisandro Dalcin NULL,
661f4259b30SLisandro Dalcin NULL,
662f4259b30SLisandro Dalcin NULL,
663f4259b30SLisandro Dalcin NULL,
664f4259b30SLisandro Dalcin /*99*/ NULL,
665f4259b30SLisandro Dalcin NULL,
666f4259b30SLisandro Dalcin NULL,
667f4259b30SLisandro Dalcin NULL,
668f4259b30SLisandro Dalcin NULL,
669f4259b30SLisandro Dalcin /*104*/ NULL,
670f4259b30SLisandro Dalcin NULL,
671f4259b30SLisandro Dalcin NULL,
672f4259b30SLisandro Dalcin NULL,
673f4259b30SLisandro Dalcin NULL,
674f4259b30SLisandro Dalcin /*109*/ NULL,
675f4259b30SLisandro Dalcin NULL,
676f4259b30SLisandro Dalcin NULL,
677f4259b30SLisandro Dalcin NULL,
678f4259b30SLisandro Dalcin NULL,
679f4259b30SLisandro Dalcin /*114*/ NULL,
680f4259b30SLisandro Dalcin NULL,
681f4259b30SLisandro Dalcin NULL,
682f4259b30SLisandro Dalcin NULL,
683421480d9SBarry Smith MatCreateSubMatricesMPI_MPIAdj,
684421480d9SBarry Smith /*119*/ NULL,
685f4259b30SLisandro Dalcin NULL,
686f4259b30SLisandro Dalcin NULL,
687f4259b30SLisandro Dalcin NULL,
688f4259b30SLisandro Dalcin NULL,
689f4259b30SLisandro Dalcin /*124*/ NULL,
690f4259b30SLisandro Dalcin NULL,
691f4259b30SLisandro Dalcin NULL,
692f4259b30SLisandro Dalcin NULL,
6938bb0f5c6SPierre Jolivet NULL,
694f4259b30SLisandro Dalcin /*129*/ NULL,
695f4259b30SLisandro Dalcin NULL,
696f4259b30SLisandro Dalcin NULL,
697f4259b30SLisandro Dalcin NULL,
698f4259b30SLisandro Dalcin NULL,
699f4259b30SLisandro Dalcin /*134*/ NULL,
700f4259b30SLisandro Dalcin NULL,
701f4259b30SLisandro Dalcin NULL,
702f4259b30SLisandro Dalcin NULL,
703f4259b30SLisandro Dalcin NULL,
704f4259b30SLisandro Dalcin /*139*/ NULL,
705f4259b30SLisandro Dalcin NULL,
706d70f29a3SPierre Jolivet NULL,
70703db1824SAlex Lindsay NULL,
708dec0b466SHong Zhang NULL};
709b97cf49bSBarry Smith
MatMPIAdjSetPreallocation_MPIAdj(Mat B,PetscInt * i,PetscInt * j,PetscInt * values)710d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMPIAdjSetPreallocation_MPIAdj(Mat B, PetscInt *i, PetscInt *j, PetscInt *values)
711d71ae5a4SJacob Faibussowitsch {
712a23d5eceSKris Buschelman Mat_MPIAdj *b = (Mat_MPIAdj *)B->data;
713e895ccc0SFande Kong PetscBool useedgeweights;
714a23d5eceSKris Buschelman
715a23d5eceSKris Buschelman PetscFunctionBegin;
7169566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(B->rmap));
7179566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(B->cmap));
7189371c9d4SSatish Balay if (values) useedgeweights = PETSC_TRUE;
7199371c9d4SSatish Balay else useedgeweights = PETSC_FALSE;
720e895ccc0SFande Kong /* Make everybody knows if they are using edge weights or not */
7215440e5dcSBarry Smith PetscCallMPI(MPIU_Allreduce(&useedgeweights, &b->useedgeweights, 1, MPI_C_BOOL, MPI_LOR, PetscObjectComm((PetscObject)B)));
72258ec18a5SLisandro Dalcin
72376bd3646SJed Brown if (PetscDefined(USE_DEBUG)) {
72476bd3646SJed Brown PetscInt ii;
72576bd3646SJed Brown
72608401ef6SPierre 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]);
727d0f46423SBarry Smith for (ii = 1; ii < B->rmap->n; ii++) {
728aed4548fSBarry 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]);
729a23d5eceSKris Buschelman }
730ad540459SPierre 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]);
73176bd3646SJed Brown }
732a23d5eceSKris Buschelman b->j = j;
733a23d5eceSKris Buschelman b->i = i;
734a23d5eceSKris Buschelman b->values = values;
735a23d5eceSKris Buschelman
736d0f46423SBarry Smith b->nz = i[B->rmap->n];
737f4259b30SLisandro Dalcin b->diag = NULL;
738a23d5eceSKris Buschelman b->symmetric = PETSC_FALSE;
739a23d5eceSKris Buschelman b->freeaij = PETSC_TRUE;
740a23d5eceSKris Buschelman
7416a09307cSBarry Smith B->ops->assemblybegin = NULL;
7426a09307cSBarry Smith B->ops->assemblyend = NULL;
7439566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
7449566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
7456a09307cSBarry Smith PetscCall(MatStashDestroy_Private(&B->stash));
7463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
747a23d5eceSKris Buschelman }
748b97cf49bSBarry Smith
MatMPIAdjCreateNonemptySubcommMat_MPIAdj(Mat A,Mat * B)749d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMPIAdjCreateNonemptySubcommMat_MPIAdj(Mat A, Mat *B)
750d71ae5a4SJacob Faibussowitsch {
7519a3665c6SJed Brown Mat_MPIAdj *a = (Mat_MPIAdj *)A->data;
7529a3665c6SJed Brown const PetscInt *ranges;
7539a3665c6SJed Brown MPI_Comm acomm, bcomm;
7549a3665c6SJed Brown MPI_Group agroup, bgroup;
755*badd099fSHan Liu PetscMPIInt i, size, nranks, *ranks;
7569a3665c6SJed Brown
7579a3665c6SJed Brown PetscFunctionBegin;
7580298fd71SBarry Smith *B = NULL;
7599566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)A, &acomm));
7609566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(acomm, &size));
7619566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRanges(A, &ranges));
7629a3665c6SJed Brown for (i = 0, nranks = 0; i < size; i++) {
7639a3665c6SJed Brown if (ranges[i + 1] - ranges[i] > 0) nranks++;
7649a3665c6SJed Brown }
7659a3665c6SJed Brown if (nranks == size) { /* All ranks have a positive number of rows, so we do not need to create a subcomm; */
7669566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)A));
7679a3665c6SJed Brown *B = A;
7683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
7699a3665c6SJed Brown }
7709a3665c6SJed Brown
7719566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nranks, &ranks));
7729a3665c6SJed Brown for (i = 0, nranks = 0; i < size; i++) {
7739a3665c6SJed Brown if (ranges[i + 1] - ranges[i] > 0) ranks[nranks++] = i;
7749a3665c6SJed Brown }
7759566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_group(acomm, &agroup));
7769566063dSJacob Faibussowitsch PetscCallMPI(MPI_Group_incl(agroup, nranks, ranks, &bgroup));
7779566063dSJacob Faibussowitsch PetscCall(PetscFree(ranks));
7789566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_create(acomm, bgroup, &bcomm));
7799566063dSJacob Faibussowitsch PetscCallMPI(MPI_Group_free(&agroup));
7809566063dSJacob Faibussowitsch PetscCallMPI(MPI_Group_free(&bgroup));
7819a3665c6SJed Brown if (bcomm != MPI_COMM_NULL) {
7829a3665c6SJed Brown PetscInt m, N;
7839a3665c6SJed Brown Mat_MPIAdj *b;
7849566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A, &m, NULL));
7859566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, NULL, &N));
7869566063dSJacob Faibussowitsch PetscCall(MatCreateMPIAdj(bcomm, m, N, a->i, a->j, a->values, B));
7879a3665c6SJed Brown b = (Mat_MPIAdj *)(*B)->data;
7889a3665c6SJed Brown b->freeaij = PETSC_FALSE;
7899566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free(&bcomm));
7909a3665c6SJed Brown }
7913ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
7929a3665c6SJed Brown }
7939a3665c6SJed Brown
MatMPIAdjToSeq_MPIAdj(Mat A,Mat * B)79466976f2fSJacob Faibussowitsch static PetscErrorCode MatMPIAdjToSeq_MPIAdj(Mat A, Mat *B)
795d71ae5a4SJacob Faibussowitsch {
796b44e7856SBarry Smith PetscInt M, N, *II, *J, NZ, nz, m, nzstart, i;
797b44e7856SBarry Smith PetscInt *Values = NULL;
798b44e7856SBarry Smith Mat_MPIAdj *adj = (Mat_MPIAdj *)A->data;
799b44e7856SBarry Smith PetscMPIInt mnz, mm, *allnz, *allm, size, *dispnz, *dispm;
800b44e7856SBarry Smith
801b44e7856SBarry Smith PetscFunctionBegin;
8029566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
8039566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &M, &N));
8049566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A, &m, NULL));
805b44e7856SBarry Smith nz = adj->nz;
80608401ef6SPierre 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]);
807462c564dSBarry Smith PetscCallMPI(MPIU_Allreduce(&nz, &NZ, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)A)));
808d4e69552SBarry Smith
8099566063dSJacob Faibussowitsch PetscCall(PetscMPIIntCast(nz, &mnz));
8109566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(size, &allnz, size, &dispnz));
8119566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&mnz, 1, MPI_INT, allnz, 1, MPI_INT, PetscObjectComm((PetscObject)A)));
8129371c9d4SSatish Balay dispnz[0] = 0;
8139371c9d4SSatish Balay for (i = 1; i < size; i++) dispnz[i] = dispnz[i - 1] + allnz[i - 1];
814b44e7856SBarry Smith if (adj->values) {
8159566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(NZ, &Values));
8169566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgatherv(adj->values, mnz, MPIU_INT, Values, allnz, dispnz, MPIU_INT, PetscObjectComm((PetscObject)A)));
817b44e7856SBarry Smith }
8189566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(NZ, &J));
8199566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgatherv(adj->j, mnz, MPIU_INT, J, allnz, dispnz, MPIU_INT, PetscObjectComm((PetscObject)A)));
8209566063dSJacob Faibussowitsch PetscCall(PetscFree2(allnz, dispnz));
8219566063dSJacob Faibussowitsch PetscCallMPI(MPI_Scan(&nz, &nzstart, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)A)));
822b44e7856SBarry Smith nzstart -= nz;
823b44e7856SBarry Smith /* shift the i[] values so they will be correct after being received */
824b44e7856SBarry Smith for (i = 0; i < m; i++) adj->i[i] += nzstart;
8259566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(M + 1, &II));
8269566063dSJacob Faibussowitsch PetscCall(PetscMPIIntCast(m, &mm));
8279566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(size, &allm, size, &dispm));
8289566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&mm, 1, MPI_INT, allm, 1, MPI_INT, PetscObjectComm((PetscObject)A)));
8299371c9d4SSatish Balay dispm[0] = 0;
8309371c9d4SSatish Balay for (i = 1; i < size; i++) dispm[i] = dispm[i - 1] + allm[i - 1];
8319566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgatherv(adj->i, mm, MPIU_INT, II, allm, dispm, MPIU_INT, PetscObjectComm((PetscObject)A)));
8329566063dSJacob Faibussowitsch PetscCall(PetscFree2(allm, dispm));
833b44e7856SBarry Smith II[M] = NZ;
834b44e7856SBarry Smith /* shift the i[] values back */
835b44e7856SBarry Smith for (i = 0; i < m; i++) adj->i[i] -= nzstart;
8369566063dSJacob Faibussowitsch PetscCall(MatCreateMPIAdj(PETSC_COMM_SELF, M, N, II, J, Values, B));
8373ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
838b44e7856SBarry Smith }
839b44e7856SBarry Smith
MatMPIAdjToSeqRankZero_MPIAdj(Mat A,Mat * B)84066976f2fSJacob Faibussowitsch static PetscErrorCode MatMPIAdjToSeqRankZero_MPIAdj(Mat A, Mat *B)
841d71ae5a4SJacob Faibussowitsch {
842252a1336SBarry Smith PetscInt M, N, *II, *J, NZ, nz, m, nzstart, i;
843252a1336SBarry Smith PetscInt *Values = NULL;
844252a1336SBarry Smith Mat_MPIAdj *adj = (Mat_MPIAdj *)A->data;
845252a1336SBarry Smith PetscMPIInt mnz, mm, *allnz = NULL, *allm, size, *dispnz, *dispm, rank;
846252a1336SBarry Smith
847252a1336SBarry Smith PetscFunctionBegin;
848252a1336SBarry Smith PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
849252a1336SBarry Smith PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)A), &rank));
850252a1336SBarry Smith PetscCall(MatGetSize(A, &M, &N));
851252a1336SBarry Smith PetscCall(MatGetLocalSize(A, &m, NULL));
852252a1336SBarry Smith nz = adj->nz;
853252a1336SBarry 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]);
854462c564dSBarry Smith PetscCallMPI(MPIU_Allreduce(&nz, &NZ, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)A)));
855252a1336SBarry Smith
856252a1336SBarry Smith PetscCall(PetscMPIIntCast(nz, &mnz));
857252a1336SBarry Smith if (!rank) PetscCall(PetscMalloc2(size, &allnz, size, &dispnz));
858252a1336SBarry Smith PetscCallMPI(MPI_Gather(&mnz, 1, MPI_INT, allnz, 1, MPI_INT, 0, PetscObjectComm((PetscObject)A)));
859252a1336SBarry Smith if (!rank) {
8609371c9d4SSatish Balay dispnz[0] = 0;
8619371c9d4SSatish Balay for (i = 1; i < size; i++) dispnz[i] = dispnz[i - 1] + allnz[i - 1];
862252a1336SBarry Smith if (adj->values) {
863252a1336SBarry Smith PetscCall(PetscMalloc1(NZ, &Values));
864252a1336SBarry Smith PetscCallMPI(MPI_Gatherv(adj->values, mnz, MPIU_INT, Values, allnz, dispnz, MPIU_INT, 0, PetscObjectComm((PetscObject)A)));
865252a1336SBarry Smith }
866252a1336SBarry Smith PetscCall(PetscMalloc1(NZ, &J));
867252a1336SBarry Smith PetscCallMPI(MPI_Gatherv(adj->j, mnz, MPIU_INT, J, allnz, dispnz, MPIU_INT, 0, PetscObjectComm((PetscObject)A)));
868252a1336SBarry Smith PetscCall(PetscFree2(allnz, dispnz));
869252a1336SBarry Smith } else {
870252a1336SBarry Smith if (adj->values) PetscCallMPI(MPI_Gatherv(adj->values, mnz, MPIU_INT, NULL, NULL, NULL, MPIU_INT, 0, PetscObjectComm((PetscObject)A)));
871252a1336SBarry Smith PetscCallMPI(MPI_Gatherv(adj->j, mnz, MPIU_INT, NULL, NULL, NULL, MPIU_INT, 0, PetscObjectComm((PetscObject)A)));
872252a1336SBarry Smith }
873252a1336SBarry Smith PetscCallMPI(MPI_Scan(&nz, &nzstart, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)A)));
874252a1336SBarry Smith nzstart -= nz;
875252a1336SBarry Smith /* shift the i[] values so they will be correct after being received */
876252a1336SBarry Smith for (i = 0; i < m; i++) adj->i[i] += nzstart;
877252a1336SBarry Smith PetscCall(PetscMPIIntCast(m, &mm));
878252a1336SBarry Smith if (!rank) {
879252a1336SBarry Smith PetscCall(PetscMalloc1(M + 1, &II));
880252a1336SBarry Smith PetscCall(PetscMalloc2(size, &allm, size, &dispm));
881252a1336SBarry Smith PetscCallMPI(MPI_Gather(&mm, 1, MPI_INT, allm, 1, MPI_INT, 0, PetscObjectComm((PetscObject)A)));
8829371c9d4SSatish Balay dispm[0] = 0;
8839371c9d4SSatish Balay for (i = 1; i < size; i++) dispm[i] = dispm[i - 1] + allm[i - 1];
884252a1336SBarry Smith PetscCallMPI(MPI_Gatherv(adj->i, mm, MPIU_INT, II, allm, dispm, MPIU_INT, 0, PetscObjectComm((PetscObject)A)));
885252a1336SBarry Smith PetscCall(PetscFree2(allm, dispm));
886252a1336SBarry Smith II[M] = NZ;
887252a1336SBarry Smith } else {
888252a1336SBarry Smith PetscCallMPI(MPI_Gather(&mm, 1, MPI_INT, NULL, 1, MPI_INT, 0, PetscObjectComm((PetscObject)A)));
889252a1336SBarry Smith PetscCallMPI(MPI_Gatherv(adj->i, mm, MPIU_INT, NULL, NULL, NULL, MPIU_INT, 0, PetscObjectComm((PetscObject)A)));
890252a1336SBarry Smith }
891252a1336SBarry Smith /* shift the i[] values back */
892252a1336SBarry Smith for (i = 0; i < m; i++) adj->i[i] -= nzstart;
893252a1336SBarry Smith if (!rank) PetscCall(MatCreateMPIAdj(PETSC_COMM_SELF, M, N, II, J, Values, B));
8943ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
895252a1336SBarry Smith }
896252a1336SBarry Smith
8979a3665c6SJed Brown /*@
89811a5261eSBarry Smith MatMPIAdjCreateNonemptySubcommMat - create the same `MATMPIADJ` matrix on a subcommunicator containing only processes owning a positive number of rows
8999a3665c6SJed Brown
9009a3665c6SJed Brown Collective
9019a3665c6SJed Brown
9024165533cSJose E. Roman Input Parameter:
9032ef1f0ffSBarry Smith . A - original `MATMPIADJ` matrix
9049a3665c6SJed Brown
9054165533cSJose E. Roman Output Parameter:
906d8a51d2aSBarry Smith . B - matrix on subcommunicator, `NULL` on MPI processes that own zero rows of `A`
9079a3665c6SJed Brown
9089a3665c6SJed Brown Level: developer
9099a3665c6SJed Brown
9109a3665c6SJed Brown Note:
9112ef1f0ffSBarry Smith The matrix `B` should be destroyed with `MatDestroy()`. The arrays are not copied, so `B` should be destroyed before `A` is destroyed.
9129a3665c6SJed Brown
913d8a51d2aSBarry Smith Developer Note:
914d8a51d2aSBarry 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.
915d8a51d2aSBarry Smith
9161cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATMPIADJ`, `MatCreateMPIAdj()`
9179a3665c6SJed Brown @*/
MatMPIAdjCreateNonemptySubcommMat(Mat A,Mat * B)918d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMPIAdjCreateNonemptySubcommMat(Mat A, Mat *B)
919d71ae5a4SJacob Faibussowitsch {
9209a3665c6SJed Brown PetscFunctionBegin;
9219a3665c6SJed Brown PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
922cac4c232SBarry Smith PetscUseMethod(A, "MatMPIAdjCreateNonemptySubcommMat_C", (Mat, Mat *), (A, B));
9233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
9249a3665c6SJed Brown }
9259a3665c6SJed Brown
9260bad9183SKris Buschelman /*MC
927fafad747SKris Buschelman MATMPIADJ - MATMPIADJ = "mpiadj" - A matrix type to be used for distributed adjacency matrices,
9280bad9183SKris Buschelman intended for use constructing orderings and partitionings.
9290bad9183SKris Buschelman
9300bad9183SKris Buschelman Level: beginner
9310bad9183SKris Buschelman
93211a5261eSBarry Smith Note:
93311a5261eSBarry Smith You can provide values to the matrix using `MatMPIAdjSetPreallocation()`, `MatCreateMPIAdj()`, or
93411a5261eSBarry Smith by calling `MatSetValues()` and `MatAssemblyBegin()` followed by `MatAssemblyEnd()`
9356a09307cSBarry Smith
9361cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatCreateMPIAdj()`, `MatMPIAdjSetPreallocation()`, `MatSetValues()`
9370bad9183SKris Buschelman M*/
MatCreate_MPIAdj(Mat B)938d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatCreate_MPIAdj(Mat B)
939d71ae5a4SJacob Faibussowitsch {
940273d9f13SBarry Smith Mat_MPIAdj *b;
941273d9f13SBarry Smith
942273d9f13SBarry Smith PetscFunctionBegin;
9434dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&b));
944b0a32e0cSBarry Smith B->data = (void *)b;
945aea10558SJacob Faibussowitsch B->ops[0] = MatOps_Values;
946273d9f13SBarry Smith B->assembled = PETSC_FALSE;
9476a09307cSBarry Smith B->preallocated = PETSC_TRUE; /* so that MatSetValues() may be used */
948273d9f13SBarry Smith
9499566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMPIAdjSetPreallocation_C", MatMPIAdjSetPreallocation_MPIAdj));
9509566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMPIAdjCreateNonemptySubcommMat_C", MatMPIAdjCreateNonemptySubcommMat_MPIAdj));
9519566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMPIAdjToSeq_C", MatMPIAdjToSeq_MPIAdj));
952252a1336SBarry Smith PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMPIAdjToSeqRankZero_C", MatMPIAdjToSeqRankZero_MPIAdj));
9539566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)B, MATMPIADJ));
9543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
955273d9f13SBarry Smith }
956273d9f13SBarry Smith
957cc4c1da9SBarry Smith /*@
958d8a51d2aSBarry Smith MatMPIAdjToSeq - Converts an parallel `MATMPIADJ` matrix to complete `MATMPIADJ` on each process (needed by sequential partitioners)
959b44e7856SBarry Smith
960d083f849SBarry Smith Logically Collective
961b44e7856SBarry Smith
962b44e7856SBarry Smith Input Parameter:
963b44e7856SBarry Smith . A - the matrix
964b44e7856SBarry Smith
965b44e7856SBarry Smith Output Parameter:
966b44e7856SBarry Smith . B - the same matrix on all processes
967b44e7856SBarry Smith
968b44e7856SBarry Smith Level: intermediate
969b44e7856SBarry Smith
9701cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATMPIADJ`, `MatCreate()`, `MatCreateMPIAdj()`, `MatSetValues()`, `MatMPIAdjToSeqRankZero()`
971b44e7856SBarry Smith @*/
MatMPIAdjToSeq(Mat A,Mat * B)972d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMPIAdjToSeq(Mat A, Mat *B)
973d71ae5a4SJacob Faibussowitsch {
974b44e7856SBarry Smith PetscFunctionBegin;
975cac4c232SBarry Smith PetscUseMethod(A, "MatMPIAdjToSeq_C", (Mat, Mat *), (A, B));
9763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
977b44e7856SBarry Smith }
978b44e7856SBarry Smith
979cc4c1da9SBarry Smith /*@
980d8a51d2aSBarry Smith MatMPIAdjToSeqRankZero - Converts an parallel `MATMPIADJ` matrix to complete `MATMPIADJ` on rank zero (needed by sequential partitioners)
981252a1336SBarry Smith
982252a1336SBarry Smith Logically Collective
983252a1336SBarry Smith
984252a1336SBarry Smith Input Parameter:
985252a1336SBarry Smith . A - the matrix
986252a1336SBarry Smith
987252a1336SBarry Smith Output Parameter:
988252a1336SBarry Smith . B - the same matrix on rank zero, not set on other ranks
989252a1336SBarry Smith
990252a1336SBarry Smith Level: intermediate
991252a1336SBarry Smith
99211a5261eSBarry Smith Note:
993252a1336SBarry Smith This routine has the advantage on systems with multiple ranks per node since only one copy of the matrix
994252a1336SBarry Smith is stored on the first node, instead of the number of ranks copies. This can allow partitioning much larger
995d5b43468SJose E. Roman parallel graph sequentially.
996252a1336SBarry Smith
9971cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATMPIADJ`, `MatCreate()`, `MatCreateMPIAdj()`, `MatSetValues()`, `MatMPIAdjToSeq()`
998252a1336SBarry Smith @*/
MatMPIAdjToSeqRankZero(Mat A,Mat * B)999d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMPIAdjToSeqRankZero(Mat A, Mat *B)
1000d71ae5a4SJacob Faibussowitsch {
1001252a1336SBarry Smith PetscFunctionBegin;
1002252a1336SBarry Smith PetscUseMethod(A, "MatMPIAdjToSeqRankZero_C", (Mat, Mat *), (A, B));
10033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1004252a1336SBarry Smith }
1005252a1336SBarry Smith
10065d83a8b1SBarry Smith /*@
1007a23d5eceSKris Buschelman MatMPIAdjSetPreallocation - Sets the array used for storing the matrix elements
1008a23d5eceSKris Buschelman
1009d083f849SBarry Smith Logically Collective
1010a23d5eceSKris Buschelman
1011a23d5eceSKris Buschelman Input Parameters:
1012fe59aa6dSJacob Faibussowitsch + B - the matrix
10132ef1f0ffSBarry Smith . i - the indices into `j` for the start of each row
1014a23d5eceSKris Buschelman . j - the column indices for each row (sorted for each row).
10152ef1f0ffSBarry Smith The indices in `i` and `j` start with zero (NOT with one).
10162ef1f0ffSBarry Smith - values - [use `NULL` if not provided] edge weights
1017a23d5eceSKris Buschelman
1018a23d5eceSKris Buschelman Level: intermediate
1019a23d5eceSKris Buschelman
1020d8a51d2aSBarry Smith Notes:
1021d8a51d2aSBarry Smith The indices in `i` and `j` start with zero (NOT with one).
1022d8a51d2aSBarry Smith
1023d8a51d2aSBarry Smith You must NOT free the `i`, `values` and `j` arrays yourself. PETSc will free them
1024d8a51d2aSBarry Smith when the matrix is destroyed; you must allocate them with `PetscMalloc()`.
1025d8a51d2aSBarry Smith
1026d8a51d2aSBarry Smith You should not include the matrix diagonal elements.
1027d8a51d2aSBarry Smith
1028d8a51d2aSBarry Smith If you already have a matrix, you can create its adjacency matrix by a call
1029d8a51d2aSBarry Smith to `MatConvert()`, specifying a type of `MATMPIADJ`.
1030d8a51d2aSBarry Smith
1031d8a51d2aSBarry Smith Possible values for `MatSetOption()` - `MAT_STRUCTURALLY_SYMMETRIC`
1032d8a51d2aSBarry Smith
1033d8a51d2aSBarry Smith Fortran Note:
1034d8a51d2aSBarry Smith From Fortran the indices and values are copied so the array space need not be provided with `PetscMalloc()`.
1035d8a51d2aSBarry Smith
1036fe59aa6dSJacob Faibussowitsch .seealso: [](ch_matrices), `Mat`, `MatCreate()`, `MatCreateMPIAdj()`, `MatSetValues()`, `MATMPIADJ`
1037a23d5eceSKris Buschelman @*/
MatMPIAdjSetPreallocation(Mat B,PetscInt * i,PetscInt * j,PetscInt * values)1038d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMPIAdjSetPreallocation(Mat B, PetscInt *i, PetscInt *j, PetscInt *values)
1039d71ae5a4SJacob Faibussowitsch {
1040273d9f13SBarry Smith PetscFunctionBegin;
1041cac4c232SBarry Smith PetscTryMethod(B, "MatMPIAdjSetPreallocation_C", (Mat, PetscInt *, PetscInt *, PetscInt *), (B, i, j, values));
10423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1043273d9f13SBarry Smith }
1044273d9f13SBarry Smith
1045b97cf49bSBarry Smith /*@C
10463eda8832SBarry Smith MatCreateMPIAdj - Creates a sparse matrix representing an adjacency list.
1047d8a51d2aSBarry Smith The matrix need not have numerical values associated with it, it is
1048b97cf49bSBarry Smith intended for ordering (to reduce bandwidth etc) and partitioning.
1049b97cf49bSBarry Smith
1050d083f849SBarry Smith Collective
1051ef5ee4d1SLois Curfman McInnes
1052b97cf49bSBarry Smith Input Parameters:
1053c2e958feSBarry Smith + comm - MPI communicator
10540752156aSBarry Smith . m - number of local rows
105558ec18a5SLisandro Dalcin . N - number of global columns
10562ef1f0ffSBarry Smith . i - the indices into `j` for the start of each row
1057329f5518SBarry Smith . j - the column indices for each row (sorted for each row).
1058d8a51d2aSBarry Smith - values - the values, optional, use `NULL` if not provided
1059b97cf49bSBarry Smith
1060b97cf49bSBarry Smith Output Parameter:
1061b97cf49bSBarry Smith . A - the matrix
1062b97cf49bSBarry Smith
106315091d37SBarry Smith Level: intermediate
106415091d37SBarry Smith
106595452b02SPatrick Sanan Notes:
1066fe59aa6dSJacob Faibussowitsch The indices in `i` and `j` start with zero (NOT with one).
1067fe59aa6dSJacob Faibussowitsch
10682ef1f0ffSBarry Smith You must NOT free the `i`, `values` and `j` arrays yourself. PETSc will free them
1069d8a51d2aSBarry Smith when the matrix is destroyed; you must allocate them with `PetscMalloc()`.
10706a09307cSBarry Smith
10716a09307cSBarry Smith You should not include the matrix diagonals.
1072b97cf49bSBarry Smith
1073e3f7e1d6Svictorle If you already have a matrix, you can create its adjacency matrix by a call
107411a5261eSBarry Smith to `MatConvert()`, specifying a type of `MATMPIADJ`.
1075ededeb1aSvictorle
107611a5261eSBarry Smith Possible values for `MatSetOption()` - `MAT_STRUCTURALLY_SYMMETRIC`
1077b97cf49bSBarry Smith
1078d8a51d2aSBarry Smith Fortran Note:
10795d83a8b1SBarry Smith From Fortran the arrays `indices` and `values` must be retained by the user until `A` is destroyed
1080d8a51d2aSBarry Smith
10811cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatCreate()`, `MatConvert()`, `MatGetOrdering()`, `MATMPIADJ`, `MatMPIAdjSetPreallocation()`
1082b97cf49bSBarry Smith @*/
MatCreateMPIAdj(MPI_Comm comm,PetscInt m,PetscInt N,PetscInt i[],PetscInt j[],PetscInt values[],Mat * A)1083c080761bSJose E. Roman PetscErrorCode MatCreateMPIAdj(MPI_Comm comm, PetscInt m, PetscInt N, PetscInt i[], PetscInt j[], PetscInt values[], Mat *A) PeNSS
1084d71ae5a4SJacob Faibussowitsch {
1085433994e6SBarry Smith PetscFunctionBegin;
10869566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, A));
10879566063dSJacob Faibussowitsch PetscCall(MatSetSizes(*A, m, PETSC_DETERMINE, PETSC_DETERMINE, N));
10889566063dSJacob Faibussowitsch PetscCall(MatSetType(*A, MATMPIADJ));
10899566063dSJacob Faibussowitsch PetscCall(MatMPIAdjSetPreallocation(*A, i, j, values));
10903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1091b97cf49bSBarry Smith }
1092