123bdbc58SBarry Smith 2d5b485abSSatish Balay /* 3d5b485abSSatish Balay Routines to compute overlapping regions of a parallel MPI matrix 4d5b485abSSatish Balay and to find submatrices that were shared across processors. 5d5b485abSSatish Balay */ 6c6db04a5SJed Brown #include <../src/mat/impls/baij/mpi/mpibaij.h> 7c6db04a5SJed Brown #include <petscbt.h> 8d5b485abSSatish Balay 9b24ad042SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Local(Mat, PetscInt, char **, PetscInt *, PetscInt **); 10b24ad042SBarry Smith static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Receive(Mat, PetscInt, PetscInt **, PetscInt **, PetscInt *); 1109573ac7SBarry Smith extern PetscErrorCode MatGetRow_MPIBAIJ(Mat, PetscInt, PetscInt *, PetscInt **, PetscScalar **); 1209573ac7SBarry Smith extern PetscErrorCode MatRestoreRow_MPIBAIJ(Mat, PetscInt, PetscInt *, PetscInt **, PetscScalar **); 13d5b485abSSatish Balay 149371c9d4SSatish Balay PetscErrorCode MatIncreaseOverlap_MPIBAIJ(Mat C, PetscInt imax, IS is[], PetscInt ov) { 15d0f46423SBarry Smith PetscInt i, N = C->cmap->N, bs = C->rmap->bs; 1636f4e84dSSatish Balay IS *is_new; 1736f4e84dSSatish Balay 183a40ed3dSBarry Smith PetscFunctionBegin; 199566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(imax, &is_new)); 20df36731bSBarry Smith /* Convert the indices into block format */ 219566063dSJacob Faibussowitsch PetscCall(ISCompressIndicesGeneral(N, C->rmap->n, bs, imax, is, is_new)); 2208401ef6SPierre Jolivet PetscCheck(ov >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative overlap specified"); 23*48a46eb9SPierre Jolivet for (i = 0; i < ov; ++i) PetscCall(MatIncreaseOverlap_MPIBAIJ_Once(C, imax, is_new)); 249566063dSJacob Faibussowitsch for (i = 0; i < imax; i++) PetscCall(ISDestroy(&is[i])); 259566063dSJacob Faibussowitsch PetscCall(ISExpandIndicesGeneral(N, N, bs, imax, is_new, is)); 269566063dSJacob Faibussowitsch for (i = 0; i < imax; i++) PetscCall(ISDestroy(&is_new[i])); 279566063dSJacob Faibussowitsch PetscCall(PetscFree(is_new)); 283a40ed3dSBarry Smith PetscFunctionReturn(0); 29d5b485abSSatish Balay } 30d5b485abSSatish Balay 31d5b485abSSatish Balay /* 32d5b485abSSatish Balay Sample message format: 33d5b485abSSatish Balay If a processor A wants processor B to process some elements corresponding 340e9b0e7eSHong Zhang to index sets is[1], is[5] 35d5b485abSSatish Balay mesg [0] = 2 (no of index sets in the mesg) 36d5b485abSSatish Balay ----------- 37d5b485abSSatish Balay mesg [1] = 1 => is[1] 38d5b485abSSatish Balay mesg [2] = sizeof(is[1]); 39d5b485abSSatish Balay ----------- 40d5b485abSSatish Balay mesg [5] = 5 => is[5] 41d5b485abSSatish Balay mesg [6] = sizeof(is[5]); 42d5b485abSSatish Balay ----------- 43d5b485abSSatish Balay mesg [7] 440e9b0e7eSHong Zhang mesg [n] data(is[1]) 45d5b485abSSatish Balay ----------- 46d5b485abSSatish Balay mesg[n+1] 47d5b485abSSatish Balay mesg[m] data(is[5]) 48d5b485abSSatish Balay ----------- 49d5b485abSSatish Balay 50d5b485abSSatish Balay Notes: 51d5b485abSSatish Balay nrqs - no of requests sent (or to be sent out) 522d4ee042Sprj- nrqr - no of requests received (which have to be or which have been processed) 53d5b485abSSatish Balay */ 549371c9d4SSatish Balay PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Once(Mat C, PetscInt imax, IS is[]) { 55df36731bSBarry Smith Mat_MPIBAIJ *c = (Mat_MPIBAIJ *)C->data; 565d0c19d7SBarry Smith const PetscInt **idx, *idx_i; 5724c049a4SHong Zhang PetscInt *n, *w3, *w4, **data, len; 58b24ad042SBarry Smith PetscMPIInt size, rank, tag1, tag2, *w2, *w1, nrqr; 59131c27b5Sprj- PetscInt Mbs, i, j, k, **rbuf, row, nrqs, msz, **outdat, **ptr; 608f9f447aSHong Zhang PetscInt *ctr, *pa, *tmp, *isz, *isz1, **xdata, **rbuf2, *d_p; 61131c27b5Sprj- PetscMPIInt *onodes1, *olengths1, *onodes2, *olengths2, proc = -1; 62f1af5d2fSBarry Smith PetscBT *table; 63749130a8SStefano Zampini MPI_Comm comm, *iscomms; 64d5b485abSSatish Balay MPI_Request *s_waits1, *r_waits1, *s_waits2, *r_waits2; 658f9f447aSHong Zhang char *t_p; 66d5b485abSSatish Balay 673a40ed3dSBarry Smith PetscFunctionBegin; 689566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)C, &comm)); 69d5b485abSSatish Balay size = c->size; 70d5b485abSSatish Balay rank = c->rank; 71df36731bSBarry Smith Mbs = c->Mbs; 72d5b485abSSatish Balay 739566063dSJacob Faibussowitsch PetscCall(PetscObjectGetNewTag((PetscObject)C, &tag1)); 749566063dSJacob Faibussowitsch PetscCall(PetscObjectGetNewTag((PetscObject)C, &tag2)); 75c7dd2924SSatish Balay 769566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(imax + 1, (PetscInt ***)&idx, imax, &n)); 7724c049a4SHong Zhang 78d5b485abSSatish Balay for (i = 0; i < imax; i++) { 799566063dSJacob Faibussowitsch PetscCall(ISGetIndices(is[i], &idx[i])); 809566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(is[i], &n[i])); 81d5b485abSSatish Balay } 82d5b485abSSatish Balay 83d5b485abSSatish Balay /* evaluate communication - mesg to who,length of mesg, and buffer space 84d5b485abSSatish Balay required. Based on this, buffers are allocated, and data copied into them*/ 859566063dSJacob Faibussowitsch PetscCall(PetscCalloc4(size, &w1, size, &w2, size, &w3, size, &w4)); 86d5b485abSSatish Balay for (i = 0; i < imax; i++) { 879566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(w4, size)); /* initialise work vector*/ 88d5b485abSSatish Balay idx_i = idx[i]; 89d5b485abSSatish Balay len = n[i]; 90d5b485abSSatish Balay for (j = 0; j < len; j++) { 91d5b485abSSatish Balay row = idx_i[j]; 9208401ef6SPierre Jolivet PetscCheck(row >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Index set cannot have negative entries"); 939566063dSJacob Faibussowitsch PetscCall(PetscLayoutFindOwner(C->rmap, row * C->rmap->bs, &proc)); 94d5b485abSSatish Balay w4[proc]++; 95d5b485abSSatish Balay } 96d5b485abSSatish Balay for (j = 0; j < size; j++) { 979371c9d4SSatish Balay if (w4[j]) { 989371c9d4SSatish Balay w1[j] += w4[j]; 999371c9d4SSatish Balay w3[j]++; 1009371c9d4SSatish Balay } 101d5b485abSSatish Balay } 102d5b485abSSatish Balay } 103d5b485abSSatish Balay 104d5b485abSSatish Balay nrqs = 0; /* no of outgoing messages */ 105d5b485abSSatish Balay msz = 0; /* total mesg length (for all proc */ 1060e9b0e7eSHong Zhang w1[rank] = 0; /* no mesg sent to itself */ 107d5b485abSSatish Balay w3[rank] = 0; 108d5b485abSSatish Balay for (i = 0; i < size; i++) { 1099371c9d4SSatish Balay if (w1[i]) { 1109371c9d4SSatish Balay w2[i] = 1; 1119371c9d4SSatish Balay nrqs++; 1129371c9d4SSatish Balay } /* there exists a message to proc i */ 113d5b485abSSatish Balay } 114d5b485abSSatish Balay /* pa - is list of processors to communicate with */ 1159566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrqs, &pa)); 116d5b485abSSatish Balay for (i = 0, j = 0; i < size; i++) { 1179371c9d4SSatish Balay if (w1[i]) { 1189371c9d4SSatish Balay pa[j] = i; 1199371c9d4SSatish Balay j++; 1209371c9d4SSatish Balay } 121d5b485abSSatish Balay } 122d5b485abSSatish Balay 123d5b485abSSatish Balay /* Each message would have a header = 1 + 2*(no of IS) + data */ 124d5b485abSSatish Balay for (i = 0; i < nrqs; i++) { 125d5b485abSSatish Balay j = pa[i]; 126d5b485abSSatish Balay w1[j] += w2[j] + 2 * w3[j]; 127d5b485abSSatish Balay msz += w1[j]; 128d5b485abSSatish Balay } 129d5b485abSSatish Balay 130c7dd2924SSatish Balay /* Determine the number of messages to expect, their lengths, from from-ids */ 1319566063dSJacob Faibussowitsch PetscCall(PetscGatherNumberOfMessages(comm, w2, w1, &nrqr)); 1329566063dSJacob Faibussowitsch PetscCall(PetscGatherMessageLengths(comm, nrqs, nrqr, w1, &onodes1, &olengths1)); 133d5b485abSSatish Balay 134c7dd2924SSatish Balay /* Now post the Irecvs corresponding to these messages */ 1359566063dSJacob Faibussowitsch PetscCall(PetscPostIrecvInt(comm, tag1, nrqr, onodes1, olengths1, &rbuf, &r_waits1)); 136d5b485abSSatish Balay 137d5b485abSSatish Balay /* Allocate Memory for outgoing messages */ 1389566063dSJacob Faibussowitsch PetscCall(PetscMalloc4(size, &outdat, size, &ptr, msz, &tmp, size, &ctr)); 1399566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(outdat, size)); 1409566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(ptr, size)); 141d5b485abSSatish Balay { 142b24ad042SBarry Smith PetscInt *iptr = tmp, ict = 0; 143d5b485abSSatish Balay for (i = 0; i < nrqs; i++) { 144d5b485abSSatish Balay j = pa[i]; 145d5b485abSSatish Balay iptr += ict; 146d5b485abSSatish Balay outdat[j] = iptr; 147d5b485abSSatish Balay ict = w1[j]; 148d5b485abSSatish Balay } 149d5b485abSSatish Balay } 150d5b485abSSatish Balay 151d5b485abSSatish Balay /* Form the outgoing messages */ 152d5b485abSSatish Balay /*plug in the headers*/ 153d5b485abSSatish Balay for (i = 0; i < nrqs; i++) { 154d5b485abSSatish Balay j = pa[i]; 155d5b485abSSatish Balay outdat[j][0] = 0; 1569566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(outdat[j] + 1, 2 * w3[j])); 157d5b485abSSatish Balay ptr[j] = outdat[j] + 2 * w3[j] + 1; 158d5b485abSSatish Balay } 159d5b485abSSatish Balay 160d5b485abSSatish Balay /* Memory for doing local proc's work*/ 161d5b485abSSatish Balay { 1629566063dSJacob Faibussowitsch PetscCall(PetscCalloc5(imax, &table, imax, &data, imax, &isz, Mbs * imax, &d_p, (Mbs / PETSC_BITS_PER_BYTE + 1) * imax, &t_p)); 163d5b485abSSatish Balay 164bbd702dbSSatish Balay for (i = 0; i < imax; i++) { 165b6410449SSatish Balay table[i] = t_p + (Mbs / PETSC_BITS_PER_BYTE + 1) * i; 166bbd702dbSSatish Balay data[i] = d_p + (Mbs)*i; 167d5b485abSSatish Balay } 168d5b485abSSatish Balay } 169d5b485abSSatish Balay 170d5b485abSSatish Balay /* Parse the IS and update local tables and the outgoing buf with the data*/ 171d5b485abSSatish Balay { 172b24ad042SBarry Smith PetscInt n_i, *data_i, isz_i, *outdat_j, ctr_j; 173f1af5d2fSBarry Smith PetscBT table_i; 174d5b485abSSatish Balay 175d5b485abSSatish Balay for (i = 0; i < imax; i++) { 1769566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(ctr, size)); 177d5b485abSSatish Balay n_i = n[i]; 178d5b485abSSatish Balay table_i = table[i]; 179d5b485abSSatish Balay idx_i = idx[i]; 180d5b485abSSatish Balay data_i = data[i]; 181d5b485abSSatish Balay isz_i = isz[i]; 182d5b485abSSatish Balay for (j = 0; j < n_i; j++) { /* parse the indices of each IS */ 183d5b485abSSatish Balay row = idx_i[j]; 1849566063dSJacob Faibussowitsch PetscCall(PetscLayoutFindOwner(C->rmap, row * C->rmap->bs, &proc)); 185d5b485abSSatish Balay if (proc != rank) { /* copy to the outgoing buffer */ 186d5b485abSSatish Balay ctr[proc]++; 187d5b485abSSatish Balay *ptr[proc] = row; 188d5b485abSSatish Balay ptr[proc]++; 189d6b45a43SBarry Smith } else { /* Update the local table */ 19026fbe8dcSKarl Rupp if (!PetscBTLookupSet(table_i, row)) data_i[isz_i++] = row; 191d5b485abSSatish Balay } 192d5b485abSSatish Balay } 193d5b485abSSatish Balay /* Update the headers for the current IS */ 194d5b485abSSatish Balay for (j = 0; j < size; j++) { /* Can Optimise this loop by using pa[] */ 195d5b485abSSatish Balay if ((ctr_j = ctr[j])) { 196d5b485abSSatish Balay outdat_j = outdat[j]; 197d5b485abSSatish Balay k = ++outdat_j[0]; 198d5b485abSSatish Balay outdat_j[2 * k] = ctr_j; 199d5b485abSSatish Balay outdat_j[2 * k - 1] = i; 200d5b485abSSatish Balay } 201d5b485abSSatish Balay } 202d5b485abSSatish Balay isz[i] = isz_i; 203d5b485abSSatish Balay } 204d5b485abSSatish Balay } 205d5b485abSSatish Balay 206d5b485abSSatish Balay /* Now post the sends */ 2079566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrqs, &s_waits1)); 208d5b485abSSatish Balay for (i = 0; i < nrqs; ++i) { 209d5b485abSSatish Balay j = pa[i]; 2109566063dSJacob Faibussowitsch PetscCallMPI(MPI_Isend(outdat[j], w1[j], MPIU_INT, j, tag1, comm, s_waits1 + i)); 211d5b485abSSatish Balay } 212d5b485abSSatish Balay 213d5b485abSSatish Balay /* No longer need the original indices*/ 214*48a46eb9SPierre Jolivet for (i = 0; i < imax; ++i) PetscCall(ISRestoreIndices(is[i], idx + i)); 2159566063dSJacob Faibussowitsch PetscCall(PetscFree2(*(PetscInt ***)&idx, n)); 216d5b485abSSatish Balay 2179566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(imax, &iscomms)); 218d5b485abSSatish Balay for (i = 0; i < imax; ++i) { 2199566063dSJacob Faibussowitsch PetscCall(PetscCommDuplicate(PetscObjectComm((PetscObject)is[i]), &iscomms[i], NULL)); 2209566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is[i])); 221d5b485abSSatish Balay } 222d5b485abSSatish Balay 223d5b485abSSatish Balay /* Do Local work*/ 2249566063dSJacob Faibussowitsch PetscCall(MatIncreaseOverlap_MPIBAIJ_Local(C, imax, table, isz, data)); 225d5b485abSSatish Balay 226d5b485abSSatish Balay /* Receive messages*/ 2279566063dSJacob Faibussowitsch PetscCallMPI(MPI_Waitall(nrqr, r_waits1, MPI_STATUSES_IGNORE)); 2289566063dSJacob Faibussowitsch PetscCallMPI(MPI_Waitall(nrqs, s_waits1, MPI_STATUSES_IGNORE)); 229d5b485abSSatish Balay 230d5b485abSSatish Balay /* Phase 1 sends are complete - deallocate buffers */ 2319566063dSJacob Faibussowitsch PetscCall(PetscFree4(outdat, ptr, tmp, ctr)); 2329566063dSJacob Faibussowitsch PetscCall(PetscFree4(w1, w2, w3, w4)); 233d5b485abSSatish Balay 2349566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrqr, &xdata)); 2359566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrqr, &isz1)); 2369566063dSJacob Faibussowitsch PetscCall(MatIncreaseOverlap_MPIBAIJ_Receive(C, nrqr, rbuf, xdata, isz1)); 237175c2bc5SJunchao Zhang if (rbuf) { 2389566063dSJacob Faibussowitsch PetscCall(PetscFree(rbuf[0])); 2399566063dSJacob Faibussowitsch PetscCall(PetscFree(rbuf)); 240175c2bc5SJunchao Zhang } 241d5b485abSSatish Balay 242d5b485abSSatish Balay /* Send the data back*/ 243d5b485abSSatish Balay /* Do a global reduction to know the buffer space req for incoming messages*/ 244d5b485abSSatish Balay { 245b24ad042SBarry Smith PetscMPIInt *rw1; 246d5b485abSSatish Balay 2479566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(size, &rw1)); 248c7dd2924SSatish Balay 249d5b485abSSatish Balay for (i = 0; i < nrqr; ++i) { 250175c2bc5SJunchao Zhang proc = onodes1[i]; 251d5b485abSSatish Balay rw1[proc] = isz1[i]; 252d5b485abSSatish Balay } 253d5b485abSSatish Balay 254c7dd2924SSatish Balay /* Determine the number of messages to expect, their lengths, from from-ids */ 2559566063dSJacob Faibussowitsch PetscCall(PetscGatherMessageLengths(comm, nrqr, nrqs, rw1, &onodes2, &olengths2)); 2569566063dSJacob Faibussowitsch PetscCall(PetscFree(rw1)); 257c7dd2924SSatish Balay } 258c7dd2924SSatish Balay /* Now post the Irecvs corresponding to these messages */ 2599566063dSJacob Faibussowitsch PetscCall(PetscPostIrecvInt(comm, tag2, nrqs, onodes2, olengths2, &rbuf2, &r_waits2)); 260d5b485abSSatish Balay 261d5b485abSSatish Balay /* Now post the sends */ 2629566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrqr, &s_waits2)); 263d5b485abSSatish Balay for (i = 0; i < nrqr; ++i) { 264175c2bc5SJunchao Zhang j = onodes1[i]; 2659566063dSJacob Faibussowitsch PetscCallMPI(MPI_Isend(xdata[i], isz1[i], MPIU_INT, j, tag2, comm, s_waits2 + i)); 266d5b485abSSatish Balay } 267d5b485abSSatish Balay 2689566063dSJacob Faibussowitsch PetscCall(PetscFree(onodes1)); 2699566063dSJacob Faibussowitsch PetscCall(PetscFree(olengths1)); 270175c2bc5SJunchao Zhang 271d5b485abSSatish Balay /* receive work done on other processors*/ 272d5b485abSSatish Balay { 273b24ad042SBarry Smith PetscMPIInt idex; 274b24ad042SBarry Smith PetscInt is_no, ct1, max, *rbuf2_i, isz_i, *data_i, jmax; 275f1af5d2fSBarry Smith PetscBT table_i; 276d5b485abSSatish Balay 277d5b485abSSatish Balay for (i = 0; i < nrqs; ++i) { 2789566063dSJacob Faibussowitsch PetscCallMPI(MPI_Waitany(nrqs, r_waits2, &idex, MPI_STATUS_IGNORE)); 279d5b485abSSatish Balay /* Process the message*/ 280999d9058SBarry Smith rbuf2_i = rbuf2[idex]; 281d5b485abSSatish Balay ct1 = 2 * rbuf2_i[0] + 1; 282999d9058SBarry Smith jmax = rbuf2[idex][0]; 283d5b485abSSatish Balay for (j = 1; j <= jmax; j++) { 284d5b485abSSatish Balay max = rbuf2_i[2 * j]; 285d5b485abSSatish Balay is_no = rbuf2_i[2 * j - 1]; 286d5b485abSSatish Balay isz_i = isz[is_no]; 287d5b485abSSatish Balay data_i = data[is_no]; 288d5b485abSSatish Balay table_i = table[is_no]; 289d5b485abSSatish Balay for (k = 0; k < max; k++, ct1++) { 290d5b485abSSatish Balay row = rbuf2_i[ct1]; 29126fbe8dcSKarl Rupp if (!PetscBTLookupSet(table_i, row)) data_i[isz_i++] = row; 292d5b485abSSatish Balay } 293d5b485abSSatish Balay isz[is_no] = isz_i; 294d5b485abSSatish Balay } 295d5b485abSSatish Balay } 2969566063dSJacob Faibussowitsch PetscCallMPI(MPI_Waitall(nrqr, s_waits2, MPI_STATUSES_IGNORE)); 297d5b485abSSatish Balay } 298d5b485abSSatish Balay 299d5b485abSSatish Balay for (i = 0; i < imax; ++i) { 3009566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(iscomms[i], isz[i], data[i], PETSC_COPY_VALUES, is + i)); 3019566063dSJacob Faibussowitsch PetscCall(PetscCommDestroy(&iscomms[i])); 302d5b485abSSatish Balay } 303d5b485abSSatish Balay 3049566063dSJacob Faibussowitsch PetscCall(PetscFree(iscomms)); 3059566063dSJacob Faibussowitsch PetscCall(PetscFree(onodes2)); 3069566063dSJacob Faibussowitsch PetscCall(PetscFree(olengths2)); 307c7dd2924SSatish Balay 3089566063dSJacob Faibussowitsch PetscCall(PetscFree(pa)); 309175c2bc5SJunchao Zhang if (rbuf2) { 3109566063dSJacob Faibussowitsch PetscCall(PetscFree(rbuf2[0])); 3119566063dSJacob Faibussowitsch PetscCall(PetscFree(rbuf2)); 312175c2bc5SJunchao Zhang } 3139566063dSJacob Faibussowitsch PetscCall(PetscFree(s_waits1)); 3149566063dSJacob Faibussowitsch PetscCall(PetscFree(r_waits1)); 3159566063dSJacob Faibussowitsch PetscCall(PetscFree(s_waits2)); 3169566063dSJacob Faibussowitsch PetscCall(PetscFree(r_waits2)); 3179566063dSJacob Faibussowitsch PetscCall(PetscFree5(table, data, isz, d_p, t_p)); 318175c2bc5SJunchao Zhang if (xdata) { 3199566063dSJacob Faibussowitsch PetscCall(PetscFree(xdata[0])); 3209566063dSJacob Faibussowitsch PetscCall(PetscFree(xdata)); 321175c2bc5SJunchao Zhang } 3229566063dSJacob Faibussowitsch PetscCall(PetscFree(isz1)); 3233a40ed3dSBarry Smith PetscFunctionReturn(0); 324d5b485abSSatish Balay } 325d5b485abSSatish Balay 326d5b485abSSatish Balay /* 327df36731bSBarry Smith MatIncreaseOverlap_MPIBAIJ_Local - Called by MatincreaseOverlap, to do 328d5b485abSSatish Balay the work on the local processor. 329d5b485abSSatish Balay 330d5b485abSSatish Balay Inputs: 331df36731bSBarry Smith C - MAT_MPIBAIJ; 332d5b485abSSatish Balay imax - total no of index sets processed at a time; 333df36731bSBarry Smith table - an array of char - size = Mbs bits. 334d5b485abSSatish Balay 335d5b485abSSatish Balay Output: 3360e9b0e7eSHong Zhang isz - array containing the count of the solution elements corresponding 337d5b485abSSatish Balay to each index set; 338d5b485abSSatish Balay data - pointer to the solutions 339d5b485abSSatish Balay */ 3409371c9d4SSatish Balay static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Local(Mat C, PetscInt imax, PetscBT *table, PetscInt *isz, PetscInt **data) { 341df36731bSBarry Smith Mat_MPIBAIJ *c = (Mat_MPIBAIJ *)C->data; 342d5b485abSSatish Balay Mat A = c->A, B = c->B; 343df36731bSBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)B->data; 344b24ad042SBarry Smith PetscInt start, end, val, max, rstart, cstart, *ai, *aj; 345b24ad042SBarry Smith PetscInt *bi, *bj, *garray, i, j, k, row, *data_i, isz_i; 346f1af5d2fSBarry Smith PetscBT table_i; 347d5b485abSSatish Balay 3483a40ed3dSBarry Smith PetscFunctionBegin; 349899cda47SBarry Smith rstart = c->rstartbs; 350899cda47SBarry Smith cstart = c->cstartbs; 351d5b485abSSatish Balay ai = a->i; 352df36731bSBarry Smith aj = a->j; 353d5b485abSSatish Balay bi = b->i; 354df36731bSBarry Smith bj = b->j; 355d5b485abSSatish Balay garray = c->garray; 356d5b485abSSatish Balay 357d5b485abSSatish Balay for (i = 0; i < imax; i++) { 358d5b485abSSatish Balay data_i = data[i]; 359d5b485abSSatish Balay table_i = table[i]; 360d5b485abSSatish Balay isz_i = isz[i]; 361d5b485abSSatish Balay for (j = 0, max = isz[i]; j < max; j++) { 362d5b485abSSatish Balay row = data_i[j] - rstart; 363d5b485abSSatish Balay start = ai[row]; 364d5b485abSSatish Balay end = ai[row + 1]; 365d5b485abSSatish Balay for (k = start; k < end; k++) { /* Amat */ 366df36731bSBarry Smith val = aj[k] + cstart; 36726fbe8dcSKarl Rupp if (!PetscBTLookupSet(table_i, val)) data_i[isz_i++] = val; 368d5b485abSSatish Balay } 369d5b485abSSatish Balay start = bi[row]; 370d5b485abSSatish Balay end = bi[row + 1]; 371d5b485abSSatish Balay for (k = start; k < end; k++) { /* Bmat */ 372df36731bSBarry Smith val = garray[bj[k]]; 37326fbe8dcSKarl Rupp if (!PetscBTLookupSet(table_i, val)) data_i[isz_i++] = val; 374d5b485abSSatish Balay } 375d5b485abSSatish Balay } 376d5b485abSSatish Balay isz[i] = isz_i; 377d5b485abSSatish Balay } 3783a40ed3dSBarry Smith PetscFunctionReturn(0); 379d5b485abSSatish Balay } 380d5b485abSSatish Balay /* 3812d4ee042Sprj- MatIncreaseOverlap_MPIBAIJ_Receive - Process the received messages, 382d5b485abSSatish Balay and return the output 383d5b485abSSatish Balay 384d5b485abSSatish Balay Input: 385d5b485abSSatish Balay C - the matrix 386d5b485abSSatish Balay nrqr - no of messages being processed. 3872d4ee042Sprj- rbuf - an array of pointers to the received requests 388d5b485abSSatish Balay 389d5b485abSSatish Balay Output: 390d5b485abSSatish Balay xdata - array of messages to be sent back 391d5b485abSSatish Balay isz1 - size of each message 392d5b485abSSatish Balay 393a8c7a070SBarry Smith For better efficiency perhaps we should malloc separately each xdata[i], 394d5b485abSSatish Balay then if a remalloc is required we need only copy the data for that one row 395a5b23f4aSJose E. Roman rather than all previous rows as it is now where a single large chunk of 396d5b485abSSatish Balay memory is used. 397d5b485abSSatish Balay 398d5b485abSSatish Balay */ 3999371c9d4SSatish Balay static PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Receive(Mat C, PetscInt nrqr, PetscInt **rbuf, PetscInt **xdata, PetscInt *isz1) { 400df36731bSBarry Smith Mat_MPIBAIJ *c = (Mat_MPIBAIJ *)C->data; 401d5b485abSSatish Balay Mat A = c->A, B = c->B; 402df36731bSBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)B->data; 403b24ad042SBarry Smith PetscInt rstart, cstart, *ai, *aj, *bi, *bj, *garray, i, j, k; 404b24ad042SBarry Smith PetscInt row, total_sz, ct, ct1, ct2, ct3, mem_estimate, oct2, l, start, end; 405d2a212eaSBarry Smith PetscInt val, max1, max2, Mbs, no_malloc = 0, *tmp, new_estimate, ctr; 406b24ad042SBarry Smith PetscInt *rbuf_i, kmax, rbuf_0; 407f1af5d2fSBarry Smith PetscBT xtable; 408d5b485abSSatish Balay 4093a40ed3dSBarry Smith PetscFunctionBegin; 410df36731bSBarry Smith Mbs = c->Mbs; 411899cda47SBarry Smith rstart = c->rstartbs; 412899cda47SBarry Smith cstart = c->cstartbs; 413d5b485abSSatish Balay ai = a->i; 414df36731bSBarry Smith aj = a->j; 415d5b485abSSatish Balay bi = b->i; 416df36731bSBarry Smith bj = b->j; 417d5b485abSSatish Balay garray = c->garray; 418d5b485abSSatish Balay 419d5b485abSSatish Balay for (i = 0, ct = 0, total_sz = 0; i < nrqr; ++i) { 420d5b485abSSatish Balay rbuf_i = rbuf[i]; 421d5b485abSSatish Balay rbuf_0 = rbuf_i[0]; 422d5b485abSSatish Balay ct += rbuf_0; 42326fbe8dcSKarl Rupp for (j = 1; j <= rbuf_0; j++) total_sz += rbuf_i[2 * j]; 424d5b485abSSatish Balay } 425d5b485abSSatish Balay 426701b8814SBarry Smith if (c->Mbs) max1 = ct * (a->nz + b->nz) / c->Mbs; 427701b8814SBarry Smith else max1 = 1; 428d5b485abSSatish Balay mem_estimate = 3 * ((total_sz > max1 ? total_sz : max1) + 1); 429175c2bc5SJunchao Zhang if (nrqr) { 4309566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(mem_estimate, &xdata[0])); 431d5b485abSSatish Balay ++no_malloc; 432175c2bc5SJunchao Zhang } 4339566063dSJacob Faibussowitsch PetscCall(PetscBTCreate(Mbs, &xtable)); 4349566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(isz1, nrqr)); 435d5b485abSSatish Balay 436d5b485abSSatish Balay ct3 = 0; 437d5b485abSSatish Balay for (i = 0; i < nrqr; i++) { /* for easch mesg from proc i */ 438d5b485abSSatish Balay rbuf_i = rbuf[i]; 439d5b485abSSatish Balay rbuf_0 = rbuf_i[0]; 440d5b485abSSatish Balay ct1 = 2 * rbuf_0 + 1; 441d5b485abSSatish Balay ct2 = ct1; 442d5b485abSSatish Balay ct3 += ct1; 443d5b485abSSatish Balay for (j = 1; j <= rbuf_0; j++) { /* for each IS from proc i*/ 4449566063dSJacob Faibussowitsch PetscCall(PetscBTMemzero(Mbs, xtable)); 445d5b485abSSatish Balay oct2 = ct2; 446d5b485abSSatish Balay kmax = rbuf_i[2 * j]; 447d5b485abSSatish Balay for (k = 0; k < kmax; k++, ct1++) { 448d5b485abSSatish Balay row = rbuf_i[ct1]; 449f1af5d2fSBarry Smith if (!PetscBTLookupSet(xtable, row)) { 450d5b485abSSatish Balay if (!(ct3 < mem_estimate)) { 451b24ad042SBarry Smith new_estimate = (PetscInt)(1.5 * mem_estimate) + 1; 4529566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(new_estimate, &tmp)); 4539566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(tmp, xdata[0], mem_estimate)); 4549566063dSJacob Faibussowitsch PetscCall(PetscFree(xdata[0])); 455d5b485abSSatish Balay xdata[0] = tmp; 4569371c9d4SSatish Balay mem_estimate = new_estimate; 4579371c9d4SSatish Balay ++no_malloc; 45826fbe8dcSKarl Rupp for (ctr = 1; ctr <= i; ctr++) xdata[ctr] = xdata[ctr - 1] + isz1[ctr - 1]; 459d5b485abSSatish Balay } 460d5b485abSSatish Balay xdata[i][ct2++] = row; 461d5b485abSSatish Balay ct3++; 462d5b485abSSatish Balay } 463d5b485abSSatish Balay } 464d5b485abSSatish Balay for (k = oct2, max2 = ct2; k < max2; k++) { 465d5b485abSSatish Balay row = xdata[i][k] - rstart; 466d5b485abSSatish Balay start = ai[row]; 467d5b485abSSatish Balay end = ai[row + 1]; 468d5b485abSSatish Balay for (l = start; l < end; l++) { 469df36731bSBarry Smith val = aj[l] + cstart; 470f1af5d2fSBarry Smith if (!PetscBTLookupSet(xtable, val)) { 471d5b485abSSatish Balay if (!(ct3 < mem_estimate)) { 472b24ad042SBarry Smith new_estimate = (PetscInt)(1.5 * mem_estimate) + 1; 4739566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(new_estimate, &tmp)); 4749566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(tmp, xdata[0], mem_estimate)); 4759566063dSJacob Faibussowitsch PetscCall(PetscFree(xdata[0])); 476d5b485abSSatish Balay xdata[0] = tmp; 4779371c9d4SSatish Balay mem_estimate = new_estimate; 4789371c9d4SSatish Balay ++no_malloc; 47926fbe8dcSKarl Rupp for (ctr = 1; ctr <= i; ctr++) xdata[ctr] = xdata[ctr - 1] + isz1[ctr - 1]; 480d5b485abSSatish Balay } 481d5b485abSSatish Balay xdata[i][ct2++] = val; 482d5b485abSSatish Balay ct3++; 483d5b485abSSatish Balay } 484d5b485abSSatish Balay } 485d5b485abSSatish Balay start = bi[row]; 486d5b485abSSatish Balay end = bi[row + 1]; 487d5b485abSSatish Balay for (l = start; l < end; l++) { 488df36731bSBarry Smith val = garray[bj[l]]; 489f1af5d2fSBarry Smith if (!PetscBTLookupSet(xtable, val)) { 490d5b485abSSatish Balay if (!(ct3 < mem_estimate)) { 491b24ad042SBarry Smith new_estimate = (PetscInt)(1.5 * mem_estimate) + 1; 4929566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(new_estimate, &tmp)); 4939566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(tmp, xdata[0], mem_estimate)); 4949566063dSJacob Faibussowitsch PetscCall(PetscFree(xdata[0])); 495d5b485abSSatish Balay xdata[0] = tmp; 4969371c9d4SSatish Balay mem_estimate = new_estimate; 4979371c9d4SSatish Balay ++no_malloc; 49826fbe8dcSKarl Rupp for (ctr = 1; ctr <= i; ctr++) xdata[ctr] = xdata[ctr - 1] + isz1[ctr - 1]; 499d5b485abSSatish Balay } 500d5b485abSSatish Balay xdata[i][ct2++] = val; 501d5b485abSSatish Balay ct3++; 502d5b485abSSatish Balay } 503d5b485abSSatish Balay } 504d5b485abSSatish Balay } 505d5b485abSSatish Balay /* Update the header*/ 506d5b485abSSatish Balay xdata[i][2 * j] = ct2 - oct2; /* Undo the vector isz1 and use only a var*/ 507d5b485abSSatish Balay xdata[i][2 * j - 1] = rbuf_i[2 * j - 1]; 508d5b485abSSatish Balay } 509d5b485abSSatish Balay xdata[i][0] = rbuf_0; 510175c2bc5SJunchao Zhang if (i + 1 < nrqr) xdata[i + 1] = xdata[i] + ct2; 511d5b485abSSatish Balay isz1[i] = ct2; /* size of each message */ 512d5b485abSSatish Balay } 5139566063dSJacob Faibussowitsch PetscCall(PetscBTDestroy(&xtable)); 5149566063dSJacob Faibussowitsch PetscCall(PetscInfo(C, "Allocated %" PetscInt_FMT " bytes, required %" PetscInt_FMT ", no of mallocs = %" PetscInt_FMT "\n", mem_estimate, ct3, no_malloc)); 5153a40ed3dSBarry Smith PetscFunctionReturn(0); 516d5b485abSSatish Balay } 517d5b485abSSatish Balay 5189371c9d4SSatish Balay PetscErrorCode MatCreateSubMatrices_MPIBAIJ(Mat C, PetscInt ismax, const IS isrow[], const IS iscol[], MatReuse scall, Mat *submat[]) { 51917df9f7cSHong Zhang IS *isrow_block, *iscol_block; 520cf4f527aSSatish Balay Mat_MPIBAIJ *c = (Mat_MPIBAIJ *)C->data; 521121971b2SHong Zhang PetscInt nmax, nstages_local, nstages, i, pos, max_no, N = C->cmap->N, bs = C->rmap->bs; 5225dc98d6aSHong Zhang Mat_SeqBAIJ *subc; 5235c39f6d9SHong Zhang Mat_SubSppt *smat; 524a2feddc5SSatish Balay 5253a40ed3dSBarry Smith PetscFunctionBegin; 52629dcf524SDmitry Karpeev /* The compression and expansion should be avoided. Doesn't point 52729dcf524SDmitry Karpeev out errors, might change the indices, hence buggey */ 5289566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(ismax + 1, &isrow_block, ismax + 1, &iscol_block)); 5299566063dSJacob Faibussowitsch PetscCall(ISCompressIndicesGeneral(N, C->rmap->n, bs, ismax, isrow, isrow_block)); 5309566063dSJacob Faibussowitsch PetscCall(ISCompressIndicesGeneral(N, C->cmap->n, bs, ismax, iscol, iscol_block)); 531cf4f527aSSatish Balay 532cf4f527aSSatish Balay /* Determine the number of stages through which submatrices are done */ 5334698041cSHong Zhang if (!C->cmap->N) nmax = 20 * 1000000 / sizeof(PetscInt); 5344698041cSHong Zhang else nmax = 20 * 1000000 / (c->Nbs * sizeof(PetscInt)); 535329f5518SBarry Smith if (!nmax) nmax = 1; 5365dc98d6aSHong Zhang 5375dc98d6aSHong Zhang if (scall == MAT_INITIAL_MATRIX) { 5384698041cSHong Zhang nstages_local = ismax / nmax + ((ismax % nmax) ? 1 : 0); /* local nstages */ 539cf4f527aSSatish Balay 540653e4784SBarry Smith /* Make sure every processor loops through the nstages */ 5411c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(&nstages_local, &nstages, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)C))); 5425dc98d6aSHong Zhang 5434698041cSHong Zhang /* Allocate memory to hold all the submatrices and dummy submatrices */ 5449566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(ismax + nstages, submat)); 5454698041cSHong Zhang } else { /* MAT_REUSE_MATRIX */ 5464698041cSHong Zhang if (ismax) { 5475dc98d6aSHong Zhang subc = (Mat_SeqBAIJ *)((*submat)[0]->data); 5485dc98d6aSHong Zhang smat = subc->submatis1; 5494698041cSHong Zhang } else { /* (*submat)[0] is a dummy matrix */ 5505c39f6d9SHong Zhang smat = (Mat_SubSppt *)(*submat)[0]->data; 5514698041cSHong Zhang } 55228b400f6SJacob Faibussowitsch PetscCheck(smat, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "MatCreateSubMatrices(...,MAT_REUSE_MATRIX,...) requires submat"); 5535dc98d6aSHong Zhang nstages = smat->nstages; 5545dc98d6aSHong Zhang } 5555dc98d6aSHong Zhang 556cf4f527aSSatish Balay for (i = 0, pos = 0; i < nstages; i++) { 557cf4f527aSSatish Balay if (pos + nmax <= ismax) max_no = nmax; 5584e195584SJunchao Zhang else if (pos >= ismax) max_no = 0; 559cf4f527aSSatish Balay else max_no = ismax - pos; 560a42ab0d6SHong Zhang 5619566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices_MPIBAIJ_local(C, max_no, isrow_block + pos, iscol_block + pos, scall, *submat + pos)); 5624e195584SJunchao Zhang if (!max_no) { 5634e195584SJunchao Zhang if (scall == MAT_INITIAL_MATRIX) { /* submat[pos] is a dummy matrix */ 5645c39f6d9SHong Zhang smat = (Mat_SubSppt *)(*submat)[pos]->data; 5654698041cSHong Zhang smat->nstages = nstages; 5664698041cSHong Zhang } 5674e195584SJunchao Zhang pos++; /* advance to next dummy matrix if any */ 5684e195584SJunchao Zhang } else pos += max_no; 569cf4f527aSSatish Balay } 57036f4e84dSSatish Balay 5715dc98d6aSHong Zhang if (scall == MAT_INITIAL_MATRIX && ismax) { 5725dc98d6aSHong Zhang /* save nstages for reuse */ 5735dc98d6aSHong Zhang subc = (Mat_SeqBAIJ *)((*submat)[0]->data); 5745dc98d6aSHong Zhang smat = subc->submatis1; 5755dc98d6aSHong Zhang smat->nstages = nstages; 5765dc98d6aSHong Zhang } 5775dc98d6aSHong Zhang 57836f4e84dSSatish Balay for (i = 0; i < ismax; i++) { 5799566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isrow_block[i])); 5809566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iscol_block[i])); 58136f4e84dSSatish Balay } 5829566063dSJacob Faibussowitsch PetscCall(PetscFree2(isrow_block, iscol_block)); 5833a40ed3dSBarry Smith PetscFunctionReturn(0); 584a2feddc5SSatish Balay } 585a2feddc5SSatish Balay 586233c2ff6SSatish Balay #if defined(PETSC_USE_CTABLE) 5879371c9d4SSatish Balay PetscErrorCode PetscGetProc(const PetscInt row, const PetscMPIInt size, const PetscInt proc_gnode[], PetscMPIInt *rank) { 588e005ede5SBarry Smith PetscInt nGlobalNd = proc_gnode[size]; 5894dc2109aSBarry Smith PetscMPIInt fproc; 590233c2ff6SSatish Balay 591233c2ff6SSatish Balay PetscFunctionBegin; 5929566063dSJacob Faibussowitsch PetscCall(PetscMPIIntCast((PetscInt)(((float)row * (float)size / (float)nGlobalNd + 0.5)), &fproc)); 59323ce1328SBarry Smith if (fproc > size) fproc = size; 594e005ede5SBarry Smith while (row < proc_gnode[fproc] || row >= proc_gnode[fproc + 1]) { 595e005ede5SBarry Smith if (row < proc_gnode[fproc]) fproc--; 596233c2ff6SSatish Balay else fproc++; 597233c2ff6SSatish Balay } 598e005ede5SBarry Smith *rank = fproc; 599233c2ff6SSatish Balay PetscFunctionReturn(0); 600233c2ff6SSatish Balay } 601233c2ff6SSatish Balay #endif 602233c2ff6SSatish Balay 603a2feddc5SSatish Balay /* -------------------------------------------------------------------------*/ 604b24ad042SBarry Smith /* This code is used for BAIJ and SBAIJ matrices (unfortunate dependency) */ 6059371c9d4SSatish Balay PetscErrorCode MatCreateSubMatrices_MPIBAIJ_local(Mat C, PetscInt ismax, const IS isrow[], const IS iscol[], MatReuse scall, Mat *submats) { 60617df9f7cSHong Zhang Mat_MPIBAIJ *c = (Mat_MPIBAIJ *)C->data; 60717df9f7cSHong Zhang Mat A = c->A; 60817df9f7cSHong Zhang Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)c->B->data, *subc; 60917df9f7cSHong Zhang const PetscInt **icol, **irow; 61017df9f7cSHong Zhang PetscInt *nrow, *ncol, start; 61117df9f7cSHong Zhang PetscMPIInt rank, size, tag0, tag2, tag3, tag4, *w1, *w2, *w3, *w4, nrqr; 612ddb3d6bcSHong Zhang PetscInt **sbuf1, **sbuf2, *sbuf2_i, i, j, k, l, ct1, ct2, **rbuf1, row, proc = -1; 61317df9f7cSHong Zhang PetscInt nrqs = 0, msz, **ptr = NULL, *req_size = NULL, *ctr = NULL, *pa, *tmp = NULL, tcol; 61417df9f7cSHong Zhang PetscInt **rbuf3 = NULL, *req_source1 = NULL, *req_source2, **sbuf_aj, **rbuf2 = NULL, max1, max2; 61517df9f7cSHong Zhang PetscInt **lens, is_no, ncols, *cols, mat_i, *mat_j, tmp2, jmax; 61617df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE) 61717df9f7cSHong Zhang PetscTable *cmap, cmap_i = NULL, *rmap, rmap_i; 61817df9f7cSHong Zhang #else 61917df9f7cSHong Zhang PetscInt **cmap, *cmap_i = NULL, **rmap, *rmap_i; 62017df9f7cSHong Zhang #endif 62196cff833SHong Zhang const PetscInt *irow_i, *icol_i; 62217df9f7cSHong Zhang PetscInt ctr_j, *sbuf1_j, *sbuf_aj_i, *rbuf1_i, kmax, *lens_i; 62317df9f7cSHong Zhang MPI_Request *s_waits1, *r_waits1, *s_waits2, *r_waits2, *r_waits3; 62417df9f7cSHong Zhang MPI_Request *r_waits4, *s_waits3, *s_waits4; 62517df9f7cSHong Zhang MPI_Comm comm; 62664f5bb2dSSatish Balay PetscScalar **rbuf4, *rbuf4_i = NULL, **sbuf_aa, *vals, *mat_a = NULL, *imat_a = NULL, *sbuf_aa_i; 62717df9f7cSHong Zhang PetscMPIInt *onodes1, *olengths1, end; 62880d31651SHong Zhang PetscInt **row2proc, *row2proc_i, *imat_ilen, *imat_j, *imat_i; 6295c39f6d9SHong Zhang Mat_SubSppt *smat_i; 63017df9f7cSHong Zhang PetscBool *issorted, colflag, iscsorted = PETSC_TRUE; 63117df9f7cSHong Zhang PetscInt *sbuf1_i, *rbuf2_i, *rbuf3_i, ilen; 632a42ab0d6SHong Zhang PetscInt bs = C->rmap->bs, bs2 = c->bs2, rstart = c->rstartbs; 633a42ab0d6SHong Zhang PetscBool ijonly = c->ijonly; /* private flag indicates only matrix data structures are requested */ 634a42ab0d6SHong Zhang PetscInt nzA, nzB, *a_i = a->i, *b_i = b->i, *a_j = a->j, *b_j = b->j, ctmp, imark, *cworkA, *cworkB; 635afb3cbd8SHong Zhang PetscScalar *vworkA = NULL, *vworkB = NULL, *a_a = a->a, *b_a = b->a; 636a42ab0d6SHong Zhang PetscInt cstart = c->cstartbs, *bmap = c->garray; 63780d31651SHong Zhang PetscBool *allrows, *allcolumns; 638a42ab0d6SHong Zhang 63917df9f7cSHong Zhang PetscFunctionBegin; 6409566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)C, &comm)); 64117df9f7cSHong Zhang size = c->size; 64217df9f7cSHong Zhang rank = c->rank; 64317df9f7cSHong Zhang 6449566063dSJacob Faibussowitsch PetscCall(PetscMalloc5(ismax, &row2proc, ismax, &cmap, ismax, &rmap, ismax + 1, &allcolumns, ismax, &allrows)); 6459566063dSJacob Faibussowitsch PetscCall(PetscMalloc5(ismax, (PetscInt ***)&irow, ismax, (PetscInt ***)&icol, ismax, &nrow, ismax, &ncol, ismax, &issorted)); 64617df9f7cSHong Zhang 64717df9f7cSHong Zhang for (i = 0; i < ismax; i++) { 6489566063dSJacob Faibussowitsch PetscCall(ISSorted(iscol[i], &issorted[i])); 64980d31651SHong Zhang if (!issorted[i]) iscsorted = issorted[i]; /* columns are not sorted! */ 6509566063dSJacob Faibussowitsch PetscCall(ISSorted(isrow[i], &issorted[i])); 65117df9f7cSHong Zhang 652ec3c739cSHong Zhang /* Check for special case: allcolumns */ 6539566063dSJacob Faibussowitsch PetscCall(ISIdentity(iscol[i], &colflag)); 6549566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(iscol[i], &ncol[i])); 6555dd0c08cSHong Zhang 6561abc2f7bSHong Zhang if (colflag && ncol[i] == c->Nbs) { 6571abc2f7bSHong Zhang allcolumns[i] = PETSC_TRUE; 6581abc2f7bSHong Zhang icol[i] = NULL; 6591abc2f7bSHong Zhang } else { 66017df9f7cSHong Zhang allcolumns[i] = PETSC_FALSE; 6619566063dSJacob Faibussowitsch PetscCall(ISGetIndices(iscol[i], &icol[i])); 6621abc2f7bSHong Zhang } 663ec3c739cSHong Zhang 664ec3c739cSHong Zhang /* Check for special case: allrows */ 6659566063dSJacob Faibussowitsch PetscCall(ISIdentity(isrow[i], &colflag)); 6669566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(isrow[i], &nrow[i])); 667ec3c739cSHong Zhang if (colflag && nrow[i] == c->Mbs) { 668ec3c739cSHong Zhang allrows[i] = PETSC_TRUE; 669ec3c739cSHong Zhang irow[i] = NULL; 670ec3c739cSHong Zhang } else { 671ec3c739cSHong Zhang allrows[i] = PETSC_FALSE; 6729566063dSJacob Faibussowitsch PetscCall(ISGetIndices(isrow[i], &irow[i])); 673ec3c739cSHong Zhang } 67417df9f7cSHong Zhang } 67517df9f7cSHong Zhang 67617df9f7cSHong Zhang if (scall == MAT_REUSE_MATRIX) { 67717df9f7cSHong Zhang /* Assumes new rows are same length as the old rows */ 67817df9f7cSHong Zhang for (i = 0; i < ismax; i++) { 67917df9f7cSHong Zhang subc = (Mat_SeqBAIJ *)(submats[i]->data); 680aed4548fSBarry Smith PetscCheck(subc->mbs == nrow[i] && subc->nbs == ncol[i], PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Cannot reuse matrix. wrong size"); 68117df9f7cSHong Zhang 68217df9f7cSHong Zhang /* Initial matrix as if empty */ 6839566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(subc->ilen, subc->mbs)); 68417df9f7cSHong Zhang 68517df9f7cSHong Zhang /* Initial matrix as if empty */ 68617df9f7cSHong Zhang submats[i]->factortype = C->factortype; 68717df9f7cSHong Zhang 68817df9f7cSHong Zhang smat_i = subc->submatis1; 68917df9f7cSHong Zhang 69017df9f7cSHong Zhang nrqs = smat_i->nrqs; 69117df9f7cSHong Zhang nrqr = smat_i->nrqr; 69217df9f7cSHong Zhang rbuf1 = smat_i->rbuf1; 69317df9f7cSHong Zhang rbuf2 = smat_i->rbuf2; 69417df9f7cSHong Zhang rbuf3 = smat_i->rbuf3; 69517df9f7cSHong Zhang req_source2 = smat_i->req_source2; 69617df9f7cSHong Zhang 69717df9f7cSHong Zhang sbuf1 = smat_i->sbuf1; 69817df9f7cSHong Zhang sbuf2 = smat_i->sbuf2; 69917df9f7cSHong Zhang ptr = smat_i->ptr; 70017df9f7cSHong Zhang tmp = smat_i->tmp; 70117df9f7cSHong Zhang ctr = smat_i->ctr; 70217df9f7cSHong Zhang 70317df9f7cSHong Zhang pa = smat_i->pa; 70417df9f7cSHong Zhang req_size = smat_i->req_size; 70517df9f7cSHong Zhang req_source1 = smat_i->req_source1; 70617df9f7cSHong Zhang 70717df9f7cSHong Zhang allcolumns[i] = smat_i->allcolumns; 708ec3c739cSHong Zhang allrows[i] = smat_i->allrows; 70917df9f7cSHong Zhang row2proc[i] = smat_i->row2proc; 71017df9f7cSHong Zhang rmap[i] = smat_i->rmap; 71117df9f7cSHong Zhang cmap[i] = smat_i->cmap; 71217df9f7cSHong Zhang } 7134698041cSHong Zhang 7144698041cSHong Zhang if (!ismax) { /* Get dummy submatrices and retrieve struct submatis1 */ 71508401ef6SPierre Jolivet PetscCheck(submats[0], PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "submats are null, cannot reuse"); 7165c39f6d9SHong Zhang smat_i = (Mat_SubSppt *)submats[0]->data; 7174698041cSHong Zhang 7184698041cSHong Zhang nrqs = smat_i->nrqs; 7194698041cSHong Zhang nrqr = smat_i->nrqr; 7204698041cSHong Zhang rbuf1 = smat_i->rbuf1; 7214698041cSHong Zhang rbuf2 = smat_i->rbuf2; 7224698041cSHong Zhang rbuf3 = smat_i->rbuf3; 7234698041cSHong Zhang req_source2 = smat_i->req_source2; 7244698041cSHong Zhang 7254698041cSHong Zhang sbuf1 = smat_i->sbuf1; 7264698041cSHong Zhang sbuf2 = smat_i->sbuf2; 7274698041cSHong Zhang ptr = smat_i->ptr; 7284698041cSHong Zhang tmp = smat_i->tmp; 7294698041cSHong Zhang ctr = smat_i->ctr; 7304698041cSHong Zhang 7314698041cSHong Zhang pa = smat_i->pa; 7324698041cSHong Zhang req_size = smat_i->req_size; 7334698041cSHong Zhang req_source1 = smat_i->req_source1; 7344698041cSHong Zhang 735c67fe082SHong Zhang allcolumns[0] = PETSC_FALSE; 7364698041cSHong Zhang } 73717df9f7cSHong Zhang } else { /* scall == MAT_INITIAL_MATRIX */ 73817df9f7cSHong Zhang /* Get some new tags to keep the communication clean */ 7399566063dSJacob Faibussowitsch PetscCall(PetscObjectGetNewTag((PetscObject)C, &tag2)); 7409566063dSJacob Faibussowitsch PetscCall(PetscObjectGetNewTag((PetscObject)C, &tag3)); 74117df9f7cSHong Zhang 74217df9f7cSHong Zhang /* evaluate communication - mesg to who, length of mesg, and buffer space 74317df9f7cSHong Zhang required. Based on this, buffers are allocated, and data copied into them*/ 7449566063dSJacob Faibussowitsch PetscCall(PetscCalloc4(size, &w1, size, &w2, size, &w3, size, &w4)); /* mesg size, initialize work vectors */ 74517df9f7cSHong Zhang 74617df9f7cSHong Zhang for (i = 0; i < ismax; i++) { 74717df9f7cSHong Zhang jmax = nrow[i]; 74817df9f7cSHong Zhang irow_i = irow[i]; 74917df9f7cSHong Zhang 7509566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(jmax, &row2proc_i)); 75117df9f7cSHong Zhang row2proc[i] = row2proc_i; 75217df9f7cSHong Zhang 75317df9f7cSHong Zhang if (issorted[i]) proc = 0; 75417df9f7cSHong Zhang for (j = 0; j < jmax; j++) { 75517df9f7cSHong Zhang if (!issorted[i]) proc = 0; 756ec3c739cSHong Zhang if (allrows[i]) row = j; 757ec3c739cSHong Zhang else row = irow_i[j]; 758ec3c739cSHong Zhang 759a42ab0d6SHong Zhang while (row >= c->rangebs[proc + 1]) proc++; 76017df9f7cSHong Zhang w4[proc]++; 76117df9f7cSHong Zhang row2proc_i[j] = proc; /* map row index to proc */ 76217df9f7cSHong Zhang } 76317df9f7cSHong Zhang for (j = 0; j < size; j++) { 7649371c9d4SSatish Balay if (w4[j]) { 7659371c9d4SSatish Balay w1[j] += w4[j]; 7669371c9d4SSatish Balay w3[j]++; 7679371c9d4SSatish Balay w4[j] = 0; 7689371c9d4SSatish Balay } 76917df9f7cSHong Zhang } 77017df9f7cSHong Zhang } 77117df9f7cSHong Zhang 77217df9f7cSHong Zhang nrqs = 0; /* no of outgoing messages */ 77317df9f7cSHong Zhang msz = 0; /* total mesg length (for all procs) */ 77417df9f7cSHong Zhang w1[rank] = 0; /* no mesg sent to self */ 77517df9f7cSHong Zhang w3[rank] = 0; 77617df9f7cSHong Zhang for (i = 0; i < size; i++) { 7779371c9d4SSatish Balay if (w1[i]) { 7789371c9d4SSatish Balay w2[i] = 1; 7799371c9d4SSatish Balay nrqs++; 7809371c9d4SSatish Balay } /* there exists a message to proc i */ 78117df9f7cSHong Zhang } 7829566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrqs, &pa)); /*(proc -array)*/ 78317df9f7cSHong Zhang for (i = 0, j = 0; i < size; i++) { 7849371c9d4SSatish Balay if (w1[i]) { 7859371c9d4SSatish Balay pa[j] = i; 7869371c9d4SSatish Balay j++; 7879371c9d4SSatish Balay } 78817df9f7cSHong Zhang } 78917df9f7cSHong Zhang 79017df9f7cSHong Zhang /* Each message would have a header = 1 + 2*(no of IS) + data */ 79117df9f7cSHong Zhang for (i = 0; i < nrqs; i++) { 79217df9f7cSHong Zhang j = pa[i]; 79317df9f7cSHong Zhang w1[j] += w2[j] + 2 * w3[j]; 79417df9f7cSHong Zhang msz += w1[j]; 79517df9f7cSHong Zhang } 7969566063dSJacob Faibussowitsch PetscCall(PetscInfo(0, "Number of outgoing messages %" PetscInt_FMT " Total message length %" PetscInt_FMT "\n", nrqs, msz)); 79717df9f7cSHong Zhang 79817df9f7cSHong Zhang /* Determine the number of messages to expect, their lengths, from from-ids */ 7999566063dSJacob Faibussowitsch PetscCall(PetscGatherNumberOfMessages(comm, w2, w1, &nrqr)); 8009566063dSJacob Faibussowitsch PetscCall(PetscGatherMessageLengths(comm, nrqs, nrqr, w1, &onodes1, &olengths1)); 80117df9f7cSHong Zhang 80217df9f7cSHong Zhang /* Now post the Irecvs corresponding to these messages */ 8039566063dSJacob Faibussowitsch PetscCall(PetscObjectGetNewTag((PetscObject)C, &tag0)); 8049566063dSJacob Faibussowitsch PetscCall(PetscPostIrecvInt(comm, tag0, nrqr, onodes1, olengths1, &rbuf1, &r_waits1)); 80517df9f7cSHong Zhang 80617df9f7cSHong Zhang /* Allocate Memory for outgoing messages */ 8079566063dSJacob Faibussowitsch PetscCall(PetscMalloc4(size, &sbuf1, size, &ptr, 2 * msz, &tmp, size, &ctr)); 8089566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(sbuf1, size)); 8099566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(ptr, size)); 81017df9f7cSHong Zhang 81117df9f7cSHong Zhang { 81217df9f7cSHong Zhang PetscInt *iptr = tmp; 81317df9f7cSHong Zhang k = 0; 81417df9f7cSHong Zhang for (i = 0; i < nrqs; i++) { 81517df9f7cSHong Zhang j = pa[i]; 81617df9f7cSHong Zhang iptr += k; 81717df9f7cSHong Zhang sbuf1[j] = iptr; 81817df9f7cSHong Zhang k = w1[j]; 81917df9f7cSHong Zhang } 82017df9f7cSHong Zhang } 82117df9f7cSHong Zhang 82217df9f7cSHong Zhang /* Form the outgoing messages. Initialize the header space */ 82317df9f7cSHong Zhang for (i = 0; i < nrqs; i++) { 82417df9f7cSHong Zhang j = pa[i]; 82517df9f7cSHong Zhang sbuf1[j][0] = 0; 8269566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(sbuf1[j] + 1, 2 * w3[j])); 82717df9f7cSHong Zhang ptr[j] = sbuf1[j] + 2 * w3[j] + 1; 82817df9f7cSHong Zhang } 82917df9f7cSHong Zhang 83017df9f7cSHong Zhang /* Parse the isrow and copy data into outbuf */ 83117df9f7cSHong Zhang for (i = 0; i < ismax; i++) { 83217df9f7cSHong Zhang row2proc_i = row2proc[i]; 8339566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(ctr, size)); 83417df9f7cSHong Zhang irow_i = irow[i]; 83517df9f7cSHong Zhang jmax = nrow[i]; 83617df9f7cSHong Zhang for (j = 0; j < jmax; j++) { /* parse the indices of each IS */ 83717df9f7cSHong Zhang proc = row2proc_i[j]; 838ec3c739cSHong Zhang if (allrows[i]) row = j; 839ec3c739cSHong Zhang else row = irow_i[j]; 840ec3c739cSHong Zhang 84117df9f7cSHong Zhang if (proc != rank) { /* copy to the outgoing buf*/ 84217df9f7cSHong Zhang ctr[proc]++; 843ec3c739cSHong Zhang *ptr[proc] = row; 84417df9f7cSHong Zhang ptr[proc]++; 84517df9f7cSHong Zhang } 84617df9f7cSHong Zhang } 84717df9f7cSHong Zhang /* Update the headers for the current IS */ 84817df9f7cSHong Zhang for (j = 0; j < size; j++) { /* Can Optimise this loop too */ 84917df9f7cSHong Zhang if ((ctr_j = ctr[j])) { 85017df9f7cSHong Zhang sbuf1_j = sbuf1[j]; 85117df9f7cSHong Zhang k = ++sbuf1_j[0]; 85217df9f7cSHong Zhang sbuf1_j[2 * k] = ctr_j; 85317df9f7cSHong Zhang sbuf1_j[2 * k - 1] = i; 85417df9f7cSHong Zhang } 85517df9f7cSHong Zhang } 85617df9f7cSHong Zhang } 85717df9f7cSHong Zhang 85817df9f7cSHong Zhang /* Now post the sends */ 8599566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrqs, &s_waits1)); 86017df9f7cSHong Zhang for (i = 0; i < nrqs; ++i) { 86117df9f7cSHong Zhang j = pa[i]; 8629566063dSJacob Faibussowitsch PetscCallMPI(MPI_Isend(sbuf1[j], w1[j], MPIU_INT, j, tag0, comm, s_waits1 + i)); 86317df9f7cSHong Zhang } 86417df9f7cSHong Zhang 86517df9f7cSHong Zhang /* Post Receives to capture the buffer size */ 8669566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrqs, &r_waits2)); 8679566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(nrqs, &req_source2, nrqs, &rbuf2, nrqs, &rbuf3)); 868175c2bc5SJunchao Zhang if (nrqs) rbuf2[0] = tmp + msz; 8699371c9d4SSatish Balay for (i = 1; i < nrqs; ++i) { rbuf2[i] = rbuf2[i - 1] + w1[pa[i - 1]]; } 87017df9f7cSHong Zhang for (i = 0; i < nrqs; ++i) { 87117df9f7cSHong Zhang j = pa[i]; 8729566063dSJacob Faibussowitsch PetscCallMPI(MPI_Irecv(rbuf2[i], w1[j], MPIU_INT, j, tag2, comm, r_waits2 + i)); 87317df9f7cSHong Zhang } 87417df9f7cSHong Zhang 87517df9f7cSHong Zhang /* Send to other procs the buf size they should allocate */ 87617df9f7cSHong Zhang /* Receive messages*/ 8779566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrqr, &s_waits2)); 8789566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(nrqr, &sbuf2, nrqr, &req_size, nrqr, &req_source1)); 87917df9f7cSHong Zhang 8809566063dSJacob Faibussowitsch PetscCallMPI(MPI_Waitall(nrqr, r_waits1, MPI_STATUSES_IGNORE)); 88117df9f7cSHong Zhang for (i = 0; i < nrqr; ++i) { 88217df9f7cSHong Zhang req_size[i] = 0; 88317df9f7cSHong Zhang rbuf1_i = rbuf1[i]; 88417df9f7cSHong Zhang start = 2 * rbuf1_i[0] + 1; 885175c2bc5SJunchao Zhang end = olengths1[i]; 8869566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(end, &sbuf2[i])); 88717df9f7cSHong Zhang sbuf2_i = sbuf2[i]; 88817df9f7cSHong Zhang for (j = start; j < end; j++) { 889ddb3d6bcSHong Zhang row = rbuf1_i[j] - rstart; 890ddb3d6bcSHong Zhang ncols = a_i[row + 1] - a_i[row] + b_i[row + 1] - b_i[row]; 89117df9f7cSHong Zhang sbuf2_i[j] = ncols; 89217df9f7cSHong Zhang req_size[i] += ncols; 89317df9f7cSHong Zhang } 894175c2bc5SJunchao Zhang req_source1[i] = onodes1[i]; 89517df9f7cSHong Zhang /* form the header */ 89617df9f7cSHong Zhang sbuf2_i[0] = req_size[i]; 89717df9f7cSHong Zhang for (j = 1; j < start; j++) sbuf2_i[j] = rbuf1_i[j]; 89817df9f7cSHong Zhang 8999566063dSJacob Faibussowitsch PetscCallMPI(MPI_Isend(sbuf2_i, end, MPIU_INT, req_source1[i], tag2, comm, s_waits2 + i)); 90017df9f7cSHong Zhang } 901ddb3d6bcSHong Zhang 9029566063dSJacob Faibussowitsch PetscCall(PetscFree(onodes1)); 9039566063dSJacob Faibussowitsch PetscCall(PetscFree(olengths1)); 904175c2bc5SJunchao Zhang 9059566063dSJacob Faibussowitsch PetscCall(PetscFree(r_waits1)); 9069566063dSJacob Faibussowitsch PetscCall(PetscFree4(w1, w2, w3, w4)); 90717df9f7cSHong Zhang 90817df9f7cSHong Zhang /* Receive messages*/ 9099566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrqs, &r_waits3)); 91017df9f7cSHong Zhang 9119566063dSJacob Faibussowitsch PetscCallMPI(MPI_Waitall(nrqs, r_waits2, MPI_STATUSES_IGNORE)); 91217df9f7cSHong Zhang for (i = 0; i < nrqs; ++i) { 9139566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(rbuf2[i][0], &rbuf3[i])); 914175c2bc5SJunchao Zhang req_source2[i] = pa[i]; 9159566063dSJacob Faibussowitsch PetscCallMPI(MPI_Irecv(rbuf3[i], rbuf2[i][0], MPIU_INT, req_source2[i], tag3, comm, r_waits3 + i)); 91617df9f7cSHong Zhang } 9179566063dSJacob Faibussowitsch PetscCall(PetscFree(r_waits2)); 91817df9f7cSHong Zhang 91917df9f7cSHong Zhang /* Wait on sends1 and sends2 */ 9209566063dSJacob Faibussowitsch PetscCallMPI(MPI_Waitall(nrqs, s_waits1, MPI_STATUSES_IGNORE)); 9219566063dSJacob Faibussowitsch PetscCallMPI(MPI_Waitall(nrqr, s_waits2, MPI_STATUSES_IGNORE)); 9229566063dSJacob Faibussowitsch PetscCall(PetscFree(s_waits1)); 9239566063dSJacob Faibussowitsch PetscCall(PetscFree(s_waits2)); 92417df9f7cSHong Zhang 92517df9f7cSHong Zhang /* Now allocate sending buffers for a->j, and send them off */ 9269566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrqr, &sbuf_aj)); 92717df9f7cSHong Zhang for (i = 0, j = 0; i < nrqr; i++) j += req_size[i]; 9289566063dSJacob Faibussowitsch if (nrqr) PetscCall(PetscMalloc1(j, &sbuf_aj[0])); 92917df9f7cSHong Zhang for (i = 1; i < nrqr; i++) sbuf_aj[i] = sbuf_aj[i - 1] + req_size[i - 1]; 93017df9f7cSHong Zhang 9319566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrqr, &s_waits3)); 93217df9f7cSHong Zhang { 93317df9f7cSHong Zhang for (i = 0; i < nrqr; i++) { 93417df9f7cSHong Zhang rbuf1_i = rbuf1[i]; 93517df9f7cSHong Zhang sbuf_aj_i = sbuf_aj[i]; 93617df9f7cSHong Zhang ct1 = 2 * rbuf1_i[0] + 1; 93717df9f7cSHong Zhang ct2 = 0; 93817df9f7cSHong Zhang for (j = 1, max1 = rbuf1_i[0]; j <= max1; j++) { 93917df9f7cSHong Zhang kmax = rbuf1[i][2 * j]; 94017df9f7cSHong Zhang for (k = 0; k < kmax; k++, ct1++) { 94117df9f7cSHong Zhang row = rbuf1_i[ct1] - rstart; 9429371c9d4SSatish Balay nzA = a_i[row + 1] - a_i[row]; 9439371c9d4SSatish Balay nzB = b_i[row + 1] - b_i[row]; 94417df9f7cSHong Zhang ncols = nzA + nzB; 9459371c9d4SSatish Balay cworkA = a_j + a_i[row]; 9469371c9d4SSatish Balay cworkB = b_j + b_i[row]; 94717df9f7cSHong Zhang 94817df9f7cSHong Zhang /* load the column indices for this row into cols */ 94917df9f7cSHong Zhang cols = sbuf_aj_i + ct2; 95017df9f7cSHong Zhang for (l = 0; l < nzB; l++) { 951a42ab0d6SHong Zhang if ((ctmp = bmap[cworkB[l]]) < cstart) cols[l] = ctmp; 952a42ab0d6SHong Zhang else break; 95317df9f7cSHong Zhang } 954a42ab0d6SHong Zhang imark = l; 955a42ab0d6SHong Zhang for (l = 0; l < nzA; l++) { cols[imark + l] = cstart + cworkA[l]; } 956a42ab0d6SHong Zhang for (l = imark; l < nzB; l++) cols[nzA + l] = bmap[cworkB[l]]; 95717df9f7cSHong Zhang ct2 += ncols; 95817df9f7cSHong Zhang } 95917df9f7cSHong Zhang } 9609566063dSJacob Faibussowitsch PetscCallMPI(MPI_Isend(sbuf_aj_i, req_size[i], MPIU_INT, req_source1[i], tag3, comm, s_waits3 + i)); 96117df9f7cSHong Zhang } 96217df9f7cSHong Zhang } 96317df9f7cSHong Zhang 96417df9f7cSHong Zhang /* create col map: global col of C -> local col of submatrices */ 96517df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE) 96617df9f7cSHong Zhang for (i = 0; i < ismax; i++) { 96717df9f7cSHong Zhang if (!allcolumns[i]) { 9689566063dSJacob Faibussowitsch PetscCall(PetscTableCreate(ncol[i], c->Nbs, &cmap[i])); 96917df9f7cSHong Zhang 97017df9f7cSHong Zhang jmax = ncol[i]; 97117df9f7cSHong Zhang icol_i = icol[i]; 97217df9f7cSHong Zhang cmap_i = cmap[i]; 973*48a46eb9SPierre Jolivet for (j = 0; j < jmax; j++) PetscCall(PetscTableAdd(cmap[i], icol_i[j] + 1, j + 1, INSERT_VALUES)); 97417df9f7cSHong Zhang } else cmap[i] = NULL; 97517df9f7cSHong Zhang } 97617df9f7cSHong Zhang #else 97717df9f7cSHong Zhang for (i = 0; i < ismax; i++) { 97817df9f7cSHong Zhang if (!allcolumns[i]) { 9799566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(c->Nbs, &cmap[i])); 98017df9f7cSHong Zhang jmax = ncol[i]; 98117df9f7cSHong Zhang icol_i = icol[i]; 98217df9f7cSHong Zhang cmap_i = cmap[i]; 983a42ab0d6SHong Zhang for (j = 0; j < jmax; j++) cmap_i[icol_i[j]] = j + 1; 98417df9f7cSHong Zhang } else cmap[i] = NULL; 98517df9f7cSHong Zhang } 98617df9f7cSHong Zhang #endif 98717df9f7cSHong Zhang 98817df9f7cSHong Zhang /* Create lens which is required for MatCreate... */ 98917df9f7cSHong Zhang for (i = 0, j = 0; i < ismax; i++) j += nrow[i]; 9909566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ismax, &lens)); 99117df9f7cSHong Zhang 992*48a46eb9SPierre Jolivet if (ismax) PetscCall(PetscCalloc1(j, &lens[0])); 99317df9f7cSHong Zhang for (i = 1; i < ismax; i++) lens[i] = lens[i - 1] + nrow[i - 1]; 99417df9f7cSHong Zhang 99517df9f7cSHong Zhang /* Update lens from local data */ 99617df9f7cSHong Zhang for (i = 0; i < ismax; i++) { 99717df9f7cSHong Zhang row2proc_i = row2proc[i]; 99817df9f7cSHong Zhang jmax = nrow[i]; 99917df9f7cSHong Zhang if (!allcolumns[i]) cmap_i = cmap[i]; 100017df9f7cSHong Zhang irow_i = irow[i]; 100117df9f7cSHong Zhang lens_i = lens[i]; 100217df9f7cSHong Zhang for (j = 0; j < jmax; j++) { 1003ec3c739cSHong Zhang if (allrows[i]) row = j; 1004ec3c739cSHong Zhang else row = irow_i[j]; /* global blocked row of C */ 1005ec3c739cSHong Zhang 100617df9f7cSHong Zhang proc = row2proc_i[j]; 100717df9f7cSHong Zhang if (proc == rank) { 1008a42ab0d6SHong Zhang /* Get indices from matA and then from matB */ 100917df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE) 1010a42ab0d6SHong Zhang PetscInt tt; 1011a42ab0d6SHong Zhang #endif 1012a42ab0d6SHong Zhang row = row - rstart; 1013a42ab0d6SHong Zhang nzA = a_i[row + 1] - a_i[row]; 1014a42ab0d6SHong Zhang nzB = b_i[row + 1] - b_i[row]; 1015a42ab0d6SHong Zhang cworkA = a_j + a_i[row]; 1016a42ab0d6SHong Zhang cworkB = b_j + b_i[row]; 1017a42ab0d6SHong Zhang 1018a42ab0d6SHong Zhang if (!allcolumns[i]) { 1019a42ab0d6SHong Zhang #if defined(PETSC_USE_CTABLE) 1020a42ab0d6SHong Zhang for (k = 0; k < nzA; k++) { 10219566063dSJacob Faibussowitsch PetscCall(PetscTableFind(cmap_i, cstart + cworkA[k] + 1, &tt)); 1022a42ab0d6SHong Zhang if (tt) lens_i[j]++; 1023a42ab0d6SHong Zhang } 1024a42ab0d6SHong Zhang for (k = 0; k < nzB; k++) { 10259566063dSJacob Faibussowitsch PetscCall(PetscTableFind(cmap_i, bmap[cworkB[k]] + 1, &tt)); 1026a42ab0d6SHong Zhang if (tt) lens_i[j]++; 1027a42ab0d6SHong Zhang } 1028a42ab0d6SHong Zhang 1029a42ab0d6SHong Zhang #else 1030a42ab0d6SHong Zhang for (k = 0; k < nzA; k++) { 1031a42ab0d6SHong Zhang if (cmap_i[cstart + cworkA[k]]) lens_i[j]++; 1032a42ab0d6SHong Zhang } 1033a42ab0d6SHong Zhang for (k = 0; k < nzB; k++) { 1034a42ab0d6SHong Zhang if (cmap_i[bmap[cworkB[k]]]) lens_i[j]++; 1035a42ab0d6SHong Zhang } 1036a42ab0d6SHong Zhang #endif 1037a42ab0d6SHong Zhang } else { /* allcolumns */ 1038a42ab0d6SHong Zhang lens_i[j] = nzA + nzB; 1039a42ab0d6SHong Zhang } 104017df9f7cSHong Zhang } 104117df9f7cSHong Zhang } 104217df9f7cSHong Zhang } 104317df9f7cSHong Zhang 104417df9f7cSHong Zhang /* Create row map: global row of C -> local row of submatrices */ 104517df9f7cSHong Zhang for (i = 0; i < ismax; i++) { 1046059f6b74SHong Zhang if (!allrows[i]) { 1047059f6b74SHong Zhang #if defined(PETSC_USE_CTABLE) 10489566063dSJacob Faibussowitsch PetscCall(PetscTableCreate(nrow[i], c->Mbs, &rmap[i])); 104917df9f7cSHong Zhang irow_i = irow[i]; 105017df9f7cSHong Zhang jmax = nrow[i]; 105117df9f7cSHong Zhang for (j = 0; j < jmax; j++) { 1052ec3c739cSHong Zhang if (allrows[i]) { 10539566063dSJacob Faibussowitsch PetscCall(PetscTableAdd(rmap[i], j + 1, j + 1, INSERT_VALUES)); 1054ec3c739cSHong Zhang } else { 10559566063dSJacob Faibussowitsch PetscCall(PetscTableAdd(rmap[i], irow_i[j] + 1, j + 1, INSERT_VALUES)); 105617df9f7cSHong Zhang } 105717df9f7cSHong Zhang } 105817df9f7cSHong Zhang #else 10599566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(c->Mbs, &rmap[i])); 106017df9f7cSHong Zhang rmap_i = rmap[i]; 106117df9f7cSHong Zhang irow_i = irow[i]; 106217df9f7cSHong Zhang jmax = nrow[i]; 106317df9f7cSHong Zhang for (j = 0; j < jmax; j++) { 1064ec3c739cSHong Zhang if (allrows[i]) rmap_i[j] = j; 1065ec3c739cSHong Zhang else rmap_i[irow_i[j]] = j; 106617df9f7cSHong Zhang } 106717df9f7cSHong Zhang #endif 1068059f6b74SHong Zhang } else rmap[i] = NULL; 1069059f6b74SHong Zhang } 107017df9f7cSHong Zhang 107117df9f7cSHong Zhang /* Update lens from offproc data */ 107217df9f7cSHong Zhang { 107317df9f7cSHong Zhang PetscInt *rbuf2_i, *rbuf3_i, *sbuf1_i; 107417df9f7cSHong Zhang 10759566063dSJacob Faibussowitsch PetscCallMPI(MPI_Waitall(nrqs, r_waits3, MPI_STATUSES_IGNORE)); 107617df9f7cSHong Zhang for (tmp2 = 0; tmp2 < nrqs; tmp2++) { 107717df9f7cSHong Zhang sbuf1_i = sbuf1[pa[tmp2]]; 107817df9f7cSHong Zhang jmax = sbuf1_i[0]; 107917df9f7cSHong Zhang ct1 = 2 * jmax + 1; 108017df9f7cSHong Zhang ct2 = 0; 108117df9f7cSHong Zhang rbuf2_i = rbuf2[tmp2]; 108217df9f7cSHong Zhang rbuf3_i = rbuf3[tmp2]; 108317df9f7cSHong Zhang for (j = 1; j <= jmax; j++) { 108417df9f7cSHong Zhang is_no = sbuf1_i[2 * j - 1]; 108517df9f7cSHong Zhang max1 = sbuf1_i[2 * j]; 108617df9f7cSHong Zhang lens_i = lens[is_no]; 108717df9f7cSHong Zhang if (!allcolumns[is_no]) cmap_i = cmap[is_no]; 108817df9f7cSHong Zhang rmap_i = rmap[is_no]; 108917df9f7cSHong Zhang for (k = 0; k < max1; k++, ct1++) { 1090059f6b74SHong Zhang if (allrows[is_no]) { 1091059f6b74SHong Zhang row = sbuf1_i[ct1]; 1092059f6b74SHong Zhang } else { 109317df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE) 10949566063dSJacob Faibussowitsch PetscCall(PetscTableFind(rmap_i, sbuf1_i[ct1] + 1, &row)); 109517df9f7cSHong Zhang row--; 109608401ef6SPierre Jolivet PetscCheck(row >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "row not found in table"); 109717df9f7cSHong Zhang #else 109817df9f7cSHong Zhang row = rmap_i[sbuf1_i[ct1]]; /* the val in the new matrix to be */ 109917df9f7cSHong Zhang #endif 1100059f6b74SHong Zhang } 110117df9f7cSHong Zhang max2 = rbuf2_i[ct1]; 110217df9f7cSHong Zhang for (l = 0; l < max2; l++, ct2++) { 110317df9f7cSHong Zhang if (!allcolumns[is_no]) { 110417df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE) 11059566063dSJacob Faibussowitsch PetscCall(PetscTableFind(cmap_i, rbuf3_i[ct2] + 1, &tcol)); 110617df9f7cSHong Zhang #else 110717df9f7cSHong Zhang tcol = cmap_i[rbuf3_i[ct2]]; 110817df9f7cSHong Zhang #endif 110917df9f7cSHong Zhang if (tcol) lens_i[row]++; 111017df9f7cSHong Zhang } else { /* allcolumns */ 111117df9f7cSHong Zhang lens_i[row]++; /* lens_i[row] += max2 ? */ 111217df9f7cSHong Zhang } 111317df9f7cSHong Zhang } 111417df9f7cSHong Zhang } 111517df9f7cSHong Zhang } 111617df9f7cSHong Zhang } 111717df9f7cSHong Zhang } 11189566063dSJacob Faibussowitsch PetscCall(PetscFree(r_waits3)); 11199566063dSJacob Faibussowitsch PetscCallMPI(MPI_Waitall(nrqr, s_waits3, MPI_STATUSES_IGNORE)); 11209566063dSJacob Faibussowitsch PetscCall(PetscFree(s_waits3)); 112117df9f7cSHong Zhang 112217df9f7cSHong Zhang /* Create the submatrices */ 112317df9f7cSHong Zhang for (i = 0; i < ismax; i++) { 1124a42ab0d6SHong Zhang PetscInt bs_tmp; 1125a42ab0d6SHong Zhang if (ijonly) bs_tmp = 1; 1126a42ab0d6SHong Zhang else bs_tmp = bs; 112717df9f7cSHong Zhang 11289566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF, submats + i)); 11299566063dSJacob Faibussowitsch PetscCall(MatSetSizes(submats[i], nrow[i] * bs_tmp, ncol[i] * bs_tmp, PETSC_DETERMINE, PETSC_DETERMINE)); 113017df9f7cSHong Zhang 11319566063dSJacob Faibussowitsch PetscCall(MatSetType(submats[i], ((PetscObject)A)->type_name)); 11329566063dSJacob Faibussowitsch PetscCall(MatSeqBAIJSetPreallocation(submats[i], bs_tmp, 0, lens[i])); 11339566063dSJacob Faibussowitsch PetscCall(MatSeqSBAIJSetPreallocation(submats[i], bs_tmp, 0, lens[i])); /* this subroutine is used by SBAIJ routines */ 113417df9f7cSHong Zhang 11355c39f6d9SHong Zhang /* create struct Mat_SubSppt and attached it to submat */ 11369566063dSJacob Faibussowitsch PetscCall(PetscNew(&smat_i)); 113717df9f7cSHong Zhang subc = (Mat_SeqBAIJ *)submats[i]->data; 113817df9f7cSHong Zhang subc->submatis1 = smat_i; 113917df9f7cSHong Zhang 114017df9f7cSHong Zhang smat_i->destroy = submats[i]->ops->destroy; 1141f68bb481SHong Zhang submats[i]->ops->destroy = MatDestroySubMatrix_SeqBAIJ; 114217df9f7cSHong Zhang submats[i]->factortype = C->factortype; 114317df9f7cSHong Zhang 114417df9f7cSHong Zhang smat_i->id = i; 114517df9f7cSHong Zhang smat_i->nrqs = nrqs; 114617df9f7cSHong Zhang smat_i->nrqr = nrqr; 114717df9f7cSHong Zhang smat_i->rbuf1 = rbuf1; 114817df9f7cSHong Zhang smat_i->rbuf2 = rbuf2; 114917df9f7cSHong Zhang smat_i->rbuf3 = rbuf3; 115017df9f7cSHong Zhang smat_i->sbuf2 = sbuf2; 115117df9f7cSHong Zhang smat_i->req_source2 = req_source2; 115217df9f7cSHong Zhang 115317df9f7cSHong Zhang smat_i->sbuf1 = sbuf1; 115417df9f7cSHong Zhang smat_i->ptr = ptr; 115517df9f7cSHong Zhang smat_i->tmp = tmp; 115617df9f7cSHong Zhang smat_i->ctr = ctr; 115717df9f7cSHong Zhang 115817df9f7cSHong Zhang smat_i->pa = pa; 115917df9f7cSHong Zhang smat_i->req_size = req_size; 116017df9f7cSHong Zhang smat_i->req_source1 = req_source1; 116117df9f7cSHong Zhang 116217df9f7cSHong Zhang smat_i->allcolumns = allcolumns[i]; 1163ec3c739cSHong Zhang smat_i->allrows = allrows[i]; 116417df9f7cSHong Zhang smat_i->singleis = PETSC_FALSE; 116517df9f7cSHong Zhang smat_i->row2proc = row2proc[i]; 116617df9f7cSHong Zhang smat_i->rmap = rmap[i]; 116717df9f7cSHong Zhang smat_i->cmap = cmap[i]; 116817df9f7cSHong Zhang } 116917df9f7cSHong Zhang 11704698041cSHong Zhang if (!ismax) { /* Create dummy submats[0] for reuse struct subc */ 11719566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF, &submats[0])); 11729566063dSJacob Faibussowitsch PetscCall(MatSetSizes(submats[0], 0, 0, PETSC_DETERMINE, PETSC_DETERMINE)); 11739566063dSJacob Faibussowitsch PetscCall(MatSetType(submats[0], MATDUMMY)); 11744698041cSHong Zhang 11755c39f6d9SHong Zhang /* create struct Mat_SubSppt and attached it to submat */ 11769566063dSJacob Faibussowitsch PetscCall(PetscNewLog(submats[0], &smat_i)); 11774698041cSHong Zhang submats[0]->data = (void *)smat_i; 11784698041cSHong Zhang 11794698041cSHong Zhang smat_i->destroy = submats[0]->ops->destroy; 1180f68bb481SHong Zhang submats[0]->ops->destroy = MatDestroySubMatrix_Dummy; 11814698041cSHong Zhang submats[0]->factortype = C->factortype; 11824698041cSHong Zhang 1183006cef31SHong Zhang smat_i->id = 0; 11844698041cSHong Zhang smat_i->nrqs = nrqs; 11854698041cSHong Zhang smat_i->nrqr = nrqr; 11864698041cSHong Zhang smat_i->rbuf1 = rbuf1; 11874698041cSHong Zhang smat_i->rbuf2 = rbuf2; 11884698041cSHong Zhang smat_i->rbuf3 = rbuf3; 11894698041cSHong Zhang smat_i->sbuf2 = sbuf2; 11904698041cSHong Zhang smat_i->req_source2 = req_source2; 11914698041cSHong Zhang 11924698041cSHong Zhang smat_i->sbuf1 = sbuf1; 11934698041cSHong Zhang smat_i->ptr = ptr; 11944698041cSHong Zhang smat_i->tmp = tmp; 11954698041cSHong Zhang smat_i->ctr = ctr; 11964698041cSHong Zhang 11974698041cSHong Zhang smat_i->pa = pa; 11984698041cSHong Zhang smat_i->req_size = req_size; 11994698041cSHong Zhang smat_i->req_source1 = req_source1; 12004698041cSHong Zhang 1201c67fe082SHong Zhang smat_i->allcolumns = PETSC_FALSE; 12024698041cSHong Zhang smat_i->singleis = PETSC_FALSE; 12034698041cSHong Zhang smat_i->row2proc = NULL; 12044698041cSHong Zhang smat_i->rmap = NULL; 12054698041cSHong Zhang smat_i->cmap = NULL; 12064698041cSHong Zhang } 12074698041cSHong Zhang 12089566063dSJacob Faibussowitsch if (ismax) PetscCall(PetscFree(lens[0])); 12099566063dSJacob Faibussowitsch PetscCall(PetscFree(lens)); 1210175c2bc5SJunchao Zhang if (sbuf_aj) { 12119566063dSJacob Faibussowitsch PetscCall(PetscFree(sbuf_aj[0])); 12129566063dSJacob Faibussowitsch PetscCall(PetscFree(sbuf_aj)); 1213175c2bc5SJunchao Zhang } 121417df9f7cSHong Zhang 121517df9f7cSHong Zhang } /* endof scall == MAT_INITIAL_MATRIX */ 121617df9f7cSHong Zhang 121717df9f7cSHong Zhang /* Post recv matrix values */ 1218bbc89d27SHong Zhang if (!ijonly) { 12199566063dSJacob Faibussowitsch PetscCall(PetscObjectGetNewTag((PetscObject)C, &tag4)); 12209566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrqs, &rbuf4)); 12219566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrqs, &r_waits4)); 122217df9f7cSHong Zhang for (i = 0; i < nrqs; ++i) { 12239566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(rbuf2[i][0] * bs2, &rbuf4[i])); 12249566063dSJacob Faibussowitsch PetscCallMPI(MPI_Irecv(rbuf4[i], rbuf2[i][0] * bs2, MPIU_SCALAR, req_source2[i], tag4, comm, r_waits4 + i)); 122517df9f7cSHong Zhang } 122617df9f7cSHong Zhang 122717df9f7cSHong Zhang /* Allocate sending buffers for a->a, and send them off */ 12289566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrqr, &sbuf_aa)); 122917df9f7cSHong Zhang for (i = 0, j = 0; i < nrqr; i++) j += req_size[i]; 12309e686edcSHong Zhang 12319566063dSJacob Faibussowitsch if (nrqr) PetscCall(PetscMalloc1(j * bs2, &sbuf_aa[0])); 1232a42ab0d6SHong Zhang for (i = 1; i < nrqr; i++) sbuf_aa[i] = sbuf_aa[i - 1] + req_size[i - 1] * bs2; 123317df9f7cSHong Zhang 12349566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrqr, &s_waits4)); 123517df9f7cSHong Zhang 123617df9f7cSHong Zhang for (i = 0; i < nrqr; i++) { 123717df9f7cSHong Zhang rbuf1_i = rbuf1[i]; 123817df9f7cSHong Zhang sbuf_aa_i = sbuf_aa[i]; 123917df9f7cSHong Zhang ct1 = 2 * rbuf1_i[0] + 1; 124017df9f7cSHong Zhang ct2 = 0; 124117df9f7cSHong Zhang for (j = 1, max1 = rbuf1_i[0]; j <= max1; j++) { 124217df9f7cSHong Zhang kmax = rbuf1_i[2 * j]; 124317df9f7cSHong Zhang for (k = 0; k < kmax; k++, ct1++) { 124417df9f7cSHong Zhang row = rbuf1_i[ct1] - rstart; 1245a42ab0d6SHong Zhang nzA = a_i[row + 1] - a_i[row]; 1246a42ab0d6SHong Zhang nzB = b_i[row + 1] - b_i[row]; 124717df9f7cSHong Zhang ncols = nzA + nzB; 124817df9f7cSHong Zhang cworkB = b_j + b_i[row]; 1249a42ab0d6SHong Zhang vworkA = a_a + a_i[row] * bs2; 1250a42ab0d6SHong Zhang vworkB = b_a + b_i[row] * bs2; 125117df9f7cSHong Zhang 125217df9f7cSHong Zhang /* load the column values for this row into vals*/ 1253a42ab0d6SHong Zhang vals = sbuf_aa_i + ct2 * bs2; 1254a42ab0d6SHong Zhang for (l = 0; l < nzB; l++) { 1255a42ab0d6SHong Zhang if ((bmap[cworkB[l]]) < cstart) { 12569566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(vals + l * bs2, vworkB + l * bs2, bs2)); 1257a42ab0d6SHong Zhang } else break; 1258a42ab0d6SHong Zhang } 1259a42ab0d6SHong Zhang imark = l; 1260*48a46eb9SPierre Jolivet for (l = 0; l < nzA; l++) PetscCall(PetscArraycpy(vals + (imark + l) * bs2, vworkA + l * bs2, bs2)); 1261*48a46eb9SPierre Jolivet for (l = imark; l < nzB; l++) PetscCall(PetscArraycpy(vals + (nzA + l) * bs2, vworkB + l * bs2, bs2)); 126217df9f7cSHong Zhang 126317df9f7cSHong Zhang ct2 += ncols; 126417df9f7cSHong Zhang } 126517df9f7cSHong Zhang } 12669566063dSJacob Faibussowitsch PetscCallMPI(MPI_Isend(sbuf_aa_i, req_size[i] * bs2, MPIU_SCALAR, req_source1[i], tag4, comm, s_waits4 + i)); 126717df9f7cSHong Zhang } 1268bbc89d27SHong Zhang } 126917df9f7cSHong Zhang 127017df9f7cSHong Zhang /* Assemble the matrices */ 127117df9f7cSHong Zhang /* First assemble the local rows */ 127217df9f7cSHong Zhang for (i = 0; i < ismax; i++) { 127317df9f7cSHong Zhang row2proc_i = row2proc[i]; 127417df9f7cSHong Zhang subc = (Mat_SeqBAIJ *)submats[i]->data; 127517df9f7cSHong Zhang imat_ilen = subc->ilen; 127617df9f7cSHong Zhang imat_j = subc->j; 127717df9f7cSHong Zhang imat_i = subc->i; 127817df9f7cSHong Zhang imat_a = subc->a; 127917df9f7cSHong Zhang 128017df9f7cSHong Zhang if (!allcolumns[i]) cmap_i = cmap[i]; 128117df9f7cSHong Zhang rmap_i = rmap[i]; 128217df9f7cSHong Zhang irow_i = irow[i]; 128317df9f7cSHong Zhang jmax = nrow[i]; 128417df9f7cSHong Zhang for (j = 0; j < jmax; j++) { 1285ec3c739cSHong Zhang if (allrows[i]) row = j; 1286ec3c739cSHong Zhang else row = irow_i[j]; 128717df9f7cSHong Zhang proc = row2proc_i[j]; 1288a42ab0d6SHong Zhang 128917df9f7cSHong Zhang if (proc == rank) { 1290a42ab0d6SHong Zhang row = row - rstart; 1291a42ab0d6SHong Zhang nzA = a_i[row + 1] - a_i[row]; 1292a42ab0d6SHong Zhang nzB = b_i[row + 1] - b_i[row]; 1293a42ab0d6SHong Zhang cworkA = a_j + a_i[row]; 1294a42ab0d6SHong Zhang cworkB = b_j + b_i[row]; 1295a42ab0d6SHong Zhang if (!ijonly) { 1296a42ab0d6SHong Zhang vworkA = a_a + a_i[row] * bs2; 1297a42ab0d6SHong Zhang vworkB = b_a + b_i[row] * bs2; 1298a42ab0d6SHong Zhang } 1299059f6b74SHong Zhang 1300059f6b74SHong Zhang if (allrows[i]) { 1301059f6b74SHong Zhang row = row + rstart; 1302059f6b74SHong Zhang } else { 1303a42ab0d6SHong Zhang #if defined(PETSC_USE_CTABLE) 13049566063dSJacob Faibussowitsch PetscCall(PetscTableFind(rmap_i, row + rstart + 1, &row)); 1305a42ab0d6SHong Zhang row--; 1306121971b2SHong Zhang 130708401ef6SPierre Jolivet PetscCheck(row >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "row not found in table"); 1308a42ab0d6SHong Zhang #else 1309a42ab0d6SHong Zhang row = rmap_i[row + rstart]; 1310a42ab0d6SHong Zhang #endif 1311059f6b74SHong Zhang } 1312a42ab0d6SHong Zhang mat_i = imat_i[row]; 1313a42ab0d6SHong Zhang if (!ijonly) mat_a = imat_a + mat_i * bs2; 1314a42ab0d6SHong Zhang mat_j = imat_j + mat_i; 131580d31651SHong Zhang ilen = imat_ilen[row]; 1316a42ab0d6SHong Zhang 1317a42ab0d6SHong Zhang /* load the column indices for this row into cols*/ 1318a42ab0d6SHong Zhang if (!allcolumns[i]) { 1319a42ab0d6SHong Zhang for (l = 0; l < nzB; l++) { 1320a42ab0d6SHong Zhang if ((ctmp = bmap[cworkB[l]]) < cstart) { 1321a42ab0d6SHong Zhang #if defined(PETSC_USE_CTABLE) 13229566063dSJacob Faibussowitsch PetscCall(PetscTableFind(cmap_i, ctmp + 1, &tcol)); 1323a42ab0d6SHong Zhang if (tcol) { 1324a42ab0d6SHong Zhang #else 1325a42ab0d6SHong Zhang if ((tcol = cmap_i[ctmp])) { 1326a42ab0d6SHong Zhang #endif 1327a42ab0d6SHong Zhang *mat_j++ = tcol - 1; 13289566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(mat_a, vworkB + l * bs2, bs2)); 1329a42ab0d6SHong Zhang mat_a += bs2; 133080d31651SHong Zhang ilen++; 1331a42ab0d6SHong Zhang } 1332a42ab0d6SHong Zhang } else break; 1333a42ab0d6SHong Zhang } 1334a42ab0d6SHong Zhang imark = l; 1335a42ab0d6SHong Zhang for (l = 0; l < nzA; l++) { 1336a42ab0d6SHong Zhang #if defined(PETSC_USE_CTABLE) 13379566063dSJacob Faibussowitsch PetscCall(PetscTableFind(cmap_i, cstart + cworkA[l] + 1, &tcol)); 1338a42ab0d6SHong Zhang if (tcol) { 1339a42ab0d6SHong Zhang #else 1340a42ab0d6SHong Zhang if ((tcol = cmap_i[cstart + cworkA[l]])) { 1341a42ab0d6SHong Zhang #endif 1342a42ab0d6SHong Zhang *mat_j++ = tcol - 1; 1343a42ab0d6SHong Zhang if (!ijonly) { 13449566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(mat_a, vworkA + l * bs2, bs2)); 1345a42ab0d6SHong Zhang mat_a += bs2; 1346a42ab0d6SHong Zhang } 134780d31651SHong Zhang ilen++; 1348a42ab0d6SHong Zhang } 1349a42ab0d6SHong Zhang } 1350a42ab0d6SHong Zhang for (l = imark; l < nzB; l++) { 1351a42ab0d6SHong Zhang #if defined(PETSC_USE_CTABLE) 13529566063dSJacob Faibussowitsch PetscCall(PetscTableFind(cmap_i, bmap[cworkB[l]] + 1, &tcol)); 1353a42ab0d6SHong Zhang if (tcol) { 1354a42ab0d6SHong Zhang #else 1355a42ab0d6SHong Zhang if ((tcol = cmap_i[bmap[cworkB[l]]])) { 1356a42ab0d6SHong Zhang #endif 1357a42ab0d6SHong Zhang *mat_j++ = tcol - 1; 1358a42ab0d6SHong Zhang if (!ijonly) { 13599566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(mat_a, vworkB + l * bs2, bs2)); 1360a42ab0d6SHong Zhang mat_a += bs2; 1361a42ab0d6SHong Zhang } 136280d31651SHong Zhang ilen++; 1363a42ab0d6SHong Zhang } 1364a42ab0d6SHong Zhang } 1365a42ab0d6SHong Zhang } else { /* allcolumns */ 1366a42ab0d6SHong Zhang for (l = 0; l < nzB; l++) { 1367a42ab0d6SHong Zhang if ((ctmp = bmap[cworkB[l]]) < cstart) { 1368a42ab0d6SHong Zhang *mat_j++ = ctmp; 13699566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(mat_a, vworkB + l * bs2, bs2)); 1370a42ab0d6SHong Zhang mat_a += bs2; 137180d31651SHong Zhang ilen++; 1372a42ab0d6SHong Zhang } else break; 1373a42ab0d6SHong Zhang } 1374a42ab0d6SHong Zhang imark = l; 1375a42ab0d6SHong Zhang for (l = 0; l < nzA; l++) { 1376a42ab0d6SHong Zhang *mat_j++ = cstart + cworkA[l]; 1377a42ab0d6SHong Zhang if (!ijonly) { 13789566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(mat_a, vworkA + l * bs2, bs2)); 1379a42ab0d6SHong Zhang mat_a += bs2; 1380a42ab0d6SHong Zhang } 138180d31651SHong Zhang ilen++; 1382a42ab0d6SHong Zhang } 1383a42ab0d6SHong Zhang for (l = imark; l < nzB; l++) { 1384a42ab0d6SHong Zhang *mat_j++ = bmap[cworkB[l]]; 1385a42ab0d6SHong Zhang if (!ijonly) { 13869566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(mat_a, vworkB + l * bs2, bs2)); 1387a42ab0d6SHong Zhang mat_a += bs2; 1388a42ab0d6SHong Zhang } 138980d31651SHong Zhang ilen++; 1390a42ab0d6SHong Zhang } 1391a42ab0d6SHong Zhang } 139280d31651SHong Zhang imat_ilen[row] = ilen; 139317df9f7cSHong Zhang } 139417df9f7cSHong Zhang } 139517df9f7cSHong Zhang } 139617df9f7cSHong Zhang 139717df9f7cSHong Zhang /* Now assemble the off proc rows */ 1398*48a46eb9SPierre Jolivet if (!ijonly) PetscCallMPI(MPI_Waitall(nrqs, r_waits4, MPI_STATUSES_IGNORE)); 139917df9f7cSHong Zhang for (tmp2 = 0; tmp2 < nrqs; tmp2++) { 140017df9f7cSHong Zhang sbuf1_i = sbuf1[pa[tmp2]]; 140117df9f7cSHong Zhang jmax = sbuf1_i[0]; 140217df9f7cSHong Zhang ct1 = 2 * jmax + 1; 140317df9f7cSHong Zhang ct2 = 0; 140417df9f7cSHong Zhang rbuf2_i = rbuf2[tmp2]; 140517df9f7cSHong Zhang rbuf3_i = rbuf3[tmp2]; 14065dd0c08cSHong Zhang if (!ijonly) rbuf4_i = rbuf4[tmp2]; 140717df9f7cSHong Zhang for (j = 1; j <= jmax; j++) { 140817df9f7cSHong Zhang is_no = sbuf1_i[2 * j - 1]; 140917df9f7cSHong Zhang rmap_i = rmap[is_no]; 141017df9f7cSHong Zhang if (!allcolumns[is_no]) cmap_i = cmap[is_no]; 141117df9f7cSHong Zhang subc = (Mat_SeqBAIJ *)submats[is_no]->data; 141217df9f7cSHong Zhang imat_ilen = subc->ilen; 141317df9f7cSHong Zhang imat_j = subc->j; 141417df9f7cSHong Zhang imat_i = subc->i; 14155dd0c08cSHong Zhang if (!ijonly) imat_a = subc->a; 141617df9f7cSHong Zhang max1 = sbuf1_i[2 * j]; 14175dd0c08cSHong Zhang for (k = 0; k < max1; k++, ct1++) { /* for each recved block row */ 141817df9f7cSHong Zhang row = sbuf1_i[ct1]; 1419059f6b74SHong Zhang 1420059f6b74SHong Zhang if (allrows[is_no]) { 1421059f6b74SHong Zhang row = sbuf1_i[ct1]; 1422059f6b74SHong Zhang } else { 142317df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE) 14249566063dSJacob Faibussowitsch PetscCall(PetscTableFind(rmap_i, row + 1, &row)); 142517df9f7cSHong Zhang row--; 142608401ef6SPierre Jolivet PetscCheck(row >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "row not found in table"); 142717df9f7cSHong Zhang #else 142817df9f7cSHong Zhang row = rmap_i[row]; 142917df9f7cSHong Zhang #endif 1430059f6b74SHong Zhang } 143117df9f7cSHong Zhang ilen = imat_ilen[row]; 143217df9f7cSHong Zhang mat_i = imat_i[row]; 14335dd0c08cSHong Zhang if (!ijonly) mat_a = imat_a + mat_i * bs2; 143417df9f7cSHong Zhang mat_j = imat_j + mat_i; 143517df9f7cSHong Zhang max2 = rbuf2_i[ct1]; 143617df9f7cSHong Zhang if (!allcolumns[is_no]) { 143717df9f7cSHong Zhang for (l = 0; l < max2; l++, ct2++) { 143817df9f7cSHong Zhang #if defined(PETSC_USE_CTABLE) 14399566063dSJacob Faibussowitsch PetscCall(PetscTableFind(cmap_i, rbuf3_i[ct2] + 1, &tcol)); 144017df9f7cSHong Zhang #else 144117df9f7cSHong Zhang tcol = cmap_i[rbuf3_i[ct2]]; 144217df9f7cSHong Zhang #endif 144317df9f7cSHong Zhang if (tcol) { 144417df9f7cSHong Zhang *mat_j++ = tcol - 1; 14455dd0c08cSHong Zhang if (!ijonly) { 14469566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(mat_a, rbuf4_i + ct2 * bs2, bs2)); 14475dd0c08cSHong Zhang mat_a += bs2; 14485dd0c08cSHong Zhang } 144917df9f7cSHong Zhang ilen++; 145017df9f7cSHong Zhang } 145117df9f7cSHong Zhang } 145217df9f7cSHong Zhang } else { /* allcolumns */ 145317df9f7cSHong Zhang for (l = 0; l < max2; l++, ct2++) { 145417df9f7cSHong Zhang *mat_j++ = rbuf3_i[ct2]; /* same global column index of C */ 14555dd0c08cSHong Zhang if (!ijonly) { 14569566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(mat_a, rbuf4_i + ct2 * bs2, bs2)); 14575dd0c08cSHong Zhang mat_a += bs2; 14585dd0c08cSHong Zhang } 145917df9f7cSHong Zhang ilen++; 146017df9f7cSHong Zhang } 146117df9f7cSHong Zhang } 146217df9f7cSHong Zhang imat_ilen[row] = ilen; 146317df9f7cSHong Zhang } 146417df9f7cSHong Zhang } 146517df9f7cSHong Zhang } 146617df9f7cSHong Zhang 146780d31651SHong Zhang if (!iscsorted) { /* sort column indices of the rows */ 146880d31651SHong Zhang MatScalar *work; 146980d31651SHong Zhang 14709566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(bs2, &work)); 147117df9f7cSHong Zhang for (i = 0; i < ismax; i++) { 147217df9f7cSHong Zhang subc = (Mat_SeqBAIJ *)submats[i]->data; 147380d31651SHong Zhang imat_ilen = subc->ilen; 147417df9f7cSHong Zhang imat_j = subc->j; 147517df9f7cSHong Zhang imat_i = subc->i; 14764b6ceb0dSHong Zhang if (!ijonly) imat_a = subc->a; 147717df9f7cSHong Zhang if (allcolumns[i]) continue; 14784b6ceb0dSHong Zhang 147917df9f7cSHong Zhang jmax = nrow[i]; 148017df9f7cSHong Zhang for (j = 0; j < jmax; j++) { 148180d31651SHong Zhang mat_i = imat_i[j]; 148280d31651SHong Zhang mat_j = imat_j + mat_i; 14834b6ceb0dSHong Zhang ilen = imat_ilen[j]; 148480d31651SHong Zhang if (ijonly) { 14859566063dSJacob Faibussowitsch PetscCall(PetscSortInt(ilen, mat_j)); 148680d31651SHong Zhang } else { 14874b6ceb0dSHong Zhang mat_a = imat_a + mat_i * bs2; 14889566063dSJacob Faibussowitsch PetscCall(PetscSortIntWithDataArray(ilen, mat_j, mat_a, bs2 * sizeof(MatScalar), work)); 148917df9f7cSHong Zhang } 149017df9f7cSHong Zhang } 149117df9f7cSHong Zhang } 14929566063dSJacob Faibussowitsch PetscCall(PetscFree(work)); 149380d31651SHong Zhang } 149417df9f7cSHong Zhang 1495bbc89d27SHong Zhang if (!ijonly) { 14969566063dSJacob Faibussowitsch PetscCall(PetscFree(r_waits4)); 14979566063dSJacob Faibussowitsch PetscCallMPI(MPI_Waitall(nrqr, s_waits4, MPI_STATUSES_IGNORE)); 14989566063dSJacob Faibussowitsch PetscCall(PetscFree(s_waits4)); 1499bbc89d27SHong Zhang } 150017df9f7cSHong Zhang 150117df9f7cSHong Zhang /* Restore the indices */ 150217df9f7cSHong Zhang for (i = 0; i < ismax; i++) { 1503*48a46eb9SPierre Jolivet if (!allrows[i]) PetscCall(ISRestoreIndices(isrow[i], irow + i)); 1504*48a46eb9SPierre Jolivet if (!allcolumns[i]) PetscCall(ISRestoreIndices(iscol[i], icol + i)); 150517df9f7cSHong Zhang } 150617df9f7cSHong Zhang 150717df9f7cSHong Zhang for (i = 0; i < ismax; i++) { 15089566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(submats[i], MAT_FINAL_ASSEMBLY)); 15099566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(submats[i], MAT_FINAL_ASSEMBLY)); 151017df9f7cSHong Zhang } 151117df9f7cSHong Zhang 15129566063dSJacob Faibussowitsch PetscCall(PetscFree5(*(PetscInt ***)&irow, *(PetscInt ***)&icol, nrow, ncol, issorted)); 15139566063dSJacob Faibussowitsch PetscCall(PetscFree5(row2proc, cmap, rmap, allcolumns, allrows)); 151417df9f7cSHong Zhang 1515bbc89d27SHong Zhang if (!ijonly) { 1516175c2bc5SJunchao Zhang if (sbuf_aa) { 15179566063dSJacob Faibussowitsch PetscCall(PetscFree(sbuf_aa[0])); 15189566063dSJacob Faibussowitsch PetscCall(PetscFree(sbuf_aa)); 1519175c2bc5SJunchao Zhang } 152017df9f7cSHong Zhang 1521*48a46eb9SPierre Jolivet for (i = 0; i < nrqs; ++i) PetscCall(PetscFree(rbuf4[i])); 15229566063dSJacob Faibussowitsch PetscCall(PetscFree(rbuf4)); 1523bbc89d27SHong Zhang } 15247a868f3eSHong Zhang c->ijonly = PETSC_FALSE; /* set back to the default */ 15253a40ed3dSBarry Smith PetscFunctionReturn(0); 1526d5b485abSSatish Balay } 1527