1d5b485abSSatish Balay /* 2d5b485abSSatish Balay Routines to compute overlapping regions of a parallel MPI matrix 3d5b485abSSatish Balay and to find submatrices that were shared across processors. 4d5b485abSSatish Balay */ 5c6db04a5SJed Brown #include <../src/mat/impls/baij/mpi/mpibaij.h> 6c6db04a5SJed Brown #include <petscbt.h> 7d5b485abSSatish Balay 8b24ad042SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Local(Mat, PetscInt, char **, PetscInt *, PetscInt **); 9b24ad042SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Receive(Mat, PetscInt, PetscInt **, PetscInt **, PetscInt *); 10d5b485abSSatish Balay 11d71ae5a4SJacob Faibussowitsch PetscErrorCode MatIncreaseOverlap_MPIBAIJ(Mat C, PetscInt imax, IS is[], PetscInt ov) 12d71ae5a4SJacob Faibussowitsch { 13bd38d1d0SPierre Jolivet PetscInt i, N = C->cmap->N, bs = C->rmap->bs, n; 14bd38d1d0SPierre Jolivet const PetscInt *idx; 1536f4e84dSSatish Balay IS *is_new; 1636f4e84dSSatish Balay 173a40ed3dSBarry Smith PetscFunctionBegin; 189566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(imax, &is_new)); 19df36731bSBarry Smith /* Convert the indices into block format */ 209566063dSJacob Faibussowitsch PetscCall(ISCompressIndicesGeneral(N, C->rmap->n, bs, imax, is, is_new)); 2108401ef6SPierre Jolivet PetscCheck(ov >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative overlap specified"); 2248a46eb9SPierre Jolivet for (i = 0; i < ov; ++i) PetscCall(MatIncreaseOverlap_MPIBAIJ_Once(C, imax, is_new)); 23bd38d1d0SPierre Jolivet for (i = 0; i < imax; i++) { 24bd38d1d0SPierre Jolivet PetscCall(ISDestroy(&is[i])); 25bd38d1d0SPierre Jolivet PetscCall(ISGetLocalSize(is_new[i], &n)); 26bd38d1d0SPierre Jolivet PetscCall(ISGetIndices(is_new[i], &idx)); 27bd38d1d0SPierre Jolivet PetscCall(ISCreateBlock(PetscObjectComm((PetscObject)is_new[i]), bs, n, idx, PETSC_COPY_VALUES, &is[i])); 28bd38d1d0SPierre Jolivet PetscCall(ISDestroy(&is_new[i])); 29bd38d1d0SPierre Jolivet } 309566063dSJacob Faibussowitsch PetscCall(PetscFree(is_new)); 313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 32d5b485abSSatish Balay } 33d5b485abSSatish Balay 34d5b485abSSatish Balay /* 35d5b485abSSatish Balay Sample message format: 36d5b485abSSatish Balay If a processor A wants processor B to process some elements corresponding 370e9b0e7eSHong Zhang to index sets is[1], is[5] 38d5b485abSSatish Balay mesg [0] = 2 (no of index sets in the mesg) 39d5b485abSSatish Balay ----------- 40d5b485abSSatish Balay mesg [1] = 1 => is[1] 41d5b485abSSatish Balay mesg [2] = sizeof(is[1]); 42d5b485abSSatish Balay ----------- 43d5b485abSSatish Balay mesg [5] = 5 => is[5] 44d5b485abSSatish Balay mesg [6] = sizeof(is[5]); 45d5b485abSSatish Balay ----------- 46d5b485abSSatish Balay mesg [7] 470e9b0e7eSHong Zhang mesg [n] data(is[1]) 48d5b485abSSatish Balay ----------- 49d5b485abSSatish Balay mesg[n+1] 50d5b485abSSatish Balay mesg[m] data(is[5]) 51d5b485abSSatish Balay ----------- 52d5b485abSSatish Balay 53d5b485abSSatish Balay Notes: 54d5b485abSSatish Balay nrqs - no of requests sent (or to be sent out) 552d4ee042Sprj- nrqr - no of requests received (which have to be or which have been processed) 56d5b485abSSatish Balay */ 57d71ae5a4SJacob Faibussowitsch PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Once(Mat C, PetscInt imax, IS is[]) 58d71ae5a4SJacob Faibussowitsch { 59df36731bSBarry Smith Mat_MPIBAIJ *c = (Mat_MPIBAIJ *)C->data; 605d0c19d7SBarry Smith const PetscInt **idx, *idx_i; 6124c049a4SHong Zhang PetscInt *n, *w3, *w4, **data, len; 626497c311SBarry Smith PetscMPIInt size, rank, tag1, tag2, *w2, *w1, nrqs, nrqr, *pa; 636497c311SBarry Smith PetscInt Mbs, **rbuf, row, msz, **outdat, **ptr; 646497c311SBarry Smith PetscInt *ctr, *tmp, *isz, *isz1, **xdata, **rbuf2, *d_p; 65131c27b5Sprj- PetscMPIInt *onodes1, *olengths1, *onodes2, *olengths2, proc = -1; 66f1af5d2fSBarry Smith PetscBT *table; 67749130a8SStefano Zampini MPI_Comm comm, *iscomms; 68d5b485abSSatish Balay MPI_Request *s_waits1, *r_waits1, *s_waits2, *r_waits2; 698f9f447aSHong Zhang char *t_p; 70d5b485abSSatish Balay 713a40ed3dSBarry Smith PetscFunctionBegin; 729566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)C, &comm)); 73d5b485abSSatish Balay size = c->size; 74d5b485abSSatish Balay rank = c->rank; 75df36731bSBarry Smith Mbs = c->Mbs; 76d5b485abSSatish Balay 779566063dSJacob Faibussowitsch PetscCall(PetscObjectGetNewTag((PetscObject)C, &tag1)); 789566063dSJacob Faibussowitsch PetscCall(PetscObjectGetNewTag((PetscObject)C, &tag2)); 79c7dd2924SSatish Balay 80864bd19bSPierre Jolivet PetscCall(PetscMalloc2(imax, (PetscInt ***)&idx, imax, &n)); 8124c049a4SHong Zhang 826497c311SBarry Smith for (PetscInt i = 0; i < imax; i++) { 839566063dSJacob Faibussowitsch PetscCall(ISGetIndices(is[i], &idx[i])); 849566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(is[i], &n[i])); 85d5b485abSSatish Balay } 86d5b485abSSatish Balay 87d5b485abSSatish Balay /* evaluate communication - mesg to who,length of mesg, and buffer space 88d5b485abSSatish Balay required. Based on this, buffers are allocated, and data copied into them*/ 899566063dSJacob Faibussowitsch PetscCall(PetscCalloc4(size, &w1, size, &w2, size, &w3, size, &w4)); 906497c311SBarry Smith for (PetscInt i = 0; i < imax; i++) { 919566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(w4, size)); /* initialise work vector*/ 92d5b485abSSatish Balay idx_i = idx[i]; 93d5b485abSSatish Balay len = n[i]; 946497c311SBarry Smith for (PetscInt j = 0; j < len; j++) { 95d5b485abSSatish Balay row = idx_i[j]; 9608401ef6SPierre Jolivet PetscCheck(row >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Index set cannot have negative entries"); 979566063dSJacob Faibussowitsch PetscCall(PetscLayoutFindOwner(C->rmap, row * C->rmap->bs, &proc)); 98d5b485abSSatish Balay w4[proc]++; 99d5b485abSSatish Balay } 1006497c311SBarry Smith for (PetscMPIInt j = 0; j < size; j++) { 1019371c9d4SSatish Balay if (w4[j]) { 1029371c9d4SSatish Balay w1[j] += w4[j]; 1039371c9d4SSatish Balay w3[j]++; 1049371c9d4SSatish Balay } 105d5b485abSSatish Balay } 106d5b485abSSatish Balay } 107d5b485abSSatish Balay 108d5b485abSSatish Balay nrqs = 0; /* no of outgoing messages */ 109d5b485abSSatish Balay msz = 0; /* total mesg length (for all proc */ 1100e9b0e7eSHong Zhang w1[rank] = 0; /* no mesg sent to itself */ 111d5b485abSSatish Balay w3[rank] = 0; 1126497c311SBarry Smith for (PetscMPIInt i = 0; i < size; i++) { 1139371c9d4SSatish Balay if (w1[i]) { 1149371c9d4SSatish Balay w2[i] = 1; 1159371c9d4SSatish Balay nrqs++; 1169371c9d4SSatish Balay } /* there exists a message to proc i */ 117d5b485abSSatish Balay } 118d5b485abSSatish Balay /* pa - is list of processors to communicate with */ 1199566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrqs, &pa)); 1206497c311SBarry Smith for (PetscMPIInt i = 0, j = 0; i < size; i++) { 1219371c9d4SSatish Balay if (w1[i]) { 1229371c9d4SSatish Balay pa[j] = i; 1239371c9d4SSatish Balay j++; 1249371c9d4SSatish Balay } 125d5b485abSSatish Balay } 126d5b485abSSatish Balay 127d5b485abSSatish Balay /* Each message would have a header = 1 + 2*(no of IS) + data */ 1286497c311SBarry Smith for (PetscMPIInt i = 0; i < nrqs; i++) { 1296497c311SBarry Smith PetscMPIInt j = pa[i]; 130d5b485abSSatish Balay w1[j] += w2[j] + 2 * w3[j]; 131d5b485abSSatish Balay msz += w1[j]; 132d5b485abSSatish Balay } 133d5b485abSSatish Balay 134c7dd2924SSatish Balay /* Determine the number of messages to expect, their lengths, from from-ids */ 1359566063dSJacob Faibussowitsch PetscCall(PetscGatherNumberOfMessages(comm, w2, w1, &nrqr)); 1369566063dSJacob Faibussowitsch PetscCall(PetscGatherMessageLengths(comm, nrqs, nrqr, w1, &onodes1, &olengths1)); 137d5b485abSSatish Balay 138c7dd2924SSatish Balay /* Now post the Irecvs corresponding to these messages */ 1399566063dSJacob Faibussowitsch PetscCall(PetscPostIrecvInt(comm, tag1, nrqr, onodes1, olengths1, &rbuf, &r_waits1)); 140d5b485abSSatish Balay 141d5b485abSSatish Balay /* Allocate Memory for outgoing messages */ 1429566063dSJacob Faibussowitsch PetscCall(PetscMalloc4(size, &outdat, size, &ptr, msz, &tmp, size, &ctr)); 1439566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(outdat, size)); 1449566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(ptr, size)); 145d5b485abSSatish Balay { 146b24ad042SBarry Smith PetscInt *iptr = tmp, ict = 0; 1476497c311SBarry Smith for (PetscMPIInt i = 0; i < nrqs; i++) { 1486497c311SBarry Smith PetscMPIInt j = pa[i]; 149d5b485abSSatish Balay iptr += ict; 150d5b485abSSatish Balay outdat[j] = iptr; 151d5b485abSSatish Balay ict = w1[j]; 152d5b485abSSatish Balay } 153d5b485abSSatish Balay } 154d5b485abSSatish Balay 155d5b485abSSatish Balay /* Form the outgoing messages */ 156d5b485abSSatish Balay /*plug in the headers*/ 1576497c311SBarry Smith for (PetscMPIInt i = 0; i < nrqs; i++) { 1586497c311SBarry Smith PetscMPIInt j = pa[i]; 159d5b485abSSatish Balay outdat[j][0] = 0; 1609566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(outdat[j] + 1, 2 * w3[j])); 161d5b485abSSatish Balay ptr[j] = outdat[j] + 2 * w3[j] + 1; 162d5b485abSSatish Balay } 163d5b485abSSatish Balay 164d5b485abSSatish Balay /* Memory for doing local proc's work*/ 165d5b485abSSatish Balay { 1669566063dSJacob Faibussowitsch PetscCall(PetscCalloc5(imax, &table, imax, &data, imax, &isz, Mbs * imax, &d_p, (Mbs / PETSC_BITS_PER_BYTE + 1) * imax, &t_p)); 167d5b485abSSatish Balay 1686497c311SBarry Smith for (PetscInt i = 0; i < imax; i++) { 169b6410449SSatish Balay table[i] = t_p + (Mbs / PETSC_BITS_PER_BYTE + 1) * i; 170bbd702dbSSatish Balay data[i] = d_p + (Mbs)*i; 171d5b485abSSatish Balay } 172d5b485abSSatish Balay } 173d5b485abSSatish Balay 174d5b485abSSatish Balay /* Parse the IS and update local tables and the outgoing buf with the data*/ 175d5b485abSSatish Balay { 1766497c311SBarry Smith PetscInt n_i, *data_i, isz_i, *outdat_j, ctr_j, k; 177f1af5d2fSBarry Smith PetscBT table_i; 178d5b485abSSatish Balay 1796497c311SBarry Smith for (PetscInt i = 0; i < imax; i++) { 1809566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(ctr, size)); 181d5b485abSSatish Balay n_i = n[i]; 182d5b485abSSatish Balay table_i = table[i]; 183d5b485abSSatish Balay idx_i = idx[i]; 184d5b485abSSatish Balay data_i = data[i]; 185d5b485abSSatish Balay isz_i = isz[i]; 1866497c311SBarry Smith for (PetscInt j = 0; j < n_i; j++) { /* parse the indices of each IS */ 187d5b485abSSatish Balay row = idx_i[j]; 1889566063dSJacob Faibussowitsch PetscCall(PetscLayoutFindOwner(C->rmap, row * C->rmap->bs, &proc)); 189d5b485abSSatish Balay if (proc != rank) { /* copy to the outgoing buffer */ 190d5b485abSSatish Balay ctr[proc]++; 191d5b485abSSatish Balay *ptr[proc] = row; 192d5b485abSSatish Balay ptr[proc]++; 193d6b45a43SBarry Smith } else { /* Update the local table */ 19426fbe8dcSKarl Rupp if (!PetscBTLookupSet(table_i, row)) data_i[isz_i++] = row; 195d5b485abSSatish Balay } 196d5b485abSSatish Balay } 197d5b485abSSatish Balay /* Update the headers for the current IS */ 1986497c311SBarry Smith for (PetscMPIInt j = 0; j < size; j++) { /* Can Optimise this loop by using pa[] */ 199d5b485abSSatish Balay if ((ctr_j = ctr[j])) { 200d5b485abSSatish Balay outdat_j = outdat[j]; 201d5b485abSSatish Balay k = ++outdat_j[0]; 202d5b485abSSatish Balay outdat_j[2 * k] = ctr_j; 203d5b485abSSatish Balay outdat_j[2 * k - 1] = i; 204d5b485abSSatish Balay } 205d5b485abSSatish Balay } 206d5b485abSSatish Balay isz[i] = isz_i; 207d5b485abSSatish Balay } 208d5b485abSSatish Balay } 209d5b485abSSatish Balay 210d5b485abSSatish Balay /* Now post the sends */ 2119566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrqs, &s_waits1)); 2126497c311SBarry Smith for (PetscMPIInt i = 0; i < nrqs; ++i) { 2136497c311SBarry Smith PetscMPIInt j = pa[i]; 2146497c311SBarry Smith PetscCallMPI(MPIU_Isend(outdat[j], w1[j], MPIU_INT, j, tag1, comm, s_waits1 + i)); 215d5b485abSSatish Balay } 216d5b485abSSatish Balay 217d5b485abSSatish Balay /* No longer need the original indices*/ 2186497c311SBarry Smith for (PetscInt i = 0; i < imax; ++i) PetscCall(ISRestoreIndices(is[i], idx + i)); 2199566063dSJacob Faibussowitsch PetscCall(PetscFree2(*(PetscInt ***)&idx, n)); 220d5b485abSSatish Balay 2219566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(imax, &iscomms)); 2226497c311SBarry Smith for (PetscInt i = 0; i < imax; ++i) { 2239566063dSJacob Faibussowitsch PetscCall(PetscCommDuplicate(PetscObjectComm((PetscObject)is[i]), &iscomms[i], NULL)); 2249566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is[i])); 225d5b485abSSatish Balay } 226d5b485abSSatish Balay 227d5b485abSSatish Balay /* Do Local work*/ 2289566063dSJacob Faibussowitsch PetscCall(MatIncreaseOverlap_MPIBAIJ_Local(C, imax, table, isz, data)); 229d5b485abSSatish Balay 230d5b485abSSatish Balay /* Receive messages*/ 2319566063dSJacob Faibussowitsch PetscCallMPI(MPI_Waitall(nrqr, r_waits1, MPI_STATUSES_IGNORE)); 2329566063dSJacob Faibussowitsch PetscCallMPI(MPI_Waitall(nrqs, s_waits1, MPI_STATUSES_IGNORE)); 233d5b485abSSatish Balay 234d5b485abSSatish Balay /* Phase 1 sends are complete - deallocate buffers */ 2359566063dSJacob Faibussowitsch PetscCall(PetscFree4(outdat, ptr, tmp, ctr)); 2369566063dSJacob Faibussowitsch PetscCall(PetscFree4(w1, w2, w3, w4)); 237d5b485abSSatish Balay 2389566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrqr, &xdata)); 2399566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrqr, &isz1)); 2409566063dSJacob Faibussowitsch PetscCall(MatIncreaseOverlap_MPIBAIJ_Receive(C, nrqr, rbuf, xdata, isz1)); 241175c2bc5SJunchao Zhang if (rbuf) { 2429566063dSJacob Faibussowitsch PetscCall(PetscFree(rbuf[0])); 2439566063dSJacob Faibussowitsch PetscCall(PetscFree(rbuf)); 244175c2bc5SJunchao Zhang } 245d5b485abSSatish Balay 246d5b485abSSatish Balay /* Send the data back*/ 247d5b485abSSatish Balay /* Do a global reduction to know the buffer space req for incoming messages*/ 248d5b485abSSatish Balay { 249b24ad042SBarry Smith PetscMPIInt *rw1; 250d5b485abSSatish Balay 2519566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(size, &rw1)); 2526497c311SBarry Smith for (PetscMPIInt i = 0; i < nrqr; ++i) { 253175c2bc5SJunchao Zhang proc = onodes1[i]; 2546497c311SBarry Smith PetscCall(PetscMPIIntCast(isz1[i], &rw1[proc])); 255d5b485abSSatish Balay } 256d5b485abSSatish Balay 257c7dd2924SSatish Balay /* Determine the number of messages to expect, their lengths, from from-ids */ 2589566063dSJacob Faibussowitsch PetscCall(PetscGatherMessageLengths(comm, nrqr, nrqs, rw1, &onodes2, &olengths2)); 2599566063dSJacob Faibussowitsch PetscCall(PetscFree(rw1)); 260c7dd2924SSatish Balay } 261c7dd2924SSatish Balay /* Now post the Irecvs corresponding to these messages */ 2629566063dSJacob Faibussowitsch PetscCall(PetscPostIrecvInt(comm, tag2, nrqs, onodes2, olengths2, &rbuf2, &r_waits2)); 263d5b485abSSatish Balay 264d5b485abSSatish Balay /* Now post the sends */ 2659566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrqr, &s_waits2)); 2666497c311SBarry Smith for (PetscMPIInt i = 0; i < nrqr; ++i) { 2676497c311SBarry Smith PetscMPIInt j = onodes1[i]; 2686497c311SBarry Smith PetscCallMPI(MPIU_Isend(xdata[i], isz1[i], MPIU_INT, j, tag2, comm, s_waits2 + i)); 269d5b485abSSatish Balay } 270d5b485abSSatish Balay 2719566063dSJacob Faibussowitsch PetscCall(PetscFree(onodes1)); 2729566063dSJacob Faibussowitsch PetscCall(PetscFree(olengths1)); 273175c2bc5SJunchao Zhang 274d5b485abSSatish Balay /* receive work done on other processors*/ 275d5b485abSSatish Balay { 276b24ad042SBarry Smith PetscMPIInt idex; 277b24ad042SBarry Smith PetscInt is_no, ct1, max, *rbuf2_i, isz_i, *data_i, jmax; 278f1af5d2fSBarry Smith PetscBT table_i; 279d5b485abSSatish Balay 2806497c311SBarry Smith for (PetscMPIInt i = 0; i < nrqs; ++i) { 2819566063dSJacob Faibussowitsch PetscCallMPI(MPI_Waitany(nrqs, r_waits2, &idex, MPI_STATUS_IGNORE)); 282d5b485abSSatish Balay /* Process the message*/ 283999d9058SBarry Smith rbuf2_i = rbuf2[idex]; 284d5b485abSSatish Balay ct1 = 2 * rbuf2_i[0] + 1; 285999d9058SBarry Smith jmax = rbuf2[idex][0]; 2866497c311SBarry Smith for (PetscInt j = 1; j <= jmax; j++) { 287d5b485abSSatish Balay max = rbuf2_i[2 * j]; 288d5b485abSSatish Balay is_no = rbuf2_i[2 * j - 1]; 289d5b485abSSatish Balay isz_i = isz[is_no]; 290d5b485abSSatish Balay data_i = data[is_no]; 291d5b485abSSatish Balay table_i = table[is_no]; 2926497c311SBarry Smith for (PetscInt k = 0; k < max; k++, ct1++) { 293d5b485abSSatish Balay row = rbuf2_i[ct1]; 29426fbe8dcSKarl Rupp if (!PetscBTLookupSet(table_i, row)) data_i[isz_i++] = row; 295d5b485abSSatish Balay } 296d5b485abSSatish Balay isz[is_no] = isz_i; 297d5b485abSSatish Balay } 298d5b485abSSatish Balay } 2999566063dSJacob Faibussowitsch PetscCallMPI(MPI_Waitall(nrqr, s_waits2, MPI_STATUSES_IGNORE)); 300d5b485abSSatish Balay } 301d5b485abSSatish Balay 3026497c311SBarry Smith for (PetscInt i = 0; i < imax; ++i) { 3039566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(iscomms[i], isz[i], data[i], PETSC_COPY_VALUES, is + i)); 3049566063dSJacob Faibussowitsch PetscCall(PetscCommDestroy(&iscomms[i])); 305d5b485abSSatish Balay } 306d5b485abSSatish Balay 3079566063dSJacob Faibussowitsch PetscCall(PetscFree(iscomms)); 3089566063dSJacob Faibussowitsch PetscCall(PetscFree(onodes2)); 3099566063dSJacob Faibussowitsch PetscCall(PetscFree(olengths2)); 310c7dd2924SSatish Balay 3119566063dSJacob Faibussowitsch PetscCall(PetscFree(pa)); 312175c2bc5SJunchao Zhang if (rbuf2) { 3139566063dSJacob Faibussowitsch PetscCall(PetscFree(rbuf2[0])); 3149566063dSJacob Faibussowitsch PetscCall(PetscFree(rbuf2)); 315175c2bc5SJunchao Zhang } 3169566063dSJacob Faibussowitsch PetscCall(PetscFree(s_waits1)); 3179566063dSJacob Faibussowitsch PetscCall(PetscFree(r_waits1)); 3189566063dSJacob Faibussowitsch PetscCall(PetscFree(s_waits2)); 3199566063dSJacob Faibussowitsch PetscCall(PetscFree(r_waits2)); 3209566063dSJacob Faibussowitsch PetscCall(PetscFree5(table, data, isz, d_p, t_p)); 321175c2bc5SJunchao Zhang if (xdata) { 3229566063dSJacob Faibussowitsch PetscCall(PetscFree(xdata[0])); 3239566063dSJacob Faibussowitsch PetscCall(PetscFree(xdata)); 324175c2bc5SJunchao Zhang } 3259566063dSJacob Faibussowitsch PetscCall(PetscFree(isz1)); 3263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 327d5b485abSSatish Balay } 328d5b485abSSatish Balay 329d5b485abSSatish Balay /* 330df36731bSBarry Smith MatIncreaseOverlap_MPIBAIJ_Local - Called by MatincreaseOverlap, to do 331d5b485abSSatish Balay the work on the local processor. 332d5b485abSSatish Balay 333d5b485abSSatish Balay Inputs: 334df36731bSBarry Smith C - MAT_MPIBAIJ; 335d5b485abSSatish Balay imax - total no of index sets processed at a time; 336df36731bSBarry Smith table - an array of char - size = Mbs bits. 337d5b485abSSatish Balay 338d5b485abSSatish Balay Output: 3390e9b0e7eSHong Zhang isz - array containing the count of the solution elements corresponding 340d5b485abSSatish Balay to each index set; 341d5b485abSSatish Balay data - pointer to the solutions 342d5b485abSSatish Balay */ 343d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Local(Mat C, PetscInt imax, PetscBT *table, PetscInt *isz, PetscInt **data) 344d71ae5a4SJacob Faibussowitsch { 345df36731bSBarry Smith Mat_MPIBAIJ *c = (Mat_MPIBAIJ *)C->data; 346d5b485abSSatish Balay Mat A = c->A, B = c->B; 347df36731bSBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)B->data; 348b24ad042SBarry Smith PetscInt start, end, val, max, rstart, cstart, *ai, *aj; 349b24ad042SBarry Smith PetscInt *bi, *bj, *garray, i, j, k, row, *data_i, isz_i; 350f1af5d2fSBarry Smith PetscBT table_i; 351d5b485abSSatish Balay 3523a40ed3dSBarry Smith PetscFunctionBegin; 353899cda47SBarry Smith rstart = c->rstartbs; 354899cda47SBarry Smith cstart = c->cstartbs; 355d5b485abSSatish Balay ai = a->i; 356df36731bSBarry Smith aj = a->j; 357d5b485abSSatish Balay bi = b->i; 358df36731bSBarry Smith bj = b->j; 359d5b485abSSatish Balay garray = c->garray; 360d5b485abSSatish Balay 361d5b485abSSatish Balay for (i = 0; i < imax; i++) { 362d5b485abSSatish Balay data_i = data[i]; 363d5b485abSSatish Balay table_i = table[i]; 364d5b485abSSatish Balay isz_i = isz[i]; 365d5b485abSSatish Balay for (j = 0, max = isz[i]; j < max; j++) { 366d5b485abSSatish Balay row = data_i[j] - rstart; 367d5b485abSSatish Balay start = ai[row]; 368d5b485abSSatish Balay end = ai[row + 1]; 369d5b485abSSatish Balay for (k = start; k < end; k++) { /* Amat */ 370df36731bSBarry Smith val = aj[k] + cstart; 37126fbe8dcSKarl Rupp if (!PetscBTLookupSet(table_i, val)) data_i[isz_i++] = val; 372d5b485abSSatish Balay } 373d5b485abSSatish Balay start = bi[row]; 374d5b485abSSatish Balay end = bi[row + 1]; 375d5b485abSSatish Balay for (k = start; k < end; k++) { /* Bmat */ 376df36731bSBarry Smith val = garray[bj[k]]; 37726fbe8dcSKarl Rupp if (!PetscBTLookupSet(table_i, val)) data_i[isz_i++] = val; 378d5b485abSSatish Balay } 379d5b485abSSatish Balay } 380d5b485abSSatish Balay isz[i] = isz_i; 381d5b485abSSatish Balay } 3823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 383d5b485abSSatish Balay } 3846497c311SBarry Smith 385d5b485abSSatish Balay /* 3862d4ee042Sprj- MatIncreaseOverlap_MPIBAIJ_Receive - Process the received messages, 387d5b485abSSatish Balay and return the output 388d5b485abSSatish Balay 389d5b485abSSatish Balay Input: 390d5b485abSSatish Balay C - the matrix 391d5b485abSSatish Balay nrqr - no of messages being processed. 3922d4ee042Sprj- rbuf - an array of pointers to the received requests 393d5b485abSSatish Balay 394d5b485abSSatish Balay Output: 395d5b485abSSatish Balay xdata - array of messages to be sent back 396d5b485abSSatish Balay isz1 - size of each message 397d5b485abSSatish Balay 398a8c7a070SBarry Smith For better efficiency perhaps we should malloc separately each xdata[i], 399d5b485abSSatish Balay then if a remalloc is required we need only copy the data for that one row 400a5b23f4aSJose E. Roman rather than all previous rows as it is now where a single large chunk of 401d5b485abSSatish Balay memory is used. 402d5b485abSSatish Balay 403d5b485abSSatish Balay */ 404d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Receive(Mat C, PetscInt nrqr, PetscInt **rbuf, PetscInt **xdata, PetscInt *isz1) 405d71ae5a4SJacob Faibussowitsch { 406df36731bSBarry Smith Mat_MPIBAIJ *c = (Mat_MPIBAIJ *)C->data; 407d5b485abSSatish Balay Mat A = c->A, B = c->B; 408df36731bSBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)B->data; 409b24ad042SBarry Smith PetscInt rstart, cstart, *ai, *aj, *bi, *bj, *garray, i, j, k; 410b24ad042SBarry Smith PetscInt row, total_sz, ct, ct1, ct2, ct3, mem_estimate, oct2, l, start, end; 411d2a212eaSBarry Smith PetscInt val, max1, max2, Mbs, no_malloc = 0, *tmp, new_estimate, ctr; 412b24ad042SBarry Smith PetscInt *rbuf_i, kmax, rbuf_0; 413f1af5d2fSBarry Smith PetscBT xtable; 414d5b485abSSatish Balay 4153a40ed3dSBarry Smith PetscFunctionBegin; 416df36731bSBarry Smith Mbs = c->Mbs; 417899cda47SBarry Smith rstart = c->rstartbs; 418899cda47SBarry Smith cstart = c->cstartbs; 419d5b485abSSatish Balay ai = a->i; 420df36731bSBarry Smith aj = a->j; 421d5b485abSSatish Balay bi = b->i; 422df36731bSBarry Smith bj = b->j; 423d5b485abSSatish Balay garray = c->garray; 424d5b485abSSatish Balay 425d5b485abSSatish Balay for (i = 0, ct = 0, total_sz = 0; i < nrqr; ++i) { 426d5b485abSSatish Balay rbuf_i = rbuf[i]; 427d5b485abSSatish Balay rbuf_0 = rbuf_i[0]; 428d5b485abSSatish Balay ct += rbuf_0; 42926fbe8dcSKarl Rupp for (j = 1; j <= rbuf_0; j++) total_sz += rbuf_i[2 * j]; 430d5b485abSSatish Balay } 431d5b485abSSatish Balay 432701b8814SBarry Smith if (c->Mbs) max1 = ct * (a->nz + b->nz) / c->Mbs; 433701b8814SBarry Smith else max1 = 1; 434d5b485abSSatish Balay mem_estimate = 3 * ((total_sz > max1 ? total_sz : max1) + 1); 435175c2bc5SJunchao Zhang if (nrqr) { 4369566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(mem_estimate, &xdata[0])); 437d5b485abSSatish Balay ++no_malloc; 438175c2bc5SJunchao Zhang } 4399566063dSJacob Faibussowitsch PetscCall(PetscBTCreate(Mbs, &xtable)); 4409566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(isz1, nrqr)); 441d5b485abSSatish Balay 442d5b485abSSatish Balay ct3 = 0; 443d5b485abSSatish Balay for (i = 0; i < nrqr; i++) { /* for easch mesg from proc i */ 444d5b485abSSatish Balay rbuf_i = rbuf[i]; 445d5b485abSSatish Balay rbuf_0 = rbuf_i[0]; 446d5b485abSSatish Balay ct1 = 2 * rbuf_0 + 1; 447d5b485abSSatish Balay ct2 = ct1; 448d5b485abSSatish Balay ct3 += ct1; 449d5b485abSSatish Balay for (j = 1; j <= rbuf_0; j++) { /* for each IS from proc i*/ 4509566063dSJacob Faibussowitsch PetscCall(PetscBTMemzero(Mbs, xtable)); 451d5b485abSSatish Balay oct2 = ct2; 452d5b485abSSatish Balay kmax = rbuf_i[2 * j]; 453d5b485abSSatish Balay for (k = 0; k < kmax; k++, ct1++) { 454d5b485abSSatish Balay row = rbuf_i[ct1]; 455f1af5d2fSBarry Smith if (!PetscBTLookupSet(xtable, row)) { 456d5b485abSSatish Balay if (!(ct3 < mem_estimate)) { 457b24ad042SBarry Smith new_estimate = (PetscInt)(1.5 * mem_estimate) + 1; 4589566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(new_estimate, &tmp)); 4599566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(tmp, xdata[0], mem_estimate)); 4609566063dSJacob Faibussowitsch PetscCall(PetscFree(xdata[0])); 461d5b485abSSatish Balay xdata[0] = tmp; 4629371c9d4SSatish Balay mem_estimate = new_estimate; 4639371c9d4SSatish Balay ++no_malloc; 46426fbe8dcSKarl Rupp for (ctr = 1; ctr <= i; ctr++) xdata[ctr] = xdata[ctr - 1] + isz1[ctr - 1]; 465d5b485abSSatish Balay } 466d5b485abSSatish Balay xdata[i][ct2++] = row; 467d5b485abSSatish Balay ct3++; 468d5b485abSSatish Balay } 469d5b485abSSatish Balay } 470d5b485abSSatish Balay for (k = oct2, max2 = ct2; k < max2; k++) { 471d5b485abSSatish Balay row = xdata[i][k] - rstart; 472d5b485abSSatish Balay start = ai[row]; 473d5b485abSSatish Balay end = ai[row + 1]; 474d5b485abSSatish Balay for (l = start; l < end; l++) { 475df36731bSBarry Smith val = aj[l] + cstart; 476f1af5d2fSBarry Smith if (!PetscBTLookupSet(xtable, val)) { 477d5b485abSSatish Balay if (!(ct3 < mem_estimate)) { 478b24ad042SBarry Smith new_estimate = (PetscInt)(1.5 * mem_estimate) + 1; 4799566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(new_estimate, &tmp)); 4809566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(tmp, xdata[0], mem_estimate)); 4819566063dSJacob Faibussowitsch PetscCall(PetscFree(xdata[0])); 482d5b485abSSatish Balay xdata[0] = tmp; 4839371c9d4SSatish Balay mem_estimate = new_estimate; 4849371c9d4SSatish Balay ++no_malloc; 48526fbe8dcSKarl Rupp for (ctr = 1; ctr <= i; ctr++) xdata[ctr] = xdata[ctr - 1] + isz1[ctr - 1]; 486d5b485abSSatish Balay } 487d5b485abSSatish Balay xdata[i][ct2++] = val; 488d5b485abSSatish Balay ct3++; 489d5b485abSSatish Balay } 490d5b485abSSatish Balay } 491d5b485abSSatish Balay start = bi[row]; 492d5b485abSSatish Balay end = bi[row + 1]; 493d5b485abSSatish Balay for (l = start; l < end; l++) { 494df36731bSBarry Smith val = garray[bj[l]]; 495f1af5d2fSBarry Smith if (!PetscBTLookupSet(xtable, val)) { 496d5b485abSSatish Balay if (!(ct3 < mem_estimate)) { 497b24ad042SBarry Smith new_estimate = (PetscInt)(1.5 * mem_estimate) + 1; 4989566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(new_estimate, &tmp)); 4999566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(tmp, xdata[0], mem_estimate)); 5009566063dSJacob Faibussowitsch PetscCall(PetscFree(xdata[0])); 501d5b485abSSatish Balay xdata[0] = tmp; 5029371c9d4SSatish Balay mem_estimate = new_estimate; 5039371c9d4SSatish Balay ++no_malloc; 50426fbe8dcSKarl Rupp for (ctr = 1; ctr <= i; ctr++) xdata[ctr] = xdata[ctr - 1] + isz1[ctr - 1]; 505d5b485abSSatish Balay } 506d5b485abSSatish Balay xdata[i][ct2++] = val; 507d5b485abSSatish Balay ct3++; 508d5b485abSSatish Balay } 509d5b485abSSatish Balay } 510d5b485abSSatish Balay } 511d5b485abSSatish Balay /* Update the header*/ 512d5b485abSSatish Balay xdata[i][2 * j] = ct2 - oct2; /* Undo the vector isz1 and use only a var*/ 513d5b485abSSatish Balay xdata[i][2 * j - 1] = rbuf_i[2 * j - 1]; 514d5b485abSSatish Balay } 515d5b485abSSatish Balay xdata[i][0] = rbuf_0; 516175c2bc5SJunchao Zhang if (i + 1 < nrqr) xdata[i + 1] = xdata[i] + ct2; 517d5b485abSSatish Balay isz1[i] = ct2; /* size of each message */ 518d5b485abSSatish Balay } 5199566063dSJacob Faibussowitsch PetscCall(PetscBTDestroy(&xtable)); 5209566063dSJacob Faibussowitsch PetscCall(PetscInfo(C, "Allocated %" PetscInt_FMT " bytes, required %" PetscInt_FMT ", no of mallocs = %" PetscInt_FMT "\n", mem_estimate, ct3, no_malloc)); 5213ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 522d5b485abSSatish Balay } 523d5b485abSSatish Balay 524d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCreateSubMatrices_MPIBAIJ(Mat C, PetscInt ismax, const IS isrow[], const IS iscol[], MatReuse scall, Mat *submat[]) 525d71ae5a4SJacob Faibussowitsch { 52617df9f7cSHong Zhang IS *isrow_block, *iscol_block; 527cf4f527aSSatish Balay Mat_MPIBAIJ *c = (Mat_MPIBAIJ *)C->data; 528*f1957bc3SPierre Jolivet PetscInt nmax, nstages, i, pos, max_no, N = C->cmap->N, bs = C->rmap->bs; 5295dc98d6aSHong Zhang Mat_SeqBAIJ *subc; 5305c39f6d9SHong Zhang Mat_SubSppt *smat; 531fdfbdca6SPierre Jolivet PetscBool sym = PETSC_FALSE, flg[2]; 532a2feddc5SSatish Balay 5333a40ed3dSBarry Smith PetscFunctionBegin; 534fdfbdca6SPierre Jolivet PetscCall(PetscObjectTypeCompare((PetscObject)C, MATMPISBAIJ, flg)); 535fdfbdca6SPierre Jolivet if (flg[0]) { 536fdfbdca6SPierre Jolivet if (isrow == iscol) sym = PETSC_TRUE; 537fdfbdca6SPierre Jolivet else { 538fdfbdca6SPierre Jolivet flg[0] = flg[1] = PETSC_TRUE; 539fdfbdca6SPierre Jolivet for (i = 0; i < ismax; i++) { 540fdfbdca6SPierre Jolivet if (isrow[i] != iscol[i]) flg[0] = PETSC_FALSE; 541fdfbdca6SPierre Jolivet PetscCall(ISGetLocalSize(iscol[0], &nmax)); 542fdfbdca6SPierre Jolivet if (nmax == C->cmap->N && flg[1]) PetscCall(ISIdentity(iscol[0], flg + 1)); 543fdfbdca6SPierre Jolivet } 544fdfbdca6SPierre Jolivet sym = (PetscBool)(flg[0] || flg[1]); 545fdfbdca6SPierre Jolivet } 546fdfbdca6SPierre Jolivet } 54729dcf524SDmitry Karpeev /* The compression and expansion should be avoided. Doesn't point 54829dcf524SDmitry Karpeev out errors, might change the indices, hence buggey */ 549864bd19bSPierre Jolivet PetscCall(PetscMalloc2(ismax, &isrow_block, ismax, &iscol_block)); 550864bd19bSPierre Jolivet PetscCall(ISCompressIndicesGeneral(C->rmap->N, C->rmap->n, bs, ismax, isrow, isrow_block)); 551864bd19bSPierre Jolivet if (isrow == iscol) { 552864bd19bSPierre Jolivet for (i = 0; i < ismax; i++) { 553864bd19bSPierre Jolivet iscol_block[i] = isrow_block[i]; 554864bd19bSPierre Jolivet PetscCall(PetscObjectReference((PetscObject)iscol_block[i])); 555864bd19bSPierre Jolivet } 556864bd19bSPierre Jolivet } else PetscCall(ISCompressIndicesGeneral(N, C->cmap->n, bs, ismax, iscol, iscol_block)); 557cf4f527aSSatish Balay 558cf4f527aSSatish Balay /* Determine the number of stages through which submatrices are done */ 5594698041cSHong Zhang if (!C->cmap->N) nmax = 20 * 1000000 / sizeof(PetscInt); 5604698041cSHong Zhang else nmax = 20 * 1000000 / (c->Nbs * sizeof(PetscInt)); 561329f5518SBarry Smith if (!nmax) nmax = 1; 5625dc98d6aSHong Zhang 5635dc98d6aSHong Zhang if (scall == MAT_INITIAL_MATRIX) { 564*f1957bc3SPierre Jolivet nstages = ismax / nmax + ((ismax % nmax) ? 1 : 0); /* local nstages */ 565cf4f527aSSatish Balay 566653e4784SBarry Smith /* Make sure every processor loops through the nstages */ 567*f1957bc3SPierre Jolivet PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &nstages, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)C))); 5685dc98d6aSHong Zhang 5694698041cSHong Zhang /* Allocate memory to hold all the submatrices and dummy submatrices */ 5709566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(ismax + nstages, submat)); 5714698041cSHong Zhang } else { /* MAT_REUSE_MATRIX */ 5724698041cSHong Zhang if (ismax) { 5735dc98d6aSHong Zhang subc = (Mat_SeqBAIJ *)((*submat)[0]->data); 5745dc98d6aSHong Zhang smat = subc->submatis1; 5754698041cSHong Zhang } else { /* (*submat)[0] is a dummy matrix */ 5765c39f6d9SHong Zhang smat = (Mat_SubSppt *)(*submat)[0]->data; 5774698041cSHong Zhang } 57828b400f6SJacob Faibussowitsch PetscCheck(smat, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "MatCreateSubMatrices(...,MAT_REUSE_MATRIX,...) requires submat"); 5795dc98d6aSHong Zhang nstages = smat->nstages; 5805dc98d6aSHong Zhang } 5815dc98d6aSHong Zhang 582cf4f527aSSatish Balay for (i = 0, pos = 0; i < nstages; i++) { 583cf4f527aSSatish Balay if (pos + nmax <= ismax) max_no = nmax; 5844e195584SJunchao Zhang else if (pos >= ismax) max_no = 0; 585cf4f527aSSatish Balay else max_no = ismax - pos; 586a42ab0d6SHong Zhang 587fdfbdca6SPierre Jolivet PetscCall(MatCreateSubMatrices_MPIBAIJ_local(C, max_no, isrow_block + pos, iscol_block + pos, scall, *submat + pos, sym)); 5884e195584SJunchao Zhang if (!max_no) { 5894e195584SJunchao Zhang if (scall == MAT_INITIAL_MATRIX) { /* submat[pos] is a dummy matrix */ 5905c39f6d9SHong Zhang smat = (Mat_SubSppt *)(*submat)[pos]->data; 5914698041cSHong Zhang smat->nstages = nstages; 5924698041cSHong Zhang } 5934e195584SJunchao Zhang pos++; /* advance to next dummy matrix if any */ 5944e195584SJunchao Zhang } else pos += max_no; 595cf4f527aSSatish Balay } 59636f4e84dSSatish Balay 5975dc98d6aSHong Zhang if (scall == MAT_INITIAL_MATRIX && ismax) { 5985dc98d6aSHong Zhang /* save nstages for reuse */ 5995dc98d6aSHong Zhang subc = (Mat_SeqBAIJ *)((*submat)[0]->data); 6005dc98d6aSHong Zhang smat = subc->submatis1; 6015dc98d6aSHong Zhang smat->nstages = nstages; 6025dc98d6aSHong Zhang } 6035dc98d6aSHong Zhang 60436f4e84dSSatish Balay for (i = 0; i < ismax; i++) { 6059566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isrow_block[i])); 6069566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iscol_block[i])); 60736f4e84dSSatish Balay } 6089566063dSJacob Faibussowitsch PetscCall(PetscFree2(isrow_block, iscol_block)); 6093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 610a2feddc5SSatish Balay } 611a2feddc5SSatish Balay 612b24ad042SBarry Smith /* This code is used for BAIJ and SBAIJ matrices (unfortunate dependency) */ 613fdfbdca6SPierre Jolivet PetscErrorCode MatCreateSubMatrices_MPIBAIJ_local(Mat C, PetscInt ismax, const IS isrow[], const IS iscol[], MatReuse scall, Mat *submats, PetscBool sym) 614d71ae5a4SJacob Faibussowitsch { 61517df9f7cSHong Zhang Mat_MPIBAIJ *c = (Mat_MPIBAIJ *)C->data; 61617df9f7cSHong Zhang Mat A = c->A; 61717df9f7cSHong Zhang Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)c->B->data, *subc; 61817df9f7cSHong Zhang const PetscInt **icol, **irow; 61917df9f7cSHong Zhang PetscInt *nrow, *ncol, start; 6206497c311SBarry Smith PetscMPIInt *pa, **row2proc, *row2proc_i, proc = -1; 6216497c311SBarry Smith PetscMPIInt rank, size, tag0, tag2, tag3, tag4, *w1, *w2, *w3, *w4, nrqr, nrqs = 0, *req_source1 = NULL, *req_source2; 6226497c311SBarry Smith PetscInt **sbuf1, **sbuf2, *sbuf2_i, ct1, ct2, **rbuf1, row; 6236497c311SBarry Smith PetscInt msz, **ptr = NULL, *req_size = NULL, *ctr = NULL, *tmp = NULL, tcol; 6246497c311SBarry Smith PetscInt **rbuf3 = NULL, **sbuf_aj, **rbuf2 = NULL, max1, max2; 62517df9f7cSHong Zhang PetscInt **lens, is_no, ncols, *cols, mat_i, *mat_j, tmp2, jmax; 62617df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE) 627eec179cfSJacob Faibussowitsch PetscHMapI *cmap, cmap_i = NULL, *rmap, rmap_i; 62817df9f7cSHong Zhang #else 62917df9f7cSHong Zhang PetscInt **cmap, *cmap_i = NULL, **rmap, *rmap_i; 63017df9f7cSHong Zhang #endif 63196cff833SHong Zhang const PetscInt *irow_i, *icol_i; 6326497c311SBarry Smith PetscInt ctr_j, *sbuf1_j, *sbuf_aj_i, *rbuf1_i, kmax, *lens_i, jcnt; 63317df9f7cSHong Zhang MPI_Request *s_waits1, *r_waits1, *s_waits2, *r_waits2, *r_waits3; 63417df9f7cSHong Zhang MPI_Request *r_waits4, *s_waits3, *s_waits4; 63517df9f7cSHong Zhang MPI_Comm comm; 63664f5bb2dSSatish Balay PetscScalar **rbuf4, *rbuf4_i = NULL, **sbuf_aa, *vals, *mat_a = NULL, *imat_a = NULL, *sbuf_aa_i; 63717df9f7cSHong Zhang PetscMPIInt *onodes1, *olengths1, end; 6386497c311SBarry Smith PetscInt *imat_ilen, *imat_j, *imat_i; 6395c39f6d9SHong Zhang Mat_SubSppt *smat_i; 64017df9f7cSHong Zhang PetscBool *issorted, colflag, iscsorted = PETSC_TRUE; 64117df9f7cSHong Zhang PetscInt *sbuf1_i, *rbuf2_i, *rbuf3_i, ilen; 642a42ab0d6SHong Zhang PetscInt bs = C->rmap->bs, bs2 = c->bs2, rstart = c->rstartbs; 643a42ab0d6SHong Zhang PetscBool ijonly = c->ijonly; /* private flag indicates only matrix data structures are requested */ 644a42ab0d6SHong Zhang PetscInt nzA, nzB, *a_i = a->i, *b_i = b->i, *a_j = a->j, *b_j = b->j, ctmp, imark, *cworkA, *cworkB; 645afb3cbd8SHong Zhang PetscScalar *vworkA = NULL, *vworkB = NULL, *a_a = a->a, *b_a = b->a; 646a42ab0d6SHong Zhang PetscInt cstart = c->cstartbs, *bmap = c->garray; 64780d31651SHong Zhang PetscBool *allrows, *allcolumns; 648a42ab0d6SHong Zhang 64917df9f7cSHong Zhang PetscFunctionBegin; 6509566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)C, &comm)); 65117df9f7cSHong Zhang size = c->size; 65217df9f7cSHong Zhang rank = c->rank; 65317df9f7cSHong Zhang 6549566063dSJacob Faibussowitsch PetscCall(PetscMalloc5(ismax, &row2proc, ismax, &cmap, ismax, &rmap, ismax + 1, &allcolumns, ismax, &allrows)); 6559566063dSJacob Faibussowitsch PetscCall(PetscMalloc5(ismax, (PetscInt ***)&irow, ismax, (PetscInt ***)&icol, ismax, &nrow, ismax, &ncol, ismax, &issorted)); 65617df9f7cSHong Zhang 6576497c311SBarry Smith for (PetscInt i = 0; i < ismax; i++) { 6589566063dSJacob Faibussowitsch PetscCall(ISSorted(iscol[i], &issorted[i])); 65980d31651SHong Zhang if (!issorted[i]) iscsorted = issorted[i]; /* columns are not sorted! */ 6609566063dSJacob Faibussowitsch PetscCall(ISSorted(isrow[i], &issorted[i])); 66117df9f7cSHong Zhang 662ec3c739cSHong Zhang /* Check for special case: allcolumns */ 6639566063dSJacob Faibussowitsch PetscCall(ISIdentity(iscol[i], &colflag)); 6649566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(iscol[i], &ncol[i])); 6655dd0c08cSHong Zhang 6661abc2f7bSHong Zhang if (colflag && ncol[i] == c->Nbs) { 6671abc2f7bSHong Zhang allcolumns[i] = PETSC_TRUE; 6681abc2f7bSHong Zhang icol[i] = NULL; 6691abc2f7bSHong Zhang } else { 67017df9f7cSHong Zhang allcolumns[i] = PETSC_FALSE; 6719566063dSJacob Faibussowitsch PetscCall(ISGetIndices(iscol[i], &icol[i])); 6721abc2f7bSHong Zhang } 673ec3c739cSHong Zhang 674ec3c739cSHong Zhang /* Check for special case: allrows */ 6759566063dSJacob Faibussowitsch PetscCall(ISIdentity(isrow[i], &colflag)); 6769566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(isrow[i], &nrow[i])); 677ec3c739cSHong Zhang if (colflag && nrow[i] == c->Mbs) { 678ec3c739cSHong Zhang allrows[i] = PETSC_TRUE; 679ec3c739cSHong Zhang irow[i] = NULL; 680ec3c739cSHong Zhang } else { 681ec3c739cSHong Zhang allrows[i] = PETSC_FALSE; 6829566063dSJacob Faibussowitsch PetscCall(ISGetIndices(isrow[i], &irow[i])); 683ec3c739cSHong Zhang } 68417df9f7cSHong Zhang } 68517df9f7cSHong Zhang 68617df9f7cSHong Zhang if (scall == MAT_REUSE_MATRIX) { 68717df9f7cSHong Zhang /* Assumes new rows are same length as the old rows */ 6886497c311SBarry Smith for (PetscInt i = 0; i < ismax; i++) { 689f4f49eeaSPierre Jolivet subc = (Mat_SeqBAIJ *)submats[i]->data; 6906497c311SBarry Smith 691aed4548fSBarry Smith PetscCheck(subc->mbs == nrow[i] && subc->nbs == ncol[i], PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Cannot reuse matrix. wrong size"); 69217df9f7cSHong Zhang 69317df9f7cSHong Zhang /* Initial matrix as if empty */ 6949566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(subc->ilen, subc->mbs)); 69517df9f7cSHong Zhang 69617df9f7cSHong Zhang /* Initial matrix as if empty */ 69717df9f7cSHong Zhang submats[i]->factortype = C->factortype; 69817df9f7cSHong Zhang 69917df9f7cSHong Zhang smat_i = subc->submatis1; 70017df9f7cSHong Zhang 70117df9f7cSHong Zhang nrqs = smat_i->nrqs; 70217df9f7cSHong Zhang nrqr = smat_i->nrqr; 70317df9f7cSHong Zhang rbuf1 = smat_i->rbuf1; 70417df9f7cSHong Zhang rbuf2 = smat_i->rbuf2; 70517df9f7cSHong Zhang rbuf3 = smat_i->rbuf3; 70617df9f7cSHong Zhang req_source2 = smat_i->req_source2; 70717df9f7cSHong Zhang 70817df9f7cSHong Zhang sbuf1 = smat_i->sbuf1; 70917df9f7cSHong Zhang sbuf2 = smat_i->sbuf2; 71017df9f7cSHong Zhang ptr = smat_i->ptr; 71117df9f7cSHong Zhang tmp = smat_i->tmp; 71217df9f7cSHong Zhang ctr = smat_i->ctr; 71317df9f7cSHong Zhang 71417df9f7cSHong Zhang pa = smat_i->pa; 71517df9f7cSHong Zhang req_size = smat_i->req_size; 71617df9f7cSHong Zhang req_source1 = smat_i->req_source1; 71717df9f7cSHong Zhang 71817df9f7cSHong Zhang allcolumns[i] = smat_i->allcolumns; 719ec3c739cSHong Zhang allrows[i] = smat_i->allrows; 72017df9f7cSHong Zhang row2proc[i] = smat_i->row2proc; 72117df9f7cSHong Zhang rmap[i] = smat_i->rmap; 72217df9f7cSHong Zhang cmap[i] = smat_i->cmap; 72317df9f7cSHong Zhang } 7244698041cSHong Zhang 7254698041cSHong Zhang if (!ismax) { /* Get dummy submatrices and retrieve struct submatis1 */ 72608401ef6SPierre Jolivet PetscCheck(submats[0], PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "submats are null, cannot reuse"); 7275c39f6d9SHong Zhang smat_i = (Mat_SubSppt *)submats[0]->data; 7284698041cSHong Zhang 7294698041cSHong Zhang nrqs = smat_i->nrqs; 7304698041cSHong Zhang nrqr = smat_i->nrqr; 7314698041cSHong Zhang rbuf1 = smat_i->rbuf1; 7324698041cSHong Zhang rbuf2 = smat_i->rbuf2; 7334698041cSHong Zhang rbuf3 = smat_i->rbuf3; 7344698041cSHong Zhang req_source2 = smat_i->req_source2; 7354698041cSHong Zhang 7364698041cSHong Zhang sbuf1 = smat_i->sbuf1; 7374698041cSHong Zhang sbuf2 = smat_i->sbuf2; 7384698041cSHong Zhang ptr = smat_i->ptr; 7394698041cSHong Zhang tmp = smat_i->tmp; 7404698041cSHong Zhang ctr = smat_i->ctr; 7414698041cSHong Zhang 7424698041cSHong Zhang pa = smat_i->pa; 7434698041cSHong Zhang req_size = smat_i->req_size; 7444698041cSHong Zhang req_source1 = smat_i->req_source1; 7454698041cSHong Zhang 746c67fe082SHong Zhang allcolumns[0] = PETSC_FALSE; 7474698041cSHong Zhang } 74817df9f7cSHong Zhang } else { /* scall == MAT_INITIAL_MATRIX */ 74917df9f7cSHong Zhang /* Get some new tags to keep the communication clean */ 7509566063dSJacob Faibussowitsch PetscCall(PetscObjectGetNewTag((PetscObject)C, &tag2)); 7519566063dSJacob Faibussowitsch PetscCall(PetscObjectGetNewTag((PetscObject)C, &tag3)); 75217df9f7cSHong Zhang 75317df9f7cSHong Zhang /* evaluate communication - mesg to who, length of mesg, and buffer space 75417df9f7cSHong Zhang required. Based on this, buffers are allocated, and data copied into them*/ 7559566063dSJacob Faibussowitsch PetscCall(PetscCalloc4(size, &w1, size, &w2, size, &w3, size, &w4)); /* mesg size, initialize work vectors */ 75617df9f7cSHong Zhang 7576497c311SBarry Smith for (PetscInt i = 0; i < ismax; i++) { 75817df9f7cSHong Zhang jmax = nrow[i]; 75917df9f7cSHong Zhang irow_i = irow[i]; 76017df9f7cSHong Zhang 7619566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(jmax, &row2proc_i)); 76217df9f7cSHong Zhang row2proc[i] = row2proc_i; 76317df9f7cSHong Zhang 76417df9f7cSHong Zhang if (issorted[i]) proc = 0; 7656497c311SBarry Smith for (PetscInt j = 0; j < jmax; j++) { 76617df9f7cSHong Zhang if (!issorted[i]) proc = 0; 767ec3c739cSHong Zhang if (allrows[i]) row = j; 768ec3c739cSHong Zhang else row = irow_i[j]; 769ec3c739cSHong Zhang 770a42ab0d6SHong Zhang while (row >= c->rangebs[proc + 1]) proc++; 77117df9f7cSHong Zhang w4[proc]++; 77217df9f7cSHong Zhang row2proc_i[j] = proc; /* map row index to proc */ 77317df9f7cSHong Zhang } 7746497c311SBarry Smith for (PetscMPIInt j = 0; j < size; j++) { 7759371c9d4SSatish Balay if (w4[j]) { 7769371c9d4SSatish Balay w1[j] += w4[j]; 7779371c9d4SSatish Balay w3[j]++; 7789371c9d4SSatish Balay w4[j] = 0; 7799371c9d4SSatish Balay } 78017df9f7cSHong Zhang } 78117df9f7cSHong Zhang } 78217df9f7cSHong Zhang 78317df9f7cSHong Zhang nrqs = 0; /* no of outgoing messages */ 78417df9f7cSHong Zhang msz = 0; /* total mesg length (for all procs) */ 78517df9f7cSHong Zhang w1[rank] = 0; /* no mesg sent to self */ 78617df9f7cSHong Zhang w3[rank] = 0; 7876497c311SBarry Smith for (PetscMPIInt i = 0; i < size; i++) { 7889371c9d4SSatish Balay if (w1[i]) { 7899371c9d4SSatish Balay w2[i] = 1; 7909371c9d4SSatish Balay nrqs++; 7919371c9d4SSatish Balay } /* there exists a message to proc i */ 79217df9f7cSHong Zhang } 7939566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrqs, &pa)); /*(proc -array)*/ 7946497c311SBarry Smith for (PetscMPIInt i = 0, j = 0; i < size; i++) { 795ac530a7eSPierre Jolivet if (w1[i]) pa[j++] = i; 79617df9f7cSHong Zhang } 79717df9f7cSHong Zhang 79817df9f7cSHong Zhang /* Each message would have a header = 1 + 2*(no of IS) + data */ 7996497c311SBarry Smith for (PetscMPIInt i = 0; i < nrqs; i++) { 8006497c311SBarry Smith PetscMPIInt j = pa[i]; 80117df9f7cSHong Zhang w1[j] += w2[j] + 2 * w3[j]; 80217df9f7cSHong Zhang msz += w1[j]; 80317df9f7cSHong Zhang } 8046497c311SBarry Smith PetscCall(PetscInfo(0, "Number of outgoing messages %d Total message length %" PetscInt_FMT "\n", nrqs, msz)); 80517df9f7cSHong Zhang 80617df9f7cSHong Zhang /* Determine the number of messages to expect, their lengths, from from-ids */ 8079566063dSJacob Faibussowitsch PetscCall(PetscGatherNumberOfMessages(comm, w2, w1, &nrqr)); 8089566063dSJacob Faibussowitsch PetscCall(PetscGatherMessageLengths(comm, nrqs, nrqr, w1, &onodes1, &olengths1)); 80917df9f7cSHong Zhang 81017df9f7cSHong Zhang /* Now post the Irecvs corresponding to these messages */ 8119566063dSJacob Faibussowitsch PetscCall(PetscObjectGetNewTag((PetscObject)C, &tag0)); 8129566063dSJacob Faibussowitsch PetscCall(PetscPostIrecvInt(comm, tag0, nrqr, onodes1, olengths1, &rbuf1, &r_waits1)); 81317df9f7cSHong Zhang 81417df9f7cSHong Zhang /* Allocate Memory for outgoing messages */ 8159566063dSJacob Faibussowitsch PetscCall(PetscMalloc4(size, &sbuf1, size, &ptr, 2 * msz, &tmp, size, &ctr)); 8169566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(sbuf1, size)); 8179566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(ptr, size)); 81817df9f7cSHong Zhang 81917df9f7cSHong Zhang { 82017df9f7cSHong Zhang PetscInt *iptr = tmp; 8216497c311SBarry Smith PetscMPIInt k = 0; 8226497c311SBarry Smith for (PetscMPIInt i = 0; i < nrqs; i++) { 8236497c311SBarry Smith PetscMPIInt j = pa[i]; 82417df9f7cSHong Zhang iptr += k; 82517df9f7cSHong Zhang sbuf1[j] = iptr; 82617df9f7cSHong Zhang k = w1[j]; 82717df9f7cSHong Zhang } 82817df9f7cSHong Zhang } 82917df9f7cSHong Zhang 83017df9f7cSHong Zhang /* Form the outgoing messages. Initialize the header space */ 8316497c311SBarry Smith for (PetscMPIInt i = 0; i < nrqs; i++) { 8326497c311SBarry Smith PetscMPIInt j = pa[i]; 8336497c311SBarry Smith 83417df9f7cSHong Zhang sbuf1[j][0] = 0; 8359566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(sbuf1[j] + 1, 2 * w3[j])); 83617df9f7cSHong Zhang ptr[j] = sbuf1[j] + 2 * w3[j] + 1; 83717df9f7cSHong Zhang } 83817df9f7cSHong Zhang 83917df9f7cSHong Zhang /* Parse the isrow and copy data into outbuf */ 8406497c311SBarry Smith for (PetscInt i = 0; i < ismax; i++) { 84117df9f7cSHong Zhang row2proc_i = row2proc[i]; 8429566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(ctr, size)); 84317df9f7cSHong Zhang irow_i = irow[i]; 84417df9f7cSHong Zhang jmax = nrow[i]; 8456497c311SBarry Smith for (PetscInt j = 0; j < jmax; j++) { /* parse the indices of each IS */ 84617df9f7cSHong Zhang proc = row2proc_i[j]; 847ec3c739cSHong Zhang if (allrows[i]) row = j; 848ec3c739cSHong Zhang else row = irow_i[j]; 849ec3c739cSHong Zhang 85017df9f7cSHong Zhang if (proc != rank) { /* copy to the outgoing buf*/ 85117df9f7cSHong Zhang ctr[proc]++; 852ec3c739cSHong Zhang *ptr[proc] = row; 85317df9f7cSHong Zhang ptr[proc]++; 85417df9f7cSHong Zhang } 85517df9f7cSHong Zhang } 85617df9f7cSHong Zhang /* Update the headers for the current IS */ 8576497c311SBarry Smith for (PetscMPIInt j = 0; j < size; j++) { /* Can Optimise this loop too */ 85817df9f7cSHong Zhang if ((ctr_j = ctr[j])) { 8596497c311SBarry Smith PetscInt k; 8606497c311SBarry Smith 86117df9f7cSHong Zhang sbuf1_j = sbuf1[j]; 86217df9f7cSHong Zhang k = ++sbuf1_j[0]; 86317df9f7cSHong Zhang sbuf1_j[2 * k] = ctr_j; 86417df9f7cSHong Zhang sbuf1_j[2 * k - 1] = i; 86517df9f7cSHong Zhang } 86617df9f7cSHong Zhang } 86717df9f7cSHong Zhang } 86817df9f7cSHong Zhang 86917df9f7cSHong Zhang /* Now post the sends */ 8709566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrqs, &s_waits1)); 8716497c311SBarry Smith for (PetscMPIInt i = 0; i < nrqs; ++i) { 8726497c311SBarry Smith PetscMPIInt j = pa[i]; 8736497c311SBarry Smith PetscCallMPI(MPIU_Isend(sbuf1[j], w1[j], MPIU_INT, j, tag0, comm, s_waits1 + i)); 87417df9f7cSHong Zhang } 87517df9f7cSHong Zhang 87617df9f7cSHong Zhang /* Post Receives to capture the buffer size */ 8779566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrqs, &r_waits2)); 8789566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(nrqs, &req_source2, nrqs, &rbuf2, nrqs, &rbuf3)); 879175c2bc5SJunchao Zhang if (nrqs) rbuf2[0] = tmp + msz; 8806497c311SBarry Smith for (PetscMPIInt i = 1; i < nrqs; ++i) rbuf2[i] = rbuf2[i - 1] + w1[pa[i - 1]]; 8816497c311SBarry Smith for (PetscMPIInt i = 0; i < nrqs; ++i) { 8826497c311SBarry Smith PetscMPIInt j = pa[i]; 8836497c311SBarry Smith PetscCallMPI(MPIU_Irecv(rbuf2[i], w1[j], MPIU_INT, j, tag2, comm, r_waits2 + i)); 88417df9f7cSHong Zhang } 88517df9f7cSHong Zhang 88617df9f7cSHong Zhang /* Send to other procs the buf size they should allocate */ 88717df9f7cSHong Zhang /* Receive messages*/ 8889566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrqr, &s_waits2)); 8899566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(nrqr, &sbuf2, nrqr, &req_size, nrqr, &req_source1)); 89017df9f7cSHong Zhang 8919566063dSJacob Faibussowitsch PetscCallMPI(MPI_Waitall(nrqr, r_waits1, MPI_STATUSES_IGNORE)); 8926497c311SBarry Smith for (PetscMPIInt i = 0; i < nrqr; ++i) { 89317df9f7cSHong Zhang req_size[i] = 0; 89417df9f7cSHong Zhang rbuf1_i = rbuf1[i]; 89517df9f7cSHong Zhang start = 2 * rbuf1_i[0] + 1; 896175c2bc5SJunchao Zhang end = olengths1[i]; 8979566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(end, &sbuf2[i])); 89817df9f7cSHong Zhang sbuf2_i = sbuf2[i]; 8996497c311SBarry Smith for (PetscInt j = start; j < end; j++) { 900ddb3d6bcSHong Zhang row = rbuf1_i[j] - rstart; 901ddb3d6bcSHong Zhang ncols = a_i[row + 1] - a_i[row] + b_i[row + 1] - b_i[row]; 90217df9f7cSHong Zhang sbuf2_i[j] = ncols; 90317df9f7cSHong Zhang req_size[i] += ncols; 90417df9f7cSHong Zhang } 905175c2bc5SJunchao Zhang req_source1[i] = onodes1[i]; 90617df9f7cSHong Zhang /* form the header */ 90717df9f7cSHong Zhang sbuf2_i[0] = req_size[i]; 9086497c311SBarry Smith for (PetscInt j = 1; j < start; j++) sbuf2_i[j] = rbuf1_i[j]; 90917df9f7cSHong Zhang 9106497c311SBarry Smith PetscCallMPI(MPIU_Isend(sbuf2_i, end, MPIU_INT, req_source1[i], tag2, comm, s_waits2 + i)); 91117df9f7cSHong Zhang } 9129566063dSJacob Faibussowitsch PetscCall(PetscFree(onodes1)); 9139566063dSJacob Faibussowitsch PetscCall(PetscFree(olengths1)); 914175c2bc5SJunchao Zhang 9159566063dSJacob Faibussowitsch PetscCall(PetscFree(r_waits1)); 9169566063dSJacob Faibussowitsch PetscCall(PetscFree4(w1, w2, w3, w4)); 91717df9f7cSHong Zhang 91817df9f7cSHong Zhang /* Receive messages*/ 9199566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrqs, &r_waits3)); 92017df9f7cSHong Zhang 9219566063dSJacob Faibussowitsch PetscCallMPI(MPI_Waitall(nrqs, r_waits2, MPI_STATUSES_IGNORE)); 9226497c311SBarry Smith for (PetscMPIInt i = 0; i < nrqs; ++i) { 9239566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(rbuf2[i][0], &rbuf3[i])); 924175c2bc5SJunchao Zhang req_source2[i] = pa[i]; 9256497c311SBarry Smith PetscCallMPI(MPIU_Irecv(rbuf3[i], rbuf2[i][0], MPIU_INT, req_source2[i], tag3, comm, r_waits3 + i)); 92617df9f7cSHong Zhang } 9279566063dSJacob Faibussowitsch PetscCall(PetscFree(r_waits2)); 92817df9f7cSHong Zhang 92917df9f7cSHong Zhang /* Wait on sends1 and sends2 */ 9309566063dSJacob Faibussowitsch PetscCallMPI(MPI_Waitall(nrqs, s_waits1, MPI_STATUSES_IGNORE)); 9319566063dSJacob Faibussowitsch PetscCallMPI(MPI_Waitall(nrqr, s_waits2, MPI_STATUSES_IGNORE)); 9329566063dSJacob Faibussowitsch PetscCall(PetscFree(s_waits1)); 9339566063dSJacob Faibussowitsch PetscCall(PetscFree(s_waits2)); 93417df9f7cSHong Zhang 93517df9f7cSHong Zhang /* Now allocate sending buffers for a->j, and send them off */ 9369566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrqr, &sbuf_aj)); 9376497c311SBarry Smith jcnt = 0; 9386497c311SBarry Smith for (PetscMPIInt i = 0; i < nrqr; i++) jcnt += req_size[i]; 9396497c311SBarry Smith if (nrqr) PetscCall(PetscMalloc1(jcnt, &sbuf_aj[0])); 9406497c311SBarry Smith for (PetscMPIInt i = 1; i < nrqr; i++) sbuf_aj[i] = sbuf_aj[i - 1] + req_size[i - 1]; 94117df9f7cSHong Zhang 9429566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrqr, &s_waits3)); 94317df9f7cSHong Zhang { 9446497c311SBarry Smith for (PetscMPIInt i = 0; i < nrqr; i++) { 94517df9f7cSHong Zhang rbuf1_i = rbuf1[i]; 94617df9f7cSHong Zhang sbuf_aj_i = sbuf_aj[i]; 94717df9f7cSHong Zhang ct1 = 2 * rbuf1_i[0] + 1; 94817df9f7cSHong Zhang ct2 = 0; 9496497c311SBarry Smith for (PetscInt j = 1, max1 = rbuf1_i[0]; j <= max1; j++) { 95017df9f7cSHong Zhang kmax = rbuf1[i][2 * j]; 9516497c311SBarry Smith for (PetscInt k = 0; k < kmax; k++, ct1++) { 9526497c311SBarry Smith PetscInt l; 95317df9f7cSHong Zhang row = rbuf1_i[ct1] - rstart; 9549371c9d4SSatish Balay nzA = a_i[row + 1] - a_i[row]; 9559371c9d4SSatish Balay nzB = b_i[row + 1] - b_i[row]; 95617df9f7cSHong Zhang ncols = nzA + nzB; 9578e3a54c0SPierre Jolivet cworkA = PetscSafePointerPlusOffset(a_j, a_i[row]); 9588e3a54c0SPierre Jolivet cworkB = PetscSafePointerPlusOffset(b_j, b_i[row]); 95917df9f7cSHong Zhang 96017df9f7cSHong Zhang /* load the column indices for this row into cols */ 96117df9f7cSHong Zhang cols = sbuf_aj_i + ct2; 96217df9f7cSHong Zhang for (l = 0; l < nzB; l++) { 963a42ab0d6SHong Zhang if ((ctmp = bmap[cworkB[l]]) < cstart) cols[l] = ctmp; 964a42ab0d6SHong Zhang else break; 96517df9f7cSHong Zhang } 966a42ab0d6SHong Zhang imark = l; 967ad540459SPierre Jolivet for (l = 0; l < nzA; l++) cols[imark + l] = cstart + cworkA[l]; 968a42ab0d6SHong Zhang for (l = imark; l < nzB; l++) cols[nzA + l] = bmap[cworkB[l]]; 96917df9f7cSHong Zhang ct2 += ncols; 97017df9f7cSHong Zhang } 97117df9f7cSHong Zhang } 9726497c311SBarry Smith PetscCallMPI(MPIU_Isend(sbuf_aj_i, req_size[i], MPIU_INT, req_source1[i], tag3, comm, s_waits3 + i)); 97317df9f7cSHong Zhang } 97417df9f7cSHong Zhang } 97517df9f7cSHong Zhang 97617df9f7cSHong Zhang /* create col map: global col of C -> local col of submatrices */ 97717df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE) 9786497c311SBarry Smith for (PetscInt i = 0; i < ismax; i++) { 97917df9f7cSHong Zhang if (!allcolumns[i]) { 980eec179cfSJacob Faibussowitsch PetscCall(PetscHMapICreateWithSize(ncol[i], cmap + i)); 98117df9f7cSHong Zhang 98217df9f7cSHong Zhang jmax = ncol[i]; 98317df9f7cSHong Zhang icol_i = icol[i]; 98417df9f7cSHong Zhang cmap_i = cmap[i]; 9856497c311SBarry Smith for (PetscInt j = 0; j < jmax; j++) PetscCall(PetscHMapISet(cmap[i], icol_i[j] + 1, j + 1)); 98617df9f7cSHong Zhang } else cmap[i] = NULL; 98717df9f7cSHong Zhang } 98817df9f7cSHong Zhang #else 9896497c311SBarry Smith for (PetscInt i = 0; i < ismax; i++) { 99017df9f7cSHong Zhang if (!allcolumns[i]) { 9919566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(c->Nbs, &cmap[i])); 99217df9f7cSHong Zhang jmax = ncol[i]; 99317df9f7cSHong Zhang icol_i = icol[i]; 99417df9f7cSHong Zhang cmap_i = cmap[i]; 9956497c311SBarry Smith for (PetscInt j = 0; j < jmax; j++) cmap_i[icol_i[j]] = j + 1; 99617df9f7cSHong Zhang } else cmap[i] = NULL; 99717df9f7cSHong Zhang } 99817df9f7cSHong Zhang #endif 99917df9f7cSHong Zhang 100017df9f7cSHong Zhang /* Create lens which is required for MatCreate... */ 10016497c311SBarry Smith jcnt = 0; 10026497c311SBarry Smith for (PetscInt i = 0; i < ismax; i++) jcnt += nrow[i]; 10039566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ismax, &lens)); 10046497c311SBarry Smith if (ismax) PetscCall(PetscCalloc1(jcnt, &lens[0])); 10056497c311SBarry Smith for (PetscInt i = 1; i < ismax; i++) lens[i] = PetscSafePointerPlusOffset(lens[i - 1], nrow[i - 1]); 100617df9f7cSHong Zhang 100717df9f7cSHong Zhang /* Update lens from local data */ 10086497c311SBarry Smith for (PetscInt i = 0; i < ismax; i++) { 100917df9f7cSHong Zhang row2proc_i = row2proc[i]; 101017df9f7cSHong Zhang jmax = nrow[i]; 101117df9f7cSHong Zhang if (!allcolumns[i]) cmap_i = cmap[i]; 101217df9f7cSHong Zhang irow_i = irow[i]; 101317df9f7cSHong Zhang lens_i = lens[i]; 10146497c311SBarry Smith for (PetscInt j = 0; j < jmax; j++) { 1015ec3c739cSHong Zhang if (allrows[i]) row = j; 1016ec3c739cSHong Zhang else row = irow_i[j]; /* global blocked row of C */ 1017ec3c739cSHong Zhang 101817df9f7cSHong Zhang proc = row2proc_i[j]; 101917df9f7cSHong Zhang if (proc == rank) { 1020a42ab0d6SHong Zhang /* Get indices from matA and then from matB */ 102117df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE) 1022a42ab0d6SHong Zhang PetscInt tt; 1023a42ab0d6SHong Zhang #endif 1024a42ab0d6SHong Zhang row = row - rstart; 1025a42ab0d6SHong Zhang nzA = a_i[row + 1] - a_i[row]; 1026a42ab0d6SHong Zhang nzB = b_i[row + 1] - b_i[row]; 1027a42ab0d6SHong Zhang cworkA = a_j + a_i[row]; 10288e3a54c0SPierre Jolivet cworkB = PetscSafePointerPlusOffset(b_j, b_i[row]); 1029a42ab0d6SHong Zhang 1030a42ab0d6SHong Zhang if (!allcolumns[i]) { 1031a42ab0d6SHong Zhang #if defined(PETSC_USE_CTABLE) 10326497c311SBarry Smith for (PetscInt k = 0; k < nzA; k++) { 1033eec179cfSJacob Faibussowitsch PetscCall(PetscHMapIGetWithDefault(cmap_i, cstart + cworkA[k] + 1, 0, &tt)); 1034a42ab0d6SHong Zhang if (tt) lens_i[j]++; 1035a42ab0d6SHong Zhang } 10366497c311SBarry Smith for (PetscInt k = 0; k < nzB; k++) { 1037eec179cfSJacob Faibussowitsch PetscCall(PetscHMapIGetWithDefault(cmap_i, bmap[cworkB[k]] + 1, 0, &tt)); 1038a42ab0d6SHong Zhang if (tt) lens_i[j]++; 1039a42ab0d6SHong Zhang } 1040a42ab0d6SHong Zhang 1041a42ab0d6SHong Zhang #else 10426497c311SBarry Smith for (PetscInt k = 0; k < nzA; k++) { 1043a42ab0d6SHong Zhang if (cmap_i[cstart + cworkA[k]]) lens_i[j]++; 1044a42ab0d6SHong Zhang } 10456497c311SBarry Smith for (PetscInt k = 0; k < nzB; k++) { 1046a42ab0d6SHong Zhang if (cmap_i[bmap[cworkB[k]]]) lens_i[j]++; 1047a42ab0d6SHong Zhang } 1048a42ab0d6SHong Zhang #endif 1049a42ab0d6SHong Zhang } else { /* allcolumns */ 1050a42ab0d6SHong Zhang lens_i[j] = nzA + nzB; 1051a42ab0d6SHong Zhang } 105217df9f7cSHong Zhang } 105317df9f7cSHong Zhang } 105417df9f7cSHong Zhang } 105517df9f7cSHong Zhang 105617df9f7cSHong Zhang /* Create row map: global row of C -> local row of submatrices */ 10576497c311SBarry Smith for (PetscInt i = 0; i < ismax; i++) { 1058059f6b74SHong Zhang if (!allrows[i]) { 1059059f6b74SHong Zhang #if defined(PETSC_USE_CTABLE) 1060eec179cfSJacob Faibussowitsch PetscCall(PetscHMapICreateWithSize(nrow[i], rmap + i)); 106117df9f7cSHong Zhang irow_i = irow[i]; 106217df9f7cSHong Zhang jmax = nrow[i]; 10636497c311SBarry Smith for (PetscInt j = 0; j < jmax; j++) { 1064ec3c739cSHong Zhang if (allrows[i]) { 1065c76ffc5fSJacob Faibussowitsch PetscCall(PetscHMapISet(rmap[i], j + 1, j + 1)); 1066ec3c739cSHong Zhang } else { 1067c76ffc5fSJacob Faibussowitsch PetscCall(PetscHMapISet(rmap[i], irow_i[j] + 1, j + 1)); 106817df9f7cSHong Zhang } 106917df9f7cSHong Zhang } 107017df9f7cSHong Zhang #else 10719566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(c->Mbs, &rmap[i])); 107217df9f7cSHong Zhang rmap_i = rmap[i]; 107317df9f7cSHong Zhang irow_i = irow[i]; 107417df9f7cSHong Zhang jmax = nrow[i]; 10756497c311SBarry Smith for (PetscInt j = 0; j < jmax; j++) { 1076ec3c739cSHong Zhang if (allrows[i]) rmap_i[j] = j; 1077ec3c739cSHong Zhang else rmap_i[irow_i[j]] = j; 107817df9f7cSHong Zhang } 107917df9f7cSHong Zhang #endif 1080059f6b74SHong Zhang } else rmap[i] = NULL; 1081059f6b74SHong Zhang } 108217df9f7cSHong Zhang 108317df9f7cSHong Zhang /* Update lens from offproc data */ 108417df9f7cSHong Zhang { 108517df9f7cSHong Zhang PetscInt *rbuf2_i, *rbuf3_i, *sbuf1_i; 108617df9f7cSHong Zhang 10879566063dSJacob Faibussowitsch PetscCallMPI(MPI_Waitall(nrqs, r_waits3, MPI_STATUSES_IGNORE)); 108817df9f7cSHong Zhang for (tmp2 = 0; tmp2 < nrqs; tmp2++) { 108917df9f7cSHong Zhang sbuf1_i = sbuf1[pa[tmp2]]; 109017df9f7cSHong Zhang jmax = sbuf1_i[0]; 109117df9f7cSHong Zhang ct1 = 2 * jmax + 1; 109217df9f7cSHong Zhang ct2 = 0; 109317df9f7cSHong Zhang rbuf2_i = rbuf2[tmp2]; 109417df9f7cSHong Zhang rbuf3_i = rbuf3[tmp2]; 10956497c311SBarry Smith for (PetscInt j = 1; j <= jmax; j++) { 109617df9f7cSHong Zhang is_no = sbuf1_i[2 * j - 1]; 109717df9f7cSHong Zhang max1 = sbuf1_i[2 * j]; 109817df9f7cSHong Zhang lens_i = lens[is_no]; 109917df9f7cSHong Zhang if (!allcolumns[is_no]) cmap_i = cmap[is_no]; 110017df9f7cSHong Zhang rmap_i = rmap[is_no]; 11016497c311SBarry Smith for (PetscInt k = 0; k < max1; k++, ct1++) { 1102059f6b74SHong Zhang if (allrows[is_no]) { 1103059f6b74SHong Zhang row = sbuf1_i[ct1]; 1104059f6b74SHong Zhang } else { 110517df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE) 1106eec179cfSJacob Faibussowitsch PetscCall(PetscHMapIGetWithDefault(rmap_i, sbuf1_i[ct1] + 1, 0, &row)); 110717df9f7cSHong Zhang row--; 110808401ef6SPierre Jolivet PetscCheck(row >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "row not found in table"); 110917df9f7cSHong Zhang #else 111017df9f7cSHong Zhang row = rmap_i[sbuf1_i[ct1]]; /* the val in the new matrix to be */ 111117df9f7cSHong Zhang #endif 1112059f6b74SHong Zhang } 111317df9f7cSHong Zhang max2 = rbuf2_i[ct1]; 11146497c311SBarry Smith for (PetscInt l = 0; l < max2; l++, ct2++) { 111517df9f7cSHong Zhang if (!allcolumns[is_no]) { 111617df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE) 1117eec179cfSJacob Faibussowitsch PetscCall(PetscHMapIGetWithDefault(cmap_i, rbuf3_i[ct2] + 1, 0, &tcol)); 111817df9f7cSHong Zhang #else 111917df9f7cSHong Zhang tcol = cmap_i[rbuf3_i[ct2]]; 112017df9f7cSHong Zhang #endif 112117df9f7cSHong Zhang if (tcol) lens_i[row]++; 112217df9f7cSHong Zhang } else { /* allcolumns */ 112317df9f7cSHong Zhang lens_i[row]++; /* lens_i[row] += max2 ? */ 112417df9f7cSHong Zhang } 112517df9f7cSHong Zhang } 112617df9f7cSHong Zhang } 112717df9f7cSHong Zhang } 112817df9f7cSHong Zhang } 112917df9f7cSHong Zhang } 11309566063dSJacob Faibussowitsch PetscCall(PetscFree(r_waits3)); 11319566063dSJacob Faibussowitsch PetscCallMPI(MPI_Waitall(nrqr, s_waits3, MPI_STATUSES_IGNORE)); 11329566063dSJacob Faibussowitsch PetscCall(PetscFree(s_waits3)); 113317df9f7cSHong Zhang 113417df9f7cSHong Zhang /* Create the submatrices */ 11356497c311SBarry Smith for (PetscInt i = 0; i < ismax; i++) { 1136a42ab0d6SHong Zhang PetscInt bs_tmp; 1137a42ab0d6SHong Zhang if (ijonly) bs_tmp = 1; 1138a42ab0d6SHong Zhang else bs_tmp = bs; 113917df9f7cSHong Zhang 11409566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF, submats + i)); 11419566063dSJacob Faibussowitsch PetscCall(MatSetSizes(submats[i], nrow[i] * bs_tmp, ncol[i] * bs_tmp, PETSC_DETERMINE, PETSC_DETERMINE)); 114217df9f7cSHong Zhang 1143fdfbdca6SPierre Jolivet PetscCall(MatSetType(submats[i], sym ? ((PetscObject)A)->type_name : MATSEQBAIJ)); 11449566063dSJacob Faibussowitsch PetscCall(MatSeqBAIJSetPreallocation(submats[i], bs_tmp, 0, lens[i])); 11459566063dSJacob Faibussowitsch PetscCall(MatSeqSBAIJSetPreallocation(submats[i], bs_tmp, 0, lens[i])); /* this subroutine is used by SBAIJ routines */ 114617df9f7cSHong Zhang 11475c39f6d9SHong Zhang /* create struct Mat_SubSppt and attached it to submat */ 11489566063dSJacob Faibussowitsch PetscCall(PetscNew(&smat_i)); 114917df9f7cSHong Zhang subc = (Mat_SeqBAIJ *)submats[i]->data; 115017df9f7cSHong Zhang subc->submatis1 = smat_i; 115117df9f7cSHong Zhang 115217df9f7cSHong Zhang smat_i->destroy = submats[i]->ops->destroy; 1153f68bb481SHong Zhang submats[i]->ops->destroy = MatDestroySubMatrix_SeqBAIJ; 115417df9f7cSHong Zhang submats[i]->factortype = C->factortype; 115517df9f7cSHong Zhang 115617df9f7cSHong Zhang smat_i->id = i; 115717df9f7cSHong Zhang smat_i->nrqs = nrqs; 115817df9f7cSHong Zhang smat_i->nrqr = nrqr; 115917df9f7cSHong Zhang smat_i->rbuf1 = rbuf1; 116017df9f7cSHong Zhang smat_i->rbuf2 = rbuf2; 116117df9f7cSHong Zhang smat_i->rbuf3 = rbuf3; 116217df9f7cSHong Zhang smat_i->sbuf2 = sbuf2; 116317df9f7cSHong Zhang smat_i->req_source2 = req_source2; 116417df9f7cSHong Zhang 116517df9f7cSHong Zhang smat_i->sbuf1 = sbuf1; 116617df9f7cSHong Zhang smat_i->ptr = ptr; 116717df9f7cSHong Zhang smat_i->tmp = tmp; 116817df9f7cSHong Zhang smat_i->ctr = ctr; 116917df9f7cSHong Zhang 117017df9f7cSHong Zhang smat_i->pa = pa; 117117df9f7cSHong Zhang smat_i->req_size = req_size; 117217df9f7cSHong Zhang smat_i->req_source1 = req_source1; 117317df9f7cSHong Zhang 117417df9f7cSHong Zhang smat_i->allcolumns = allcolumns[i]; 1175ec3c739cSHong Zhang smat_i->allrows = allrows[i]; 117617df9f7cSHong Zhang smat_i->singleis = PETSC_FALSE; 117717df9f7cSHong Zhang smat_i->row2proc = row2proc[i]; 117817df9f7cSHong Zhang smat_i->rmap = rmap[i]; 117917df9f7cSHong Zhang smat_i->cmap = cmap[i]; 118017df9f7cSHong Zhang } 118117df9f7cSHong Zhang 11824698041cSHong Zhang if (!ismax) { /* Create dummy submats[0] for reuse struct subc */ 11839566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF, &submats[0])); 11849566063dSJacob Faibussowitsch PetscCall(MatSetSizes(submats[0], 0, 0, PETSC_DETERMINE, PETSC_DETERMINE)); 11859566063dSJacob Faibussowitsch PetscCall(MatSetType(submats[0], MATDUMMY)); 11864698041cSHong Zhang 11875c39f6d9SHong Zhang /* create struct Mat_SubSppt and attached it to submat */ 11884dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&smat_i)); 11894698041cSHong Zhang submats[0]->data = (void *)smat_i; 11904698041cSHong Zhang 11914698041cSHong Zhang smat_i->destroy = submats[0]->ops->destroy; 1192f68bb481SHong Zhang submats[0]->ops->destroy = MatDestroySubMatrix_Dummy; 11934698041cSHong Zhang submats[0]->factortype = C->factortype; 11944698041cSHong Zhang 1195006cef31SHong Zhang smat_i->id = 0; 11964698041cSHong Zhang smat_i->nrqs = nrqs; 11974698041cSHong Zhang smat_i->nrqr = nrqr; 11984698041cSHong Zhang smat_i->rbuf1 = rbuf1; 11994698041cSHong Zhang smat_i->rbuf2 = rbuf2; 12004698041cSHong Zhang smat_i->rbuf3 = rbuf3; 12014698041cSHong Zhang smat_i->sbuf2 = sbuf2; 12024698041cSHong Zhang smat_i->req_source2 = req_source2; 12034698041cSHong Zhang 12044698041cSHong Zhang smat_i->sbuf1 = sbuf1; 12054698041cSHong Zhang smat_i->ptr = ptr; 12064698041cSHong Zhang smat_i->tmp = tmp; 12074698041cSHong Zhang smat_i->ctr = ctr; 12084698041cSHong Zhang 12094698041cSHong Zhang smat_i->pa = pa; 12104698041cSHong Zhang smat_i->req_size = req_size; 12114698041cSHong Zhang smat_i->req_source1 = req_source1; 12124698041cSHong Zhang 1213c67fe082SHong Zhang smat_i->allcolumns = PETSC_FALSE; 12144698041cSHong Zhang smat_i->singleis = PETSC_FALSE; 12154698041cSHong Zhang smat_i->row2proc = NULL; 12164698041cSHong Zhang smat_i->rmap = NULL; 12174698041cSHong Zhang smat_i->cmap = NULL; 12184698041cSHong Zhang } 12194698041cSHong Zhang 12209566063dSJacob Faibussowitsch if (ismax) PetscCall(PetscFree(lens[0])); 12219566063dSJacob Faibussowitsch PetscCall(PetscFree(lens)); 1222175c2bc5SJunchao Zhang if (sbuf_aj) { 12239566063dSJacob Faibussowitsch PetscCall(PetscFree(sbuf_aj[0])); 12249566063dSJacob Faibussowitsch PetscCall(PetscFree(sbuf_aj)); 1225175c2bc5SJunchao Zhang } 122617df9f7cSHong Zhang 122717df9f7cSHong Zhang } /* endof scall == MAT_INITIAL_MATRIX */ 122817df9f7cSHong Zhang 122917df9f7cSHong Zhang /* Post recv matrix values */ 1230bbc89d27SHong Zhang if (!ijonly) { 12319566063dSJacob Faibussowitsch PetscCall(PetscObjectGetNewTag((PetscObject)C, &tag4)); 12329566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrqs, &rbuf4)); 12339566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrqs, &r_waits4)); 12346497c311SBarry Smith for (PetscMPIInt i = 0; i < nrqs; ++i) { 12359566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(rbuf2[i][0] * bs2, &rbuf4[i])); 12366497c311SBarry Smith PetscCallMPI(MPIU_Irecv(rbuf4[i], rbuf2[i][0] * bs2, MPIU_SCALAR, req_source2[i], tag4, comm, r_waits4 + i)); 123717df9f7cSHong Zhang } 123817df9f7cSHong Zhang 123917df9f7cSHong Zhang /* Allocate sending buffers for a->a, and send them off */ 12409566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrqr, &sbuf_aa)); 12416497c311SBarry Smith jcnt = 0; 12426497c311SBarry Smith for (PetscMPIInt i = 0; i < nrqr; i++) jcnt += req_size[i]; 12436497c311SBarry Smith if (nrqr) PetscCall(PetscMalloc1(jcnt * bs2, &sbuf_aa[0])); 12446497c311SBarry Smith for (PetscMPIInt i = 1; i < nrqr; i++) sbuf_aa[i] = sbuf_aa[i - 1] + req_size[i - 1] * bs2; 124517df9f7cSHong Zhang 12469566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrqr, &s_waits4)); 124717df9f7cSHong Zhang 12486497c311SBarry Smith for (PetscMPIInt i = 0; i < nrqr; i++) { 124917df9f7cSHong Zhang rbuf1_i = rbuf1[i]; 125017df9f7cSHong Zhang sbuf_aa_i = sbuf_aa[i]; 125117df9f7cSHong Zhang ct1 = 2 * rbuf1_i[0] + 1; 125217df9f7cSHong Zhang ct2 = 0; 12536497c311SBarry Smith for (PetscInt j = 1, max1 = rbuf1_i[0]; j <= max1; j++) { 125417df9f7cSHong Zhang kmax = rbuf1_i[2 * j]; 12556497c311SBarry Smith for (PetscInt k = 0; k < kmax; k++, ct1++) { 12566497c311SBarry Smith PetscInt l; 12576497c311SBarry Smith 125817df9f7cSHong Zhang row = rbuf1_i[ct1] - rstart; 1259a42ab0d6SHong Zhang nzA = a_i[row + 1] - a_i[row]; 1260a42ab0d6SHong Zhang nzB = b_i[row + 1] - b_i[row]; 126117df9f7cSHong Zhang ncols = nzA + nzB; 12628e3a54c0SPierre Jolivet cworkB = PetscSafePointerPlusOffset(b_j, b_i[row]); 12638e3a54c0SPierre Jolivet vworkA = PetscSafePointerPlusOffset(a_a, a_i[row] * bs2); 12648e3a54c0SPierre Jolivet vworkB = PetscSafePointerPlusOffset(b_a, b_i[row] * bs2); 126517df9f7cSHong Zhang 126617df9f7cSHong Zhang /* load the column values for this row into vals*/ 1267a42ab0d6SHong Zhang vals = sbuf_aa_i + ct2 * bs2; 1268a42ab0d6SHong Zhang for (l = 0; l < nzB; l++) { 1269ac530a7eSPierre Jolivet if ((bmap[cworkB[l]]) < cstart) PetscCall(PetscArraycpy(vals + l * bs2, vworkB + l * bs2, bs2)); 1270ac530a7eSPierre Jolivet else break; 1271a42ab0d6SHong Zhang } 1272a42ab0d6SHong Zhang imark = l; 127348a46eb9SPierre Jolivet for (l = 0; l < nzA; l++) PetscCall(PetscArraycpy(vals + (imark + l) * bs2, vworkA + l * bs2, bs2)); 127448a46eb9SPierre Jolivet for (l = imark; l < nzB; l++) PetscCall(PetscArraycpy(vals + (nzA + l) * bs2, vworkB + l * bs2, bs2)); 127517df9f7cSHong Zhang 127617df9f7cSHong Zhang ct2 += ncols; 127717df9f7cSHong Zhang } 127817df9f7cSHong Zhang } 12796497c311SBarry Smith PetscCallMPI(MPIU_Isend(sbuf_aa_i, req_size[i] * bs2, MPIU_SCALAR, req_source1[i], tag4, comm, s_waits4 + i)); 128017df9f7cSHong Zhang } 1281bbc89d27SHong Zhang } 128217df9f7cSHong Zhang 128317df9f7cSHong Zhang /* Assemble the matrices */ 128417df9f7cSHong Zhang /* First assemble the local rows */ 12856497c311SBarry Smith for (PetscInt i = 0; i < ismax; i++) { 128617df9f7cSHong Zhang row2proc_i = row2proc[i]; 128717df9f7cSHong Zhang subc = (Mat_SeqBAIJ *)submats[i]->data; 128817df9f7cSHong Zhang imat_ilen = subc->ilen; 128917df9f7cSHong Zhang imat_j = subc->j; 129017df9f7cSHong Zhang imat_i = subc->i; 129117df9f7cSHong Zhang imat_a = subc->a; 129217df9f7cSHong Zhang 129317df9f7cSHong Zhang if (!allcolumns[i]) cmap_i = cmap[i]; 129417df9f7cSHong Zhang rmap_i = rmap[i]; 129517df9f7cSHong Zhang irow_i = irow[i]; 129617df9f7cSHong Zhang jmax = nrow[i]; 12976497c311SBarry Smith for (PetscInt j = 0; j < jmax; j++) { 1298ec3c739cSHong Zhang if (allrows[i]) row = j; 1299ec3c739cSHong Zhang else row = irow_i[j]; 130017df9f7cSHong Zhang proc = row2proc_i[j]; 1301a42ab0d6SHong Zhang 130217df9f7cSHong Zhang if (proc == rank) { 1303a42ab0d6SHong Zhang row = row - rstart; 1304a42ab0d6SHong Zhang nzA = a_i[row + 1] - a_i[row]; 1305a42ab0d6SHong Zhang nzB = b_i[row + 1] - b_i[row]; 1306a42ab0d6SHong Zhang cworkA = a_j + a_i[row]; 13078e3a54c0SPierre Jolivet cworkB = PetscSafePointerPlusOffset(b_j, b_i[row]); 1308a42ab0d6SHong Zhang if (!ijonly) { 1309a42ab0d6SHong Zhang vworkA = a_a + a_i[row] * bs2; 13108e3a54c0SPierre Jolivet vworkB = PetscSafePointerPlusOffset(b_a, b_i[row] * bs2); 1311a42ab0d6SHong Zhang } 1312059f6b74SHong Zhang 1313059f6b74SHong Zhang if (allrows[i]) { 1314059f6b74SHong Zhang row = row + rstart; 1315059f6b74SHong Zhang } else { 1316a42ab0d6SHong Zhang #if defined(PETSC_USE_CTABLE) 1317eec179cfSJacob Faibussowitsch PetscCall(PetscHMapIGetWithDefault(rmap_i, row + rstart + 1, 0, &row)); 1318a42ab0d6SHong Zhang row--; 1319121971b2SHong Zhang 132008401ef6SPierre Jolivet PetscCheck(row >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "row not found in table"); 1321a42ab0d6SHong Zhang #else 1322a42ab0d6SHong Zhang row = rmap_i[row + rstart]; 1323a42ab0d6SHong Zhang #endif 1324059f6b74SHong Zhang } 1325a42ab0d6SHong Zhang mat_i = imat_i[row]; 13268e3a54c0SPierre Jolivet if (!ijonly) mat_a = PetscSafePointerPlusOffset(imat_a, mat_i * bs2); 13278e3a54c0SPierre Jolivet mat_j = PetscSafePointerPlusOffset(imat_j, mat_i); 132880d31651SHong Zhang ilen = imat_ilen[row]; 1329a42ab0d6SHong Zhang 1330a42ab0d6SHong Zhang /* load the column indices for this row into cols*/ 1331a42ab0d6SHong Zhang if (!allcolumns[i]) { 13326497c311SBarry Smith PetscInt l; 13336497c311SBarry Smith 1334a42ab0d6SHong Zhang for (l = 0; l < nzB; l++) { 1335a42ab0d6SHong Zhang if ((ctmp = bmap[cworkB[l]]) < cstart) { 1336a42ab0d6SHong Zhang #if defined(PETSC_USE_CTABLE) 1337eec179cfSJacob Faibussowitsch PetscCall(PetscHMapIGetWithDefault(cmap_i, ctmp + 1, 0, &tcol)); 1338a42ab0d6SHong Zhang if (tcol) { 1339a42ab0d6SHong Zhang #else 1340a42ab0d6SHong Zhang if ((tcol = cmap_i[ctmp])) { 1341a42ab0d6SHong Zhang #endif 1342a42ab0d6SHong Zhang *mat_j++ = tcol - 1; 13439566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(mat_a, vworkB + l * bs2, bs2)); 1344a42ab0d6SHong Zhang mat_a += bs2; 134580d31651SHong Zhang ilen++; 1346a42ab0d6SHong Zhang } 1347a42ab0d6SHong Zhang } else break; 1348a42ab0d6SHong Zhang } 1349a42ab0d6SHong Zhang imark = l; 13506497c311SBarry Smith for (PetscInt l = 0; l < nzA; l++) { 1351a42ab0d6SHong Zhang #if defined(PETSC_USE_CTABLE) 1352eec179cfSJacob Faibussowitsch PetscCall(PetscHMapIGetWithDefault(cmap_i, cstart + cworkA[l] + 1, 0, &tcol)); 1353a42ab0d6SHong Zhang if (tcol) { 1354a42ab0d6SHong Zhang #else 1355a42ab0d6SHong Zhang if ((tcol = cmap_i[cstart + cworkA[l]])) { 1356a42ab0d6SHong Zhang #endif 1357a42ab0d6SHong Zhang *mat_j++ = tcol - 1; 1358a42ab0d6SHong Zhang if (!ijonly) { 13599566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(mat_a, vworkA + l * bs2, bs2)); 1360a42ab0d6SHong Zhang mat_a += bs2; 1361a42ab0d6SHong Zhang } 136280d31651SHong Zhang ilen++; 1363a42ab0d6SHong Zhang } 1364a42ab0d6SHong Zhang } 1365a42ab0d6SHong Zhang for (l = imark; l < nzB; l++) { 1366a42ab0d6SHong Zhang #if defined(PETSC_USE_CTABLE) 1367eec179cfSJacob Faibussowitsch PetscCall(PetscHMapIGetWithDefault(cmap_i, bmap[cworkB[l]] + 1, 0, &tcol)); 1368a42ab0d6SHong Zhang if (tcol) { 1369a42ab0d6SHong Zhang #else 1370a42ab0d6SHong Zhang if ((tcol = cmap_i[bmap[cworkB[l]]])) { 1371a42ab0d6SHong Zhang #endif 1372a42ab0d6SHong Zhang *mat_j++ = tcol - 1; 1373a42ab0d6SHong Zhang if (!ijonly) { 13749566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(mat_a, vworkB + l * bs2, bs2)); 1375a42ab0d6SHong Zhang mat_a += bs2; 1376a42ab0d6SHong Zhang } 137780d31651SHong Zhang ilen++; 1378a42ab0d6SHong Zhang } 1379a42ab0d6SHong Zhang } 1380a42ab0d6SHong Zhang } else { /* allcolumns */ 13816497c311SBarry Smith PetscInt l; 1382a42ab0d6SHong Zhang for (l = 0; l < nzB; l++) { 1383a42ab0d6SHong Zhang if ((ctmp = bmap[cworkB[l]]) < cstart) { 1384a42ab0d6SHong Zhang *mat_j++ = ctmp; 13859566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(mat_a, vworkB + l * bs2, bs2)); 1386a42ab0d6SHong Zhang mat_a += bs2; 138780d31651SHong Zhang ilen++; 1388a42ab0d6SHong Zhang } else break; 1389a42ab0d6SHong Zhang } 1390a42ab0d6SHong Zhang imark = l; 1391a42ab0d6SHong Zhang for (l = 0; l < nzA; l++) { 1392a42ab0d6SHong Zhang *mat_j++ = cstart + cworkA[l]; 1393a42ab0d6SHong Zhang if (!ijonly) { 13949566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(mat_a, vworkA + l * bs2, bs2)); 1395a42ab0d6SHong Zhang mat_a += bs2; 1396a42ab0d6SHong Zhang } 139780d31651SHong Zhang ilen++; 1398a42ab0d6SHong Zhang } 1399a42ab0d6SHong Zhang for (l = imark; l < nzB; l++) { 1400a42ab0d6SHong Zhang *mat_j++ = bmap[cworkB[l]]; 1401a42ab0d6SHong Zhang if (!ijonly) { 14029566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(mat_a, vworkB + l * bs2, bs2)); 1403a42ab0d6SHong Zhang mat_a += bs2; 1404a42ab0d6SHong Zhang } 140580d31651SHong Zhang ilen++; 1406a42ab0d6SHong Zhang } 1407a42ab0d6SHong Zhang } 140880d31651SHong Zhang imat_ilen[row] = ilen; 140917df9f7cSHong Zhang } 141017df9f7cSHong Zhang } 141117df9f7cSHong Zhang } 141217df9f7cSHong Zhang 141317df9f7cSHong Zhang /* Now assemble the off proc rows */ 141448a46eb9SPierre Jolivet if (!ijonly) PetscCallMPI(MPI_Waitall(nrqs, r_waits4, MPI_STATUSES_IGNORE)); 141517df9f7cSHong Zhang for (tmp2 = 0; tmp2 < nrqs; tmp2++) { 141617df9f7cSHong Zhang sbuf1_i = sbuf1[pa[tmp2]]; 141717df9f7cSHong Zhang jmax = sbuf1_i[0]; 141817df9f7cSHong Zhang ct1 = 2 * jmax + 1; 141917df9f7cSHong Zhang ct2 = 0; 142017df9f7cSHong Zhang rbuf2_i = rbuf2[tmp2]; 142117df9f7cSHong Zhang rbuf3_i = rbuf3[tmp2]; 14225dd0c08cSHong Zhang if (!ijonly) rbuf4_i = rbuf4[tmp2]; 14236497c311SBarry Smith for (PetscInt j = 1; j <= jmax; j++) { 142417df9f7cSHong Zhang is_no = sbuf1_i[2 * j - 1]; 142517df9f7cSHong Zhang rmap_i = rmap[is_no]; 142617df9f7cSHong Zhang if (!allcolumns[is_no]) cmap_i = cmap[is_no]; 142717df9f7cSHong Zhang subc = (Mat_SeqBAIJ *)submats[is_no]->data; 142817df9f7cSHong Zhang imat_ilen = subc->ilen; 142917df9f7cSHong Zhang imat_j = subc->j; 143017df9f7cSHong Zhang imat_i = subc->i; 14315dd0c08cSHong Zhang if (!ijonly) imat_a = subc->a; 143217df9f7cSHong Zhang max1 = sbuf1_i[2 * j]; 14336497c311SBarry Smith for (PetscInt k = 0; k < max1; k++, ct1++) { /* for each recved block row */ 143417df9f7cSHong Zhang row = sbuf1_i[ct1]; 1435059f6b74SHong Zhang 1436059f6b74SHong Zhang if (allrows[is_no]) { 1437059f6b74SHong Zhang row = sbuf1_i[ct1]; 1438059f6b74SHong Zhang } else { 143917df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE) 1440eec179cfSJacob Faibussowitsch PetscCall(PetscHMapIGetWithDefault(rmap_i, row + 1, 0, &row)); 144117df9f7cSHong Zhang row--; 144208401ef6SPierre Jolivet PetscCheck(row >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "row not found in table"); 144317df9f7cSHong Zhang #else 144417df9f7cSHong Zhang row = rmap_i[row]; 144517df9f7cSHong Zhang #endif 1446059f6b74SHong Zhang } 144717df9f7cSHong Zhang ilen = imat_ilen[row]; 144817df9f7cSHong Zhang mat_i = imat_i[row]; 14495dd0c08cSHong Zhang if (!ijonly) mat_a = imat_a + mat_i * bs2; 145017df9f7cSHong Zhang mat_j = imat_j + mat_i; 145117df9f7cSHong Zhang max2 = rbuf2_i[ct1]; 145217df9f7cSHong Zhang if (!allcolumns[is_no]) { 14536497c311SBarry Smith for (PetscInt l = 0; l < max2; l++, ct2++) { 145417df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE) 1455eec179cfSJacob Faibussowitsch PetscCall(PetscHMapIGetWithDefault(cmap_i, rbuf3_i[ct2] + 1, 0, &tcol)); 145617df9f7cSHong Zhang #else 145717df9f7cSHong Zhang tcol = cmap_i[rbuf3_i[ct2]]; 145817df9f7cSHong Zhang #endif 145917df9f7cSHong Zhang if (tcol) { 146017df9f7cSHong Zhang *mat_j++ = tcol - 1; 14615dd0c08cSHong Zhang if (!ijonly) { 14629566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(mat_a, rbuf4_i + ct2 * bs2, bs2)); 14635dd0c08cSHong Zhang mat_a += bs2; 14645dd0c08cSHong Zhang } 146517df9f7cSHong Zhang ilen++; 146617df9f7cSHong Zhang } 146717df9f7cSHong Zhang } 146817df9f7cSHong Zhang } else { /* allcolumns */ 14696497c311SBarry Smith for (PetscInt l = 0; l < max2; l++, ct2++) { 147017df9f7cSHong Zhang *mat_j++ = rbuf3_i[ct2]; /* same global column index of C */ 14715dd0c08cSHong Zhang if (!ijonly) { 14729566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(mat_a, rbuf4_i + ct2 * bs2, bs2)); 14735dd0c08cSHong Zhang mat_a += bs2; 14745dd0c08cSHong Zhang } 147517df9f7cSHong Zhang ilen++; 147617df9f7cSHong Zhang } 147717df9f7cSHong Zhang } 147817df9f7cSHong Zhang imat_ilen[row] = ilen; 147917df9f7cSHong Zhang } 148017df9f7cSHong Zhang } 148117df9f7cSHong Zhang } 148217df9f7cSHong Zhang 148380d31651SHong Zhang if (!iscsorted) { /* sort column indices of the rows */ 148480d31651SHong Zhang MatScalar *work; 148580d31651SHong Zhang 14869566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(bs2, &work)); 14876497c311SBarry Smith for (PetscInt i = 0; i < ismax; i++) { 148817df9f7cSHong Zhang subc = (Mat_SeqBAIJ *)submats[i]->data; 148980d31651SHong Zhang imat_ilen = subc->ilen; 149017df9f7cSHong Zhang imat_j = subc->j; 149117df9f7cSHong Zhang imat_i = subc->i; 14924b6ceb0dSHong Zhang if (!ijonly) imat_a = subc->a; 149317df9f7cSHong Zhang if (allcolumns[i]) continue; 14944b6ceb0dSHong Zhang 149517df9f7cSHong Zhang jmax = nrow[i]; 14966497c311SBarry Smith for (PetscInt j = 0; j < jmax; j++) { 149780d31651SHong Zhang mat_i = imat_i[j]; 149880d31651SHong Zhang mat_j = imat_j + mat_i; 14994b6ceb0dSHong Zhang ilen = imat_ilen[j]; 150080d31651SHong Zhang if (ijonly) { 15019566063dSJacob Faibussowitsch PetscCall(PetscSortInt(ilen, mat_j)); 150280d31651SHong Zhang } else { 15034b6ceb0dSHong Zhang mat_a = imat_a + mat_i * bs2; 15049566063dSJacob Faibussowitsch PetscCall(PetscSortIntWithDataArray(ilen, mat_j, mat_a, bs2 * sizeof(MatScalar), work)); 150517df9f7cSHong Zhang } 150617df9f7cSHong Zhang } 150717df9f7cSHong Zhang } 15089566063dSJacob Faibussowitsch PetscCall(PetscFree(work)); 150980d31651SHong Zhang } 151017df9f7cSHong Zhang 1511bbc89d27SHong Zhang if (!ijonly) { 15129566063dSJacob Faibussowitsch PetscCall(PetscFree(r_waits4)); 15139566063dSJacob Faibussowitsch PetscCallMPI(MPI_Waitall(nrqr, s_waits4, MPI_STATUSES_IGNORE)); 15149566063dSJacob Faibussowitsch PetscCall(PetscFree(s_waits4)); 1515bbc89d27SHong Zhang } 151617df9f7cSHong Zhang 151717df9f7cSHong Zhang /* Restore the indices */ 15186497c311SBarry Smith for (PetscInt i = 0; i < ismax; i++) { 151948a46eb9SPierre Jolivet if (!allrows[i]) PetscCall(ISRestoreIndices(isrow[i], irow + i)); 152048a46eb9SPierre Jolivet if (!allcolumns[i]) PetscCall(ISRestoreIndices(iscol[i], icol + i)); 152117df9f7cSHong Zhang } 152217df9f7cSHong Zhang 15236497c311SBarry Smith for (PetscInt i = 0; i < ismax; i++) { 15249566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(submats[i], MAT_FINAL_ASSEMBLY)); 15259566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(submats[i], MAT_FINAL_ASSEMBLY)); 152617df9f7cSHong Zhang } 152717df9f7cSHong Zhang 15289566063dSJacob Faibussowitsch PetscCall(PetscFree5(*(PetscInt ***)&irow, *(PetscInt ***)&icol, nrow, ncol, issorted)); 15299566063dSJacob Faibussowitsch PetscCall(PetscFree5(row2proc, cmap, rmap, allcolumns, allrows)); 153017df9f7cSHong Zhang 1531bbc89d27SHong Zhang if (!ijonly) { 1532175c2bc5SJunchao Zhang if (sbuf_aa) { 15339566063dSJacob Faibussowitsch PetscCall(PetscFree(sbuf_aa[0])); 15349566063dSJacob Faibussowitsch PetscCall(PetscFree(sbuf_aa)); 1535175c2bc5SJunchao Zhang } 153617df9f7cSHong Zhang 15376497c311SBarry Smith for (PetscMPIInt i = 0; i < nrqs; ++i) PetscCall(PetscFree(rbuf4[i])); 15389566063dSJacob Faibussowitsch PetscCall(PetscFree(rbuf4)); 1539bbc89d27SHong Zhang } 15407a868f3eSHong Zhang c->ijonly = PETSC_FALSE; /* set back to the default */ 15413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1542d5b485abSSatish Balay } 1543